Go to the documentation of this file.
54 wordList bTypes(
U.boundaryField().size());
59 if (isA<fixedValueFvPatchVectorField>(pf))
61 bTypes[patchI] = fixedValueFvPatchScalarField::typeName;
65 bTypes[patchI] = zeroGradientFvPatchScalarField::typeName;
82 mesh_.time().timeName(),
118 mesh_.time().timeName(),
130 const icoModel& model = mesh_.lookupObject<icoModel>
135 return model.nuEff();
139 const cmpModel& model = mesh_.lookupObject<cmpModel>
144 return model.muEff();
155 mesh_.time().timeName(),
175 const bool loadFromFiles
179 mesh_(refCast<const fvMesh>(obr)),
186 resetOnStartUp_(
false),
217 log_.readIfPresent(
"log",
dict);
219 if (log_)
Info<<
type() <<
" " << name_ <<
" output:" <<
nl;
232 dict.
lookup(
"resetOnStartUp") >> resetOnStartUp_;
244 if (log_)
Info<<
type() <<
" " << name_ <<
" output:" <<
nl;
255 word schemeVar =
T.name();
261 word divScheme(
"div(phi," + schemeVar +
")");
262 word laplacianScheme(
"laplacian(" + DT.name() +
"," + schemeVar +
")");
265 scalar relaxCoeff = 0.0;
266 if (mesh_.relaxEquation(schemeVar))
268 relaxCoeff = mesh_.equationRelaxationFactor(schemeVar);
277 for (
label i = 0; i <= nCorr_; i++)
290 fvOptions_.constrain(
TEqn);
298 for (
label i = 0; i <= nCorr_; i++)
311 fvOptions_.constrain(
TEqn);
319 <<
"Incompatible dimensions for phi: " <<
phi.dimensions() <<
nl
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
const fvMesh & mesh_
Reference to the mesh database.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
virtual void write()
Calculate the scalarTransport and write.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
A class for handling words, derived from string.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
bool read(const char *, int32_t &)
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
static const word propertiesName
Default name of the turbulence properties dictionary.
Calculate the divergence of the given field.
wordList boundaryTypes() const
Return the boundary types for the scalar field.
word UName_
Name of velocity field (optional)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Calculate the matrix for the divergence of the given field and flux.
Registry of regIOobjects.
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
scalarTransport(const scalarTransport &)
Disallow default bitwise copy construct.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
virtual void end()
Execute at the final time-loop, currently does nothing.
virtual void execute()
Execute, currently does nothing.
volScalarField & transportedField()
Return reference to registered transported field.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Calculate the matrix for the laplacian of the field.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
A list of keyword definitions, which are a keyword followed by any number of values (e....
A scalar instance of fvMatrix.
tmp< volScalarField > DT(const surfaceScalarField &phi) const
Return the diffusivity field.
Generic dimensioned Type class.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Calculate the matrix for implicit and explicit sources.
fvScalarMatrix TEqn(fvm::div(phi, T) - fvm::laplacian(alphaEff, T)==radiation->ST(rhoCpRef, T)+fvOptions(T))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
virtual void read(const dictionary &)
Read the scalarTransport data.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
virtual void timeSet()
Called when time was set at the end of the Time::operator++.
Templated abstract base class for single-phase incompressible turbulence models.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
const dimensionSet dimVolume(pow3(dimLength))
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Generic GeometricField class.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
defineTypeNameAndDebug(combustionModel, 0)
virtual ~scalarTransport()
Destructor.
word name(const complex &)
Return a string representation of a complex.
Calulate the matrix for the first temporal derivative.