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-2021 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 {
40  defineTypeNameAndDebug(interRegionHeatTransfer, 0);
42  (
43  fvModel,
44  interRegionHeatTransfer,
45  dictionary
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  heatTransferModel_ = heatTransferModel::New(coeffs(), *this);
63  }
64 }
65 
66 
69 {
70  return
71  refCast<const interRegionHeatTransfer>(nbrModel()).heatTransferModel_;
72 }
73 
74 
75 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
76 
78 (
79  const word& name,
80  const word& modelType,
81  const dictionary& dict,
82  const fvMesh& mesh
83 )
84 :
85  interRegionModel(name, modelType, dict, mesh),
86  semiImplicit_(false),
87  TName_(word::null),
88  TNbrName_(word::null),
89  heatTransferModel_(nullptr)
90 {
91  readCoeffs();
92 }
93 
94 
95 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
96 
98 {}
99 
100 
101 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102 
104 {
105  const basicThermo& thermo =
106  mesh().lookupObject<basicThermo>(basicThermo::dictName);
107 
108  return wordList(1, thermo.he().name());
109 }
110 
111 
113 (
114  fvMatrix<scalar>& eqn,
115  const word& fieldName
116 ) const
117 {
118  const volScalarField& he = eqn.psi();
119 
120  const volScalarField& T =
121  mesh().lookupObject<volScalarField>(TName_);
122  const volScalarField& Tnbr =
123  nbrMesh().lookupObject<volScalarField>(TNbrName_);
124 
125  tmp<volScalarField> tTnbrMapped =
126  volScalarField::New(TName_ + "nbrMapped", T);
127  interpolate(Tnbr, tTnbrMapped->primitiveFieldRef());
128  volScalarField& TnbrMapped = tTnbrMapped.ref();
129 
130  // Get the heat transfer coefficient field
131  tmp<volScalarField> tHtcAoV;
132  if (master())
133  {
134  tmp<volScalarField> mask =
136  (
137  "mask",
138  mesh(),
140  );
141  tmp<volScalarField> oneNbr =
143  (
144  "one",
145  nbrMesh(),
147  );
148  interpolate(oneNbr(), mask.ref().primitiveFieldRef());
149  tHtcAoV =
150  mask
151  *heatTransferModel_->htc()
152  *heatTransferModel_->AoV();
153  }
154  else
155  {
156  tmp<volScalarField> tHtcNbr =
157  nbrHeatTransferModel().htc()
158  *nbrHeatTransferModel().AoV();
159  tHtcAoV =
161  (
162  tHtcNbr().name(),
163  mesh(),
164  dimensionedScalar(tHtcNbr().dimensions(), 0)
165  );
166  interpolate(tHtcNbr(), tHtcAoV.ref().primitiveFieldRef());
167  }
168  const volScalarField& htcAoV = tHtcAoV();
169 
170  if (semiImplicit_)
171  {
172  if (he.dimensions() == dimEnergy/dimMass)
173  {
174  const basicThermo& thermo =
175  mesh().lookupObject<basicThermo>(basicThermo::dictName);
176 
177  const volScalarField htcAoVByCpv(htcAoV/thermo.Cpv());
178 
179  eqn +=
180  htcAoV*(TnbrMapped - T)
181  + htcAoVByCpv*he - fvm::Sp(htcAoVByCpv, he);
182  }
183  else if (he.dimensions() == dimTemperature)
184  {
185  eqn += htcAoV*TnbrMapped - fvm::Sp(htcAoV, he);
186  }
187  }
188  else
189  {
190  eqn += htcAoV*(TnbrMapped - T);
191  }
192 }
193 
194 
196 (
197  const volScalarField& rho,
198  fvMatrix<scalar>& eqn,
199  const word& fieldName
200 ) const
201 {
202  addSup(eqn, fieldName);
203 }
204 
205 
207 {
208  if (master())
209  {
210  heatTransferModel_->correct();
211  }
212 }
213 
214 
216 {
217  if (interRegionModel::read(dict))
218  {
219  readCoeffs();
220  return true;
221  }
222  else
223  {
224  return false;
225  }
226 }
227 
228 
229 // ************************************************************************* //
virtual void correct()
Correct the model.
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
const word & name() const
Return name.
Definition: IOobject.H:303
fluidReactionThermo & thermo
Definition: createFields.H:28
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:77
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Base class for heat transfer coefficient modelling used in heat transfer fvModels. Area per unit volume [1/m] (AoV) must be provided as a value in the coefficients dictionary or as a field in constant.
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:282
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
const dimensionSet dimless
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
virtual tmp< volScalarField > Cpv() const =0
Heat capacity at constant pressure/volume [J/kg/K].
Macros for easy insertion into run-time selection tables.
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
virtual ~interRegionHeatTransfer()
Destructor.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
interRegionHeatTransfer(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from dictionary.
const dimensionSet & dimensions() const
Return dimensions.
Base class for inter-region exchange.
const heatTransferModel & nbrHeatTransferModel() const
Get the neighbour heat transfer model.
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
virtual void addSup(fvMatrix< scalar > &eqn, const word &fieldName) const
Source term to energy equation.
static const word null
An empty word.
Definition: word.H:77
Volume integrate volField creating a volField.
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
virtual bool read(const dictionary &dict)
Read dictionary.
static autoPtr< heatTransferModel > New(const dictionary &dict, const fvMesh &mesh)
Select from dictionary and mesh.
static const word dictName
Name of the thermophysical properties dictionary.
Definition: basicThermo.H:129
thermo he()
const dimensionSet dimEnergy
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimMass
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Internal & ref()
Return a reference to the dimensioned internal field.
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
List< word > wordList
A List of words.
Definition: fileName.H:54
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual bool read(const dictionary &dict)
Read dictionary.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A class for managing temporary objects.
Definition: PtrList.H:53
Calculate the matrix for implicit and explicit sources.
const dimensionSet dimTemperature
Namespace for OpenFOAM.