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 | Copyright (C) 2015 OpenCFD Ltd.
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 
28 #include "surfaceFields.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34  namespace fv
35  {
36  defineTypeNameAndDebug(tabulatedNTUHeatTransfer, 0);
38  (
39  option,
40  tabulatedNTUHeatTransfer,
41  dictionary
42  );
43  }
44 
45  template<>
46  const char* Foam::NamedEnum
47  <
48  Foam::fv::tabulatedNTUHeatTransfer::geometryModeType,
49  2
50  >::names[] =
51  {
52  "calculated",
53  "user"
54  };
55 }
56 
58 Foam::fv::tabulatedNTUHeatTransfer::geometryModelNames_;
59 
60 
61 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
62 
64 Foam::fv::tabulatedNTUHeatTransfer::ntuTable()
65 {
66  if (!ntuTable_.valid())
67  {
68  ntuTable_.reset(new interpolation2DTable<scalar>(coeffs_));
69  }
70 
71  return ntuTable_();
72 }
73 
74 
76 (
77  const fvMesh& mesh
78 ) const
79 {
80  if (!mesh.foundObject<basicThermo>("thermophysicalProperties"))
81  {
83  << " on mesh " << mesh.name()
84  << " could not find thermophysicalProperties "
85  << exit(FatalError);
86  }
87 
88  return mesh.lookupObject<basicThermo>("thermophysicalProperties");
89 }
90 
91 
92 void Foam::fv::tabulatedNTUHeatTransfer::initialiseGeometry()
93 {
94  if (Ain_ < 0)
95  {
96  geometryMode_ =
97  geometryModelNames_.read(coeffs_.lookup("geometryMode"));
98 
99  Info<< "Region " << mesh_.name() << " " << type() << " " << name_ << " "
100  << geometryModelNames_[geometryMode_] << " geometry:" << nl;
101 
102  switch (geometryMode_)
103  {
104  case gmCalculated:
105  {
106  const fvMesh& nbrMesh =
107  mesh_.time().lookupObject<fvMesh>(nbrRegionName());
108 
109  word inletPatchName(coeffs_.lookup("inletPatch"));
110  word inletPatchNbrName(coeffs_.lookup("inletPatchNbr"));
111 
112  Info<< " Inlet patch : " << inletPatchName << nl
113  << " Inlet patch neighbour : " << inletPatchNbrName
114  << nl;
115 
116  label patchI = mesh_.boundary().findPatchID(inletPatchName);
117  label patchINbr =
118  nbrMesh.boundary().findPatchID(inletPatchNbrName);
119 
120  scalar alpha(readScalar(coeffs_.lookup("inletBlockageRatio")));
121 
122  if ((alpha < 0) || (alpha > 1))
123  {
125  << "Inlet patch blockage ratio must be between 0 and 1"
126  << ". Current value: " << alpha
127  << abort(FatalError);
128  }
129 
130  scalar alphaNbr
131  (
132  readScalar(coeffs_.lookup("inletBlockageRatioNbr"))
133  );
134 
135  if ((alphaNbr < 0) || (alphaNbr > 1))
136  {
138  << "Inlet patch neighbour blockage ratio must be "
139  << "between 0 and 1. Current value: " << alphaNbr
140  << abort(FatalError);
141  }
142 
143  Info<< " Inlet blockage ratio : " << alpha << nl
144  << " Inlet blockage ratio neighbour : " << alphaNbr
145  << nl;
146 
147  Ain_ =
148  (scalar(1) - alpha)
149  *gSum(mesh_.magSf().boundaryField()[patchI]);
150 
151  AinNbr_ =
152  (scalar(1) - alphaNbr)
153  *gSum(nbrMesh.magSf().boundaryField()[patchINbr]);
154 
155  scalar beta(readScalar(coeffs_.lookup("coreBlockageRatio")));
156 
157  if ((beta < 0) || (beta > 1))
158  {
160  << "Core volume blockage ratio must be between 0 and 1"
161  << ". Current value: " << beta
162  << abort(FatalError);
163  }
164 
165  Info<< " Core volume blockage ratio : " << beta << nl;
166 
167  Vcore_ = (scalar(1) - beta)*meshInterp().V();
168 
169  break;
170  }
171  case gmUser:
172  {
173  coeffs_.lookup("Ain") >> Ain_;
174  coeffs_.lookup("AinNbr") >> AinNbr_;
175 
176  if (!coeffs_.readIfPresent("Vcore", Vcore_))
177  {
178  Vcore_ = meshInterp().V();
179  }
180 
181  break;
182  }
183  default:
184  {
186  << "Unhandled enumeration " << geometryMode_
187  << abort(FatalError);
188  }
189  }
190 
191  Info<< " Inlet area local : " << Ain_ << nl
192  << " Inlet area neighbour : " << AinNbr_ << nl
193  << " Core volume : " << Vcore_ << nl
194  << endl;
195  }
196 }
197 
198 
199 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
200 
201 Foam::fv::tabulatedNTUHeatTransfer::tabulatedNTUHeatTransfer
202 (
203  const word& name,
204  const word& modelType,
205  const dictionary& dict,
206  const fvMesh& mesh
207 )
208 :
209  interRegionHeatTransferModel(name, modelType, dict, mesh),
210  UName_(coeffs_.lookupOrDefault<word>("UName", "U")),
211  UNbrName_(coeffs_.lookupOrDefault<word>("UNbrName", "U")),
212  rhoName_(coeffs_.lookupOrDefault<word>("rhoName", "rho")),
213  rhoNbrName_(coeffs_.lookupOrDefault<word>("rhoNbrName", "rho")),
214  ntuTable_(),
215  geometryMode_(gmCalculated),
216  Ain_(-1),
217  AinNbr_(-1),
218  Vcore_(-1)
219 {}
220 
221 
222 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
223 
224 Foam::fv::tabulatedNTUHeatTransfer::~tabulatedNTUHeatTransfer()
225 {}
226 
227 
228 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
229 
230 void Foam::fv::tabulatedNTUHeatTransfer::calculateHtc()
231 {
232  initialiseGeometry();
233 
234  const fvMesh& nbrMesh = mesh_.time().lookupObject<fvMesh>(nbrRegionName());
235 
236  const basicThermo& thermo = this->thermo(mesh_);
237  const basicThermo& thermoNbr = this->thermo(nbrMesh);
238  const volScalarField Cp(thermo.Cp());
239  const volScalarField CpNbr(thermoNbr.Cp());
240 
241  // Calculate scaled mass flow for primary region
242  const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
243  const volScalarField& rho = mesh_.lookupObject<volScalarField>(rhoName_);
244  const scalarField mDot(mag(U)*rho*Ain_);
245 
246  // Calculate scaled mass flow for neighbour region
247  const volVectorField& UNbr =
248  nbrMesh.lookupObject<volVectorField>(UNbrName_);
249  const scalarField UMagNbr(mag(UNbr));
250  const scalarField UMagNbrMapped(interpolate(UMagNbr));
251  const scalarField& rhoNbr =
252  nbrMesh.lookupObject<volScalarField>(rhoNbrName_).internalField();
253  const scalarField rhoNbrMapped(interpolate(rhoNbr));
254  const scalarField mDotNbr(UMagNbrMapped*rhoNbrMapped*AinNbr_);
255 
256 
257  scalarField& htcc = htc_.internalField();
258  const interpolation2DTable<Foam::scalar>& ntuTable = this->ntuTable();
259 
260  forAll(htcc, cellI)
261  {
262  scalar Cpc = Cp[cellI];
263  scalar CpcNbr = CpNbr[cellI];
264  scalar mDotc = mDot[cellI];
265  scalar mDotcNbr = mDotNbr[cellI];
266  scalar Cmin = min(Cpc*mDotc, CpcNbr*mDotcNbr);
267  scalar ntu = ntuTable(mDotc, mDotcNbr);
268 
269  htcc[cellI] = Cmin*ntu/Vcore_;
270  }
271 }
272 
273 
274 bool Foam::fv::tabulatedNTUHeatTransfer::read(const dictionary& dict)
275 {
276  if (option::read(dict))
277  {
278  coeffs_.readIfPresent("UName", UName_);
279  coeffs_.readIfPresent("UNbrName", UNbrName_);
280  coeffs_.readIfPresent("rhoName", rhoName_);
281  coeffs_.readIfPresent("rhoNbrName", rhoNbrName_);
282 
283  // Force geometry re-initialisation
284  Ain_ = -1;
285  initialiseGeometry();
286 
287  return true;
288  }
289  else
290  {
291  return false;
292  }
293 }
294 
295 
296 // ************************************************************************* //
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
tabulatedNTUHeatTransfer.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::MULES::interpolate
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Definition: IMULESTemplates.C:40
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:52
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:564
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
thermo
rhoThermo & thermo
Definition: setRegionFluidFields.H:3
U
U
Definition: pEqn.H:46
Foam::interpolation2DTable
2D table interpolation. The data must be in ascending order in both dimensions x and y.
Definition: interpolation2DTable.H:52
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:59
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
readScalar
#define readScalar
Definition: doubleScalar.C:38
rho
rho
Definition: pEqn.H:3
internalField
conserve internalField()+
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
fv
labelList fv(nPoints)
Foam::fv::defineTypeNameAndDebug
defineTypeNameAndDebug(option, 0)
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::fv::addToRunTimeSelectionTable
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
scalarField
volScalarField scalarField(fieldObject, mesh)
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::NamedEnum
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52