Loading...
Searching...
No Matches
effectivenessTable.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-2023 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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 "effectivenessTable.H"
29#include "basicThermo.H"
30#include "surfaceInterpolate.H"
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace heatExchangerModels
38{
41 (
45 );
46}
47}
48
49
50// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51
52void Foam::heatExchangerModels::effectivenessTable::writeFileHeader
53(
54 Ostream& os
55) const
56{
57 writeFile::writeHeader(os, "Heat exchanger source");
58 writeFile::writeCommented(os, "Time");
59 writeFile::writeTabbed(os, "Net mass flux [kg/s]");
60 writeFile::writeTabbed(os, "Total heat exchange [W]");
61 writeFile::writeTabbed(os, "Secondary inlet T [K]");
62 writeFile::writeTabbed(os, "Reference T [K]");
63 writeFile::writeTabbed(os, "Effectiveness [-]");
64
65 if (secondaryCpPtr_)
66 {
67 writeFile::writeTabbed(os, "Secondary outlet T [K]");
68 }
70 os << endl;
71}
72
73
74// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
75
77(
78 const fvMesh& mesh,
79 const word& name,
80 const dictionary& coeffs
81)
82:
84 userPrimaryInletT_(false),
85 targetQdotActive_(false),
86 secondaryCpPtr_
87 (
88 Function1<scalar>::NewIfPresent
89 (
90 "secondaryCp",
91 coeffs,
92 word::null,
93 &mesh
94 )
95 ),
96 eTable_(),
97 targetQdotCalcInterval_(5),
98 secondaryMassFlowRate_(0),
99 secondaryInletT_(0),
100 primaryInletT_(0),
101 targetQdot_(0),
102 targetQdotRelax_(0.5),
103 sumPhi_(0),
104 Qt_(0),
105 Tref_(0),
106 effectiveness_(0)
108 writeFileHeader(file());
109}
110
111
112// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
113
115{
117
119}
120
121
124(
125 const labelList& cells
126)
127{
128 const auto& thermo = mesh_.lookupObject<basicThermo>(basicThermo::dictName);
129 const auto& phi = mesh_.lookupObject<surfaceScalarField>(phiName_);
130 const auto& T = mesh_.lookupObject<volScalarField>(TName_);
131
134
135 sumPhi_ = 0;
136 scalar sumMagPhi = 0;
137 scalar CpfMean = 0;
138 scalar primaryInletTfMean = 0;
139 forAll(faceId_, i)
140 {
141 const label facei = faceId_[i];
142 if (facePatchId_[i] != -1)
143 {
144 const label patchi = facePatchId_[i];
145 const scalar phii = phi.boundaryField()[patchi][facei]*faceSign_[i];
146 const scalar magPhii = mag(phii);
147
148 sumPhi_ += phii;
149 sumMagPhi += magPhii;
150
151 const scalar Cpfi = Cpf.boundaryField()[patchi][facei];
152 const scalar Tfi = Tf.boundaryField()[patchi][facei];
153
154 CpfMean += Cpfi*magPhii;
155 primaryInletTfMean += Tfi*magPhii;
156 }
157 else
158 {
159 const scalar phii = phi[facei]*faceSign_[i];
160 const scalar magPhii = mag(phii);
161
162 sumPhi_ += phii;
163 sumMagPhi += magPhii;
164
165 CpfMean += Cpf[facei]*magPhii;
166 primaryInletTfMean += Tf[facei]*magPhii;
167 }
168 }
169 reduce(CpfMean, sumOp<scalar>());
170 reduce(sumPhi_, sumOp<scalar>());
171 reduce(sumMagPhi, sumOp<scalar>());
172
173 CpfMean /= sumMagPhi + ROOTVSMALL;
174
175 scalar primaryInletT = primaryInletT_;
176 if (!userPrimaryInletT_)
177 {
178 reduce(primaryInletTfMean, sumOp<scalar>());
179 primaryInletT = primaryInletTfMean/(sumMagPhi + ROOTVSMALL);
180 }
181
182 effectiveness_ = eTable_()(mag(sumPhi_), secondaryMassFlowRate_);
183
184 const scalar alpha = effectiveness_*CpfMean*mag(sumPhi_);
185
186 Qt_ = alpha*(secondaryInletT_ - primaryInletT);
187
188 if
189 (
190 targetQdotActive_
191 && (mesh_.time().timeIndex() % targetQdotCalcInterval_ == 0)
192 )
193 {
194 const scalar dT = (targetQdot_ - Qt_)/(alpha + ROOTVSMALL);
195 secondaryInletT_ += targetQdotRelax_*dT;
196 }
197
198 // start with a copy
199 scalarField deltaTCells(T, cells);
200 Tref_ = 0;
201 if (Qt_ > 0)
202 {
203 Tref_ = gMax(deltaTCells);
204 for (scalar& delta : deltaTCells)
205 {
206 delta = Foam::max(Tref_ - delta, scalar(0));
207 }
208 }
209 else
210 {
211 Tref_ = gMin(deltaTCells);
212 for (scalar& delta : deltaTCells)
213 {
214 delta = Foam::max(delta - Tref_, scalar(0));
215 }
216 }
217
218 const auto& U = mesh_.lookupObject<volVectorField>(UName_);
219 const scalarField& V = mesh_.V();
220 scalar sumWeight = 0;
221 forAll(cells, i)
222 {
223 const label celli = cells[i];
224 sumWeight += V[celli]*mag(U[celli])*deltaTCells[i];
226 reduce(sumWeight, sumOp<scalar>());
227
228 return Qt_*deltaTCells/(sumWeight + ROOTVSMALL);
229}
230
231
233{
234 if (!writeFile::read(dict))
235 {
236 return false;
237 }
238
239 Info<< incrIndent << indent << "- using model: " << type() << endl;
240
241 coeffs_.readEntry("secondaryMassFlowRate", secondaryMassFlowRate_);
242 coeffs_.readEntry("secondaryInletT", secondaryInletT_);
243
244 if (coeffs_.readIfPresent("primaryInletT", primaryInletT_))
245 {
246 userPrimaryInletT_ = true;
247
248 Info<< indent
249 << "- using user-specified primary flow inlet temperature: "
250 << primaryInletT_ << endl;
251 }
252 else
253 {
254 Info<< indent
255 << "- using flux-weighted primary flow inlet temperature"
256 << endl;
257 }
258
259 if (coeffs_.readIfPresent("targetQdot", targetQdot_))
260 {
261 targetQdotActive_ = true;
262
263 Info<< indent << "- using target heat rejection: " << targetQdot_ << nl;
264
265 coeffs_.readIfPresent
266 (
267 "targetQdotCalcInterval",
268 targetQdotCalcInterval_
269 );
270
271 Info<< indent << "- updating secondary inlet temperature every "
272 << targetQdotCalcInterval_ << " iterations" << nl;
273
274 coeffs_.readIfPresent("targetQdotRelax", targetQdotRelax_);
275
276 Info<< indent << "- temperature relaxation: "
277 << targetQdotRelax_ << endl;
278 }
279
280 UName_ = coeffs_.getOrDefault<word>("U", "U");
281 TName_ = coeffs_.getOrDefault<word>("T", "T");
282 phiName_ = coeffs_.getOrDefault<word>("phi", "phi");
283 coeffs_.readEntry("faceZone", faceZoneName_);
286
287 return true;
288}
289
290
292(
293 const bool log
294)
295{
296 if (log)
297 {
298 Info<< nl
299 << type() << ": " << name_ << nl << incrIndent
300 << indent << "Net mass flux [kg/s] : " << sumPhi_ << nl
301 << indent << "Total heat exchange [W] : " << Qt_ << nl
302 << indent << "Secondary inlet T [K] : " << secondaryInletT_<< nl
303 << indent << "Reference T [K] : " << Tref_ << nl
304 << indent << "Effectiveness [-] : " << effectiveness_
305 << decrIndent;
306 }
307
308 if (Pstream::master())
309 {
310 Ostream& os = file();
311 writeCurrentTime(os);
312
313 os << tab << sumPhi_
314 << tab << Qt_
315 << tab << secondaryInletT_
316 << tab << Tref_
317 << tab << effectiveness_;
318
319 if (secondaryCpPtr_)
320 {
321 // Secondary Cp as a function of the starting secondary temperature
322 const scalar secondaryCp = secondaryCpPtr_->value(secondaryInletT_);
323 const scalar secondaryOutletT =
324 secondaryInletT_ - Qt_/(secondaryMassFlowRate_*secondaryCp);
325
326 if (log)
327 {
328 Info << nl << incrIndent << indent
329 << "Secondary outlet T [K] : " << secondaryOutletT
330 << decrIndent;
331 }
332
333 os << tab << secondaryOutletT;
334 }
335 os << endl;
336 }
337
338 Info<< nl << endl;
339}
340
341
342// ************************************************************************* //
scalar delta
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition Function1.H:92
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
Abstract base-class for fluid and solid thermodynamic properties.
Definition basicThermo.H:62
static const word dictName
The dictionary name ("thermophysicalProperties").
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
virtual OFstream & file()
Return access to the file (if only 1).
Definition writeFile.C:270
virtual void writeCurrentTime(Ostream &os) const
Write the current time to stream.
Definition writeFile.C:354
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Base class for heat exchanger models to handle various characteristics for the heatExchangerSource fv...
virtual void initialise()
Initialise data members of the model.
const word & name_
Reference to the name of the fvOption source.
const fvMesh & mesh_
Reference to the mesh.
word UName_
Name of operand velocity field.
heatExchangerModel(const fvMesh &mesh, const word &name, const dictionary &coeffs)
Construct from components.
word phiName_
Name of operand flux field.
const dictionary & coeffs_
Dictionary containing coefficients specific to the chosen model.
word faceZoneName_
Name of the faceZone at the heat exchanger inlet.
labelList faceId_
Local list of face IDs.
word TName_
Name of operand temperature field.
labelList facePatchId_
Local list of patch IDs per face.
labelList faceSign_
List of +1/-1 representing face flip map (1 use as is, -1 negate).
virtual const word & U() const
Return const reference to the name of velocity field.
A heat exchanger model where heat exchange is calculated via an effectiveness table.
virtual void initialise()
Initialise data members of the model.
virtual bool read(const dictionary &dict)
Read top-level dictionary.
effectivenessTable(const fvMesh &mesh, const word &name, const dictionary &coeffs)
Construct from components.
virtual void write(const bool log)
Write data to stream and files.
virtual tmp< scalarField > energyDensity(const labelList &cells)
Return energy density per unit length [J/m3/m].
2D table interpolation. The data must be in ascending order in both dimensions x and y.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
const volScalarField & T
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
const cellShapeList & cells
const expr V(m.psi().mesh().V())
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.
A namespace for various heat exchanger model implementations.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< label > labelList
A List of labels.
Definition List.H:62
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition Ostream.H:490
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Ostream & indent(Ostream &os)
Indent stream.
Definition Ostream.H:481
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
Type gMin(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition Ostream.H:499
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Type gMax(const FieldField< Field, Type > &f)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
constexpr char tab
The tab '\t' character(0x09).
Definition Ostream.H:49
volScalarField & alpha
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299