tabulatedNTUHeatTransfer.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) 2015-2018 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
30 #include "surfaceFields.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  namespace fv
37  {
38  defineTypeNameAndDebug(tabulatedNTUHeatTransfer, 0);
40  (
41  option,
42  tabulatedNTUHeatTransfer,
43  dictionary
44  );
45  }
46 }
47 
48 const Foam::Enum
49 <
51 >
53 ({
54  { geometryModeType::gmCalculated, "calculated" },
55  { geometryModeType::gmUser, "user" },
56 });
57 
58 
59 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
60 
62 Foam::fv::tabulatedNTUHeatTransfer::ntuTable()
63 {
64  if (!ntuTable_.valid())
65  {
66  ntuTable_.reset(new interpolation2DTable<scalar>(coeffs_));
67  }
68 
69  return *ntuTable_;
70 }
71 
72 
73 const Foam::basicThermo& Foam::fv::tabulatedNTUHeatTransfer::thermo
74 (
75  const fvMesh& mesh
76 ) const
77 {
78  if (!mesh.foundObject<basicThermo>(basicThermo::dictName))
79  {
81  << " on mesh " << mesh.name()
82  << " could not find " << basicThermo::dictName
83  << exit(FatalError);
84  }
85 
86  return mesh.lookupObject<basicThermo>(basicThermo::dictName);
87 }
88 
89 
90 void Foam::fv::tabulatedNTUHeatTransfer::initialiseGeometry()
91 {
92  if (Ain_ < 0)
93  {
94  geometryMode_ = geometryModelNames_.get("geometryMode", coeffs_);
95 
96  Info<< "Region " << mesh_.name() << " " << type() << " " << name_
97  << " " << geometryModelNames_[geometryMode_] << " geometry:" << nl;
98 
99  switch (geometryMode_)
100  {
101  case gmCalculated:
102  {
103  const fvMesh& nbrMesh =
104  mesh_.time().lookupObject<fvMesh>(nbrRegionName());
105 
106  word inletPatchName(coeffs_.get<word>("inletPatch"));
107  word inletPatchNbrName(coeffs_.get<word>("inletPatchNbr"));
108 
109  Info<< " Inlet patch : " << inletPatchName << nl
110  << " Inlet patch neighbour : " << inletPatchNbrName
111  << nl;
112 
113  label patchI = mesh_.boundary().findPatchID(inletPatchName);
114  label patchINbr =
115  nbrMesh.boundary().findPatchID(inletPatchNbrName);
116 
117  scalar alpha(coeffs_.get<scalar>("inletBlockageRatio"));
118 
119  if (alpha < 0 || alpha > 1)
120  {
122  << "Inlet patch blockage ratio must be between 0 and 1"
123  << ". Current value: " << alpha
124  << abort(FatalError);
125  }
126 
127  scalar alphaNbr(coeffs_.get<scalar>("inletBlockageRatioNbr"));
128 
129  if (alphaNbr < 0 || alphaNbr > 1)
130  {
132  << "Inlet patch neighbour blockage ratio must be "
133  << "between 0 and 1. Current value: " << alphaNbr
134  << abort(FatalError);
135  }
136 
137  Info<< " Inlet blockage ratio : " << alpha << nl
138  << " Inlet blockage ratio neighbour : " << alphaNbr
139  << nl;
140 
141  Ain_ =
142  (scalar(1) - alpha)
143  *gSum(mesh_.magSf().boundaryField()[patchI]);
144 
145  AinNbr_ =
146  (scalar(1) - alphaNbr)
147  *gSum(nbrMesh.magSf().boundaryField()[patchINbr]);
148 
149  scalar beta(coeffs_.get<scalar>("coreBlockageRatio"));
150 
151  if (beta < 0 || beta > 1)
152  {
154  << "Core volume blockage ratio must be between 0 and 1"
155  << ". Current value: " << beta
156  << abort(FatalError);
157  }
158 
159  Info<< " Core volume blockage ratio : " << beta << nl;
160 
161  Vcore_ = (scalar(1) - beta)*meshInterp().V();
162 
163  break;
164  }
165  case gmUser:
166  {
167  coeffs_.readEntry("Ain", Ain_);
168  coeffs_.readEntry("AinNbr", AinNbr_);
169 
170  if (!coeffs_.readIfPresent("Vcore", Vcore_))
171  {
172  Vcore_ = meshInterp().V();
173  }
174 
175  break;
176  }
177  default:
178  {
180  << "Unhandled enumeration " << geometryMode_
181  << abort(FatalError);
182  }
183  }
184 
185  Info<< " Inlet area local : " << Ain_ << nl
186  << " Inlet area neighbour : " << AinNbr_ << nl
187  << " Core volume : " << Vcore_ << nl
188  << endl;
189  }
190 }
191 
192 
193 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
194 
196 (
197  const word& name,
198  const word& modelType,
199  const dictionary& dict,
200  const fvMesh& mesh
201 )
202 :
204  UName_(coeffs_.lookupOrDefault<word>("U", "U")),
205  UNbrName_(coeffs_.lookupOrDefault<word>("UNbr", "U")),
206  rhoName_(coeffs_.lookupOrDefault<word>("rho", "rho")),
207  rhoNbrName_(coeffs_.lookupOrDefault<word>("rhoNbr", "rho")),
208  ntuTable_(),
209  geometryMode_(gmCalculated),
210  Ain_(-1),
211  AinNbr_(-1),
212  Vcore_(-1)
213 {}
214 
215 
216 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
217 
219 {}
220 
221 
222 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
223 
225 {
226  initialiseGeometry();
227 
228  const fvMesh& nbrMesh = mesh_.time().lookupObject<fvMesh>(nbrRegionName());
229 
230  const basicThermo& thermo = this->thermo(mesh_);
231  const basicThermo& thermoNbr = this->thermo(nbrMesh);
232  const volScalarField Cp(thermo.Cp());
233  const volScalarField CpNbr(thermoNbr.Cp());
234 
235  // Calculate scaled mass flow for primary region
236  const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
237  const volScalarField& rho = mesh_.lookupObject<volScalarField>(rhoName_);
238  const scalarField mDot(mag(U)*rho*Ain_);
239 
240  // Calculate scaled mass flow for neighbour region
241  const volVectorField& UNbr =
242  nbrMesh.lookupObject<volVectorField>(UNbrName_);
243  const scalarField UMagNbr(mag(UNbr));
244  const scalarField UMagNbrMapped(interpolate(UMagNbr));
245  const scalarField& rhoNbr =
246  nbrMesh.lookupObject<volScalarField>(rhoNbrName_).internalField();
247  const scalarField rhoNbrMapped(interpolate(rhoNbr));
248  const scalarField mDotNbr(UMagNbrMapped*rhoNbrMapped*AinNbr_);
249 
250 
251  scalarField& htcc = htc_.primitiveFieldRef();
252  const interpolation2DTable<Foam::scalar>& ntuTable = this->ntuTable();
253 
254  forAll(htcc, cellI)
255  {
256  scalar Cpc = Cp[cellI];
257  scalar CpcNbr = CpNbr[cellI];
258  scalar mDotc = mDot[cellI];
259  scalar mDotcNbr = mDotNbr[cellI];
260  scalar Cmin = min(Cpc*mDotc, CpcNbr*mDotcNbr);
261  scalar ntu = ntuTable(mDotc, mDotcNbr);
262 
263  htcc[cellI] = Cmin*ntu/Vcore_;
264  }
265 }
266 
267 
269 {
270  if (option::read(dict))
271  {
272  coeffs_.readIfPresent("U", UName_);
273  coeffs_.readIfPresent("UNbr", UNbrName_);
274  coeffs_.readIfPresent("rho", rhoName_);
275  coeffs_.readIfPresent("rhoNbr", rhoNbrName_);
276 
277  // Force geometry re-initialisation
278  Ain_ = -1;
279  initialiseGeometry();
280 
281  return true;
282  }
283 
284  return false;
285 }
286 
287 
288 // ************************************************************************* //
Foam::fv::tabulatedNTUHeatTransfer::geometryModeType
geometryModeType
Definition: tabulatedNTUHeatTransfer.H:133
Foam::fv::tabulatedNTUHeatTransfer::geometryModelNames_
static const Enum< geometryModeType > geometryModelNames_
Definition: tabulatedNTUHeatTransfer.H:139
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:51
Foam::fv::tabulatedNTUHeatTransfer::read
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: tabulatedNTUHeatTransfer.C:268
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
tabulatedNTUHeatTransfer.H
Foam::basicThermo::Cp
virtual tmp< volScalarField > Cp() const =0
Heat capacity at constant pressure [J/kg/K].
Foam::fv::interRegionHeatTransferModel
Base class for inter region heat exchange. The derived classes must provide the heat transfer coeffis...
Definition: interRegionHeatTransferModel.H:59
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
thermo
psiReactionThermo & thermo
Definition: createFields.H:28
Foam::fv::tabulatedNTUHeatTransfer::calculateHtc
virtual void calculateHtc()
Calculate the heat transfer coefficient.
Definition: tabulatedNTUHeatTransfer.C:224
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
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
surfaceFields.H
Foam::surfaceFields.
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
rho
rho
Definition: readInitialConditions.H:96
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::interpolation2DTable
2D table interpolation. The data must be in ascending order in both dimensions x and y.
Definition: interpolation2DTable.H:53
Foam::interpolate
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::Field< scalar >
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
Foam::fv::option::coeffs_
dictionary coeffs_
Dictionary containing source coefficients.
Definition: fvOption.H:88
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
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::fv::tabulatedNTUHeatTransfer::~tabulatedNTUHeatTransfer
virtual ~tabulatedNTUHeatTransfer()
Destructor.
Definition: tabulatedNTUHeatTransfer.C:218
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
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::fv::option::read
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvOptionIO.C:55
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
fv
labelList fv(nPoints)
Foam::fv::defineTypeNameAndDebug
defineTypeNameAndDebug(option, 0)
U
U
Definition: pEqn.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::fv::addToRunTimeSelectionTable
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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
Foam::fv::tabulatedNTUHeatTransfer::tabulatedNTUHeatTransfer
tabulatedNTUHeatTransfer(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from dictionary.
Definition: tabulatedNTUHeatTransfer.C:196
Cp
const volScalarField & Cp
Definition: EEqn.H:7
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
Foam::GeometricField< scalar, fvPatchField, volMesh >