reactionDrivenPhaseChange.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "multicomponentThermo.H"
29 
30 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace fv
35 {
38  (
39  fvModel,
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 void Foam::fv::reactionDrivenPhaseChange::readCoeffs(const dictionary& dict)
50 {
52 }
53 
54 
55 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
56 
58 (
59  const word& name,
60  const word& modelType,
61  const fvMesh& mesh,
62  const dictionary& dict
63 )
64 :
66  (
67  name,
68  modelType,
69  mesh,
70  dict,
71  readSpecies(coeffs(modelType, dict), true)
72  ),
73  fluid_(mesh().lookupObject<phaseSystem>(phaseSystem::propertiesName)),
74  phase1_(fluid_.phases()[phaseNames().first()]),
75  phase2_(fluid_.phases()[phaseNames().second()])
76 {
77  readCoeffs(coeffs(dict));
78 
79  const ThermoRefPair<multicomponentThermo> mcThermos =
80  thermos().thermos<multicomponentThermo>();
81 
82  if (!mcThermos.either())
83  {
85  << "Model " << name << " of type " << modelType
86  << " applied to two pure phases " << phase1_.name()
87  << " and " << phase2_.name() << exit(FatalError);
88  }
89 }
90 
91 
92 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93 
96 {
97  const volScalarField::Internal& alpha1 = phase1_, alpha2 = phase2_;
98 
99  const labelPair specieis = this->specieis(mDoti);
100 
101  const ThermoRefPair<multicomponentThermo> mcThermos =
102  thermos().thermos<multicomponentThermo>();
103 
106  (
107  name() + ":mDot_" + species()[mDoti],
108  mesh(),
110  );
111 
112  if (mcThermos.valid().first())
113  {
114  tResult.ref() += alpha1*phase1_.R(specieis.first());
115  }
116 
117  if (mcThermos.valid().second())
118  {
119  tResult.ref() -= alpha2*phase2_.R(specieis.second());
120  }
121 
122  return tResult;
123 }
124 
125 
127 {
128  if (phaseChange::read(dict))
129  {
130  readCoeffs(coeffs(dict));
131  return true;
132  }
133  else
134  {
135  return false;
136  }
137 }
138 
139 
140 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static tmp< DimensionedField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Mesh &mesh, const dimensionSet &, const PrimitiveField< Type > &)
Return a temporary field constructed from name, mesh,.
const Type & second() const
Return second.
Definition: PairI.H:121
const Type & first() const
Return first.
Definition: PairI.H:107
Class containing a pair of thermo references. Handles down-casting to more specific thermo types by c...
Definition: ThermoRefPair.H:51
const Pair< bool > & valid() const
Access the validity flags.
bool either() const
Return if either validity flag is set.
ThermoRefPair< OtherThermoType > thermos() const
Cast to a different thermo type, with error checking.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
Finite volume model abstract base class.
Definition: fvModel.H:60
static const dictionary & coeffs(const word &modelType, const dictionary &)
Return the coefficients sub-dictionary for a given model type.
Definition: fvModelI.H:31
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:57
Base class for phase change models.
Definition: phaseChange.H:61
void reReadSpecies(const dictionary &dict) const
Re-read the names of the transferring species.
Definition: phaseChange.C:140
virtual tmp< DimensionedField< scalar, volMesh > > mDot() const
Return the total phase change rate.
Definition: phaseChange.C:498
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: phaseChange.C:624
const ThermoRefPair< basicThermo > & thermos() const
Return the thermo references.
Definition: phaseChangeI.H:32
Model for mass-diffusion rate limited phase change between two phases.
virtual bool read(const dictionary &dict)
Read source dictionary.
reactionDrivenPhaseChange(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
Base-class for multi-component thermodynamic properties.
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:145
Class to represent a system of phases.
Definition: phaseSystem.H:74
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define WarningInFunction
Report a warning using Foam::Warning.
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
labelList second(const UList< labelPair > &p)
Definition: patchToPatch.C:49
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
const dimensionSet dimTime
const dimensionSet dimDensity
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
labelList fv(nPoints)
dictionary dict