referenceTemperature.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) 2022 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 
28 #include "referenceTemperature.H"
29 #include "fvMesh.H"
30 #include "basicThermo.H"
31 #include "surfaceInterpolate.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace heatExchangerModels
39 {
40  defineTypeNameAndDebug(referenceTemperature, 0);
42  (
43  heatExchangerModel,
44  referenceTemperature,
45  dictionary
46  );
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 Foam::scalar
54 Foam::heatExchangerModels::referenceTemperature::primaryNetMassFlux() const
55 {
57 
58  scalar sumPhi = 0;
59 
60  forAll(faceId_, i)
61  {
62  const label facei = faceId_[i];
63  if (facePatchId_[i] != -1)
64  {
65  const label patchi = facePatchId_[i];
66  sumPhi += phi.boundaryField()[patchi][facei]*faceSign_[i];
67  }
68  else
69  {
70  sumPhi += phi[facei]*faceSign_[i];
71  }
72  }
73  reduce(sumPhi, sumOp<scalar>());
74 
75  return sumPhi;
76 }
77 
78 
79 Foam::scalar
80 Foam::heatExchangerModels::referenceTemperature::primaryInletTemperature() const
81 {
82  const auto& phi = mesh_.lookupObject<surfaceScalarField>(phiName_);
83  const auto& T = mesh_.lookupObject<volScalarField>(TName_);
84 
86 
87  scalar sumMagPhi = 0;
88  scalar primaryInletTfMean = 0;
89 
90  forAll(faceId_, i)
91  {
92  const label facei = faceId_[i];
93  if (facePatchId_[i] != -1)
94  {
95  const label patchi = facePatchId_[i];
96  const scalar phii = phi.boundaryField()[patchi][facei]*faceSign_[i];
97  const scalar magPhii = mag(phii);
98 
99  sumMagPhi += magPhii;
100 
101  const scalar Tfi = Tf.boundaryField()[patchi][facei];
102 
103  primaryInletTfMean += Tfi*magPhii;
104  }
105  else
106  {
107  const scalar phii = phi[facei]*faceSign_[i];
108  const scalar magPhii = mag(phii);
109 
110  sumMagPhi += magPhii;
111 
112  primaryInletTfMean += Tf[facei]*magPhii;
113  }
114  }
115  reduce(sumMagPhi, sumOp<scalar>());
116  reduce(primaryInletTfMean, sumOp<scalar>());
117 
118  return primaryInletTfMean/(sumMagPhi + ROOTVSMALL);
119 }
120 
121 
123 Foam::heatExchangerModels::referenceTemperature::temperatureDifferences
124 (
125  const labelList& cells
126 ) const
127 {
128  const auto& T = mesh_.lookupObject<volScalarField>(TName_);
129  const scalarField TCells(T, cells);
130  scalarField deltaTCells(cells.size(), Zero);
131 
132  if (Qt_ > 0)
133  {
134  forAll(deltaTCells, i)
135  {
136  deltaTCells[i] = max(Tref_ - TCells[i], scalar(0));
137  }
138  }
139  else
140  {
141  forAll(deltaTCells, i)
142  {
143  deltaTCells[i] = max(TCells[i] - Tref_, scalar(0));
144  }
145  }
146 
147  return deltaTCells;
148 }
149 
150 
151 Foam::scalar
152 Foam::heatExchangerModels::referenceTemperature::weight
153 (
154  const labelList& cells,
155  const scalarField& deltaTCells
156 ) const
157 {
158  scalar sumWeight = 0;
159 
160  const auto& U = mesh_.lookupObject<volVectorField>(UName_);
161  const scalarField& V = mesh_.V();
162 
163  forAll(cells, i)
164  {
165  const label celli = cells[i];
166  sumWeight += V[celli]*mag(U[celli])*deltaTCells[i];
167  }
168  reduce(sumWeight, sumOp<scalar>());
169 
170  return sumWeight;
171 }
172 
173 
174 void Foam::heatExchangerModels::referenceTemperature::writeFileHeader
175 (
176  Ostream& os
177 ) const
178 {
179  writeFile::writeHeader(os, "Effectiveness heat exchanger source");
180  writeFile::writeCommented(os, "Time");
181  writeFile::writeTabbed(os, "Net mass flux [kg/s]");
182  writeFile::writeTabbed(os, "Total heat exchange [W]");
183  writeFile::writeTabbed(os, "Reference T [K]");
184  os << endl;
185 }
186 
187 
188 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
189 
191 (
192  const fvMesh& mesh,
193  const word& name,
194  const dictionary& coeffs
195 )
196 :
197  heatExchangerModel(mesh, name, coeffs),
198  targetQdotPtr_
199  (
200  Function1<scalar>::New
201  (
202  "targetQdot",
203  coeffs,
204  &mesh_
205  )
206  ),
207  TrefTablePtr_(nullptr),
208  sumPhi_(0),
209  Qt_(0),
210  Tref_(0)
211 {
212  writeFileHeader(file());
213 }
214 
215 
216 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
217 
219 {
221 }
222 
223 
226 (
227  const labelList& cells
228 )
229 {
230  sumPhi_ = primaryNetMassFlux();
231 
232  Qt_ = targetQdotPtr_->value(mesh_.time().value());
233 
234  if (TrefTablePtr_)
235  {
236  const scalar primaryInletT = primaryInletTemperature();
237  Tref_ = TrefTablePtr_()(mag(sumPhi_), primaryInletT);
238  }
239 
240  const scalarField deltaTCells(temperatureDifferences(cells));
241 
242  const scalar sumWeight = weight(cells, deltaTCells);
243 
244  return Qt_*deltaTCells/(sumWeight + ROOTVSMALL);
245 }
246 
247 
249 (
250  const dictionary& dict
251 )
252 {
253  if (!writeFile::read(dict))
254  {
255  return false;
256  }
257 
258  Info<< incrIndent << indent << "- using model: " << type() << endl;
259 
260  if (coeffs_.readIfPresent("Tref", Tref_))
261  {
262  Info<< indent << "- using constant reference temperature: " << Tref_
263  << endl;
264  }
265  else
266  {
267  TrefTablePtr_.reset(new interpolation2DTable<scalar>(coeffs_));
268 
269  Info<< indent << "- using reference temperature table"
270  << endl;
271  }
272 
273  UName_ = coeffs_.getOrDefault<word>("U", "U");
274  TName_ = coeffs_.getOrDefault<word>("T", "T");
275  phiName_ = coeffs_.getOrDefault<word>("phi", "phi");
276  coeffs_.readEntry("faceZone", faceZoneName_);
278  Info<< decrIndent;
279 
280  return true;
281 }
282 
283 
285 {
286  if (log)
287  {
288  Info<< nl
289  << type() << ": " << name_ << nl << incrIndent
290  << indent << "Net mass flux [kg/s] : " << sumPhi_ << nl
291  << indent << "Total heat exchange [W] : " << Qt_ << nl
292  << indent << "Reference T [K] : " << Tref_
293  << decrIndent;
294  }
295 
296  if (Pstream::master())
297  {
298  Ostream& os = file();
299  writeCurrentTime(os);
300 
301  os << tab << sumPhi_
302  << tab << Qt_
303  << tab << Tref_
304  << endl;
305  }
306 
307  if (log) Info<< nl << endl;
308 }
309 
310 
311 // ************************************************************************* //
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
labelList faceId_
Local list of face IDs.
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...
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:493
static void writeHeader(Ostream &os, const word &fieldName)
addToRunTimeSelectionTable(heatExchangerModel, effectivenessTable, dictionary)
dimensionedScalar log(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
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
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
virtual void initialise()
Initialise data members of the model.
constexpr char tab
The tab &#39;\t&#39; character(0x09)
Definition: Ostream.H:49
labelList faceSign_
List of +1/-1 representing face flip map (1 use as is, -1 negate)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
Macros for easy insertion into run-time selection tables.
Base class for heat exchanger models to handle various characteristics for the heatExchangerSource fv...
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
word phiName_
Name of operand flux field.
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
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const cellShapeList & cells
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const fvMesh & mesh_
Reference to the mesh.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
OBJstream os(runTime.globalPath()/outputName)
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:511
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
virtual bool read(const dictionary &dict)
Read top-level dictionary.
labelList facePatchId_
Local list of patch IDs per face.
U
Definition: pEqn.H:72
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1094
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual void write(const bool log)
Write data to stream and files.
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual void initialise()
Initialise data members of the model.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:502
referenceTemperature(const fvMesh &mesh, const word &name, const dictionary &coeffs)
Construct from components.
void reduce(T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) using linear/tree communication schedule.
virtual tmp< scalarField > energyDensity(const labelList &cells)
Return energy density per unit length [J/m3/m].
Namespace for OpenFOAM.
defineTypeNameAndDebug(effectivenessTable, 0)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127