coefficientPhaseChange.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) 2021-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 
26 #include "coefficientPhaseChange.H"
27 #include "fvcGrad.H"
28 #include "multicomponentThermo.H"
30 
31 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace fv
36 {
39 }
40 }
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 void Foam::fv::coefficientPhaseChange::readCoeffs(const dictionary& dict)
45 {
47 
48  C_.read(dict);
49 }
50 
51 
53 Foam::fv::coefficientPhaseChange::timesY1
54 (
55  tmp<volScalarField::Internal> mDot
56 ) const
57 {
58  const ThermoRefPair<multicomponentThermo> mcThermos =
59  thermos().thermos<multicomponentThermo>();
60 
61  if (!mcThermos.valid().first() || species().empty())
62  {
63  return mDot;
64  }
65 
66  if (species().size() == 1)
67  {
68  return mcThermos.first().Y()[specieis().first()]*mDot;
69  }
70 
71  tmp<volScalarField::Internal> tY1 =
73  (
74  typedName("Y1"),
75  mcThermos.first().Y()[specieis(0).first()]
76  );
77 
78  for (label mDoti = 1; mDoti < species().size(); ++ mDoti)
79  {
80  tY1.ref() += mcThermos.first().Y()[specieis(mDoti).first()];
81  }
82 
83  return tY1*mDot;
84 }
85 
86 
88 Foam::fv::coefficientPhaseChange::mDotByAlpha1Y1() const
89 {
90  return C_*mag(fvc::grad(alpha1_))()();
91 }
92 
93 
94 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
95 
97 (
98  const word& name,
99  const word& modelType,
100  const fvMesh& mesh,
101  const dictionary& dict
102 )
103 :
105  (
106  name,
107  modelType,
108  mesh,
109  dict,
110  readSpecies(coeffs(modelType, dict), false)
111  ),
112  C_("C", dimMass/dimArea/dimTime, NaN),
113  alpha1_(mesh().lookupObject<volScalarField>(alphaNames().first()))
114 {
115  readCoeffs(coeffs(dict));
116 }
117 
118 
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120 
123 {
124  return timesY1(alpha1_*mDotByAlpha1Y1());
125 }
126 
127 
130 {
131  const ThermoRefPair<multicomponentThermo> mcThermos =
132  thermos().thermos<multicomponentThermo>();
133 
134  tmp<volScalarField::Internal> tmDot = alpha1_*mDotByAlpha1Y1();
135 
136  if (mcThermos.valid().first())
137  {
138  const labelPair specieis = this->specieis(mDoti);
139 
140  tmDot.ref() *= mcThermos.first().Y()[specieis.first()];
141  }
142 
143  return tmDot;
144 }
145 
146 
148 (
149  const volScalarField& alpha,
150  const volScalarField& rho,
151  fvMatrix<scalar>& eqn
152 ) const
153 {
154  const label i = index(alphaNames(), eqn.psi().name());
155 
156  if (i != -1)
157  {
158  const volScalarField::Internal mDotByAlpha1(timesY1(mDotByAlpha1Y1()));
159 
160  if (i == 0)
161  {
162  eqn -= fvm::Sp(mDotByAlpha1, eqn.psi());
163  }
164  else
165  {
166  eqn +=
167  mDotByAlpha1*alpha1_
168  - correction(fvm::Sp(mDotByAlpha1, eqn.psi()));
169  }
170  }
171  else
172  {
174  }
175 }
176 
177 
179 (
180  const volScalarField& alpha,
181  const volScalarField& rho,
182  const volScalarField& Yi,
183  fvMatrix<scalar>& eqn
184 ) const
185 {
186  const label i = index(phaseNames(), eqn.psi().group());
187 
188  const ThermoRefPair<multicomponentThermo>& mcThermos =
189  thermos().thermos<multicomponentThermo>();
190 
191  if
192  (
193  i == 0
194  && mcThermos.valid().first()
195  && mcThermos.first().containsSpecie(Yi.member())
196  && species().found(Yi.member())
197  )
198  {
199  eqn -= fvm::Sp(alpha1_*mDotByAlpha1Y1(), Yi);
200  }
201  else
202  {
203  phaseChange::addSup(alpha, rho, Yi, eqn);
204  }
205 }
206 
207 
209 {
210  if (phaseChange::read(dict))
211  {
212  readCoeffs(coeffs(dict));
213  return true;
214  }
215  else
216  {
217  return false;
218  }
219 }
220 
221 
222 // ************************************************************************* //
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,.
Generic GeometricField class.
static word member(const word &name)
Return member (name without the extension)
Definition: IOobject.C:146
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.
ThermoRefPair< OtherThermoType > thermos() const
Cast to a different thermo type, with error checking.
const ThermoType & first() const
Access the first thermo.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
void read(const dictionary &, const unitConversion &defaultUnits=NullObjectRef< unitConversion >())
Update the value of dimensioned<Type>
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
VolField< Type > & psi()
Definition: fvMatrix.H:289
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
This simple model generates a phase change between two phases calculated from the following expressio...
virtual tmp< DimensionedField< scalar, volMesh > > mDot() const
Return the total phase change rate.
virtual bool read(const dictionary &dict)
Read source dictionary.
coefficientPhaseChange(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
void addSup(const volScalarField &alpha, const volScalarField &rho, fvMatrix< scalar > &eqn) const
Override the compressible continuity equation to add.
virtual void addSup(fvMatrix< scalar > &eqn) const
Add a source term to a field-less proxy equation.
Definition: massTransfer.C:225
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 bool read(const dictionary &dict)
Read source dictionary.
Definition: phaseChange.C:624
void addSup(const volScalarField &alpha, const volScalarField &rho, const volScalarField &heOrYi, fvMatrix< scalar > &eqn) const
Override the energy equation to add the phase change heat, or.
Definition: phaseChange.C:551
Base-class for multi-component thermodynamic properties.
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)
Calculate the gradient of the given field.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
Namespace for OpenFOAM.
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
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
const dimensionSet dimTime
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:181
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimMass
const dimensionSet dimArea
labelList fv(nPoints)
dictionary dict