functionObjects::age: add optinoal diffusion term and convergence control
to enable the calculation of the residence time for a fluid; mainly used in HVAC analysis. E.g. residence time of air inside a ventilated room, see the new tutorial roomResidenceTime. Contributed by Tobias Holzmann
This commit is contained in:
@ -24,8 +24,9 @@ License
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "age.H"
|
||||
#include "fvmDdt.H"
|
||||
#include "fvmDiv.H"
|
||||
#include "fvmLaplacian.H"
|
||||
#include "turbulenceModel.H"
|
||||
#include "inletOutletFvPatchField.H"
|
||||
#include "wallFvPatch.H"
|
||||
#include "zeroGradientFvPatchField.H"
|
||||
@ -65,6 +66,26 @@ Foam::wordList Foam::functionObjects::age::patchTypes() const
|
||||
}
|
||||
|
||||
|
||||
bool Foam::functionObjects::age::converged
|
||||
(
|
||||
const int nCorr,
|
||||
const scalar initialResidual
|
||||
) const
|
||||
{
|
||||
if (initialResidual < tolerance_)
|
||||
{
|
||||
Info<< "Field " << typeName
|
||||
<< " converged in " << nCorr << " correctors\n" << endl;
|
||||
|
||||
return true;
|
||||
}
|
||||
else
|
||||
{
|
||||
return false;
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
Foam::functionObjects::age::age
|
||||
@ -74,12 +95,7 @@ Foam::functionObjects::age::age
|
||||
const dictionary& dict
|
||||
)
|
||||
:
|
||||
fvMeshFunctionObject(name, runTime, dict),
|
||||
phiName_(),
|
||||
rhoName_(),
|
||||
nCorr_(0),
|
||||
schemesField_()
|
||||
|
||||
fvMeshFunctionObject(name, runTime, dict)
|
||||
{
|
||||
read(dict);
|
||||
}
|
||||
@ -97,10 +113,10 @@ bool Foam::functionObjects::age::read(const dictionary& dict)
|
||||
{
|
||||
phiName_ = dict.lookupOrDefault<word>("phi", "phi");
|
||||
rhoName_ = dict.lookupOrDefault<word>("rho", "rho");
|
||||
|
||||
dict.readIfPresent("nCorr", nCorr_);
|
||||
|
||||
nCorr_ = dict.lookupOrDefault<int>("nCorr", 5);
|
||||
schemesField_ = dict.lookupOrDefault<word>("schemesField", typeName);
|
||||
diffusion_ = dict.lookupOrDefault<Switch>("diffusion", false);
|
||||
tolerance_ = dict.lookupOrDefault<scalar>("tolerance", 1e-5);
|
||||
|
||||
return true;
|
||||
}
|
||||
@ -114,7 +130,7 @@ bool Foam::functionObjects::age::execute()
|
||||
|
||||
bool Foam::functionObjects::age::write()
|
||||
{
|
||||
volScalarField t
|
||||
volScalarField age
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
@ -152,36 +168,83 @@ bool Foam::functionObjects::age::write()
|
||||
const volScalarField& rho =
|
||||
mesh_.lookupObject<volScalarField>(rhoName_);
|
||||
|
||||
for (label i = 0; i <= nCorr_; ++ i)
|
||||
tmp<volScalarField> tmuEff;
|
||||
word laplacianScheme;
|
||||
|
||||
if (diffusion_)
|
||||
{
|
||||
fvScalarMatrix tEqn
|
||||
tmuEff =
|
||||
mesh_.lookupObject<turbulenceModel>
|
||||
(
|
||||
turbulenceModel::propertiesName
|
||||
).muEff();
|
||||
|
||||
laplacianScheme =
|
||||
"laplacian(" + tmuEff().name() + ',' + schemesField_ + ")";
|
||||
}
|
||||
|
||||
for (int i=0; i<=nCorr_; i++)
|
||||
{
|
||||
fvScalarMatrix ageEqn
|
||||
(
|
||||
fvm::div(phi, t, divScheme) == rho
|
||||
fvm::div(phi, age, divScheme) == rho
|
||||
);
|
||||
|
||||
tEqn.relax(relaxCoeff);
|
||||
if (diffusion_)
|
||||
{
|
||||
ageEqn -= fvm::laplacian(tmuEff(), age, laplacianScheme);
|
||||
}
|
||||
|
||||
tEqn.solve(schemesField_);
|
||||
ageEqn.relax(relaxCoeff);
|
||||
|
||||
if (converged(i, ageEqn.solve(schemesField_).initialResidual()))
|
||||
{
|
||||
break;
|
||||
};
|
||||
}
|
||||
}
|
||||
else
|
||||
{
|
||||
for (label i = 0; i <= nCorr_; ++ i)
|
||||
tmp<volScalarField> tnuEff;
|
||||
word laplacianScheme;
|
||||
|
||||
if (diffusion_)
|
||||
{
|
||||
fvScalarMatrix tEqn
|
||||
tnuEff =
|
||||
mesh_.lookupObject<turbulenceModel>
|
||||
(
|
||||
turbulenceModel::propertiesName
|
||||
).nuEff();
|
||||
|
||||
laplacianScheme =
|
||||
"laplacian(" + tnuEff().name() + ',' + schemesField_ + ")";
|
||||
}
|
||||
|
||||
for (int i=0; i<=nCorr_; i++)
|
||||
{
|
||||
fvScalarMatrix ageEqn
|
||||
(
|
||||
fvm::div(phi, t, divScheme) == dimensionedScalar(1)
|
||||
fvm::div(phi, age, divScheme) == dimensionedScalar(1)
|
||||
);
|
||||
|
||||
tEqn.relax(relaxCoeff);
|
||||
if (diffusion_)
|
||||
{
|
||||
ageEqn -= fvm::laplacian(tnuEff(), age, laplacianScheme);
|
||||
}
|
||||
|
||||
tEqn.solve(schemesField_);
|
||||
ageEqn.relax(relaxCoeff);
|
||||
|
||||
if (converged(i, ageEqn.solve(schemesField_).initialResidual()))
|
||||
{
|
||||
break;
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
Info<< "Min/max age:" << min(t).value() << ' ' << max(t).value() << endl;
|
||||
Info<< "Min/max age:" << min(age).value()
|
||||
<< ' ' << max(age).value() << endl;
|
||||
|
||||
t.write();
|
||||
age.write();
|
||||
|
||||
return true;
|
||||
}
|
||||
|
||||
@ -40,12 +40,14 @@ Description
|
||||
|
||||
Usage
|
||||
\table
|
||||
Property | Description | Required | Default value
|
||||
phi | The name of the flux field | no | phi
|
||||
rho | The name of the density field | no | rho
|
||||
nCorr | The number of correctors | no | 0
|
||||
Property | Description | Required | Default value
|
||||
phi | The name of the flux field | no | phi
|
||||
rho | The name of the density field | no | rho
|
||||
nCorr | The maximum number of correctors | no | 5
|
||||
schemesField | The name of the field from which schemes are taken | \\
|
||||
no | age
|
||||
diffusion | Switch to turn on/off the diffusion term | no | off
|
||||
tolerance | Solver residual control | no | 1e-5
|
||||
\endtable
|
||||
|
||||
\verbatim
|
||||
@ -57,9 +59,6 @@ Usage
|
||||
writeControl writeTime;
|
||||
writeInterval 1;
|
||||
|
||||
phi phi;
|
||||
rho rho;
|
||||
nCorr 10;
|
||||
schemesField k;
|
||||
}
|
||||
\endverbatim
|
||||
@ -99,17 +98,26 @@ class age
|
||||
word rhoName_;
|
||||
|
||||
//- Number of corrections
|
||||
label nCorr_;
|
||||
int nCorr_;
|
||||
|
||||
//- Name of field from which schemes are taken
|
||||
word schemesField_;
|
||||
|
||||
//- Switch to turn on/off the diffusion term
|
||||
Switch diffusion_;
|
||||
|
||||
//- Convergence tolerance
|
||||
scalar tolerance_;
|
||||
|
||||
|
||||
// Private Member Functions
|
||||
|
||||
//- The list of patch types for the age field
|
||||
wordList patchTypes() const;
|
||||
|
||||
//- Check convergence
|
||||
bool converged(const int nCorr, const scalar initialResidual) const;
|
||||
|
||||
|
||||
public:
|
||||
|
||||
|
||||
Reference in New Issue
Block a user