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-2021 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 {
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 :
118  interRegionOption
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  IOobject::scopedName(type(), "htc"),
135  mesh.time().timeName(),
136  mesh,
137  IOobject::NO_READ,
138  IOobject::NO_WRITE
139  ),
140  mesh,
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  coeffs_.readEntry("semiImplicit", semiImplicit_);
152  }
153 }
154 
155 
156 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
157 
159 (
160  fvMatrix<scalar>& eqn,
161  const label fieldi
162 )
163 {
164  setNbrModel();
165 
166  correct();
167 
168  const volScalarField& he = eqn.psi();
169 
170  const auto& T = mesh_.lookupObject<volScalarField>(TName_);
171 
172  auto tTmapped = volScalarField::New
173  (
174  IOobject::scopedName(type(), "Tmapped"),
175  T
176  );
177  auto& Tmapped = tTmapped.ref();
178 
179  const auto& nbrMesh = mesh_.time().lookupObject<fvMesh>(nbrRegionName_);
180 
181  const auto& Tnbr = nbrMesh.lookupObject<volScalarField>(TNbrName_);
182 
183  interpolate(Tnbr, Tmapped.primitiveFieldRef());
184 
185  if (debug)
186  {
187  Info<< "Volumetric integral of htc: "
188  << fvc::domainIntegrate(htc_).value()
189  << endl;
190 
191  if (mesh_.time().writeTime())
192  {
193  Tmapped.write();
194  htc_.write();
195  }
196  }
197 
198  if (semiImplicit_)
199  {
200  if (he.dimensions() == dimEnergy/dimMass)
201  {
202  const auto* thermoPtr =
203  mesh_.findObject<basicThermo>(basicThermo::dictName);
204 
205  if (thermoPtr)
206  {
207  const basicThermo& thermo = *thermoPtr;
208 
209  volScalarField htcByCpv(htc_/thermo.Cpv());
210 
211  eqn += htc_*(Tmapped - T) + htcByCpv*he - fvm::Sp(htcByCpv, he);
212 
213  if (debug)
214  {
215  const dimensionedScalar energy =
216  fvc::domainIntegrate(htc_*(he/thermo.Cp() - Tmapped));
217 
218  Info<< "Energy exchange from region " << nbrMesh.name()
219  << " To " << mesh_.name() << " : " << energy.value()
220  << endl;
221  }
222  }
223  else
224  {
226  << " on mesh " << mesh_.name()
227  << " could not find object basicThermo."
228  << " The available objects are: "
229  << mesh_.names()
230  << exit(FatalError);
231  }
232  }
233  else if (he.dimensions() == dimTemperature)
234  {
235  eqn += htc_*Tmapped - fvm::Sp(htc_, he);
236  }
237  }
238  else
239  {
240  eqn += htc_*(Tmapped - T);
241  }
242 }
243 
244 
246 (
247  const volScalarField& rho,
248  fvMatrix<scalar>& eqn,
249  const label fieldi
250 )
251 {
252  addSup(eqn, fieldi);
253 }
254 
255 
257 {
259  {
260  return true;
261  }
262 
263  return false;
264 }
265 
266 
267 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
volScalarField & he
Definition: YEEqn.H:52
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
dictionary dict
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Intermediate class for handling inter-region heat exchanges.
void setNbrModel()
Set the neighbour interRegionHeatTransferModel.
const word zeroGradientType
A zeroGradient patch field type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition: fvOption.H:157
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:142
interRegionHeatTransferModel(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from dictionary.
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
word nbrModelName_
Name of the model in the neighbour mesh.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Template invariant parts for fvPatchField.
Definition: fvPatchField.H:77
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
Source term to energy equation.
fv::options & fvOptions
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
bool firstIter_
Flag to determine the first iteration.
virtual bool read(const dictionary &dict)
Read dictionary.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
psiReactionThermo & thermo
Definition: createFields.H:28
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
word timeName
Definition: getTimeIndex.H:3
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:485
dynamicFvMesh & mesh
bool active() const noexcept
True if source is active.
Definition: fvOptionI.H:42
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
void correct()
Correct to calculate the inter-region heat transfer coefficient.
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
Volume integrate volField creating a volField.
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
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:64
const word & name() const noexcept
Return const access to the source name.
Definition: fvOptionI.H:24
word nbrRegionName_
Name of the neighbour region to map.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition: IOobject.C:456
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
bool semiImplicit_
Flag to activate semi-implicit coupling.
int debug
Static debugging option.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/chef2/andy/OpenFOAM/release/v2312/OpenFOAM-v2312/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=clamp((rho - rholSat)/(rhovSat - rholSat), zero_one{});alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:63
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimEnergy
interRegionHeatTransferModel * nbrModel_
Pointer to neighbour interRegionHeatTransferModel.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:41
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
Information stream (stdout output on master, null elsewhere)
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
virtual bool read(const dictionary &dict)
Read source dictionary.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
dictionary coeffs_
Dictionary containing source coefficients.
Definition: fvOption.H:152
Calculate the finiteVolume matrix for implicit and explicit sources.
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127