59void Foam::functionObjects::BilgerMixtureFraction::calcBilgerMixtureFraction()
81 auto& f_Bilger = *resultPtr;
83 auto&
Y = thermo_.
Y();
85 f_Bilger = -o2RequiredOx_;
90 *(nAtomsC_[i] + nAtomsS_[i] + 0.25*nAtomsH_[i] - 0.5*nAtomsO_[i])
93 f_Bilger /= o2RequiredFuelOx_;
94 f_Bilger.clamp_range(zero_one{});
98bool Foam::functionObjects::BilgerMixtureFraction::readComposition
104 auto&
Y = thermo_.Y();
111 subDict.getCheckOrDefault<scalar>
119 if (
sum(comp) < SMALL)
122 <<
"No composition is given" <<
nl
123 <<
"Valid species are:" <<
nl
130 const word fractionBasisType
132 subDict.getOrDefault<word>(
"fractionBasis",
"mass")
135 if (fractionBasisType ==
"mass")
140 else if (fractionBasisType ==
"mole")
147 comp[i] *= thermo_.W(i);
155 <<
"The given fractionBasis type is invalid" <<
nl
156 <<
"Valid fractionBasis types are" <<
nl
157 <<
" \"mass\" (default)" <<
nl
168Foam::scalar Foam::functionObjects::BilgerMixtureFraction::o2Required
177 comp[i]/thermo_.W(i)*(nAtomsC_[i] + nAtomsS_[i] + 0.25*nAtomsH_[i]);
184Foam::scalar Foam::functionObjects::BilgerMixtureFraction::o2Present
192 o2pres += comp[i]/thermo_.W(i)*nAtomsO_[i];
209 phaseName_(
dict.getOrDefault<
word>(
"phase",
word::null)),
215 IOobject::groupName(
"f_Bilger", phaseName_)
225 nSpecies_(thermo_.
Y().size()),
227 o2RequiredFuelOx_(0),
228 nAtomsC_(nSpecies_, 0),
229 nAtomsS_(nSpecies_, 0),
230 nAtomsH_(nSpecies_, 0),
231 nAtomsO_(nSpecies_, 0),
232 Yoxidiser_(nSpecies_, 0),
237 calcBilgerMixtureFraction();
260 nSpecies_ = thermo_.Y().size();
265 <<
"Number of input species is zero"
269 nAtomsC_.resize(nSpecies_, 0);
270 nAtomsS_.resize(nSpecies_, 0);
271 nAtomsH_.resize(nSpecies_, 0);
272 nAtomsO_.resize(nSpecies_, 0);
274 auto&
Y = thermo_.Y();
277 typedef BasicChemistryModel<psiReactionThermo> psiChemistryModelType;
278 typedef BasicChemistryModel<rhoReactionThermo> rhoChemistryModelType;
280 const auto* psiChemPtr =
281 mesh_.cfindObject<psiChemistryModelType>(
"chemistryProperties");
283 const auto* rhoChemPtr =
284 mesh_.cfindObject<rhoChemistryModelType>(
"chemistryProperties");
286 autoPtr<speciesCompositionTable> speciesCompPtr;
290 speciesCompPtr.reset((*psiChemPtr).thermo().specieComposition());
294 speciesCompPtr.reset((*rhoChemPtr).thermo().specieComposition());
299 <<
"BasicChemistryModel not found"
305 const List<specieElement>& curSpecieComposition =
306 (speciesCompPtr.ref())[speciesTab[i]];
308 forAll(curSpecieComposition, j)
310 const word&
e = curSpecieComposition[j].name();
311 const label nAtoms = curSpecieComposition[j].nAtoms();
315 nAtomsC_[i] = nAtoms;
319 nAtomsS_[i] = nAtoms;
323 nAtomsH_[i] = nAtoms;
327 nAtomsO_[i] = nAtoms;
332 if (
sum(nAtomsO_) == 0)
335 <<
"No specie contains oxygen"
339 Yoxidiser_.
resize(nSpecies_, 0);
340 Yfuel_.resize(nSpecies_, 0);
344 !readComposition(
dict.subDict(
"oxidiser"), Yoxidiser_)
345 || !readComposition(
dict.subDict(
"fuel"), Yfuel_)
351 o2RequiredOx_ = o2Required(Yoxidiser_) - o2Present(Yoxidiser_);
353 if (o2RequiredOx_ > 0)
356 <<
"Oxidiser composition contains not enough oxygen" <<
endl
357 <<
"Mixed up fuel and oxidiser compositions?"
361 const scalar o2RequiredFuel = o2Required(Yfuel_) - o2Present(Yfuel_);
363 if (o2RequiredFuel < 0)
366 <<
"Fuel composition contains too much oxygen" <<
endl
367 <<
"Mixed up fuel and oxidiser compositions?"
371 o2RequiredFuelOx_ = o2RequiredFuel - o2RequiredOx_;
373 if (
mag(o2RequiredFuelOx_) < SMALL)
376 <<
"Fuel and oxidiser have the same composition"
386 calcBilgerMixtureFraction();
401 <<
" writing field " << resultName_ <<
endl;
403 return writeObject(resultName_);
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Basic chemistry model templated on thermodynamics.
void clamp_range(const dimensioned< MinMax< Type > > &range)
Clamp field values (in-place) to the specified range.
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void resize(const label len)
Adjust allocated size of list.
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.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
Specialization of basicMultiComponentMixture for a mixture consisting of a number for molecular speci...
virtual scalar W(const label speciei) const =0
Molecular weight of the given specie [kg/kmol].
Abstract base-class for fluid and solid thermodynamic properties.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
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.
Calculates the Bilger mixture-fraction field (i.e. ) based on the elemental composition of the mixtur...
BilgerMixtureFraction(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
virtual bool clear()
Clear the Bilger mixture-fraction field from registry.
virtual bool execute()
Calculate the Bilger mixture-fraction field.
virtual bool write()
Write Bilger mixture-fraction field.
virtual bool read(const dictionary &)
Read the BilgerMixtureFraction data.
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.
const ObjectType & lookupObject(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
bool writeObject(const word &fieldName)
Write field if present in the (sub) objectRegistry.
bool clearObject(const word &fieldName)
Clear field from the (sub) objectRegistry if present.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
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...
bool store()
Register object with its registry and transfer ownership to the registry.
static constexpr scalarRange ge0() noexcept
A greater-equals zero bound.
A class for handling words, derived from Foam::string.
static const word null
An empty word.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
PtrList< volScalarField > & Y
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const word dictName("faMeshDefinition")
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
hashedWordList speciesTable
A table of species as a hashedWordList.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
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...
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.
errorManipArg< error, int > exit(error &err, const int errNo=1)
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.
Type definitions for thermo-physics models.