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-2018 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  << nl << 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 
135 Foam::fv::effectivenessHeatExchangerSource::effectivenessHeatExchangerSource
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 field from which the temperature
164  // is obtained
165 
166  const basicThermo& thermo =
167  mesh_.lookupObject<basicThermo>(basicThermo::dictName);
168 
169  fieldNames_.setSize(1, thermo.he().name());
170 
171  applied_.setSize(1, false);
172 
173  eTable_.reset(new interpolation2DTable<scalar>(coeffs_));
174 
175  initialise();
176 }
177 
178 
179 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
180 
182 (
183  const volScalarField& rho,
184  fvMatrix<scalar>& eqn,
185  const label
186 )
187 {
188  const basicThermo& thermo =
189  mesh_.lookupObject<basicThermo>(basicThermo::dictName);
190 
191  const surfaceScalarField Cpf(fvc::interpolate(thermo.Cp()));
192 
193  const surfaceScalarField& phi =
194  mesh_.lookupObject<surfaceScalarField>(phiName_);
195 
196  const volScalarField& T = mesh_.lookupObject<volScalarField>(TName_);
198 
199  scalar sumPhi = 0;
200  scalar sumMagPhi = 0;
201  scalar CpfMean = 0;
202  scalar primaryInletTfMean = 0;
203  forAll(faceId_, i)
204  {
205  label facei = faceId_[i];
206  if (facePatchId_[i] != -1)
207  {
208  label patchi = facePatchId_[i];
209  scalar phii = phi.boundaryField()[patchi][facei]*faceSign_[i];
210 
211  sumPhi += phii;
212 
213  scalar Cpfi = Cpf.boundaryField()[patchi][facei];
214  scalar Tfi = Tf.boundaryField()[patchi][facei];
215  scalar magPhii = mag(phii);
216 
217  sumMagPhi += magPhii;
218  CpfMean += Cpfi*magPhii;
219  primaryInletTfMean += Tfi*magPhii;
220  }
221  else
222  {
223  scalar phii = phi[facei]*faceSign_[i];
224  scalar magPhii = mag(phii);
225 
226  sumPhi += phii;
227  sumMagPhi += magPhii;
228  CpfMean += Cpf[facei]*magPhii;
229  primaryInletTfMean += Tf[facei]*magPhii;
230  }
231  }
232  reduce(CpfMean, sumOp<scalar>());
233  reduce(sumPhi, sumOp<scalar>());
234  reduce(sumMagPhi, sumOp<scalar>());
235  CpfMean /= sumMagPhi + ROOTVSMALL;
236 
237  scalar primaryInletT = primaryInletT_;
238  if (!userPrimaryInletT_)
239  {
240  reduce(primaryInletTfMean, sumOp<scalar>());
241  primaryInletT = primaryInletTfMean/(sumMagPhi + ROOTVSMALL);
242  }
243 
244  const scalar alpha =
245  eTable_()(mag(sumPhi), secondaryMassFlowRate_)
246  *CpfMean*mag(sumPhi);
247 
248  const scalar Qt = alpha*(secondaryInletT_ - primaryInletT);
249 
250  if
251  (
252  targetQdotActive_
253  && (mesh_.time().timeIndex() % targetQdotCalcInterval_ == 0)
254  )
255  {
256  scalar dT = (targetQdot_ - Qt)/(alpha + ROOTVSMALL);
257  secondaryInletT_ += targetQdotRelax_*dT;
258  }
259 
260  const scalarField TCells(T, cells_);
261  scalar Tref = 0;
262  scalarField deltaTCells(cells_.size(), Zero);
263  if (Qt > 0)
264  {
265  Tref = gMax(TCells);
266  forAll(deltaTCells, i)
267  {
268  deltaTCells[i] = max(Tref - TCells[i], 0.0);
269  }
270  }
271  else
272  {
273  Tref = gMin(TCells);
274  forAll(deltaTCells, i)
275  {
276  deltaTCells[i] = max(TCells[i] - Tref, 0.0);
277  }
278  }
279 
280 
281  const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
282  const scalarField& V = mesh_.V();
283  scalar sumWeight = 0;
284  forAll(cells_, i)
285  {
286  label celli = cells_[i];
287  sumWeight += V[celli]*mag(U[celli])*deltaTCells[i];
288  }
289  reduce(sumWeight, sumOp<scalar>());
290 
291  if (this->V() > VSMALL && mag(Qt) > VSMALL)
292  {
293  scalarField& heSource = eqn.source();
294 
295  forAll(cells_, i)
296  {
297  label celli = cells_[i];
298  heSource[celli] -=
299  Qt*V[celli]*mag(U[celli])*deltaTCells[i]
300  /(sumWeight + ROOTVSMALL);
301  }
302  }
303 
304  Info<< type() << ": " << name() << nl << incrIndent
305  << indent << "Net mass flux [Kg/s] : " << sumPhi << nl
306  << indent << "Total heat exchange [W] : " << Qt << nl
307  << indent << "Secondary inlet T [K] : " << secondaryInletT_ << nl
308  << indent << "Tref [K] : " << Tref << nl
309  << indent << "Effectiveness : "
310  << eTable_()(mag(sumPhi), secondaryMassFlowRate_) << decrIndent
311  << nl << endl;
312 }
313 
314 
316 {
318  {
319  UName_ = coeffs_.lookupOrDefault<word>("U", "U");
320  TName_ = coeffs_.lookupOrDefault<word>("T", "T");
321  phiName_ = coeffs_.lookupOrDefault<word>("phi", "phi");
322  coeffs_.readEntry("faceZone", faceZoneName_);
323 
324  coeffs_.readEntry("secondaryMassFlowRate", secondaryMassFlowRate_);
325  coeffs_.readEntry("secondaryInletT", secondaryInletT_);
326 
327  if (coeffs_.readIfPresent("primaryInletT", primaryInletT_))
328  {
329  userPrimaryInletT_ = true;
330  Info<< type() << " " << this->name() << ": " << indent << nl
331  << "employing user-specified primary flow inlet temperature: "
332  << primaryInletT_ << endl;
333  }
334  else
335  {
336  Info<< type() << " " << this->name() << ": " << indent << nl
337  << "employing flux-weighted primary flow inlet temperature"
338  << endl;
339  }
340 
341  if (coeffs_.readIfPresent("targetQdot", targetQdot_))
342  {
343  targetQdotActive_ = true;
344  Info<< indent << "employing target heat rejection of "
345  << targetQdot_ << nl;
346 
347  coeffs_.readIfPresent
348  (
349  "targetQdotCalcInterval",
350  targetQdotCalcInterval_
351  );
352 
353  Info<< indent << "updating secondary inlet temperature every "
354  << targetQdotCalcInterval_ << " iterations" << nl;
355 
356  coeffs_.readIfPresent("targetQdotRelax", targetQdotRelax_);
357 
358  Info<< indent << "temperature relaxation: "
359  << targetQdotRelax_ << endl;
360 
361  }
362 
363  return true;
364  }
365 
366  return false;
367 }
368 
369 
370 // ************************************************************************* //
basicThermo.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fv::cellSetOption
Cell-set options abtract base class. Provides a base set of controls, e.g.:
Definition: cellSetOption.H:72
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
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:315
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)
Scalar.
Definition: effectivenessHeatExchangerSource.H:325
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:435
fvMatrix.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::incrIndent
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:314
rho
rho
Definition: readInitialConditions.H:96
Foam::fv::option::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:82
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:477
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:268
coupledPolyPatch.H
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
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:274
Foam::polyBoundaryMesh::whichPatch
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Definition: polyBoundaryMesh.C:805
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:265
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:84
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::decrIndent
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:321
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 given a name, return -1 if not found.
Definition: ZoneMesh.C:484
Foam::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:307
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)
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
effectivenessHeatExchangerSource.H
Foam::BitOps::count
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:74
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:295
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:271
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::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592