effectivenessHeatExchanger.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2024 OpenFOAM Foundation
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 
27 #include "fvMatrix.H"
28 #include "basicThermo.H"
29 #include "surfaceInterpolate.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace fv
37 {
41  (
42  fvModel,
44  dictionary,
45  effectivenessHeatExchangerSource,
46  "effectivenessHeatExchangerSource"
47  );
48 }
49 }
50 
51 
52 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53 
54 void Foam::fv::effectivenessHeatExchanger::readCoeffs()
55 {
56  secondaryMassFlowRate_ =
57  coeffs().lookup<scalar>("secondaryMassFlowRate", dimMass/dimTime);
58  secondaryInletT_ =
59  coeffs().lookup<scalar>("secondaryInletT", dimTemperature);
60  primaryInletT_ =
61  coeffs().lookup<scalar>("primaryInletT", dimTemperature);
62 
63  eTable_.reset
64  (
66  (
67  "effectiveness",
70  dimless,
71  coeffs()
72  ).ptr()
73  );
74 
75  UName_ = coeffs().lookupOrDefault<word>("U", "U");
76  TName_ = coeffs().lookupOrDefault<word>("T", "T");
77  phiName_ = coeffs().lookupOrDefault<word>("phi", "phi");
78 
79  faceZoneName_ = coeffs().lookup<word>("faceZone");
80 }
81 
82 
83 void Foam::fv::effectivenessHeatExchanger::setZone()
84 {
85  zoneIndex_ = mesh().faceZones().findIndex(faceZoneName_);
86  if (zoneIndex_ < 0)
87  {
89  << type() << " " << this->name() << ": "
90  << " Unknown face zone name: " << faceZoneName_
91  << ". Valid face zones are: " << mesh().faceZones().toc()
92  << nl << exit(FatalError);
93  }
94 
95  const faceZone& fZone = mesh().faceZones()[zoneIndex_];
96 
97  faceId_.setSize(fZone.size());
98  facePatchId_.setSize(fZone.size());
99  faceSign_.setSize(fZone.size());
100 
101  label count = 0;
102  forAll(fZone, i)
103  {
104  const label facei = fZone[i];
105  label faceId = -1;
106  label facePatchId = -1;
107  if (mesh().isInternalFace(facei))
108  {
109  faceId = facei;
110  facePatchId = -1;
111  }
112  else
113  {
114  facePatchId = mesh().boundaryMesh().whichPatch(facei);
115  const polyPatch& pp = mesh().boundaryMesh()[facePatchId];
116  if (isA<coupledPolyPatch>(pp))
117  {
118  if (refCast<const coupledPolyPatch>(pp).owner())
119  {
120  faceId = pp.whichFace(facei);
121  }
122  else
123  {
124  faceId = -1;
125  }
126  }
127  else if (!isA<emptyPolyPatch>(pp))
128  {
129  faceId = facei - pp.start();
130  }
131  else
132  {
133  faceId = -1;
134  facePatchId = -1;
135  }
136  }
137 
138  if (faceId >= 0)
139  {
140  if (fZone.flipMap()[i])
141  {
142  faceSign_[count] = -1;
143  }
144  else
145  {
146  faceSign_[count] = 1;
147  }
148  faceId_[count] = faceId;
149  facePatchId_[count] = facePatchId;
150  count++;
151  }
152  }
153  faceId_.setSize(count);
154  facePatchId_.setSize(count);
155  faceSign_.setSize(count);
156 
157  calculateTotalArea(faceZoneArea_);
158 }
159 
160 
161 void Foam::fv::effectivenessHeatExchanger::calculateTotalArea
162 (
163  scalar& area
164 ) const
165 {
166  area = 0;
167  forAll(faceId_, i)
168  {
169  const label facei = faceId_[i];
170  if (facePatchId_[i] != -1)
171  {
172  label patchi = facePatchId_[i];
173  area += mesh().magSf().boundaryField()[patchi][facei];
174  }
175  else
176  {
177  area += mesh().magSf()[facei];
178  }
179  }
180  reduce(area, sumOp<scalar>());
181 }
182 
183 
184 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
185 
187 (
188  const word& name,
189  const word& modelType,
190  const fvMesh& mesh,
191  const dictionary& dict
192 )
193 :
194  fvModel(name, modelType, mesh, dict),
195  set_(mesh, coeffs()),
196  secondaryMassFlowRate_(NaN),
197  secondaryInletT_(NaN),
198  primaryInletT_(NaN),
199  eTable_(nullptr),
200  UName_(word::null),
201  TName_(word::null),
202  phiName_(word::null),
203  faceZoneName_(word::null),
204  zoneIndex_(-1),
205  faceId_(),
206  facePatchId_(),
207  faceSign_(),
208  faceZoneArea_(NaN)
209 {
210  readCoeffs();
211  setZone();
212 }
213 
214 
215 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
216 
218 {
219  const basicThermo& thermo =
220  mesh().lookupObject<basicThermo>(physicalProperties::typeName);
221 
222  return wordList(1, thermo.he().name());
223 }
224 
225 
227 (
228  const volScalarField& rho,
229  const volScalarField& he,
230  fvMatrix<scalar>& eqn
231 ) const
232 {
233  const basicThermo& thermo =
234  mesh().lookupObject<basicThermo>(physicalProperties::typeName);
235 
236  const surfaceScalarField Cpf(fvc::interpolate(thermo.Cp()));
237 
238  const surfaceScalarField& phi =
239  mesh().lookupObject<surfaceScalarField>(phiName_);
240 
241  scalar totalphi = 0;
242  scalar CpfMean = 0;
243  forAll(faceId_, i)
244  {
245  label facei = faceId_[i];
246  if (facePatchId_[i] != -1)
247  {
248  label patchi = facePatchId_[i];
249  totalphi += phi.boundaryField()[patchi][facei]*faceSign_[i];
250 
251  CpfMean +=
252  Cpf.boundaryField()[patchi][facei]
253  *mesh().magSf().boundaryField()[patchi][facei];
254  }
255  else
256  {
257  totalphi += phi[facei]*faceSign_[i];
258  CpfMean += Cpf[facei]*mesh().magSf()[facei];
259  }
260  }
261  reduce(CpfMean, sumOp<scalar>());
262  reduce(totalphi, sumOp<scalar>());
263 
264  const scalar Qt =
265  eTable_->value(mag(totalphi), secondaryMassFlowRate_)
266  *(secondaryInletT_ - primaryInletT_)
267  *(CpfMean/faceZoneArea_)*mag(totalphi);
268 
269  const labelUList cells = set_.cells();
270 
271  const volScalarField& T = mesh().lookupObject<volScalarField>(TName_);
272  const scalarField TCells(T, cells);
273  scalar Tref = 0;
274  if (Qt > 0)
275  {
276  Tref = max(TCells);
277  reduce(Tref, maxOp<scalar>());
278  }
279  else
280  {
281  Tref = min(TCells);
282  reduce(Tref, minOp<scalar>());
283  }
284 
285  scalarField deltaTCells(cells.size(), 0);
286  forAll(deltaTCells, i)
287  {
288  if (Qt > 0)
289  {
290  deltaTCells[i] = max(Tref - TCells[i], 0.0);
291  }
292  else
293  {
294  deltaTCells[i] = max(TCells[i] - Tref, 0.0);
295  }
296  }
297 
298  const volVectorField& U = mesh().lookupObject<volVectorField>(UName_);
299  const scalarField& V = mesh().V();
300  scalar sumWeight = 0;
301 
302  forAll(cells, i)
303  {
304  sumWeight += V[cells[i]]*mag(U[cells[i]])*deltaTCells[i];
305  }
306  reduce(sumWeight, sumOp<scalar>());
307 
308  if (mag(Qt) > vSmall)
309  {
310  scalarField& heSource = eqn.source();
311 
312  forAll(cells, i)
313  {
314  heSource[cells[i]] -=
315  Qt*V[cells[i]]*mag(U[cells[i]])*deltaTCells[i]/sumWeight;
316  }
317  }
318 
319  if (debug && Pstream::master())
320  {
321  Info<< indent << "Net mass flux [Kg/s] = " << totalphi << nl;
322  Info<< indent << "Total energy exchange [W] = " << Qt << nl;
323  Info<< indent << "Tref [K] = " << Tref << nl;
324  Info<< indent << "Effectiveness : "
325  << eTable_->value(mag(totalphi), secondaryMassFlowRate_) << endl;
326  }
327 }
328 
329 
331 {
332  set_.movePoints();
333  return true;
334 }
335 
336 
338 (
339  const polyTopoChangeMap& map
340 )
341 {
342  set_.topoChange(map);
343 }
344 
345 
347 {
348  set_.mapMesh(map);
349 }
350 
351 
353 (
354  const polyDistributionMap& map
355 )
356 {
357  set_.distribute(map);
358 }
359 
360 
362 {
363  if (fvModel::read(dict))
364  {
365  set_.read(coeffs());
366  readCoeffs();
367  setZone();
368  return true;
369  }
370  else
371  {
372  return false;
373  }
374 }
375 
376 
377 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
static autoPtr< Function2< Type > > New(const word &name, const Function2s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function2New.C:32
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:78
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Field< Type > & source()
Definition: fvMatrix.H:307
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:59
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:199
Heat exchanger model, based on an effectiveness.
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
effectivenessHeatExchanger(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read dictionary.
virtual void addSup(const volScalarField &rho, const volScalarField &he, fvMatrix< scalar > &eqn) const
Explicit and implicit source for compressible equation.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
const cellShapeList & cells
label faceId(-1)
U
Definition: pEqn.H:72
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
addBackwardCompatibleToRunTimeSelectionTable(fvConstraint, fixedTemperature, dictionary, fixedTemperatureConstraint, "fixedTemperatureConstraint")
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
List< word > wordList
A List of words.
Definition: fileName.H:54
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
messageStream Info
const dimensionSet dimTemperature
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimMass
error FatalError
label count(const ListType &l, typename ListType::const_reference x)
Count the number of occurrences of a value in a list.
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:227
static const char nl
Definition: Ostream.H:266
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
thermo he()
labelList fv(nPoints)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31