44namespace functionObjects
62 if (!foundObject<volScalarField>(fieldName_))
76 store(fieldName_, tfldPtr);
78 if (phaseName_ !=
"none")
84 return lookupObjectRef<volScalarField>(fieldName_);
94 const word Dname(
"D" +
s.name());
103 mesh_.time().timeName(),
113 if (nutName_ !=
"none")
124 findObject<incompressible::turbulenceModel>
134 alphaD_ *
turb->nu() + alphaDt_ *
turb->nut()
142 findObject<compressible::turbulenceModel>
152 alphaD_ *
turb->mu() + alphaDt_ *
turb->mut()
163 mesh_.time().timeName(),
184 fieldName_(
dict.getOrDefault<
word>(
"field",
"s")),
185 phiName_(
dict.getOrDefault<
word>(
"phi",
"phi")),
186 rhoName_(
dict.getOrDefault<
word>(
"rho",
"rho")),
187 nutName_(
dict.getOrDefault<
word>(
"nut",
"none")),
188 phaseName_(
dict.getOrDefault<
word>(
"phase",
"none")),
189 phasePhiCompressedName_
191 dict.getOrDefault<
word>(
"phasePhiCompressed",
"alphaPhiUn")
196 resetOnStartUp_(false),
197 schemesField_(
"unknown-schemesField"),
199 bounded01_(
dict.getOrDefault(
"bounded01", true))
226 dict.readIfPresent(
"phi", phiName_);
227 dict.readIfPresent(
"rho", rhoName_);
228 dict.readIfPresent(
"nut", nutName_);
229 dict.readIfPresent(
"phase", phaseName_);
230 dict.readIfPresent(
"bounded01", bounded01_);
232 schemesField_ =
dict.getOrDefault(
"schemesField", fieldName_);
233 constantD_ =
dict.readIfPresent(
"D", D_);
234 alphaD_ =
dict.getOrDefault<scalar>(
"alphaD", 1);
235 alphaDt_ =
dict.getOrDefault<scalar>(
"alphaDt", 1);
237 dict.readIfPresent(
"nCorr", nCorr_);
238 dict.readIfPresent(
"resetOnStartUp", resetOnStartUp_);
240 if (
dict.found(
"fvOptions"))
242 fvOptions_.reset(
dict.subDict(
"fvOptions"));
261 word divScheme(
"div(phi," + schemesField_ +
")");
262 word laplacianScheme(
"laplacian(" +
D.name() +
"," + schemesField_ +
")");
265 scalar relaxCoeff = 0.0;
266 if (mesh_.relaxEquation(schemesField_))
268 relaxCoeff = mesh_.equationRelaxationFactor(schemesField_);
272 if (phaseName_ !=
"none")
283 D.dimensions().reset(limitedPhiAlpha.dimensions()/
dimLength);
287 for (label i = 0; i <= nCorr_; i++)
292 +
fvm::div(limitedPhiAlpha,
s, divScheme)
298 sEqn.
relax(relaxCoeff);
299 fvOptions_.constrain(sEqn);
300 sEqn.
solve(mesh_.solverDict(schemesField_));
302 tTPhiUD = sEqn.
flux();
322 for (label i = 0; i <= nCorr_; i++)
334 sEqn.
relax(relaxCoeff);
336 fvOptions_.constrain(sEqn);
338 sEqn.
solve(mesh_.solverDict(schemesField_));
343 for (label i = 0; i <= nCorr_; i++)
354 sEqn.
relax(relaxCoeff);
356 fvOptions_.constrain(sEqn);
358 sEqn.
solve(mesh_.solverDict(schemesField_));
364 <<
"Incompatible dimensions for phi: " <<
phi.dimensions() <<
nl
378 <<
tab << fieldName_ <<
nl
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
compressible::turbulenceModel & turb
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
static word timeName(const scalar t, const int precision=precision_)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Abstract base-class for Time/database function objects.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
bool store(word &fieldName, const tmp< ObjectType > &tfield, bool cacheable=false)
Store the field in the (sub) objectRegistry under the given name.
Evolves a passive scalar transport equation.
virtual ~scalarTransport()
Destructor.
virtual bool execute()
Calculate the scalarTransport.
virtual bool write()
Do nothing.
virtual bool read(const dictionary &)
Read the scalarTransport data.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
const Time & time() const
Return the top-level database.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
A class representing the concept of a field of 1 used to avoid unnecessary manipulations for objects ...
void setFluxRequired(const word &name) const
Get flux-required for given name, or default.
A class for managing temporary objects.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for the laplacian of the field.
Calculate the finiteVolume matrix for implicit and explicit sources.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar pos(const dimensionedScalar &ds)
const dimensionSet dimless
Dimensionless.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Ostream & endl(Ostream &os)
Add newline and flush stream.
static constexpr const zero Zero
Global zero (0)
const dimensionSet dimVolume(pow3(dimLength))
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
errorManipArg< error, int > exit(error &err, const int errNo=1)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
constexpr char nl
The newline '\n' character (0x0a)
constexpr char tab
The tab '\t' character(0x09)
const dimensionedScalar & D