effectivenessHeatExchangerSource.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-2023 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 {
40  (
41  fvModel,
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 void Foam::fv::effectivenessHeatExchangerSource::readCoeffs()
52 {
53  secondaryMassFlowRate_ = coeffs().lookup<scalar>("secondaryMassFlowRate");
54  secondaryInletT_ = coeffs().lookup<scalar>("secondaryInletT");
55  primaryInletT_ = coeffs().lookup<scalar>("primaryInletT");
56 
57  eTable_.reset(Function2<scalar>::New("effectiveness", coeffs()).ptr());
58 
59  UName_ = coeffs().lookupOrDefault<word>("U", "U");
60  TName_ = coeffs().lookupOrDefault<word>("T", "T");
61  phiName_ = coeffs().lookupOrDefault<word>("phi", "phi");
62 
63  faceZoneName_ = coeffs().lookup<word>("faceZone");
64 }
65 
66 
67 void Foam::fv::effectivenessHeatExchangerSource::setZone()
68 {
69  zoneID_ = mesh().faceZones().findZoneID(faceZoneName_);
70  if (zoneID_ < 0)
71  {
73  << type() << " " << this->name() << ": "
74  << " Unknown face zone name: " << faceZoneName_
75  << ". Valid face zones are: " << mesh().faceZones().names()
76  << nl << exit(FatalError);
77  }
78 
79  const faceZone& fZone = mesh().faceZones()[zoneID_];
80 
81  faceId_.setSize(fZone.size());
82  facePatchId_.setSize(fZone.size());
83  faceSign_.setSize(fZone.size());
84 
85  label count = 0;
86  forAll(fZone, i)
87  {
88  const label facei = fZone[i];
89  label faceId = -1;
90  label facePatchId = -1;
91  if (mesh().isInternalFace(facei))
92  {
93  faceId = facei;
94  facePatchId = -1;
95  }
96  else
97  {
98  facePatchId = mesh().boundaryMesh().whichPatch(facei);
99  const polyPatch& pp = mesh().boundaryMesh()[facePatchId];
100  if (isA<coupledPolyPatch>(pp))
101  {
102  if (refCast<const coupledPolyPatch>(pp).owner())
103  {
104  faceId = pp.whichFace(facei);
105  }
106  else
107  {
108  faceId = -1;
109  }
110  }
111  else if (!isA<emptyPolyPatch>(pp))
112  {
113  faceId = facei - pp.start();
114  }
115  else
116  {
117  faceId = -1;
118  facePatchId = -1;
119  }
120  }
121 
122  if (faceId >= 0)
123  {
124  if (fZone.flipMap()[i])
125  {
126  faceSign_[count] = -1;
127  }
128  else
129  {
130  faceSign_[count] = 1;
131  }
132  faceId_[count] = faceId;
133  facePatchId_[count] = facePatchId;
134  count++;
135  }
136  }
137  faceId_.setSize(count);
138  facePatchId_.setSize(count);
139  faceSign_.setSize(count);
140 
141  calculateTotalArea(faceZoneArea_);
142 }
143 
144 
145 void Foam::fv::effectivenessHeatExchangerSource::calculateTotalArea
146 (
147  scalar& area
148 ) const
149 {
150  area = 0;
151  forAll(faceId_, i)
152  {
153  const label facei = faceId_[i];
154  if (facePatchId_[i] != -1)
155  {
156  label patchi = facePatchId_[i];
157  area += mesh().magSf().boundaryField()[patchi][facei];
158  }
159  else
160  {
161  area += mesh().magSf()[facei];
162  }
163  }
164  reduce(area, sumOp<scalar>());
165 }
166 
167 
168 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
169 
171 (
172  const word& name,
173  const word& modelType,
174  const fvMesh& mesh,
175  const dictionary& dict
176 )
177 :
178  fvModel(name, modelType, mesh, dict),
179  set_(mesh, coeffs()),
180  secondaryMassFlowRate_(NaN),
181  secondaryInletT_(NaN),
182  primaryInletT_(NaN),
183  eTable_(nullptr),
184  UName_(word::null),
185  TName_(word::null),
186  phiName_(word::null),
187  faceZoneName_(word::null),
188  zoneID_(-1),
189  faceId_(),
190  facePatchId_(),
191  faceSign_(),
192  faceZoneArea_(NaN)
193 {
194  readCoeffs();
195  setZone();
196 }
197 
198 
199 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
200 
202 {
203  const basicThermo& thermo =
204  mesh().lookupObject<basicThermo>(physicalProperties::typeName);
205 
206  return wordList(1, thermo.he().name());
207 }
208 
209 
211 (
212  const volScalarField& rho,
213  fvMatrix<scalar>& eqn,
214  const word& fieldName
215 ) const
216 {
217  const basicThermo& thermo =
218  mesh().lookupObject<basicThermo>(physicalProperties::typeName);
219 
220  const surfaceScalarField Cpf(fvc::interpolate(thermo.Cp()));
221 
222  const surfaceScalarField& phi =
223  mesh().lookupObject<surfaceScalarField>(phiName_);
224 
225  scalar totalphi = 0;
226  scalar CpfMean = 0;
227  forAll(faceId_, i)
228  {
229  label facei = faceId_[i];
230  if (facePatchId_[i] != -1)
231  {
232  label patchi = facePatchId_[i];
233  totalphi += phi.boundaryField()[patchi][facei]*faceSign_[i];
234 
235  CpfMean +=
236  Cpf.boundaryField()[patchi][facei]
237  *mesh().magSf().boundaryField()[patchi][facei];
238  }
239  else
240  {
241  totalphi += phi[facei]*faceSign_[i];
242  CpfMean += Cpf[facei]*mesh().magSf()[facei];
243  }
244  }
245  reduce(CpfMean, sumOp<scalar>());
246  reduce(totalphi, sumOp<scalar>());
247 
248  const scalar Qt =
249  eTable_->value(mag(totalphi), secondaryMassFlowRate_)
250  *(secondaryInletT_ - primaryInletT_)
251  *(CpfMean/faceZoneArea_)*mag(totalphi);
252 
253  const labelUList cells = set_.cells();
254 
255  const volScalarField& T = mesh().lookupObject<volScalarField>(TName_);
256  const scalarField TCells(T, cells);
257  scalar Tref = 0;
258  if (Qt > 0)
259  {
260  Tref = max(TCells);
261  reduce(Tref, maxOp<scalar>());
262  }
263  else
264  {
265  Tref = min(TCells);
266  reduce(Tref, minOp<scalar>());
267  }
268 
269  scalarField deltaTCells(cells.size(), 0);
270  forAll(deltaTCells, i)
271  {
272  if (Qt > 0)
273  {
274  deltaTCells[i] = max(Tref - TCells[i], 0.0);
275  }
276  else
277  {
278  deltaTCells[i] = max(TCells[i] - Tref, 0.0);
279  }
280  }
281 
282  const volVectorField& U = mesh().lookupObject<volVectorField>(UName_);
283  const scalarField& V = mesh().V();
284  scalar sumWeight = 0;
285 
286  forAll(cells, i)
287  {
288  sumWeight += V[cells[i]]*mag(U[cells[i]])*deltaTCells[i];
289  }
290  reduce(sumWeight, sumOp<scalar>());
291 
292  if (mag(Qt) > vSmall)
293  {
294  scalarField& heSource = eqn.source();
295 
296  forAll(cells, i)
297  {
298  heSource[cells[i]] -=
299  Qt*V[cells[i]]*mag(U[cells[i]])*deltaTCells[i]/sumWeight;
300  }
301  }
302 
303  if (debug && Pstream::master())
304  {
305  Info<< indent << "Net mass flux [Kg/s] = " << totalphi << nl;
306  Info<< indent << "Total energy exchange [W] = " << Qt << nl;
307  Info<< indent << "Tref [K] = " << Tref << nl;
308  Info<< indent << "Effectiveness : "
309  << eTable_->value(mag(totalphi), secondaryMassFlowRate_) << endl;
310  }
311 }
312 
313 
315 {
316  set_.movePoints();
317  return true;
318 }
319 
320 
322 (
323  const polyTopoChangeMap& map
324 )
325 {
326  set_.topoChange(map);
327 }
328 
329 
331 {
332  set_.mapMesh(map);
333 }
334 
335 
337 (
338  const polyDistributionMap& map
339 )
340 {
341  set_.distribute(map);
342 }
343 
344 
346 {
347  if (fvModel::read(dict))
348  {
349  set_.read(coeffs());
350  readCoeffs();
351  setZone();
352  return true;
353  }
354  else
355  {
356  return false;
357  }
358 }
359 
360 
361 // ************************************************************************* //
#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 dictionary &dict)
Selector.
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:160
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
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:101
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:40
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:187
Heat exchanger source model, in which the heat exchanger is defined as a selection of cells.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
virtual void addSup(const volScalarField &rho, fvMatrix< scalar > &eqn, const word &fieldName) const
Explicit and implicit source for compressible equation.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read dictionary.
effectivenessHeatExchangerSource(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
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:306
label patchi
const cellShapeList & cells
label faceId(-1)
U
Definition: pEqn.H:72
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
static Type NaN()
Return a primitive with all components set to NaN.
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:251
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
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)
error FatalError
label count(const ListType &l, typename ListType::const_reference x)
Count the number of occurrences of a value in a list.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
static const char nl
Definition: Ostream.H:260
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
labelList fv(nPoints)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31