singleComponentPhaseChange.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::singleComponentPhaseChange::readCoeffs()
45 {
46  if
47  (
49  && specie_ != coeffs().lookup<word>("specie")
50  )
51  {
53  << "Cannot change the specie of a " << typeName << " model "
54  << "at run time" << exit(FatalIOError);
55  }
56 
57  energySemiImplicit_ =
58  coeffs().lookup<bool>("energySemiImplicit");
59 }
60 
61 
62 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
63 
65 (
66  const word& name,
67  const word& modelType,
68  const fvMesh& mesh,
69  const dictionary& dict,
70  const Pair<bool>& fluidThermosRequired,
71  const Pair<bool>& specieThermosRequired
72 )
73 :
75  (
76  name,
77  modelType,
78  mesh,
79  dict,
80  fluidThermosRequired,
81  specieThermosRequired
82  ),
83  specie_
84  (
85  specieThermos().valid().first() || specieThermos().valid().second()
86  ? coeffs().lookup<word>("specie")
87  : word::null
88  ),
89  specieis_
90  (
91  specieThermos().valid().first()
92  ? specieThermos().first().species()[specie_]
93  : -1,
94  specieThermos().valid().second()
95  ? specieThermos().second().species()[specie_]
96  : -1
97  ),
98  energySemiImplicit_(false)
99 {
100  readCoeffs();
101 }
102 
103 
104 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
105 
107 (
108  const volScalarField& alpha,
109  const volScalarField& rho,
110  const volScalarField& heOrYi,
111  fvMatrix<scalar>& eqn
112 ) const
113 {
114  const label i = index(phaseNames(), alpha.group());
115  const label s = sign(phaseNames(), alpha.group());
116 
117  // Energy equation
118  if (index(heNames(), heOrYi.name()) != -1)
119  {
120  const volScalarField& p = this->p();
121  const volScalarField Tchange(vifToVf(this->Tchange()));
122 
123  const volScalarField::Internal mDot(this->mDot());
124 
125  // Direct transfer of energy due to mass transfer
127  specieThermos().valid()[i]
128  ? vfToVif(specieThermos()[i].hsi(specieis_[i], p, Tchange))
129  : vfToVif(thermos()[i].hs(p, Tchange));
130  if (energySemiImplicit_)
131  {
132  eqn += -fvm::SuSp(-s*mDot, heOrYi) + s*mDot*(hs - heOrYi);
133  }
134  else
135  {
136  eqn += s*mDot*hs;
137  }
138 
139  // Absolute enthalpies at the interface
141  for (label j = 0; j < 2; ++ j)
142  {
143  has[j] =
144  specieThermos().valid()[j]
145  ? vfToVif(specieThermos()[j].hai(specieis_[j], p, Tchange))
146  : vfToVif(thermos()[j].ha(p, Tchange));
147  }
148 
149  // Latent heat of phase change
150  eqn -=
151  (i == 0 ? 1 - Lfraction() : Lfraction())
152  *mDot*(has.second() - has.first());
153  }
154  // Mass fraction equation
155  else if
156  (
157  specieThermos().valid()[i]
158  && specieThermos()[i].containsSpecie(heOrYi.member())
159  )
160  {
161  // The transferring specie
162  if (heOrYi.member() == specie())
163  {
164  eqn += s*mDot();
165  }
166  // A non-transferring specie. Do nothing.
167  else
168  {}
169  }
170  // Something else. Fall back.
171  else
172  {
173  massTransfer::addSup(alpha, rho, heOrYi, eqn);
174  }
175 }
176 
177 
179 {
180  if (phaseChange::read(dict))
181  {
182  readCoeffs();
183  return true;
184  }
185  else
186  {
187  return false;
188  }
189 }
190 
191 
192 // ************************************************************************* //
scalar hs(const scalar p, const scalar T) const
Definition: EtoHthermo.H:11
scalar ha(const scalar p, const scalar T) const
Definition: EtoHthermo.H:20
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...
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.
Definition: phaseChange.H:59
const ThermoRefPair< multicomponentThermo > & specieThermos() const
Return the specie thermo references.
Definition: phaseChangeI.H:45
Base class for phase change models in which only a single component changes phase....
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.
singleComponentPhaseChange(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict, const Pair< bool > &fluidThermosRequired, const Pair< bool > &specieThermosRequired)
Construct from explicit source name and mesh.
Base class of the thermophysical property types.
Definition: specie.H:66
A class for managing temporary objects.
Definition: tmp.H:55
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))
bool valid(const PtrList< ModelType > &l)
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
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
labelList second(const UList< labelPair > &p)
Definition: patchToPatch.C:49
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
IOerror FatalIOError
labelList fv(nPoints)
dictionary dict
volScalarField & p