30#include "surfaceInterpolate.H"
52void Foam::heatExchangerModels::effectivenessTable::writeFileHeader
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 [-]");
67 writeFile::writeTabbed(
os,
"Secondary outlet T [K]");
84 userPrimaryInletT_(false),
85 targetQdotActive_(false),
97 targetQdotCalcInterval_(5),
98 secondaryMassFlowRate_(0),
102 targetQdotRelax_(0.5),
108 writeFileHeader(
file());
136 scalar sumMagPhi = 0;
138 scalar primaryInletTfMean = 0;
141 const label facei = faceId_[i];
142 if (facePatchId_[i] != -1)
144 const label patchi = facePatchId_[i];
145 const scalar phii =
phi.boundaryField()[patchi][facei]*faceSign_[i];
146 const scalar magPhii =
mag(phii);
149 sumMagPhi += magPhii;
151 const scalar Cpfi = Cpf.boundaryField()[patchi][facei];
152 const scalar Tfi = Tf.boundaryField()[patchi][facei];
154 CpfMean += Cpfi*magPhii;
155 primaryInletTfMean += Tfi*magPhii;
159 const scalar phii =
phi[facei]*faceSign_[i];
160 const scalar magPhii =
mag(phii);
163 sumMagPhi += magPhii;
165 CpfMean += Cpf[facei]*magPhii;
166 primaryInletTfMean += Tf[facei]*magPhii;
169 reduce(CpfMean, sumOp<scalar>());
170 reduce(sumPhi_, sumOp<scalar>());
171 reduce(sumMagPhi, sumOp<scalar>());
173 CpfMean /= sumMagPhi + ROOTVSMALL;
175 scalar primaryInletT = primaryInletT_;
176 if (!userPrimaryInletT_)
178 reduce(primaryInletTfMean, sumOp<scalar>());
179 primaryInletT = primaryInletTfMean/(sumMagPhi + ROOTVSMALL);
182 effectiveness_ = eTable_()(
mag(sumPhi_), secondaryMassFlowRate_);
184 const scalar
alpha = effectiveness_*CpfMean*
mag(sumPhi_);
186 Qt_ =
alpha*(secondaryInletT_ - primaryInletT);
191 && (mesh_.time().timeIndex() % targetQdotCalcInterval_ == 0)
194 const scalar dT = (targetQdot_ - Qt_)/(
alpha + ROOTVSMALL);
195 secondaryInletT_ += targetQdotRelax_*dT;
203 Tref_ =
gMax(deltaTCells);
204 for (scalar&
delta : deltaTCells)
211 Tref_ =
gMin(deltaTCells);
212 for (scalar&
delta : deltaTCells)
220 scalar sumWeight = 0;
223 const label celli =
cells[i];
224 sumWeight +=
V[celli]*
mag(
U[celli])*deltaTCells[i];
228 return Qt_*deltaTCells/(sumWeight + ROOTVSMALL);
234 if (!writeFile::read(
dict))
241 coeffs_.readEntry(
"secondaryMassFlowRate", secondaryMassFlowRate_);
242 coeffs_.readEntry(
"secondaryInletT", secondaryInletT_);
244 if (coeffs_.readIfPresent(
"primaryInletT", primaryInletT_))
246 userPrimaryInletT_ =
true;
249 <<
"- using user-specified primary flow inlet temperature: "
250 << primaryInletT_ <<
endl;
255 <<
"- using flux-weighted primary flow inlet temperature"
259 if (coeffs_.readIfPresent(
"targetQdot", targetQdot_))
261 targetQdotActive_ =
true;
263 Info<<
indent <<
"- using target heat rejection: " << targetQdot_ <<
nl;
265 coeffs_.readIfPresent
267 "targetQdotCalcInterval",
268 targetQdotCalcInterval_
271 Info<<
indent <<
"- updating secondary inlet temperature every "
272 << targetQdotCalcInterval_ <<
" iterations" <<
nl;
274 coeffs_.readIfPresent(
"targetQdotRelax", targetQdotRelax_);
277 << targetQdotRelax_ <<
endl;
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_);
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_
311 writeCurrentTime(
os);
315 <<
tab << secondaryInletT_
317 <<
tab << effectiveness_;
322 const scalar secondaryCp = secondaryCpPtr_->value(secondaryInletT_);
323 const scalar secondaryOutletT =
324 secondaryInletT_ - Qt_/(secondaryMassFlowRate_*secondaryCp);
329 <<
"Secondary outlet T [K] : " << secondaryOutletT
333 os <<
tab << secondaryOutletT;
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...
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,...
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Abstract base-class for fluid and solid thermodynamic properties.
static const word dictName
The dictionary name ("thermophysicalProperties").
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
virtual OFstream & file()
Return access to the file (if only 1).
virtual void writeCurrentTime(Ostream &os) const
Write the current time to stream.
Mesh data needed to do the Finite Volume discretisation.
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.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
OBJstream os(runTime.globalPath()/outputName)
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.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< label > labelList
A List of labels.
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.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Ostream & indent(Ostream &os)
Indent stream.
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.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
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).
constexpr char tab
The tab '\t' character(0x09).
#define forAll(list, i)
Loop across all elements in list.