Go to the documentation of this file.
45 namespace functionObjects
59 inletOutletFvPatchField<scalar>::typeName
66 result[patchi] = zeroGradientFvPatchField<scalar>::typeName;
74 bool Foam::functionObjects::age::converged
77 const scalar initialResidual
80 if (initialResidual < tolerance_)
82 Info<<
"Field " << typeName
83 <<
" converged in " <<
nCorr <<
" correctors"
93 template<
class GeoField>
95 Foam::functionObjects::age::newField
105 scopedName(baseName),
139 phiName_ =
dict.getOrDefault<
word>(
"phi",
"phi");
140 rhoName_ =
dict.getOrDefault<
word>(
"rho",
"rho");
141 schemesField_ =
dict.getOrDefault<
word>(
"schemesField", typeName);
142 tolerance_ =
dict.getOrDefault<scalar>(
"tolerance", 1
e-5);
143 nCorr_ =
dict.getOrDefault<
int>(
"nCorr", 5);
144 diffusion_ =
dict.getOrDefault<
bool>(
"diffusion",
false);
160 mesh_.time().timeName(),
172 const word divScheme(
"div(phi," + schemesField_ +
")");
175 scalar relaxCoeff = 0;
176 if (mesh_.relaxEquation(schemesField_))
178 relaxCoeff = mesh_.equationRelaxationFactor(schemesField_);
197 word laplacianScheme;
208 "laplacian(" + tmuEff().name() +
',' + schemesField_ +
")";
211 for (
int i = 0; i <= nCorr_; ++i)
223 ageEqn.
relax(relaxCoeff);
227 if (converged(i, ageEqn.
solve().initialResidual()))
238 word laplacianScheme;
249 "laplacian(" + tnuEff().name() +
',' + schemesField_ +
")";
252 for (
int i = 0; i <= nCorr_; ++i)
265 ageEqn.
relax(relaxCoeff);
269 if (converged(i, ageEqn.
solve().initialResidual()))
278 Info<<
"Min/max age:"
279 <<
min(
age).value() <<
' '
284 word fieldName = typeName;
286 return store(fieldName, tage);
292 return writeObject(typeName);
static autoPtr< T > New(Args &&... args)
Construct autoPtr of T with forwarding arguments.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A class for handling words, derived from Foam::string.
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
bool read(const char *buf, int32_t &val)
Same as readInt32.
static const word propertiesName
Default name of the turbulence properties dictionary.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Calculate the matrix for the divergence of the given field and flux.
#define forAll(list, i)
Loop across all elements in list.
wordList patchTypes(nPatches)
List< word > wordList
A List of words.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
virtual bool read(const dictionary &)
Read the data.
messageStream Info
Information stream (stdout output on master, null elsewhere)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
virtual bool write()
Write.
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Calculate the matrix for the laplacian of the field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
age(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
virtual bool read(const dictionary &dict)
Read optional controls.
Macros for easy insertion into run-time selection tables.
Generic dimensioned Type class.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Calculates and writes out the time taken for a particle to travel from an inlet to the location....
virtual bool execute()
Execute.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
const fvMesh & mesh_
Reference to the fvMesh.
const dimensionedScalar e
Elementary charge.
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
defineTypeNameAndDebug(ObukhovLength, 0)
const polyBoundaryMesh & patches
A special matrix type and solver, designed for finite volume solutions of scalar equations....
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Templated abstract base class for single-phase incompressible turbulence models.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)