twoResistanceHeatTransfer.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::twoResistanceHeatTransfer
26 
27 Description
28  Model for heat transfer between two phases. Two heat transfer coefficients
29  are used to calculate the heat fluxes that result from the temperature
30  differences between the phases and the shared interface.
31 
32 Usage
33  Example usage:
34  \verbatim
35  heatTransfer
36  {
37  type twoResistanceHeatTransfer;
38  libs ("libmultiphaseEulerFvModels.so");
39 
40  phases (air water);
41 
42  blending
43  {
44  type linear;
45  minPartlyContinuousAlpha.air 0;
46  minFullyContinuousAlpha.air 1;
47  minPartlyContinuousAlpha.water 0;
48  minFullyContinuousAlpha.water 1;
49  }
50 
51  models
52  {
53  air_dispersedIn_water_inThe_air
54  {
55  type spherical;
56  residualAlpha 1e-4;
57  }
58  air_dispersedIn_water_inThe_water
59  {
60  type RanzMarshall;
61  residualAlpha 1e-4;
62  }
63  water_dispersedIn_air_inThe_air
64  {
65  type RanzMarshall;
66  residualAlpha 1e-4;
67  }
68  water_dispersedIn_air_inThe_water
69  {
70  type spherical;
71  residualAlpha 1e-4;
72  }
73  }
74  }
75  \endverbatim
76 
77 SourceFiles
78  twoResistanceHeatTransfer.C
79 
80 \*---------------------------------------------------------------------------*/
81 
82 #ifndef twoResistanceHeatTransfer_H
83 #define twoResistanceHeatTransfer_H
84 
85 #include "phaseSystem.H"
86 #include "heatTransferModel.H"
87 
88 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
89 
90 namespace Foam
91 {
92 namespace fv
93 {
94 
95 /*---------------------------------------------------------------------------*\
96  Class twoResistanceHeatTransfer Declaration
97 \*---------------------------------------------------------------------------*/
98 
100 :
101  public fvModel
102 {
103 private:
104 
105  // Private Type Definitions
106 
107  //- Table of heat transfer models
108  typedef
110  <
114  > modelsTable;
115 
116  //- Table of heat transfer coefficients
117  typedef
118  HashTable
119  <
123  > KsTable;
124 
125 
126  // Private Data
127 
128  //- Phase system
129  const phaseSystem& fluid_;
130 
131  //- Blended-sided heat transfer coefficient models
132  modelsTable models_;
133 
134  //- Names of the affected energy fields
135  wordList heNames_;
136 
137  //- Counter for the evaluations of the energy equation sources
138  mutable label energyEquationIndex_;
139 
140  //- Cached heat transfer coefficients
141  mutable KsTable Ks_;
142 
143 
144  // Private Member Functions
145 
146  //- Non-virtual read
147  void readCoeffs(const dictionary& dict);
148 
149  //- Get the heat transfer coefficients
150  const KsTable& Ks() const;
151 
152 
153 public:
154 
155  //- Runtime type information
156  TypeName("twoResistanceHeatTransfer");
157 
158 
159  // Constructors
160 
161  //- Construct from explicit source name and mesh
163  (
164  const word& name,
165  const word& modelType,
166  const fvMesh& mesh,
167  const dictionary& dict
168  );
169 
170 
171  // Member Functions
172 
173  // Access
174 
175  //- Return the heat transfer coefficients between a pair of phases
177  (
178  const phaseModel& phase1,
179  const phaseModel& phase2
180  ) const;
181 
182  //- Return the heat transfer coefficients between a pair of phases
183  // with a specified residual volume fraction
185  (
186  const phaseModel& phase1,
187  const phaseModel& phase2,
188  const scalar residualAlpha
189  ) const;
190 
191 
192  // Checks
193 
194  //- Return the list of fields for which the fvModel adds source term
195  // to the transport equation
196  virtual wordList addSupFields() const;
197 
198 
199  // Sources
200 
201  //- Add a source term to the phase-compressible energy equation
202  void addSup
203  (
204  const volScalarField& alpha,
205  const volScalarField& rho,
206  const volScalarField& he,
207  fvMatrix<scalar>& eqn
208  ) const;
209 
210 
211  // Mesh changes
212 
213  //- Update for mesh motion
214  virtual bool movePoints();
215 
216  //- Update topology using the given map
217  virtual void topoChange(const polyTopoChangeMap&);
218 
219  //- Update from another mesh using the given map
220  virtual void mapMesh(const polyMeshMap&);
221 
222  //- Redistribute or update using the given distribution map
223  virtual void distribute(const polyDistributionMap&);
224 
225 
226  //- Correct the fvModel
227  virtual void correct();
228 
229 
230  // IO
231 
232  //- Read source dictionary
233  virtual bool read(const dictionary& dict);
234 };
235 
236 
237 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
238 
239 } // End namespace fv
240 } // End namespace Foam
241 
242 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
243 
244 #endif
245 
246 // ************************************************************************* //
Generic GeometricField class.
A HashTable specialisation for hashing pointers.
Definition: HashPtrTable.H:67
An STL-conforming hash table.
Definition: HashTable.H:127
An ordered pair of two objects of type <Type> with first() and second() elements.
Definition: Pair.H:66
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
Finite volume model abstract base class.
Definition: fvModel.H:60
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 between two phases. Two heat transfer coefficients are used to calculate the ...
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
void addSup(const volScalarField &alpha, const volScalarField &rho, const volScalarField &he, fvMatrix< scalar > &eqn) const
Add a source term to the phase-compressible energy equation.
virtual void correct()
Correct the fvModel.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read source dictionary.
TypeName("twoResistanceHeatTransfer")
Runtime type information.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
twoResistanceHeatTransfer(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
Word-pair based class used for keying interface models in hash tables.
Class to represent a system of phases.
Definition: phaseSystem.H:74
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
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
thermo he()
labelList fv(nPoints)
dictionary dict