heatTransferLimitedPhaseChange.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) 2024-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::heatTransferLimitedPhaseChange
26 
27 Description
28  Model for heat transfer rate limited phase change between two phases.
29 
30  The interface between the two phases is assumed to be at a saturated
31  condition. This allows the temperature of the interface to be evaluated
32  from a user-supplied saturation curve. This temperature then defines the
33  heat flux being transferred to the interface from the surrounding fluid.
34  The imbalance in the heat fluxes on either side of the interface is then
35  divided by the latent heat of phase change in order to get the rate at
36  which mass is being changed from one phase to the other.
37 
38  This model only supports pure phases. A two-resistance heat transfer model
39  must also be in operation between the two changing phases.
40 
41 Usage
42  Example usage:
43  \verbatim
44  phaseChange
45  {
46  type heatTransferLimitedPhaseChange;
47  libs ("libmultiphaseEulerFvModels.so");
48 
49  phases (steam water);
50 
51  energySemiImplicit yes;
52  pressureImplicit no;
53 
54  saturationTemperature
55  {
56  type constant;
57  value 372.76;
58  }
59  }
60  \endverbatim
61 
62 SourceFiles
63  heatTransferLimitedPhaseChange.C
64 
65 \*---------------------------------------------------------------------------*/
66 
67 #ifndef heatTransferLimitedPhaseChange_H
68 #define heatTransferLimitedPhaseChange_H
69 
70 #include "phaseChange.H"
71 #include "phaseSystem.H"
73 
74 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
75 
76 namespace Foam
77 {
78 
79 // Forward declaration of classes
80 namespace solvers
81 {
82  class multiphaseEuler;
83 }
84 
85 namespace fv
86 {
87 
88 /*---------------------------------------------------------------------------*\
89  Class heatTransferLimitedPhaseChange Declaration
90 \*---------------------------------------------------------------------------*/
91 
93 :
94  public phaseChange
95 {
96 private:
97 
98  // Private Data
99 
100  //- Solver
101  const solvers::multiphaseEuler& solver_;
102 
103  //- Phase system
104  const phaseSystem& fluid_;
105 
106  //- Phase 1
107  const phaseModel& phase1_;
108 
109  //- Phase 2
110  const phaseModel& phase2_;
111 
112  //- The saturation curve
113  autoPtr<saturationTemperatureModel> saturationModelPtr_;
114 
115  //- Should this phase change be linearised in the pressure equation?
116  bool pressureImplicit_;
117 
118  //- Counter for the evaluations of the pressure equation sources
119  mutable label pressureEquationIndex_;
120 
121  //- The phase change rate
122  mutable volScalarField::Internal mDot_;
123 
124  //- The derivative of the phase change rate w.r.t. the pressure. Only
125  // valid if pressure implicit is selected.
126  mutable autoPtr<volScalarField::Internal> dmDotdpPtr_;
127 
128 
129  // Private Member Functions
130 
131  //- Non-virtual read
132  void readCoeffs(const dictionary& dict);
133 
134  //- Correct the phase change rate
135  void correctMDot() const;
136 
137 
138 public:
139 
140  //- Runtime type information
141  TypeName("heatTransferLimitedPhaseChange");
142 
143 
144  // Constructors
145 
146  //- Construct from explicit source name and mesh
148  (
149  const word& name,
150  const word& modelType,
151  const fvMesh& mesh,
152  const dictionary& dict
153  );
154 
155 
156  // Member Functions
157 
158  // Evaluation
159 
160  //- Return the fraction of the latent heat that is transferred into
161  // the second phase
163 
164 
165  // Sources
166 
167  //- Return the mass transfer rate
169 
170  //- Use phaseChange's source functions
171  using phaseChange::addSup;
172 
173  //- Override the pressure equation to add the mass transfer rate
174  // linearised in the pressure
175  void addSup
176  (
177  const volScalarField& alpha,
178  const volScalarField& rho,
179  fvMatrix<scalar>& eqn
180  ) const;
181 
182 
183  //- Correct the fvModel
184  virtual void correct();
185 
186 
187  // IO
188 
189  //- Read source dictionary
190  virtual bool read(const dictionary& dict);
191 };
192 
193 
194 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195 
196 } // End namespace fv
197 } // End namespace Foam
198 
199 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200 
201 #endif
202 
203 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.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
Model for heat transfer rate limited phase change between two phases.
virtual tmp< DimensionedField< scalar, volMesh > > Lfraction() const
Return the fraction of the latent heat that is transferred into.
TypeName("heatTransferLimitedPhaseChange")
Runtime type information.
virtual tmp< DimensionedField< scalar, volMesh > > mDot() const
Return the 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
Use phaseChange's source functions.
Definition: phaseChange.C:551
heatTransferLimitedPhaseChange(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
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 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
Class to represent a system of phases.
Definition: phaseSystem.H:74
Solver module for a system of any number of compressible fluid phases with a common pressure,...
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
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