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