phaseChange.H
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 Class
25  Foam::fv::phaseChange
26 
27 Description
28  Base class for phase change models
29 
30 SourceFiles
31  phaseChange.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef phaseChange_H
36 #define phaseChange_H
37 
38 #include "hashedWordList.H"
39 #include "massTransfer.H"
40 #include "ThermoRefPair.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 class basicThermo;
48 class fluidThermo;
49 class multicomponentThermo;
50 class fluidMulticomponentThermo;
51 
52 namespace fv
53 {
54 
55 /*---------------------------------------------------------------------------*\
56  Class phaseChange Declaration
57 \*---------------------------------------------------------------------------*/
58 
59 class phaseChange
60 :
61  public massTransfer
62 {
63 private:
64 
65  // Private Data
66 
67  //- The thermo references
68  const ThermoRefPair<basicThermo> thermos_;
69 
70  //- Names of the energy fields
71  const Pair<word> heNames_;
72 
73  //- The names of the transferring species
74  hashedWordList species_;
75 
76  //- The indices of the transferring species in the two phases, or -1 if
77  // the phase is not multicomponent
78  List<labelPair> specieis_;
79 
80  //- Whether or not to linearise the energy source
81  bool energySemiImplicit_;
82 
83 
84  // Private Member Functions
85 
86  //- Non-virtual read
87  void readCoeffs(const dictionary& dict);
88 
89  //- Initialise the indices of the transferring species in the two
90  // phases, or -1 if the phase is not multicomponent
91  const List<labelPair> initSpecieis() const;
92 
93 
94 protected:
95 
96  // Protected Member Functions
97 
98  //- Read the names of the transferring specie
99  wordList readSpecie(const dictionary& dict, const bool required) const;
100 
101  //- Read the names of the transferring species
102  wordList readSpecies(const dictionary& dict, const bool required) const;
103 
104  //- Re-read the names of the transferring specie
105  void reReadSpecie(const dictionary& dict) const;
106 
107  //- Re-read the names of the transferring species
108  void reReadSpecies(const dictionary& dict) const;
109 
110  //- Set the names of the transferring species
111  void setSpecies
112  (
113  const word& name,
114  const word& modelType,
115  const wordList& species
116  );
117 
118  //- Set the names of the transferring species
119  void setSpecies(const wordList& species);
120 
121  //- Re-set the names of the transferring species
122  void reSetSpecies(const wordList& species);
123 
124  //- Access the pressure field
125  const volScalarField& p() const;
126 
127  //- Add a boundary field to the given internal field
129  (
131  );
132 
133  //- Remove the boundary field from the given geometric field
135  (
136  const tmp<volScalarField>& tvf
137  );
138 
139 
140 public:
141 
142  //- Runtime type information
143  TypeName("phaseChange");
144 
145 
146  // Constructors
147 
148  //- Construct from explicit source name and mesh
150  (
151  const word& name,
152  const word& modelType,
153  const fvMesh& mesh,
154  const dictionary& dict,
155  const wordList& species
156  );
157 
158 
159  // Member Functions
160 
161  // Access
162 
163  //- Return the thermo references
164  inline const ThermoRefPair<basicThermo>& thermos() const;
165 
166  //- Return the fluid thermo references
168  fluidThermos(const bool, const bool) const;
169 
170  //- Return the multicomponent thermo references
172  multicomponentThermos(const bool, const bool) const;
173 
174  //- Return the fluid multicomponent thermo references
176  fluidMulticomponentThermos(const bool, const bool) const;
177 
178  //- Return the names of the energy fields
179  inline const Pair<word>& heNames() const;
180 
181  //- Return the names of the transferring species. Empty if neither
182  // thermo is multicomponent.
183  inline const hashedWordList& species() const;
184 
185  //- Return the indices of the transferring specie in the two
186  // phases, or -1 if the phase is not multicomponent
187  const labelPair& specieis(const label mDoti = -1) const;
188 
189 
190  // Evaluation
191 
192  //- Return the temperature at which the phases are considered to be
193  // changing. By default this is considered to be the temperature
194  // of the "source" phase (i.e., the phase for which mDot is
195  // negative), but this can be overridden to account for heat
196  // transfer modelling or similar.
198 
199  //- Return the fraction of the latent heat that is transferred into
200  // the second phase. By default this is weighted by the phase
201  // thermal conductivities, but this can be overridden to account
202  // for heat transfer modelling or similar.
204 
205  //- Return the latent heat
207  (
208  const label mDoti = -1
209  ) const;
210 
211  //- Return the latent heat for a given changing temperature
213  (
215  const label mDoti = -1
216  ) const;
217 
218  //- Return the latent heat for a patch and a given changing
219  // temperature
221  (
222  const label patchi,
223  const scalarField& Tchange,
224  const label mDoti = -1
225  ) const;
226 
227 
228  // Sources
229 
230  //- Return the total phase change rate
232 
233  //- Return the mass transfer rate of a specie
235  (
236  const label mDoti
237  ) const;
238 
239  //- Override the energy equation to add the phase change heat, or
240  // the species equations to add the relevant mass sources
241  void addSup
242  (
243  const volScalarField& alpha,
244  const volScalarField& rho,
245  const volScalarField& heOrYi,
246  fvMatrix<scalar>& eqn
247  ) const;
248 
249 
250  //- Read source dictionary
251  virtual bool read(const dictionary& dict);
252 };
253 
254 
255 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
256 
257 } // End namespace fv
258 } // End namespace Foam
259 
260 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
261 
262 #include "phaseChangeI.H"
263 
264 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
265 
266 #endif
267 
268 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
Class containing a pair of thermo references. Handles down-casting to more specific thermo types by c...
Definition: ThermoRefPair.H:51
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
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:96
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:69
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:57
Base class for mass transfers between phases.
Definition: massTransfer.H:55
tmp< volScalarField::Internal > rho(const label i) const
Return the density.
Definition: massTransfer.C:92
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
const labelPair & specieis(const label mDoti=-1) const
Return the indices of the transferring specie in the two.
Definition: phaseChange.C:366
wordList readSpecie(const dictionary &dict, const bool required) const
Read the names of the transferring specie.
Definition: phaseChange.C:74
const ThermoRefPair< multicomponentThermo > multicomponentThermos(const bool, const bool) const
Return the multicomponent thermo references.
Definition: phaseChange.C:334
virtual tmp< DimensionedField< scalar, volMesh > > Lfraction() const
Return the fraction of the latent heat that is transferred into.
Definition: phaseChange.C:423
const ThermoRefPair< fluidThermo > fluidThermos(const bool, const bool) const
Return the fluid thermo references.
Definition: phaseChange.C:317
void reReadSpecie(const dictionary &dict) const
Re-read the names of the transferring specie.
Definition: phaseChange.C:129
wordList readSpecies(const dictionary &dict, const bool required) const
Read the names of the transferring species.
Definition: phaseChange.C:96
static tmp< DimensionedField< scalar, volMesh > > vfToVif(const tmp< volScalarField > &tvf)
Remove the boundary field from the given geometric field.
Definition: phaseChange.C:277
virtual tmp< DimensionedField< scalar, volMesh > > mDot() const
Return the total phase change rate.
Definition: phaseChange.C:498
TypeName("phaseChange")
Runtime type information.
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
void setSpecies(const word &name, const word &modelType, const wordList &species)
Set the names of the transferring species.
Definition: phaseChange.C:152
const ThermoRefPair< fluidMulticomponentThermo > fluidMulticomponentThermos(const bool, const bool) const
Return the fluid multicomponent thermo references.
Definition: phaseChange.C:351
void reSetSpecies(const wordList &species)
Re-set the names of the transferring species.
Definition: phaseChange.C:228
tmp< DimensionedField< scalar, volMesh > > L(const label mDoti=-1) const
Return the latent heat.
Definition: phaseChange.C:433
const hashedWordList & species() const
Return the names of the transferring species. Empty if neither.
Definition: phaseChangeI.H:44
const volScalarField & p() const
Access the pressure field.
Definition: phaseChange.C:239
virtual tmp< DimensionedField< scalar, volMesh > > Tchange() const
Return the temperature at which the phases are considered to be.
Definition: phaseChange.C:393
phaseChange(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict, const wordList &species)
Construct from explicit source name and mesh.
Definition: phaseChange.C:292
const ThermoRefPair< basicThermo > & thermos() const
Return the thermo references.
Definition: phaseChangeI.H:32
const Pair< word > & heNames() const
Return the names of the energy fields.
Definition: phaseChangeI.H:38
static tmp< volScalarField > vifToVf(const tmp< DimensionedField< scalar, volMesh >> &tvif)
Add a boundary field to the given internal field.
Definition: phaseChange.C:254
A wordList with hashed indices for faster lookup by name.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
label patchi
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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
labelList fv(nPoints)
dictionary dict