42 #define PhiDimensionErrorInFunction(phi) \
43 FatalErrorInFunction \
44 << "Incompatible dimensions for " << phi.name() << ": " \
45 << phi.dimensions() << nl \
46 << "Dimensions should be " << dimMass/dimTime << " or " \
47 << dimVolume/dimTime << exit(FatalError)
53 namespace functionObjects
81 PhiPatchFieldTypes[
patchi] =
82 p.boundaryField()[
patchi].fixesValue()
113 Foam::functionObjects::phaseScalarTransport::alphaPhi()
126 tmp<surfaceScalarField> tAlphaPhi
140 const word laplacianScheme =
"laplacian(" + pName_ +
")";
146 auto writeDDt = [&](
const label i)
159 Info<<
type() <<
": Writing " << DDtAlpha.name() <<
endl;
162 if (debug && time_.writeTime())
168 nonOrthogonalSolutionControl& control =
169 mesh_.lookupObjectRef<nonOrthogonalSolutionControl>
171 solutionControl::typeName
177 while (control.correctNonOrthogonal())
186 PhiEqn.solve(pName_);
188 if (control.finalNonOrthogonalIter())
190 alphaPhi += PhiEqn.flux();
199 while (control.correctNonOrthogonal())
208 PhiEqn.solve(pName_);
210 if (control.finalNonOrthogonalIter())
212 alphaPhi += PhiEqn.flux();
222 if (debug && time_.writeTime())
232 Foam::functionObjects::phaseScalarTransport::D
237 const word Dname(
"D" + s_.name());
250 const word& nameNoPhase = momentumTransportModel::typeName;
282 fieldName_(
dict.lookup(
"field")),
301 <<
"Field \"" << fieldName_ <<
"\" does not have a phase extension "
302 <<
"in its name. If it is associated with \"phaseA\" then it "
303 <<
"should be named \"" << fieldName_ <<
".phaseA\"."
323 solveAlphaPhi_ =
dict.lookupOrDefault<
bool>(
"solveAlphaPhi",
false);
331 const word defaultAlphaPhiName =
336 :
dict.lookupOrDefault<
word>(
"alphaPhi", defaultAlphaPhiName);
337 phiName_ =
dict.lookupOrDefault<
word>(
"phi",
"phi");
344 pName_ =
dict.lookupOrDefault<
word>(
"p",
"p");
345 schemesField_ =
dict.lookupOrDefault<
word>(
"schemesField", fieldName_);
356 dict.lookup(
"D") >> D_;
361 dict.lookupBackwardsCompatible<scalar>({
"alphal",
"alphaD"});
363 dict.lookupBackwardsCompatible<scalar>({
"alphat",
"alphaDt"});
367 nCorr_ =
dict.lookupOrDefaultBackwardsCompatible<
label>
369 {
"nCorrectors",
"nCorr"},
372 residualAlpha_ =
dict.lookupOrDefault<scalar>(
"residualAlpha", rootSmall);
373 writeAlphaField_ =
dict.lookupOrDefault<
bool>(
"writeAlphaField",
true);
381 return wordList{alphaName_, alphaPhiName_, phiName_, pName_};
397 const scalar relaxCoeff =
398 mesh_.solution().relaxEquation(schemesField_)
399 ? mesh_.solution().equationRelaxationFactor(schemesField_)
409 for (
int i=0; i<=nCorr_; i++)
418 "div(" + alphaPhi.
name() +
"," + schemesField_ +
")"
435 "laplacian(" +
D.
name() +
"," + schemesField_ +
")"
439 fieldEqn.
relax(relaxCoeff);
441 fieldEqn.
solve(schemesField_);
450 for (
int i=0; i<=nCorr_; i++)
459 "div(" + alphaPhi.
name() +
"," + schemesField_ +
")"
476 "laplacian(" +
D.
name() +
"," + schemesField_ +
")"
480 fieldEqn.
relax(relaxCoeff);
482 fieldEqn.
solve(schemesField_);
501 if (writeAlphaField_)
510 "alpha" + fieldName_.capitalise(),
#define forAll(list, i)
Loop across all elements in list.
Macros for easy insertion into run-time selection tables.
const dimensionSet & dimensions() const
Return dimensions.
static const char *const typeName
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
word group() const
Return group (extension part of name)
word member() const
Return member (name without the extension)
const word & name() const
Return name.
static word groupName(Name name, const word &group)
Enum read(Istream &) const
Read a word from Istream and return the corresponding.
virtual Ostream & write(const char)=0
Write character.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
label size() const
Return the number of elements in the UPtrList.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
const word & name() const
Return const reference to name.
Abstract base-class for Time/database functionObjects.
const Time & time_
Reference to time.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
virtual bool read(const dictionary &)
Read optional controls.
Evolves a passive scalar transport equation within one phase of a multiphase simulation....
virtual wordList fields() const
Return the list of fields required.
virtual ~phaseScalarTransport()
Destructor.
virtual bool execute()
Solve for the evolution of the field.
virtual bool write()
Do nothing. The field is registered and written automatically.
virtual bool read(const dictionary &)
Read the settings from the given dictionary.
phaseScalarTransport(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
static const NamedEnum< diffusivityType, 3 > diffusivityTypeNames_
Diffusivity type names.
Finite volume constraints.
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
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).
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
const fvSchemes & schemes() const
Return the fvSchemes.
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
void setFluxRequired(const word &name) const
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
A class for managing temporary objects.
A class for handling words, derived from string.
static const word null
An empty word.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Calculate the first temporal derivative.
Calculate the divergence of the given field.
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.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
compressibleMomentumTransportModel momentumTransportModel
const char *const group
Group name for atomic constants.
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
static const coefficient D("D", dimTemperature, 257.14)
errorManipArg< error, int > exit(error &err, const int errNo=1)
List< word > wordList
A List of words.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & endl(Ostream &os)
Add newline and flush stream.
fvMatrix< scalar > fvScalarMatrix
const dimensionSet dimKinematicViscosity
SurfaceField< scalar > surfaceScalarField
const dimensionSet dimLength
const dimensionSet dimTime
const dimensionSet dimVolume
VolField< scalar > volScalarField
word typedName(Name name)
Return the name of the object within the given type.
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimMass
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
#define PhiDimensionErrorInFunction(phi)
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))