effectivenessHeatExchangerSource.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) 2013-2015 OpenFOAM Foundation
9  Copyright (C) 2016-2020 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 "fvMesh.H"
31 #include "fvMatrix.H"
33 #include "basicThermo.H"
34 #include "coupledPolyPatch.H"
35 #include "surfaceInterpolate.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 namespace fv
42 {
43  defineTypeNameAndDebug(effectivenessHeatExchangerSource, 0);
45  (
46  option,
47  effectivenessHeatExchangerSource,
48  dictionary
49  );
50 }
51 }
52 
53 
54 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55 
56 void Foam::fv::effectivenessHeatExchangerSource::initialise()
57 {
59 
60  if (zoneID < 0)
61  {
63  << type() << " " << this->name() << ": "
64  << " Unknown face zone name: " << faceZoneName_
65  << ". Valid face zones are: " << mesh_.faceZones().names()
66  << exit(FatalError);
67  }
68 
69  const faceZone& fZone = mesh_.faceZones()[zoneID];
70 
71  faceId_.setSize(fZone.size());
72  facePatchId_.setSize(fZone.size());
73  faceSign_.setSize(fZone.size());
74 
75  label count = 0;
76  forAll(fZone, i)
77  {
78  label facei = fZone[i];
79  label faceId = -1;
80  label facePatchId = -1;
81  if (mesh_.isInternalFace(facei))
82  {
83  faceId = facei;
84  facePatchId = -1;
85  }
86  else
87  {
88  facePatchId = mesh_.boundaryMesh().whichPatch(facei);
89  const polyPatch& pp = mesh_.boundaryMesh()[facePatchId];
90  if (isA<coupledPolyPatch>(pp))
91  {
92  if (refCast<const coupledPolyPatch>(pp).owner())
93  {
94  faceId = pp.whichFace(facei);
95  }
96  else
97  {
98  faceId = -1;
99  }
100  }
101  else if (!isA<emptyPolyPatch>(pp))
102  {
103  faceId = facei - pp.start();
104  }
105  else
106  {
107  faceId = -1;
108  facePatchId = -1;
109  }
110  }
111 
112  if (faceId >= 0)
113  {
114  if (fZone.flipMap()[i])
115  {
116  faceSign_[count] = -1;
117  }
118  else
119  {
120  faceSign_[count] = 1;
121  }
122  faceId_[count] = faceId;
123  facePatchId_[count] = facePatchId;
124  count++;
125  }
126  }
130 }
131 
132 
133 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
134 
136 (
137  const word& name,
138  const word& modelType,
139  const dictionary& dict,
140  const fvMesh& mesh
141 )
142 :
143  cellSetOption(name, modelType, dict, mesh),
144  secondaryMassFlowRate_(0),
145  secondaryInletT_(0),
146  primaryInletT_(0),
147  userPrimaryInletT_(false),
148  targetQdotActive_(false),
149  targetQdot_(0),
150  targetQdotCalcInterval_(5),
151  targetQdotRelax_(0.5),
152  eTable_(),
153  UName_("U"),
154  TName_("T"),
155  phiName_("phi"),
156  faceZoneName_("unknown-faceZone"),
157  faceId_(),
158  facePatchId_(),
159  faceSign_()
160 {
161  read(dict);
162 
163  // Set the field name to that of the energy
164  // field from which the temperature is obtained
165 
166  const auto& thermo = mesh_.lookupObject<basicThermo>(basicThermo::dictName);
167 
168  fieldNames_.setSize(1, thermo.he().name());
169 
170  applied_.setSize(1, false);
171 
172  eTable_.reset(new interpolation2DTable<scalar>(coeffs_));
173 
174  initialise();
175 }
176 
177 
178 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
179 
181 (
182  const volScalarField& rho,
183  fvMatrix<scalar>& eqn,
184  const label
185 )
186 {
187  const auto& thermo = mesh_.lookupObject<basicThermo>(basicThermo::dictName);
188 
189  const surfaceScalarField Cpf(fvc::interpolate(thermo.Cp()));
190 
191  const auto& phi = mesh_.lookupObject<surfaceScalarField>(phiName_);
192 
193  const auto& T = mesh_.lookupObject<volScalarField>(TName_);
195 
196  scalar sumPhi = 0;
197  scalar sumMagPhi = 0;
198  scalar CpfMean = 0;
199  scalar primaryInletTfMean = 0;
200  forAll(faceId_, i)
201  {
202  label facei = faceId_[i];
203  if (facePatchId_[i] != -1)
204  {
205  label patchi = facePatchId_[i];
206  scalar phii = phi.boundaryField()[patchi][facei]*faceSign_[i];
207 
208  sumPhi += phii;
209 
210  scalar Cpfi = Cpf.boundaryField()[patchi][facei];
211  scalar Tfi = Tf.boundaryField()[patchi][facei];
212  scalar magPhii = mag(phii);
213 
214  sumMagPhi += magPhii;
215  CpfMean += Cpfi*magPhii;
216  primaryInletTfMean += Tfi*magPhii;
217  }
218  else
219  {
220  scalar phii = phi[facei]*faceSign_[i];
221  scalar magPhii = mag(phii);
222 
223  sumPhi += phii;
224  sumMagPhi += magPhii;
225  CpfMean += Cpf[facei]*magPhii;
226  primaryInletTfMean += Tf[facei]*magPhii;
227  }
228  }
229  reduce(CpfMean, sumOp<scalar>());
230  reduce(sumPhi, sumOp<scalar>());
231  reduce(sumMagPhi, sumOp<scalar>());
232  CpfMean /= sumMagPhi + ROOTVSMALL;
233 
234  scalar primaryInletT = primaryInletT_;
235  if (!userPrimaryInletT_)
236  {
237  reduce(primaryInletTfMean, sumOp<scalar>());
238  primaryInletT = primaryInletTfMean/(sumMagPhi + ROOTVSMALL);
239  }
240 
241  const scalar alpha =
242  eTable_()(mag(sumPhi), secondaryMassFlowRate_)
243  *CpfMean*mag(sumPhi);
244 
245  const scalar Qt = alpha*(secondaryInletT_ - primaryInletT);
246 
247  if
248  (
249  targetQdotActive_
250  && (mesh_.time().timeIndex() % targetQdotCalcInterval_ == 0)
251  )
252  {
253  scalar dT = (targetQdot_ - Qt)/(alpha + ROOTVSMALL);
254  secondaryInletT_ += targetQdotRelax_*dT;
255  }
256 
257  const scalarField TCells(T, cells_);
258  scalar Tref = 0;
259  scalarField deltaTCells(cells_.size(), Zero);
260  if (Qt > 0)
261  {
262  Tref = gMax(TCells);
263  forAll(deltaTCells, i)
264  {
265  deltaTCells[i] = max(Tref - TCells[i], 0.0);
266  }
267  }
268  else
269  {
270  Tref = gMin(TCells);
271  forAll(deltaTCells, i)
272  {
273  deltaTCells[i] = max(TCells[i] - Tref, 0.0);
274  }
275  }
276 
277 
278  const auto& U = mesh_.lookupObject<volVectorField>(UName_);
279  const scalarField& V = mesh_.V();
280  scalar sumWeight = 0;
281  forAll(cells_, i)
282  {
283  label celli = cells_[i];
284  sumWeight += V[celli]*mag(U[celli])*deltaTCells[i];
285  }
286  reduce(sumWeight, sumOp<scalar>());
287 
288  if (this->V() > VSMALL && mag(Qt) > VSMALL)
289  {
290  scalarField& heSource = eqn.source();
291 
292  forAll(cells_, i)
293  {
294  label celli = cells_[i];
295  heSource[celli] -=
296  Qt*V[celli]*mag(U[celli])*deltaTCells[i]
297  /(sumWeight + ROOTVSMALL);
298  }
299  }
300 
301  Info<< type() << ": " << name() << nl << incrIndent
302  << indent << "Net mass flux [Kg/s] : " << sumPhi << nl
303  << indent << "Total heat exchange [W] : " << Qt << nl
304  << indent << "Secondary inlet T [K] : " << secondaryInletT_ << nl
305  << indent << "Tref [K] : " << Tref << nl
306  << indent << "Effectiveness : "
307  << eTable_()(mag(sumPhi), secondaryMassFlowRate_) << decrIndent
308  << nl << endl;
309 }
310 
311 
313 {
315  {
316  UName_ = coeffs_.getOrDefault<word>("U", "U");
317  TName_ = coeffs_.getOrDefault<word>("T", "T");
318  phiName_ = coeffs_.getOrDefault<word>("phi", "phi");
319  coeffs_.readEntry("faceZone", faceZoneName_);
320 
321  coeffs_.readEntry("secondaryMassFlowRate", secondaryMassFlowRate_);
322  coeffs_.readEntry("secondaryInletT", secondaryInletT_);
323 
324  if (coeffs_.readIfPresent("primaryInletT", primaryInletT_))
325  {
326  userPrimaryInletT_ = true;
327  Info<< type() << " " << this->name() << ": " << indent << nl
328  << "employing user-specified primary flow inlet temperature: "
329  << primaryInletT_ << endl;
330  }
331  else
332  {
333  Info<< type() << " " << this->name() << ": " << indent << nl
334  << "employing flux-weighted primary flow inlet temperature"
335  << endl;
336  }
337 
338  if (coeffs_.readIfPresent("targetQdot", targetQdot_))
339  {
340  targetQdotActive_ = true;
341  Info<< indent << "employing target heat rejection of "
342  << targetQdot_ << nl;
343 
344  coeffs_.readIfPresent
345  (
346  "targetQdotCalcInterval",
347  targetQdotCalcInterval_
348  );
349 
350  Info<< indent << "updating secondary inlet temperature every "
351  << targetQdotCalcInterval_ << " iterations" << nl;
352 
353  coeffs_.readIfPresent("targetQdotRelax", targetQdotRelax_);
354 
355  Info<< indent << "temperature relaxation: "
356  << targetQdotRelax_ << endl;
357  }
358 
359  return true;
360  }
361 
362  return false;
363 }
364 
365 
366 // ************************************************************************* //
basicThermo.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fv::cellSetOption
Intermediate abstract class for handling cell-set options for the derived fvOptions.
Definition: cellSetOption.H:163
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::fv::option::name
const word & name() const
Return const access to the source name.
Definition: fvOptionI.H:30
Foam::fv::effectivenessHeatExchangerSource::read
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: effectivenessHeatExchangerSource.C:312
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::fv::effectivenessHeatExchangerSource::addSup
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
Add explicit/implicit contribution to momentum equation.
Definition: effectivenessHeatExchangerSource.H:464
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::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
fvMatrix.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::incrIndent
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:327
rho
rho
Definition: readInitialConditions.H:88
Foam::fv::option::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:139
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:486
Foam::interpolation2DTable< scalar >
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::fv::effectivenessHeatExchangerSource::faceId_
labelList faceId_
Local list of face IDs.
Definition: effectivenessHeatExchangerSource.H:409
coupledPolyPatch.H
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
faceId
label faceId(-1)
Foam::fv::cellSetOption::read
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: cellSetOption.C:255
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::fv::effectivenessHeatExchangerSource::faceSign_
labelList faceSign_
List of +1/-1 representing face flip map (1 use as is, -1 negate)
Definition: effectivenessHeatExchangerSource.H:415
Foam::polyBoundaryMesh::whichPatch
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Definition: polyBoundaryMesh.C:810
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
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::effectivenessHeatExchangerSource::faceZoneName_
word faceZoneName_
Name of the faceZone at the heat exchange inlet.
Definition: effectivenessHeatExchangerSource.H:406
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.
zoneID
const labelIOList & zoneID
Definition: interpolatedFaces.H:22
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::decrIndent
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:334
Foam::fvc::interpolate
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.
Foam::ZoneMesh::findZoneID
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition: ZoneMesh.C:484
Foam::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:320
Foam::ZoneMesh::names
wordList names() const
A list of the zone names.
Definition: ZoneMesh.C:340
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
fv
labelList fv(nPoints)
U
U
Definition: pEqn.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::nl
constexpr char nl
Definition: Ostream.H:385
effectivenessHeatExchangerSource.H
Foam::BitOps::count
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:77
Foam::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
Definition: primitiveMeshI.H:102
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::fvMatrix::source
Field< Type > & source()
Definition: fvMatrix.H:297
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::fv::effectivenessHeatExchangerSource::facePatchId_
labelList facePatchId_
Local list of patch ID per face.
Definition: effectivenessHeatExchangerSource.H:412
Foam::fv::defineTypeNameAndDebug
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
Foam::fv::addToRunTimeSelectionTable
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::fv::effectivenessHeatExchangerSource::effectivenessHeatExchangerSource
effectivenessHeatExchangerSource(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
Definition: effectivenessHeatExchangerSource.C:136
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592