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