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-2024 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
76 {
77  const volScalarField& he = eqn.psi();
78 
79  const volScalarField& T =
80  mesh().lookupObject<volScalarField>
81  (
82  IOobject::groupName(TName_, phaseName_)
83  );
84 
85  tmp<volScalarField> mask =
86  volScalarField::New("mask", mesh(), dimensionedScalar(dimless, 0));
87  UIndirectList<scalar>(mask.ref().primitiveFieldRef(), set_.cells()) = 1;
88  const volScalarField htcAv
89  (
90  alpha*mask*heatTransferCoefficientModel_->htc()*heatTransferAv_->Av()
91  );
92 
93  if (semiImplicit_)
94  {
95  if (he.dimensions() == dimEnergy/dimMass)
96  {
97  const basicThermo& thermo =
98  mesh().lookupObject<basicThermo>
99  (
100  IOobject::groupName(physicalProperties::typeName, phaseName_)
101  );
102 
103  const volScalarField htcAvByCpv(htcAv/thermo.Cpv());
104 
105  eqn += htcAv*(Ta_ - T) + htcAvByCpv*he - fvm::Sp(htcAvByCpv, he);
106  }
107  else if (he.dimensions() == dimTemperature)
108  {
109  eqn += htcAv*Ta_ - fvm::Sp(htcAv, he);
110  }
111  }
112  else
113  {
114  eqn += htcAv*(Ta_ - T);
115  }
116 }
117 
118 
119 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
120 
122 (
123  const word& name,
124  const word& modelType,
125  const fvMesh& mesh,
126  const dictionary& dict
127 )
128 :
129  fvModel(name, modelType, mesh, dict),
130  set_(mesh, coeffs()),
131  phaseName_(word::null),
132  semiImplicit_(false),
133  TName_(word::null),
134  Ta_("Ta", dimTemperature, NaN),
135  heatTransferAv_(nullptr),
136  heatTransferCoefficientModel_(nullptr)
137 {
138  readCoeffs();
139 }
140 
141 
142 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
143 
145 {}
146 
147 
148 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
149 
151 {
152  const basicThermo& thermo =
153  mesh().lookupObject<basicThermo>
154  (
155  IOobject::groupName(physicalProperties::typeName, phaseName_)
156  );
157 
158  return wordList(1, thermo.he().name());
159 }
160 
161 
163 (
164  const volScalarField& he,
165  fvMatrix<scalar>& eqn
166 ) const
167 {
168  add(geometricOneField(), eqn);
169 }
170 
171 
173 (
174  const volScalarField& rho,
175  const volScalarField& he,
176  fvMatrix<scalar>& eqn
177 ) const
178 {
179  add(geometricOneField(), eqn);
180 }
181 
182 
184 (
185  const volScalarField& alpha,
186  const volScalarField& rho,
187  const volScalarField& he,
188  fvMatrix<scalar>& eqn
189 ) const
190 {
191  add(alpha, eqn);
192 }
193 
194 
196 {
197  heatTransferCoefficientModel_->correct();
198 }
199 
200 
202 {
203  set_.movePoints();
204  return true;
205 }
206 
207 
209 {
210  set_.topoChange(map);
211 }
212 
213 
215 {
216  set_.mapMesh(map);
217 }
218 
219 
221 {
222  set_.distribute(map);
223 }
224 
225 
227 {
228  if (fvModel::read(dict))
229  {
230  readCoeffs();
231  return true;
232  }
233  else
234  {
235  return false;
236  }
237 }
238 
239 
240 // ************************************************************************* //
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,.
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: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
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:199
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 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:201
virtual ~heatTransfer()
Destructor.
Definition: heatTransfer.C:144
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
Definition: heatTransfer.C:150
virtual void correct()
Correct the model.
Definition: heatTransfer.C:195
virtual void addSup(const volScalarField &he, fvMatrix< scalar > &eqn) const
Source term to energy equation.
Definition: heatTransfer.C:163
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: heatTransfer.C:208
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: heatTransfer.C:220
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: heatTransfer.C:226
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: heatTransfer.C:214
heatTransfer(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from dictionary.
Definition: heatTransfer.C:122
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 > &)
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:54
const dimensionSet dimEnergy
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
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:64
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