interRegionHeatTransfer.C
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 \*---------------------------------------------------------------------------*/
25 
27 #include "basicThermo.H"
28 #include "fvmSup.H"
30 #include "fvcVolumeIntegrate.H"
31 #include "fvModels.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace fv
39 {
42  (
43  fvModel,
46  );
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 void Foam::fv::interRegionHeatTransfer::readCoeffs()
54 {
55  semiImplicit_ = coeffs().lookup<bool>("semiImplicit");
56 
57  TName_ = coeffs().lookupOrDefault<word>("T", "T");
58  TNbrName_ = coeffs().lookupOrDefault<word>("TNbr", "T");
59 
60  if (master())
61  {
62  heatTransferAv_.reset(new heatTransferAv(coeffs(), mesh()));
63 
64  heatTransferCoefficientModel_ =
66  }
67 }
68 
69 
71 Foam::fv::interRegionHeatTransfer::nbrHeatTransfer() const
72 {
73  return refCast<const interRegionHeatTransfer>(nbrModel());
74 }
75 
76 
77 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
78 
80 (
81  const word& name,
82  const word& modelType,
83  const fvMesh& mesh,
84  const dictionary& dict
85 )
86 :
87  interRegionModel(name, modelType, mesh, dict),
88  semiImplicit_(false),
89  TName_(word::null),
90  TNbrName_(word::null),
91  heatTransferAv_(nullptr),
92  heatTransferCoefficientModel_(nullptr)
93 {
94  readCoeffs();
95 }
96 
97 
98 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
99 
101 {}
102 
103 
104 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
105 
107 {
108  const basicThermo& thermo =
109  mesh().lookupObject<basicThermo>(physicalProperties::typeName);
110 
111  return wordList(1, thermo.he().name());
112 }
113 
114 
116 (
117  const volScalarField& he,
118  fvMatrix<scalar>& eqn
119 ) const
120 {
121  const volScalarField& T =
122  mesh().lookupObject<volScalarField>(TName_);
123 
124  tmp<volScalarField> tTnbr = volScalarField::New(TNbrName_, T);
126  (
127  nbrMesh().lookupObject<volScalarField>(TNbrName_),
128  tTnbr->primitiveFieldRef()
129  );
130  const volScalarField& Tnbr = tTnbr();
131 
132  // Get the heat transfer coefficient field
133  tmp<volScalarField> tHtcAv;
134  if (master())
135  {
136  tmp<volScalarField> mask =
138  (
139  "mask",
140  mesh(),
142  );
143  tmp<volScalarField> oneNbr =
145  (
146  "one",
147  nbrMesh(),
149  );
150  interpolate(oneNbr(), mask.ref().primitiveFieldRef());
151  tHtcAv =
152  mask
153  *heatTransferCoefficientModel_->htc()
154  *heatTransferAv_->Av();
155  }
156  else
157  {
158  tmp<volScalarField> tHtcNbr =
159  nbrHeatTransfer().heatTransferCoefficientModel_->htc()
160  *nbrHeatTransfer().heatTransferAv_->Av();
161  tHtcAv =
163  (
164  tHtcNbr().name(),
165  mesh(),
166  dimensionedScalar(tHtcNbr().dimensions(), 0)
167  );
168  interpolate(tHtcNbr(), tHtcAv.ref().primitiveFieldRef());
169  }
170  const volScalarField& htcAv = tHtcAv();
171 
172  if (semiImplicit_)
173  {
174  if (he.dimensions() == dimEnergy/dimMass)
175  {
176  const basicThermo& thermo =
177  mesh().lookupObject<basicThermo>(physicalProperties::typeName);
178 
179  const volScalarField htcAvByCpv(htcAv/thermo.Cpv());
180 
181  eqn +=
182  htcAv*(Tnbr - T)
183  + htcAvByCpv*he - fvm::Sp(htcAvByCpv, he);
184  }
185  else if (he.dimensions() == dimTemperature)
186  {
187  eqn += htcAv*Tnbr - fvm::Sp(htcAv, he);
188  }
189  }
190  else
191  {
192  eqn += htcAv*(Tnbr - T);
193  }
194 }
195 
196 
198 (
199  const volScalarField& rho,
200  const volScalarField& he,
201  fvMatrix<scalar>& eqn
202 ) const
203 {
204  addSup(he, eqn);
205 }
206 
207 
209 {
210  if (master())
211  {
212  heatTransferCoefficientModel_->correct();
213  }
214 }
215 
216 
218 {
220  return true;
221 }
222 
223 
225 {
227 }
228 
229 
231 {
233 }
234 
235 
237 (
238  const polyDistributionMap&
239 )
240 {
242 }
243 
244 
246 {
248  {
249  readCoeffs();
250  return true;
251  }
252  else
253  {
254  return false;
255  }
256 }
257 
258 
259 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:78
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
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
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:59
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:53
static autoPtr< heatTransferCoefficientModel > New(const dictionary &dict, const fvMesh &mesh)
Select from dictionary and mesh.
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.
Base class for inter-region exchange.
bool master() const
Return whether the master region.
virtual bool read(const dictionary &dict)
Read dictionary.
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 managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
Volume integrate volField creating a volField.
Calculate the matrix for implicit and explicit sources.
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:54
const dimensionSet dimEnergy
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
const dimensionSet dimTemperature
const dimensionSet dimMass
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
thermo he()
labelList fv(nPoints)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31