70 store(fieldName_, tfldPtr);
77Foam::tmp<Foam::volScalarField>
78Foam::functionObjects::energyTransport::kappaEff()
const
84 const turbType* turbPtr = findObject<turbType>
109Foam::tmp<Foam::volScalarField>
110Foam::functionObjects::energyTransport::rho()
const
122 trho.ref() = lookupObject<volScalarField>(rhoName_);
128Foam::tmp<Foam::volScalarField>
129Foam::functionObjects::energyTransport::Cp()
const
133 tmp<volScalarField>
tCp(phases_[0]*Cps_[0]);
135 for (label i = 1; i < phases_.size(); i++)
137 tCp.ref() += phases_[i]*Cps_[i];
152Foam::tmp<Foam::volScalarField>
153Foam::functionObjects::energyTransport::kappa()
const
157 tmp<volScalarField> tkappa(phases_[0]*kappas_[0]);
159 for (label i = 1; i < phases_.size(); i++)
161 tkappa.ref() += phases_[i]*kappas_[i];
201 multiphaseThermo_(
dict.subOrEmptyDict(
"phaseThermos")),
212 fieldName_(
dict.getOrDefault<
word>(
"field",
"T")),
213 schemesField_(
"unknown-schemesField"),
214 phiName_(
dict.getOrDefault<
word>(
"phi",
"phi")),
215 rhoName_(
dict.getOrDefault<
word>(
"rho",
"rho")),
222 if (!multiphaseThermo_.
empty())
224 Cps_.setSize(multiphaseThermo_.
size());
225 kappas_.setSize(Cps_.size());
226 phaseNames_.
setSize(Cps_.size());
231 const word& key = iter().keyword();
238 <<
"Found non-dictionary entry " << iter()
239 <<
" in top-level dictionary " << multiphaseThermo_
245 phaseNames_[
phasei] = key;
272 phases_.setSize(phaseNames_.
size());
287 if (Cp_.value() == 0.0 || kappa_.value() == 0.0)
290 <<
" Multiphase thermo dictionary not found and Cp/kappa "
291 <<
" for single phase are zero. Please entry either"
300 s.correctBoundaryConditions();
313 dict.readIfPresent(
"phi", phiName_);
314 dict.readIfPresent(
"rho", rhoName_);
316 schemesField_ =
dict.getOrDefault(
"schemesField", fieldName_);
318 dict.readIfPresent(
"nCorr", nCorr_);
319 dict.readIfPresent(
"tolerance", tol_);
321 if (
dict.found(
"fvOptions"))
323 fvOptions_.reset(
dict.subDict(
"fvOptions"));
342 word divScheme(
"div(phi," + schemesField_ +
")");
345 "laplacian(kappaEff," + schemesField_ +
")"
349 scalar relaxCoeff = 0;
350 mesh_.relaxEquation(schemesField_, relaxCoeff);
353 bool converged =
false;
361 for (
int i = 0; i <= nCorr_; ++i)
370 fvOptions_(rhoCp_,
s)
373 sEqn.
relax(relaxCoeff);
375 fvOptions_.constrain(sEqn);
378 converged = (sEqn.
solve(schemesField_).initialResidual() < tol_);
379 if (converged)
break;
396 for (
int i = 0; i <= nCorr_; ++i)
404 fvOptions_(trhoCp.ref(),
s)
407 sEqn.relax(relaxCoeff);
409 fvOptions_.constrain(sEqn);
412 converged = (sEqn.solve(schemesField_).initialResidual() < tol_);
413 if (converged)
break;
419 <<
"Incompatible dimensions for phi: " <<
phi.dimensions() <<
nl
427 <<
s.name() <<
" is converged." <<
nl
428 <<
tab <<
"initial-residual tolerance: " << tol_ <<
nl
429 <<
tab <<
"outer iteration: " << iter <<
nl;
441 <<
tab << fieldName_ <<
nl
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
bool empty() const noexcept
True if the list is empty.
label size() const noexcept
The number of elements in list.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
@ REGISTER
Request registration (bool: true).
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
void setSize(label n)
Alias for resize().
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
void size(const label n)
Older name for setAddressableSize.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
const dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary pointer if present (and it is a dictionary) otherwise return nullptr...
Abstract base-class for Time/database function objects.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
Computes the simplified energy transport equation in single-phase or two-phase flow,...
virtual bool execute()
Execute the function-object operations.
virtual bool write()
Write the function object output.
energyTransport(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
virtual bool read(const dictionary &)
Read the function-object dictionary.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
bool foundObject(const word &fieldName) const
Find object (eg, a field) in the (sub) objectRegistry.
bool store(word &fieldName, const tmp< ObjectType > &tfield, bool cacheable=false)
Store the field in the (sub) objectRegistry under the given name.
ObjectType & lookupObjectRef(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
const Time & time() const
Return time database.
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
const Time & time() const
Return the top-level database.
Type * getObjectPtr(const word &name, const bool recursive=false) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
const tmp< volScalarField > & tCp
const volScalarField & Cp
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
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))
const surfaceScalarField rhoCpPhi(fvc::interpolate(fluid.Cp()) *rhoPhi)
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
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)
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
IncompressibleTurbulenceModel< transportModel > turbulenceModel
const dimensionSet dimless
Dimensionless.
const dimensionSet dimEnergy
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
fvMatrix< scalar > fvScalarMatrix
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.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
const dimensionSet dimVolume(pow3(dimLength))
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
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).
tmp< volScalarField > trho
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.