interRegionHeatTransferModel.C
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  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "basicThermo.H"
31 #include "fvmSup.H"
33 #include "fvcVolumeIntegrate.H"
34 #include "fvOptionList.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace fv
41 {
42  defineTypeNameAndDebug(interRegionHeatTransferModel, 0);
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * Protected member functions * * * * * * * * * * * //
48 
50 {
51  if (!firstIter_)
52  {
53  return;
54  }
55 
56  const auto& nbrMesh = mesh_.time().lookupObject<fvMesh>(nbrRegionName_);
57 
58  const auto& fvOptions = nbrMesh.lookupObject<optionList>("fvOptions");
59 
60  bool nbrModelFound = false;
61 
62  forAll(fvOptions, i)
63  {
64  if (fvOptions[i].name() == nbrModelName_)
65  {
67  (
68  refCast<const interRegionHeatTransferModel>(fvOptions[i])
69  );
70  nbrModelFound = true;
71  break;
72  }
73  }
74 
75  if (!nbrModelFound)
76  {
78  << "Neighbour model not found" << nbrModelName_
79  << " in region " << nbrMesh.name() << nl
80  << exit(FatalError);
81  }
82 
83  firstIter_ = false;
84 
85  // Set nbr model's nbr model to avoid construction order problems
87 }
88 
89 
91 {
92  if (master_)
93  {
94  if (mesh_.time().timeIndex() != timeIndex_)
95  {
96  calculateHtc();
97  timeIndex_ = mesh_.time().timeIndex();
98  }
99  }
100  else
101  {
102  nbrModel().correct();
103  interpolate(nbrModel().htc(), htc_);
104  }
105 }
106 
107 
108 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
109 
111 (
112  const word& name,
113  const word& modelType,
114  const dictionary& dict,
115  const fvMesh& mesh
116 )
117 :
119  (
120  name,
121  modelType,
122  dict,
123  mesh
124  ),
125  nbrModelName_(coeffs_.get<word>("nbrModel")),
126  nbrModel_(nullptr),
127  firstIter_(true),
128  semiImplicit_(false),
129  timeIndex_(-1),
130  htc_
131  (
132  IOobject
133  (
134  type() + ":htc",
135  mesh.time().timeName(),
136  mesh,
139  ),
140  mesh,
142  zeroGradientFvPatchScalarField::typeName
143  ),
144  TName_(coeffs_.getOrDefault<word>("T", "T")),
145  TNbrName_(coeffs_.getOrDefault<word>("TNbr", "T"))
146 {
147  if (active())
148  {
149  coeffs_.readEntry("fields", fieldNames_);
150  applied_.setSize(fieldNames_.size(), false);
151 
152  coeffs_.readEntry("semiImplicit", semiImplicit_);
153  }
154 }
155 
156 
157 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
158 
160 (
161  fvMatrix<scalar>& eqn,
162  const label fieldi
163 )
164 {
165  setNbrModel();
166 
167  correct();
168 
169  const volScalarField& he = eqn.psi();
170 
171  const auto& T = mesh_.lookupObject<volScalarField>(TName_);
172 
173  auto tTmapped = tmp<volScalarField>::New
174  (
175  IOobject
176  (
177  type() + ":Tmapped",
178  mesh_.time().timeName(),
179  mesh_,
182  ),
183  T
184  );
185 
186  auto& Tmapped = tTmapped.ref();
187 
188  const auto& nbrMesh = mesh_.time().lookupObject<fvMesh>(nbrRegionName_);
189 
190  const auto& Tnbr = nbrMesh.lookupObject<volScalarField>(TNbrName_);
191 
192  interpolate(Tnbr, Tmapped.primitiveFieldRef());
193 
194  if (debug)
195  {
196  Info<< "Volumetric integral of htc: "
197  << fvc::domainIntegrate(htc_).value()
198  << endl;
199 
200  if (mesh_.time().writeTime())
201  {
202  Tmapped.write();
203  htc_.write();
204  }
205  }
206 
207  if (semiImplicit_)
208  {
209  if (he.dimensions() == dimEnergy/dimMass)
210  {
211  const auto* thermoPtr =
212  mesh_.findObject<basicThermo>(basicThermo::dictName);
213 
214  if (thermoPtr)
215  {
216  const basicThermo& thermo = *thermoPtr;
217 
218  volScalarField htcByCpv(htc_/thermo.Cpv());
219 
220  eqn += htc_*(Tmapped - T) + htcByCpv*he - fvm::Sp(htcByCpv, he);
221 
222  if (debug)
223  {
224  const dimensionedScalar energy =
225  fvc::domainIntegrate(htc_*(he/thermo.Cp() - Tmapped));
226 
227  Info<< "Energy exchange from region " << nbrMesh.name()
228  << " To " << mesh_.name() << " : " << energy.value()
229  << endl;
230  }
231  }
232  else
233  {
235  << " on mesh " << mesh_.name()
236  << " could not find object basicThermo."
237  << " The available objects are: "
238  << mesh_.names()
239  << exit(FatalError);
240  }
241  }
242  else if (he.dimensions() == dimTemperature)
243  {
244  eqn += htc_*Tmapped - fvm::Sp(htc_, he);
245  }
246  }
247  else
248  {
249  eqn += htc_*(Tmapped - T);
250  }
251 }
252 
253 
255 (
256  const volScalarField& rho,
257  fvMatrix<scalar>& eqn,
258  const label fieldi
259 )
260 {
261  addSup(eqn, fieldi);
262 }
263 
264 
266 {
268  {
269  return true;
270  }
271 
272  return false;
273 }
274 
275 
276 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::fv::interRegionOption::nbrRegionName_
word nbrRegionName_
Name of the neighbour region to map.
Definition: interRegionOption.H:119
basicThermo.H
Foam::fv::interRegionHeatTransferModel::correct
void correct()
Correct to calculate the inter-region heat transfer coefficient.
Definition: interRegionHeatTransferModel.C:90
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::fv::option::name
const word & name() const
Return const access to the source name.
Definition: fvOptionI.H:30
Foam::fv::interRegionHeatTransferModel
Intermediate class for handling inter-region heat exchanges.
Definition: interRegionHeatTransferModel.H:142
Foam::fvc::domainIntegrate
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcVolumeIntegrate.C:88
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::basicThermo
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:63
interRegionHeatTransferModel.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
rho
rho
Definition: readInitialConditions.H:88
Foam::fvMatrix::psi
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:287
Foam::fv::interRegionHeatTransferModel::nbrModel_
interRegionHeatTransferModel * nbrModel_
Pointer to neighbour interRegionHeatTransferModel.
Definition: interRegionHeatTransferModel.H:154
Foam::fv::option::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:139
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::interpolate
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
fvOptionList.H
Foam::fv::interRegionOption
Intermediate class for handling inter-region exchanges.
Definition: interRegionOption.H:107
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
correct
fvOptions correct(rho)
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::fv::interRegionHeatTransferModel::nbrModelName_
word nbrModelName_
Name of the model in the neighbour mesh.
Definition: interRegionHeatTransferModel.H:151
Foam::fv::interRegionHeatTransferModel::interRegionHeatTransferModel
interRegionHeatTransferModel(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from dictionary.
Definition: interRegionHeatTransferModel.C:111
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::dictionary::dictName
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionary.H:458
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned< scalar >
Foam::Ostream::write
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::fv::interRegionHeatTransferModel::addSup
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
Source term to energy equation.
Definition: interRegionHeatTransferModel.C:160
Foam::fv::interRegionHeatTransferModel::firstIter_
bool firstIter_
Flag to determine the first iteration.
Definition: interRegionHeatTransferModel.H:157
fvmSup.H
Calculate the matrix for implicit and explicit sources.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::fv::interRegionHeatTransferModel::setNbrModel
void setNbrModel()
Set the neighbour interRegionHeatTransferModel.
Definition: interRegionHeatTransferModel.C:49
fv
labelList fv(nPoints)
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
he
volScalarField & he
Definition: YEEqn.H:52
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::GeometricField::ref
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Definition: GeometricField.C:749
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
fvcVolumeIntegrate.H
Volume integrate volField creating a volField.
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
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 source dictionary.
Definition: interRegionHeatTransferModel.C:265
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:55
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:275
Foam::fv::defineTypeNameAndDebug
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::fv::optionList
List of finite volume options.
Definition: fvOptionList.H:69
zeroGradientFvPatchFields.H
Foam::fv::interRegionOption::read
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: interRegionOption.C:127