multicomponentPhaseChange.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-2024 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 "basicThermo.H"
28 #include "multicomponentThermo.H"
30 
31 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace fv
36 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 void Foam::fv::multicomponentPhaseChange::readCoeffs()
45 {
46  if (species_ != coeffs().lookup<wordList>("species"))
47  {
49  << "Cannot change the species of a " << typeName << " model "
50  << "at run time" << exit(FatalIOError);
51  }
52 
53  energySemiImplicit_ =
54  coeffs().lookup<bool>("energySemiImplicit");
55 }
56 
57 
58 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59 
61 (
62  const word& name,
63  const word& modelType,
64  const fvMesh& mesh,
65  const dictionary& dict,
66  const Pair<bool>& fluidThermosRequired
67 )
68 :
70  (
71  name,
72  modelType,
73  mesh,
74  dict,
75  fluidThermosRequired,
76  {true, true}
77  ),
78  species_(coeffs().lookup<wordList>("species")),
79  energySemiImplicit_(false)
80 {
81  readCoeffs();
82 }
83 
84 
85 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
86 
89 {
92  (
93  "mDot",
94  mesh(),
96  );
97 
98  forAll(species(), mDoti)
99  {
100  tmDot.ref() += mDot(mDoti);
101  }
102 
103  return tmDot;
104 }
105 
106 
108 (
109  const volScalarField& alpha,
110  const volScalarField& rho,
111  const volScalarField& heOrYi,
112  fvMatrix<scalar>& eqn
113 ) const
114 {
115  const label i = index(phaseNames(), alpha.group());
116  const label s = sign(phaseNames(), alpha.group());
117 
118  // Energy equation
119  if (index(heNames(), heOrYi.name()) != -1)
120  {
121  const volScalarField& p = this->p();
122  const volScalarField Tchange(vifToVf(this->Tchange()));
123 
124  forAll(species(), mDoti)
125  {
126  const label speciei =
127  specieThermos()[i].species()[species()[mDoti]];
128 
129  const volScalarField::Internal mDot(this->mDot(mDoti));
130 
131  // Direct transfer of energy due to mass transfer
133  vfToVif(specieThermos()[i].hsi(speciei, p, Tchange));
134  if (energySemiImplicit_)
135  {
136  eqn += -fvm::SuSp(-s*mDot, heOrYi) + s*mDot*(hs - heOrYi);
137  }
138  else
139  {
140  eqn += s*mDot*hs;
141  }
142 
143  // Absolute enthalpies at the interface
145  for (label j = 0; j < 2; ++ j)
146  {
147  has[j] = vfToVif(specieThermos()[j].hai(speciei, p, Tchange));
148  }
149 
150  // Latent heat of phase change
151  eqn -=
152  (i == 0 ? 1 - Lfraction() : Lfraction())
153  *mDot*(has.second() - has.first());
154  }
155  }
156  // Mass fraction equation
157  else if
158  (
159  specieThermos().valid()[i]
160  && specieThermos()[i].containsSpecie(heOrYi.member())
161  )
162  {
163  // A transferring specie
164  if (species().found(heOrYi.member()))
165  {
166  eqn += s*mDot(species()[heOrYi.member()]);
167  }
168  // A non-transferring specie. Do nothing.
169  else
170  {}
171  }
172  // Something else. Fall back.
173  else
174  {
176  }
177 }
178 
179 
181 {
182  if (phaseChange::read(dict))
183  {
184  readCoeffs();
185  return true;
186  }
187  else
188  {
189  return false;
190  }
191 }
192 
193 
194 // ************************************************************************* //
scalar hs(const scalar p, const scalar T) const
Definition: EtoHthermo.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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 > > New(const word &name, const Mesh &mesh, const dimensionSet &, const Field< 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:149
const word & name() const
Return name.
Definition: IOobject.H:310
const Type & second() const
Return second.
Definition: Pair.H:110
const Type & first() const
Return first.
Definition: Pair.H:98
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:59
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: massTransfer.C:242
virtual void addSup(fvMatrix< scalar > &eqn) const
Add a source term to a field-less proxy equation.
Definition: massTransfer.C:218
Base class for phase change models in which multiple components may change phase. This can only be ap...
multicomponentPhaseChange(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict, const Pair< bool > &fluidThermosRequired)
Construct from explicit source name and mesh.
virtual tmp< DimensionedField< scalar, volMesh > > mDot() const
Return the total mass transfer rate.
virtual bool read(const dictionary &dict)
Read source dictionary.
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.
Base class for phase change models.
Definition: phaseChange.H:59
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
defineTypeNameAndDebug(bound, 0)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
dimensionedScalar sign(const dimensionedScalar &ds)
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimTime
const dimensionSet dimDensity
IOerror FatalIOError
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
labelList fv(nPoints)
dictionary dict
volScalarField & p