heatTransfer.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) 2021-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 
26 #include "heatTransfer.H"
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::heatTransfer::readCoeffs()
54 {
55  phaseName_ = coeffs().lookupOrDefault<word>("phase", word::null);
56 
57  semiImplicit_ = coeffs().lookup<bool>("semiImplicit");
58 
59  TName_ = coeffs().lookupOrDefault<word>("T", "T");
60 
61  Ta_ = dimensionedScalar("Ta", dimTemperature, coeffs());
62 
63  heatTransferAv_.reset(new heatTransferAv(coeffs(), mesh()));
64 
65  heatTransferCoefficientModel_ =
67 }
68 
69 
70 template<class AlphaFieldType>
71 void Foam::fv::heatTransfer::add
72 (
73  const AlphaFieldType& alpha,
74  fvMatrix<scalar>& eqn,
75  const word& fieldName
76 ) const
77 {
78  const volScalarField& he = eqn.psi();
79 
80  const volScalarField& T =
81  mesh().lookupObject<volScalarField>
82  (
83  IOobject::groupName(TName_, phaseName_)
84  );
85 
86  tmp<volScalarField> mask =
87  volScalarField::New("mask", mesh(), dimensionedScalar(dimless, 0));
88  UIndirectList<scalar>(mask.ref().primitiveFieldRef(), set_.cells()) = 1;
89  const volScalarField htcAv
90  (
91  alpha*mask*heatTransferCoefficientModel_->htc()*heatTransferAv_->Av()
92  );
93 
94  if (semiImplicit_)
95  {
96  if (he.dimensions() == dimEnergy/dimMass)
97  {
98  const basicThermo& thermo =
99  mesh().lookupObject<basicThermo>
100  (
101  IOobject::groupName(physicalProperties::typeName, phaseName_)
102  );
103 
104  const volScalarField htcAvByCpv(htcAv/thermo.Cpv());
105 
106  eqn += htcAv*(Ta_ - T) + htcAvByCpv*he - fvm::Sp(htcAvByCpv, he);
107  }
108  else if (he.dimensions() == dimTemperature)
109  {
110  eqn += htcAv*Ta_ - fvm::Sp(htcAv, he);
111  }
112  }
113  else
114  {
115  eqn += htcAv*(Ta_ - T);
116  }
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
121 
123 (
124  const word& name,
125  const word& modelType,
126  const fvMesh& mesh,
127  const dictionary& dict
128 )
129 :
130  fvModel(name, modelType, mesh, dict),
131  set_(mesh, coeffs()),
132  phaseName_(word::null),
133  semiImplicit_(false),
134  TName_(word::null),
135  Ta_("Ta", dimTemperature, NaN),
136  heatTransferAv_(nullptr),
137  heatTransferCoefficientModel_(nullptr)
138 {
139  readCoeffs();
140 }
141 
142 
143 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
144 
146 {}
147 
148 
149 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150 
152 {
153  const basicThermo& thermo =
154  mesh().lookupObject<basicThermo>
155  (
156  IOobject::groupName(physicalProperties::typeName, phaseName_)
157  );
158 
159  return wordList(1, thermo.he().name());
160 }
161 
162 
164 (
165  fvMatrix<scalar>& eqn,
166  const word& fieldName
167 ) const
168 {
169  add(geometricOneField(), eqn, fieldName);
170 }
171 
172 
174 (
175  const volScalarField& rho,
176  fvMatrix<scalar>& eqn,
177  const word& fieldName
178 ) const
179 {
180  add(geometricOneField(), eqn, fieldName);
181 }
182 
183 
185 (
186  const volScalarField& alpha,
187  const volScalarField& rho,
188  fvMatrix<scalar>& eqn,
189  const word& fieldName
190 ) const
191 {
192  add(alpha, eqn, fieldName);
193 }
194 
195 
197 {
198  heatTransferCoefficientModel_->correct();
199 }
200 
201 
203 {
204  set_.movePoints();
205  return true;
206 }
207 
208 
210 {
211  set_.topoChange(map);
212 }
213 
214 
216 {
217  set_.mapMesh(map);
218 }
219 
220 
222 {
223  set_.distribute(map);
224 }
225 
226 
228 {
229  if (fvModel::read(dict))
230  {
231  readCoeffs();
232  return true;
233  }
234  else
235  {
236  return false;
237  }
238 }
239 
240 
241 // ************************************************************************* //
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 >> &)
Return a temporary field constructed from name,.
static word groupName(Name name, const word &group)
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:160
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
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
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:40
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:187
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:34
static autoPtr< heatTransferCoefficientModel > New(const dictionary &dict, const fvMesh &mesh)
Select from dictionary and mesh.
Model for heat exchange. Requires specification of an ambient temperature with which to exchange heat...
Definition: heatTransfer.H:94
virtual bool movePoints()
Update for mesh motion.
Definition: heatTransfer.C:202
virtual ~heatTransfer()
Destructor.
Definition: heatTransfer.C:145
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
Definition: heatTransfer.C:151
virtual void correct()
Correct the model.
Definition: heatTransfer.C:196
virtual void addSup(fvMatrix< scalar > &eqn, const word &fieldName) const
Source term to energy equation.
Definition: heatTransfer.C:164
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: heatTransfer.C:209
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: heatTransfer.C:221
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: heatTransfer.C:227
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: heatTransfer.C:215
heatTransfer(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from dictionary.
Definition: heatTransfer.C:123
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
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
static const word null
An empty word.
Definition: word.H:77
Volume integrate volField creating a volField.
Calculate the matrix for implicit and explicit sources.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
static Type NaN()
Return a primitive with all components set to NaN.
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
List< word > wordList
A List of words.
Definition: fileName.H:54
const dimensionSet dimEnergy
const dimensionSet dimless
const dimensionSet dimTemperature
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
const dimensionSet dimMass
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
thermo he()
labelList fv(nPoints)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31