Skip to content
Snippets Groups Projects
Commit 43b3fa2c authored by Henry Weller's avatar Henry Weller
Browse files

functionObjects::scalarTransport: simplified, standardized, rationalized

tutorials/incompressible/pisoFoam/les/pitzDaily: Added scalarTransport
functionObject to demonstrate the new functionality
parent 421d3ffd
Branches
Tags
No related merge requests found
Showing
with 125 additions and 130 deletions
...@@ -6,12 +6,14 @@ ...@@ -6,12 +6,14 @@
| \\/ M anipulation | | | \\/ M anipulation | |
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/
T scalarTransport
{ {
#includeEtc "caseDicts/postProcessing/scalarTransport/scalarTransport.cfg" type scalarTransport;
libs ("libsolverFunctionObjects.so");
userDT true; field s;
DT 1e-09; schemesField s;
D 1e-09;
} }
// ************************************************************************* // // ************************************************************************* //
...@@ -7,16 +7,6 @@ ...@@ -7,16 +7,6 @@
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/
type scalarTransport; type scalarTransport;
libs ("libutilityFunctionObjects.so"); libs ("libsolverFunctionObjects.so");
writeControl timeStep;
writeInterval 1;
write true;
log false;
resetOnStartUp false;
autoSchemes true;
fvOptions {};
// ************************************************************************* // // ************************************************************************* //
...@@ -20,6 +20,6 @@ writeInterval 1e-2; ...@@ -20,6 +20,6 @@ writeInterval 1e-2;
// transportProperties settings // transportProperties settings
DT DT [ 0 2 -1 0 0 0 0 ] 1e-9; DT DT [ 0 2 -1 0 0 0 0 ] 1e-9;
#includeEtc "caseDicts/postProcessing/scalarTransport/scalarTransportDict.cfg" #includeEtc "caseDicts/solvers/scalarTransport/scalarTransportDict.cfg"
// ************************************************************************* // // ************************************************************************* //
...@@ -25,13 +25,8 @@ License ...@@ -25,13 +25,8 @@ License
#include "scalarTransport.H" #include "scalarTransport.H"
#include "surfaceFields.H" #include "surfaceFields.H"
#include "dictionary.H"
#include "fixedValueFvPatchFields.H"
#include "zeroGradientFvPatchFields.H"
#include "fvScalarMatrix.H"
#include "fvmDdt.H" #include "fvmDdt.H"
#include "fvmDiv.H" #include "fvmDiv.H"
#include "fvcDiv.H"
#include "fvmLaplacian.H" #include "fvmLaplacian.H"
#include "fvmSup.H" #include "fvmSup.H"
#include "turbulentTransportModel.H" #include "turbulentTransportModel.H"
...@@ -58,30 +53,7 @@ namespace functionObjects ...@@ -58,30 +53,7 @@ namespace functionObjects
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * // // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
Foam::wordList Foam::functionObjects::scalarTransport::boundaryTypes() const Foam::tmp<Foam::volScalarField> Foam::functionObjects::scalarTransport::D
{
const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
wordList bTypes(U.boundaryField().size());
forAll(bTypes, patchi)
{
const fvPatchField<vector>& pf = U.boundaryField()[patchi];
if (isA<fixedValueFvPatchVectorField>(pf))
{
bTypes[patchi] = fixedValueFvPatchScalarField::typeName;
}
else
{
bTypes[patchi] = zeroGradientFvPatchScalarField::typeName;
}
}
return bTypes;
}
Foam::tmp<Foam::volScalarField> Foam::functionObjects::scalarTransport::DT
( (
const surfaceScalarField& phi const surfaceScalarField& phi
) const ) const
...@@ -89,7 +61,9 @@ Foam::tmp<Foam::volScalarField> Foam::functionObjects::scalarTransport::DT ...@@ -89,7 +61,9 @@ Foam::tmp<Foam::volScalarField> Foam::functionObjects::scalarTransport::DT
typedef incompressible::turbulenceModel icoModel; typedef incompressible::turbulenceModel icoModel;
typedef compressible::turbulenceModel cmpModel; typedef compressible::turbulenceModel cmpModel;
if (userDT_) word Dname("D" + s_.name());
if (constantD_)
{ {
return tmp<volScalarField> return tmp<volScalarField>
( (
...@@ -97,14 +71,14 @@ Foam::tmp<Foam::volScalarField> Foam::functionObjects::scalarTransport::DT ...@@ -97,14 +71,14 @@ Foam::tmp<Foam::volScalarField> Foam::functionObjects::scalarTransport::DT
( (
IOobject IOobject
( (
"DT", Dname,
mesh_.time().timeName(), mesh_.time().timeName(),
mesh_.time(), mesh_.time(),
IOobject::NO_READ, IOobject::NO_READ,
IOobject::NO_WRITE IOobject::NO_WRITE
), ),
mesh_, mesh_,
dimensionedScalar("DT", phi.dimensions()/dimLength, DT_) dimensionedScalar(Dname, phi.dimensions()/dimLength, D_)
) )
); );
} }
...@@ -134,14 +108,14 @@ Foam::tmp<Foam::volScalarField> Foam::functionObjects::scalarTransport::DT ...@@ -134,14 +108,14 @@ Foam::tmp<Foam::volScalarField> Foam::functionObjects::scalarTransport::DT
( (
IOobject IOobject
( (
"DT", Dname,
mesh_.time().timeName(), mesh_.time().timeName(),
mesh_.time(), mesh_.time(),
IOobject::NO_READ, IOobject::NO_READ,
IOobject::NO_WRITE IOobject::NO_WRITE
), ),
mesh_, mesh_,
dimensionedScalar("DT", phi.dimensions()/dimLength, 0.0) dimensionedScalar(Dname, phi.dimensions()/dimLength, 0.0)
) )
); );
} }
...@@ -158,30 +132,24 @@ Foam::functionObjects::scalarTransport::scalarTransport ...@@ -158,30 +132,24 @@ Foam::functionObjects::scalarTransport::scalarTransport
) )
: :
fvMeshFunctionObject(name, runTime, dict), fvMeshFunctionObject(name, runTime, dict),
DT_(0.0), fieldName_(dict.lookupOrDefault<word>("field", "s")),
D_(0),
nCorr_(0), nCorr_(0),
fvOptions_(mesh_), fvOptions_(mesh_),
T_ s_
( (
IOobject IOobject
( (
name, fieldName_,
mesh_.time().timeName(), mesh_.time().timeName(),
mesh_, mesh_,
IOobject::READ_IF_PRESENT, IOobject::MUST_READ,
IOobject::AUTO_WRITE IOobject::AUTO_WRITE
), ),
mesh_, mesh_
dimensionedScalar("zero", dimless, 0.0),
boundaryTypes()
) )
{ {
read(dict); read(dict);
if (resetOnStartUp_)
{
T_ == dimensionedScalar("zero", dimless, 0.0);
}
} }
...@@ -198,22 +166,21 @@ bool Foam::functionObjects::scalarTransport::read(const dictionary& dict) ...@@ -198,22 +166,21 @@ bool Foam::functionObjects::scalarTransport::read(const dictionary& dict)
fvMeshFunctionObject::read(dict); fvMeshFunctionObject::read(dict);
phiName_ = dict.lookupOrDefault<word>("phi", "phi"); phiName_ = dict.lookupOrDefault<word>("phi", "phi");
UName_ = dict.lookupOrDefault<word>("U", "U");
rhoName_ = dict.lookupOrDefault<word>("rho", "rho"); rhoName_ = dict.lookupOrDefault<word>("rho", "rho");
schemesField_ = dict.lookupOrDefault<word>("schemesField", fieldName_);
userDT_ = false; constantD_ = false;
if (dict.readIfPresent("DT", DT_)) if (dict.readIfPresent("D", D_))
{ {
userDT_ = true; constantD_ = true;
} }
dict.lookup("resetOnStartUp") >> resetOnStartUp_;
dict.readIfPresent("nCorr", nCorr_); dict.readIfPresent("nCorr", nCorr_);
dict.lookup("autoSchemes") >> autoSchemes_; if (dict.found("fvOptions"))
{
fvOptions_.reset(dict.subDict("fvOptions")); fvOptions_.reset(dict.subDict("fvOptions"));
}
return true; return true;
} }
...@@ -226,24 +193,17 @@ bool Foam::functionObjects::scalarTransport::execute(const bool postProcess) ...@@ -226,24 +193,17 @@ bool Foam::functionObjects::scalarTransport::execute(const bool postProcess)
const surfaceScalarField& phi = const surfaceScalarField& phi =
mesh_.lookupObject<surfaceScalarField>(phiName_); mesh_.lookupObject<surfaceScalarField>(phiName_);
// calculate the diffusivity // Calculate the diffusivity
volScalarField DT(this->DT(phi)); volScalarField D(this->D(phi));
// set schemes
word schemeVar = T_.name();
if (autoSchemes_)
{
schemeVar = UName_;
}
word divScheme("div(phi," + schemeVar + ")"); word divScheme("div(phi," + schemesField_ + ")");
word laplacianScheme("laplacian(" + DT.name() + "," + schemeVar + ")"); word laplacianScheme("laplacian(" + D.name() + "," + schemesField_ + ")");
// set under-relaxation coeff // Set under-relaxation coeff
scalar relaxCoeff = 0.0; scalar relaxCoeff = 0.0;
if (mesh_.relaxEquation(schemeVar)) if (mesh_.relaxEquation(schemesField_))
{ {
relaxCoeff = mesh_.equationRelaxationFactor(schemeVar); relaxCoeff = mesh_.equationRelaxationFactor(schemesField_);
} }
if (phi.dimensions() == dimMass/dimTime) if (phi.dimensions() == dimMass/dimTime)
...@@ -251,44 +211,42 @@ bool Foam::functionObjects::scalarTransport::execute(const bool postProcess) ...@@ -251,44 +211,42 @@ bool Foam::functionObjects::scalarTransport::execute(const bool postProcess)
const volScalarField& rho = const volScalarField& rho =
mesh_.lookupObject<volScalarField>(rhoName_); mesh_.lookupObject<volScalarField>(rhoName_);
// solve
for (label i = 0; i <= nCorr_; i++) for (label i = 0; i <= nCorr_; i++)
{ {
fvScalarMatrix TEqn fvScalarMatrix sEqn
( (
fvm::ddt(rho, T_) fvm::ddt(rho, s_)
+ fvm::div(phi, T_, divScheme) + fvm::div(phi, s_, divScheme)
- fvm::laplacian(DT, T_, laplacianScheme) - fvm::laplacian(D, s_, laplacianScheme)
== ==
fvOptions_(rho, T_) fvOptions_(rho, s_)
); );
TEqn.relax(relaxCoeff); sEqn.relax(relaxCoeff);
fvOptions_.constrain(TEqn); fvOptions_.constrain(sEqn);
TEqn.solve(mesh_.solverDict(schemeVar)); sEqn.solve(mesh_.solverDict(schemesField_));
} }
} }
else if (phi.dimensions() == dimVolume/dimTime) else if (phi.dimensions() == dimVolume/dimTime)
{ {
// solve
for (label i = 0; i <= nCorr_; i++) for (label i = 0; i <= nCorr_; i++)
{ {
fvScalarMatrix TEqn fvScalarMatrix sEqn
( (
fvm::ddt(T_) fvm::ddt(s_)
+ fvm::div(phi, T_, divScheme) + fvm::div(phi, s_, divScheme)
- fvm::laplacian(DT, T_, laplacianScheme) - fvm::laplacian(D, s_, laplacianScheme)
== ==
fvOptions_(T_) fvOptions_(s_)
); );
TEqn.relax(relaxCoeff); sEqn.relax(relaxCoeff);
fvOptions_.constrain(TEqn); fvOptions_.constrain(sEqn);
TEqn.solve(mesh_.solverDict(schemeVar)); sEqn.solve(mesh_.solverDict(schemesField_));
} }
} }
else else
......
...@@ -28,16 +28,13 @@ Group ...@@ -28,16 +28,13 @@ Group
grpSolversFunctionObjects grpSolversFunctionObjects
Description Description
This function object evolves a passive scalar transport equation. The This function object evolves a passive scalar transport equation.
field in ininitially zero, to which sources are added. The field name
is assigned the name of the function object. Boundary conditions are
automatically applied, based on the velocity boundary conditions.
- the field can be zeroed on start-up using the resetOnStartUp flag - To specify the field name set the 'field' entry
- to employ the same numerical schemes as the flow velocity, use the - To employ the same numerical schemes as another field set
autoSchemes flag the 'schemesField' entry,
- the diffusivity can be set manually using the DT entry, or retrieved - The diffusivity can be set manually using the 'D' entry, or retrieved
from the turbulence model (if applicable) from the turbulence model (if applicable).
SeeAlso SeeAlso
Foam::functionObjects::fvMeshFunctionObject Foam::functionObjects::fvMeshFunctionObject
...@@ -52,7 +49,6 @@ SourceFiles ...@@ -52,7 +49,6 @@ SourceFiles
#include "fvMeshFunctionObject.H" #include "fvMeshFunctionObject.H"
#include "volFields.H" #include "volFields.H"
#include "surfaceFieldsFwd.H"
#include "fvOptionList.H" #include "fvOptionList.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
...@@ -72,44 +68,38 @@ class scalarTransport ...@@ -72,44 +68,38 @@ class scalarTransport
{ {
// Private data // Private data
//- Name of field to process
word fieldName_;
//- Name of flux field (optional) //- Name of flux field (optional)
word phiName_; word phiName_;
//- Name of velocity field (optional)
word UName_;
//- Name of density field (optional) //- Name of density field (optional)
word rhoName_; word rhoName_;
//- Diffusion coefficient (optional) //- Diffusion coefficient (optional)
scalar DT_; scalar D_;
//- Flag to indicate whether user DT_ is used //- Flag to indicate whether a constant, uniform D_ is specified
bool userDT_; bool constantD_;
//- Flag to reset scalar field on start-up
bool resetOnStartUp_;
//- Number of corrector iterations (optional) //- Number of corrector iterations (optional)
label nCorr_; label nCorr_;
//- Flag to employ schemes for velocity for scalar transport //- Name of field whose schemes are used (optional)
bool autoSchemes_; word schemesField_;
//- Run-time selectable finite volume options, e.g. sources, constraints //- Run-time selectable finite volume options, e.g. sources, constraints
fv::optionList fvOptions_; fv::optionList fvOptions_;
//- The scalar field //- The scalar field
volScalarField T_; volScalarField s_;
// Private Member Functions // Private Member Functions
//- Return the boundary types for the scalar field
wordList boundaryTypes() const;
//- Return the diffusivity field //- Return the diffusivity field
tmp<volScalarField> DT(const surfaceScalarField& phi) const; tmp<volScalarField> D(const surfaceScalarField& phi) const;
//- Disallow default bitwise copy construct //- Disallow default bitwise copy construct
scalarTransport(const scalarTransport&); scalarTransport(const scalarTransport&);
......
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: dev |
| \\ / A nd | Web: www.OpenFOAM.org |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class volScalarField;
object s;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
dimensions [0 0 0 0 0 0 0];
internalField uniform 0;
boundaryField
{
inlet
{
type fixedValue;
value uniform 1;
}
outlet
{
type inletOutlet;
inletValue uniform 0;
value uniform 0;
}
upperWall
{
type zeroGradient;
}
lowerWall
{
type zeroGradient;
}
frontAndBack
{
type empty;
}
}
// ************************************************************************* //
...@@ -128,6 +128,8 @@ functions ...@@ -128,6 +128,8 @@ functions
} }
); );
} }
#includeFunc scalarTransport
} }
// ************************************************************************* // // ************************************************************************* //
...@@ -30,6 +30,7 @@ divSchemes ...@@ -30,6 +30,7 @@ divSchemes
default none; default none;
div(phi,U) Gauss LUST grad(U); div(phi,U) Gauss LUST grad(U);
div(phi,k) Gauss limitedLinear 1; div(phi,k) Gauss limitedLinear 1;
div(phi,s) bounded Gauss limitedLinear 1;
div((nuEff*dev2(T(grad(U))))) Gauss linear; div((nuEff*dev2(T(grad(U))))) Gauss linear;
} }
......
...@@ -39,7 +39,7 @@ solvers ...@@ -39,7 +39,7 @@ solvers
relTol 0; relTol 0;
} }
"(U|k|B|nuTilda)" "(U|k|B|nuTilda|s)"
{ {
solver smoothSolver; solver smoothSolver;
smoother GaussSeidel; smoother GaussSeidel;
......
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment