interRegionHeatTransferModel.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::fv::interRegionHeatTransferModel
28 
29 Group
30  grpFvOptionsSources
31 
32 Description
33  Base class for inter region heat exchange. The derived classes must
34  provide the heat transfer coeffisine (htc) which is used as follows
35  in the energy equation.
36 
37  \f[
38  -htc*Tmapped + Sp(htc, T)
39  \f]
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #ifndef interRegionHeatTransferModel_H
44 #define interRegionHeatTransferModel_H
45 
46 #include "interRegionOption.H"
47 #include "volFields.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 namespace fv
54 {
55 
56 /*---------------------------------------------------------------------------*\
57  Class interRegionHeatTransferModel Declaration
58 \*---------------------------------------------------------------------------*/
59 
61 :
62  public interRegionOption
63 {
64 protected:
65 
66  // Protected data
67 
68  //- Name of the model in the neighbour mesh
70 
71  //- Pointer to neighbour interRegionHeatTransferModel
73 
74  //- First iteration
75  bool firstIter_;
76 
77  //- Flag to activate semi-implicit coupling
78  bool semiImplicit_;
79 
80  //- Time index - used for updating htc
82 
83  //- Heat transfer coefficient [W/m2/k] times area/volume [1/m]
85 
86  //- Name of temperature field; default = "T"
87  word TName_;
88 
89  //- Name of neighbour temperature field; default = "T"
91 
92 
93  // Protected member functions
94 
95  //- Set the neighbour interRegionHeatTransferModel
96  void setNbrModel();
97 
98  //- Inherit correct from interRegionOption
100 
101  //- Correct to calculate the inter-region heat transfer coefficient
102  void correct();
103 
104  //- Interpolate field with nbrModel specified
105  template<class Type>
107  (
109  const Field<Type>& field
110  ) const;
111 
112  //- Interpolate field without nbrModel specified
113  template<class Type>
115  (
116  const Field<Type>& field
117  ) const;
118 
119  //- Interpolate field with nbrModel specified
120  template<class Type>
121  void interpolate
122  (
124  const Field<Type>& field,
125  Field<Type>& result
126  ) const;
127 
128  //- Interpolate field without nbrModel specified
129  template<class Type>
130  void interpolate
131  (
132  const Field<Type>& field,
133  Field<Type>& result
134  ) const;
135 
136 
137 public:
138 
139  //- Runtime type information
140  TypeName("interRegionHeatTransferModel");
141 
142 
143  // Constructors
144 
145  //- Construct from dictionary
147  (
148  const word& name,
149  const word& modelType,
150  const dictionary& dict,
151  const fvMesh& mesh
152  );
153 
154 
155  //- Destructor
156  virtual ~interRegionHeatTransferModel() = default;
157 
158 
159  // Member Functions
160 
161  // Access
162 
163  //- Return const access to the neighbour region name
164  inline const word& nbrRegionName() const;
165 
166  //- Return const access to the mapToMap pointer
167  inline const meshToMesh& meshInterp() const;
168 
169  //- Return the heat transfer coefficient
170  inline const volScalarField& htc() const;
171 
172  //- Return const access to the neighbour model
173  inline const interRegionHeatTransferModel& nbrModel() const;
174 
175  //- Return access to the neighbour model
177 
178  //- Source term to energy equation
179  virtual void addSup
180  (
181  fvMatrix<scalar>& eqn,
182  const label fieldi
183  );
184 
185  //- Source term to compressible energy equation
186  virtual void addSup
187  (
188  const volScalarField& rho,
189  fvMatrix<scalar>& eqn,
190  const label fieldi
191  );
192 
193  //- Calculate heat transfer coefficient
194  virtual void calculateHtc() = 0;
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 
212 
213 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
214 
215 #ifdef NoRepository
217 #endif
218 
219 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220 
221 #endif
222 
223 // ************************************************************************* //
Foam::fv::option::correct
virtual void correct(volScalarField &field)
Definition: fvOption.C:307
volFields.H
Foam::fv::interRegionHeatTransferModel::nbrRegionName
const word & nbrRegionName() const
Return const access to the neighbour region name.
Definition: interRegionHeatTransferModelI.H:31
interRegionHeatTransferModelI.H
Foam::fv::interRegionHeatTransferModel::correct
void correct()
Correct to calculate the inter-region heat transfer coefficient.
Definition: interRegionHeatTransferModel.C:89
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fv::interRegionHeatTransferModel::timeIndex_
label timeIndex_
Time index - used for updating htc.
Definition: interRegionHeatTransferModel.H:80
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::fv::option::name
const word & name() const
Return const access to the source name.
Definition: fvOptionI.H:30
Foam::fv::interRegionHeatTransferModel::TName_
word TName_
Name of temperature field; default = "T".
Definition: interRegionHeatTransferModel.H:86
Foam::fv::interRegionHeatTransferModel
Base class for inter region heat exchange. The derived classes must provide the heat transfer coeffis...
Definition: interRegionHeatTransferModel.H:59
rho
rho
Definition: readInitialConditions.H:96
Foam::fv::interRegionHeatTransferModel::nbrModel_
interRegionHeatTransferModel * nbrModel_
Pointer to neighbour interRegionHeatTransferModel.
Definition: interRegionHeatTransferModel.H:71
Foam::label
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:62
Foam::meshToMesh
Class to calculate the cell-addressing between two overlapping meshes.
Definition: meshToMesh.H:64
Foam::fv::interRegionHeatTransferModel::TypeName
TypeName("interRegionHeatTransferModel")
Runtime type information.
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::fv::interRegionOption
Base class for inter-region exchange.
Definition: interRegionOption.H:53
interRegionOption.H
Foam::fv::interRegionHeatTransferModel::calculateHtc
virtual void calculateHtc()=0
Calculate heat transfer coefficient.
Foam::fv::interRegionHeatTransferModel::nbrModel
const interRegionHeatTransferModel & nbrModel() const
Return const access to the neighbour model.
Definition: interRegionHeatTransferModelI.H:59
Foam::fv::interRegionHeatTransferModel::nbrModelName_
word nbrModelName_
Name of the model in the neighbour mesh.
Definition: interRegionHeatTransferModel.H:68
Foam::fv::interRegionHeatTransferModel::interpolate
tmp< Field< Type > > interpolate(const interRegionHeatTransferModel &nbrModel, const Field< Type > &field) const
Interpolate field with nbrModel specified.
Foam::fv::interRegionHeatTransferModel::interRegionHeatTransferModel
interRegionHeatTransferModel(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from dictionary.
Definition: interRegionHeatTransferModel.C:110
Foam::fv::interRegionHeatTransferModel::htc
const volScalarField & htc() const
Return the heat transfer coefficient.
Definition: interRegionHeatTransferModelI.H:52
field
rDeltaTY field()
Foam::fv::interRegionHeatTransferModel::htc_
volScalarField htc_
Heat transfer coefficient [W/m2/k] times area/volume [1/m].
Definition: interRegionHeatTransferModel.H:83
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::fv::interRegionHeatTransferModel::~interRegionHeatTransferModel
virtual ~interRegionHeatTransferModel()=default
Destructor.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::fv::interRegionHeatTransferModel::addSup
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
Source term to energy equation.
Definition: interRegionHeatTransferModel.C:159
Foam::fv::interRegionHeatTransferModel::firstIter_
bool firstIter_
First iteration.
Definition: interRegionHeatTransferModel.H:74
Foam::fv::interRegionHeatTransferModel::setNbrModel
void setNbrModel()
Set the neighbour interRegionHeatTransferModel.
Definition: interRegionHeatTransferModel.C:48
Foam::fv::interRegionHeatTransferModel::TNbrName_
word TNbrName_
Name of neighbour temperature field; default = "T".
Definition: interRegionHeatTransferModel.H:89
fv
labelList fv(nPoints)
interRegionHeatTransferModelTemplates.C
Foam::fv::option::mesh
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvOptionI.H:36
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::fv::interRegionHeatTransferModel::read
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: interRegionHeatTransferModel.C:265
Foam::fv::interRegionHeatTransferModel::meshInterp
const meshToMesh & meshInterp() const
Return const access to the mapToMap pointer.
Definition: interRegionHeatTransferModelI.H:38
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::fv::interRegionHeatTransferModel::semiImplicit_
bool semiImplicit_
Flag to activate semi-implicit coupling.
Definition: interRegionHeatTransferModel.H:77