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-2023 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  nbrRegion other;
50 
51  interpolationMethod intersection;
52  master true;
53 
54  semiImplicit no;
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  interRegionHeatTransfer.C
69 
70 \*---------------------------------------------------------------------------*/
71 
72 #ifndef interRegionHeatTransfer_H
73 #define interRegionHeatTransfer_H
74 
75 #include "interRegionModel.H"
77 #include "heatTransferAv.H"
78 
79 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
80 
81 namespace Foam
82 {
83 namespace fv
84 {
85 
86 /*---------------------------------------------------------------------------*\
87  Class interRegionHeatTransfer Declaration
88 \*---------------------------------------------------------------------------*/
89 
91 :
92  public interRegionModel
93 {
94  // Private data
95 
96  //- Flag to activate semi-implicit coupling
97  bool semiImplicit_;
98 
99  //- Name of temperature field; default = "T"
100  word TName_;
101 
102  //- Name of neighbour temperature field; default = "T"
103  word TNbrName_;
104 
105  //- The heat transfer area per unit volume
106  autoPtr<heatTransferAv> heatTransferAv_;
107 
108  //- The heat transfer model
109  autoPtr<heatTransferCoefficientModel> heatTransferCoefficientModel_;
110 
111 
112  // Private member functions
113 
114  //- Non-virtual read
115  void readCoeffs();
116 
117  //- Get the neighbour heat transfer
118  const interRegionHeatTransfer& nbrHeatTransfer() const;
119 
120 
121 public:
122 
123  //- Runtime type information
124  TypeName("interRegionHeatTransfer");
125 
126 
127  // Constructors
128 
129  //- Construct from dictionary
131  (
132  const word& name,
133  const word& modelType,
134  const fvMesh& mesh,
135  const dictionary& dict
136  );
137 
138 
139  //- Destructor
140  virtual ~interRegionHeatTransfer();
141 
142 
143  // Member Functions
144 
145  // Access
146 
147  //- Return the heat transfer coefficient
148  inline const volScalarField& htc() const;
149 
150 
151  // Checks
152 
153  //- Return the list of fields for which the fvModel adds source term
154  // to the transport equation
155  virtual wordList addSupFields() const;
156 
157 
158  // Sources
159 
160  //- Source term to energy equation
161  virtual void addSup
162  (
163  fvMatrix<scalar>& eqn,
164  const word& fieldName
165  ) const;
166 
167  //- Source term to compressible energy equation
168  virtual void addSup
169  (
170  const volScalarField& rho,
171  fvMatrix<scalar>& eqn,
172  const word& fieldName
173  ) const;
174 
175 
176  // Correction
177 
178  //- Correct the model
179  virtual void correct();
180 
181 
182  // Mesh changes
183 
184  //- Update for mesh motion
185  virtual bool movePoints();
186 
187  //- Update topology using the given map
188  virtual void topoChange(const polyTopoChangeMap&);
189 
190  //- Update from another mesh using the given map
191  virtual void mapMesh(const polyMeshMap&);
192 
193  //- Redistribute or update using the given distribution map
194  virtual void distribute(const polyDistributionMap&);
195 
196 
197  // IO
198 
199  //- Read dictionary
200  virtual bool read(const dictionary& dict);
201 };
202 
203 
204 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
205 
206 } // End namespace fv
207 } // End namespace Foam
208 
209 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210 
211 #endif
212 
213 // ************************************************************************* //
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:160
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:101
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:34
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:28
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(fvMatrix< scalar > &eqn, const word &fieldName) 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.
labelList fv(nPoints)
dictionary dict