heatTransfer.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::heatTransfer
26 
27 Description
28  Model for heat exchange. Requires specification of an ambient temperature
29  with which to exchange heat, and a model for the heat transfer coefficient
30  (htc) and the area per unit volume (Av). These are then used to apply the
31  following source to the energy equation:
32 
33  \f[
34  -htc*Av*(T_a - T)
35  \f]
36 
37  If the semiImplicit option is set, then this becomes:
38 
39  \f[
40  -htc*Av*(T_a - T) + htc*Av/Cp*h - Sp(htc*Av/Cp, h);
41  \f]
42 
43 Usage
44  Example usage:
45  \verbatim
46  heatTransfer
47  {
48  type heatTransfer;
49 
50  cellZone c0;
51 
52  semiImplicit no;
53 
54  Ta 300;
55 
56  Av 200;
57 
58  heatTransferCoefficientModel constant;
59 
60  htc 10;
61  }
62  \endverbatim
63 
64 See also
65  Foam::fv::heatTransferCoefficientModel
66 
67 SourceFiles
68  heatTransfer.C
69 
70 \*---------------------------------------------------------------------------*/
71 
72 #ifndef heatTransfer_H
73 #define heatTransfer_H
74 
75 #include "fvModel.H"
76 #include "fvCellZone.H"
78 #include "heatTransferAv.H"
79 
80 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
81 
82 namespace Foam
83 {
84 namespace fv
85 {
86 
87 /*---------------------------------------------------------------------------*\
88  Class heatTransfer Declaration
89 \*---------------------------------------------------------------------------*/
90 
91 class heatTransfer
92 :
93  public fvModel
94 {
95  // Private data
96 
97  //- The cellzone the model applies to
98  fvCellZone zone_;
99 
100  //- Name of the phase
101  word phaseName_;
102 
103  //- Flag to activate semi-implicit coupling
104  bool semiImplicit_;
105 
106  //- Name of temperature field; default = "T"
107  word TName_;
108 
109  //- Ambient temperature
110  dimensionedScalar Ta_;
111 
112  //- The heat transfer area per unit volume
113  autoPtr<heatTransferAv> heatTransferAv_;
114 
115  //- The heat transfer model
116  autoPtr<heatTransferCoefficientModel> heatTransferCoefficientModel_;
117 
118 
119  // Private member functions
120 
121  //- Non-virtual read
122  void readCoeffs(const dictionary& dict);
123 
124  //- Source term to energy equation
125  template<class AlphaFieldType>
126  inline void add
127  (
128  const AlphaFieldType& alpha,
129  fvMatrix<scalar>& eqn
130  ) const;
131 
132 
133 public:
134 
135  //- Runtime type information
136  TypeName("heatTransfer");
137 
138 
139  // Constructors
140 
141  //- Construct from dictionary
143  (
144  const word& name,
145  const word& modelType,
146  const fvMesh& mesh,
147  const dictionary& dict
148  );
149 
150 
151  //- Destructor
152  virtual ~heatTransfer();
153 
154 
155  // Member Functions
156 
157  // Checks
158 
159  //- Return the list of fields for which the fvModel adds source term
160  // to the transport equation
161  virtual wordList addSupFields() const;
162 
163 
164  // Sources
165 
166  //- Source term to energy equation
167  virtual void addSup
168  (
169  const volScalarField& he,
170  fvMatrix<scalar>& eqn
171  ) const;
172 
173  //- Source term to compressible energy equation
174  virtual void addSup
175  (
176  const volScalarField& rho,
177  const volScalarField& he,
178  fvMatrix<scalar>& eqn
179  ) const;
180 
181  //- Source term to phase energy equation
182  virtual void addSup
183  (
184  const volScalarField& alpha,
185  const volScalarField& rho,
186  const volScalarField& he,
187  fvMatrix<scalar>& eqn
188  ) const;
189 
190 
191  // Correction
192 
193  //- Correct the model
194  virtual void correct();
195 
196 
197  // Mesh changes
198 
199  //- Update for mesh motion
200  virtual bool movePoints();
201 
202  //- Update topology using the given map
203  virtual void topoChange(const polyTopoChangeMap&);
204 
205  //- Update from another mesh using the given map
206  virtual void mapMesh(const polyMeshMap&);
207 
208  //- Redistribute or update using the given distribution map
209  virtual void distribute(const polyDistributionMap&);
210 
211 
212  // IO
213 
214  //- Read dictionary
215  virtual bool read(const dictionary& dict);
216 };
217 
218 
219 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220 
221 } // End namespace fv
222 } // End namespace Foam
223 
224 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
225 
226 #endif
227 
228 // ************************************************************************* //
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
cellZone selection or generation class with caching of zone volume
Definition: fvCellZone.H:94
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 exchange. Requires specification of an ambient temperature with which to exchange heat...
Definition: heatTransfer.H:93
virtual bool movePoints()
Update for mesh motion.
Definition: heatTransfer.C:201
virtual ~heatTransfer()
Destructor.
Definition: heatTransfer.C:144
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
Definition: heatTransfer.C:150
TypeName("heatTransfer")
Runtime type information.
virtual void correct()
Correct the model.
Definition: heatTransfer.C:195
virtual void addSup(const volScalarField &he, fvMatrix< scalar > &eqn) const
Source term to energy equation.
Definition: heatTransfer.C:163
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: heatTransfer.C:208
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: heatTransfer.C:220
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: heatTransfer.C:226
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: heatTransfer.C:214
heatTransfer(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from dictionary.
Definition: heatTransfer.C:122
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.
thermo he()
labelList fv(nPoints)
dictionary dict