oneResistanceHeatTransfer.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::oneResistanceHeatTransfer
26 
27 Description
28  Model for heat transfer between two phases. One heat transfer coefficient
29  is used to calculate the heat flux from the temperature difference between
30  the two phases.
31 
32 Usage
33  Example usage:
34  \verbatim
35  heatTransfer
36  {
37  type oneResistanceHeatTransfer;
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
54  {
55  type RanzMarshall;
56  residualAlpha 1e-4;
57  }
58  water_dispersedIn_air
59  {
60  type RanzMarshall;
61  residualAlpha 1e-4;
62  }
63  }
64  }
65  \endverbatim
66 
67 SourceFiles
68  oneResistanceHeatTransfer.C
69 
70 \*---------------------------------------------------------------------------*/
71 
72 #ifndef oneResistanceHeatTransfer_H
73 #define oneResistanceHeatTransfer_H
74 
75 #include "phaseSystem.H"
76 #include "heatTransferModel.H"
77 
78 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
79 
80 namespace Foam
81 {
82 namespace fv
83 {
84 
85 /*---------------------------------------------------------------------------*\
86  Class oneResistanceHeatTransfer Declaration
87 \*---------------------------------------------------------------------------*/
88 
90 :
91  public fvModel
92 {
93 private:
94 
95  // Private Type Definitions
96 
97  //- Table of heat transfer models
98  typedef
100  <
104  > modelsTable;
105 
106  //- Table of heat transfer coefficients
107  typedef
108  HashTable
109  <
113  > KsTable;
114 
115 
116  // Private Data
117 
118  //- Phase system
119  const phaseSystem& fluid_;
120 
121  //- Blended heat transfer coefficient models
122  modelsTable models_;
123 
124  //- Names of the affected energy fields
125  wordList heNames_;
126 
127  //- Counter for the evaluations of the energy equation sources
128  mutable label energyEquationIndex_;
129 
130  //- Cached heat transfer coefficients
131  mutable KsTable Ks_;
132 
133 
134  // Private Member Functions
135 
136  //- Non-virtual read
137  void readCoeffs(const dictionary& dict);
138 
139  //- Get the heat transfer coefficients
140  const KsTable& Ks() const;
141 
142 
143 public:
144 
145  //- Runtime type information
146  TypeName("oneResistanceHeatTransfer");
147 
148 
149  // Constructors
150 
151  //- Construct from explicit source name and mesh
153  (
154  const word& name,
155  const word& modelType,
156  const fvMesh& mesh,
157  const dictionary& dict
158  );
159 
160 
161  // Member Functions
162 
163  // Access
164 
165  //- Return the heat transfer coefficient between a pair of phases
167  (
168  const phaseModel& phase1,
169  const phaseModel& phase2
170  ) const;
171 
172 
173  // Checks
174 
175  //- Return the list of fields for which the fvModel adds source term
176  // to the transport equation
177  virtual wordList addSupFields() const;
178 
179 
180  // Sources
181 
182  //- Add a source term to the phase-compressible energy equation
183  void addSup
184  (
185  const volScalarField& alpha,
186  const volScalarField& rho,
187  const volScalarField& he,
188  fvMatrix<scalar>& eqn
189  ) const;
190 
191 
192  // Mesh changes
193 
194  //- Update for mesh motion
195  virtual bool movePoints();
196 
197  //- Update topology using the given map
198  virtual void topoChange(const polyTopoChangeMap&);
199 
200  //- Update from another mesh using the given map
201  virtual void mapMesh(const polyMeshMap&);
202 
203  //- Redistribute or update using the given distribution map
204  virtual void distribute(const polyDistributionMap&);
205 
206 
207  //- Correct the fvModel
208  virtual void correct();
209 
210 
211  // IO
212 
213  //- Read source dictionary
214  virtual bool read(const dictionary& dict);
215 };
216 
217 
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 
220 } // End namespace fv
221 } // End namespace Foam
222 
223 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
224 
225 #endif
226 
227 // ************************************************************************* //
Generic GeometricField class.
A HashTable specialisation for hashing pointers.
Definition: HashPtrTable.H:67
An STL-conforming hash table.
Definition: HashTable.H:127
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. One heat transfer coefficient is used to calculate the he...
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.
tmp< volScalarField > K(const phaseModel &phase1, const phaseModel &phase2) const
Return the heat transfer coefficient between a pair of phases.
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.
oneResistanceHeatTransfer(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
TypeName("oneResistanceHeatTransfer")
Runtime type information.
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 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
thermo he()
labelList fv(nPoints)
dictionary dict