31#include "surfaceInterpolate.H"
54Foam::heatExchangerModels::referenceTemperature::primaryNetMassFlux()
const
66 sumPhi +=
phi.boundaryField()[patchi][facei]*
faceSign_[i];
73 reduce(sumPhi, sumOp<scalar>());
80Foam::heatExchangerModels::referenceTemperature::primaryInletTemperature()
const
88 scalar primaryInletTfMean = 0;
92 const label facei = faceId_[i];
93 if (facePatchId_[i] != -1)
95 const label patchi = facePatchId_[i];
96 const scalar phii =
phi.boundaryField()[patchi][facei]*faceSign_[i];
97 const scalar magPhii =
mag(phii);
101 const scalar Tfi = Tf.boundaryField()[patchi][facei];
103 primaryInletTfMean += Tfi*magPhii;
107 const scalar phii =
phi[facei]*faceSign_[i];
108 const scalar magPhii =
mag(phii);
110 sumMagPhi += magPhii;
112 primaryInletTfMean += Tf[facei]*magPhii;
115 reduce(sumMagPhi, sumOp<scalar>());
116 reduce(primaryInletTfMean, sumOp<scalar>());
118 return primaryInletTfMean/(sumMagPhi + ROOTVSMALL);
123Foam::heatExchangerModels::referenceTemperature::temperatureDifferences
136 deltaTCells[i] =
max(Tref_ - TCells[i], scalar(0));
143 deltaTCells[i] =
max(TCells[i] - Tref_, scalar(0));
152Foam::heatExchangerModels::referenceTemperature::weight
158 scalar sumWeight = 0;
165 const label celli =
cells[i];
166 sumWeight +=
V[celli]*
mag(
U[celli])*deltaTCells[i];
168 reduce(sumWeight, sumOp<scalar>());
174void Foam::heatExchangerModels::referenceTemperature::writeFileHeader
179 writeFile::writeHeader(
os,
"Effectiveness heat exchanger source");
180 writeFile::writeCommented(
os,
"Time");
181 writeFile::writeTabbed(
os,
"Net mass flux [kg/s]");
182 writeFile::writeTabbed(
os,
"Total heat exchange [W]");
183 writeFile::writeTabbed(
os,
"Reference T [K]");
207 TrefTablePtr_(nullptr),
212 writeFileHeader(
file());
230 sumPhi_ = primaryNetMassFlux();
232 Qt_ = targetQdotPtr_->value(mesh_.time().value());
236 const scalar primaryInletT = primaryInletTemperature();
237 Tref_ = TrefTablePtr_()(
mag(sumPhi_), primaryInletT);
242 const scalar sumWeight = weight(
cells, deltaTCells);
244 return Qt_*deltaTCells/(sumWeight + ROOTVSMALL);
253 if (!writeFile::read(
dict))
260 if (coeffs_.readIfPresent(
"Tref", Tref_))
262 Info<<
indent <<
"- using constant reference temperature: " << Tref_
267 TrefTablePtr_.reset(
new interpolation2DTable<scalar>(coeffs_));
269 Info<<
indent <<
"- using reference temperature table"
273 UName_ = coeffs_.getOrDefault<word>(
"U",
"U");
274 TName_ = coeffs_.getOrDefault<word>(
"T",
"T");
275 phiName_ = coeffs_.getOrDefault<word>(
"phi",
"phi");
276 coeffs_.readEntry(
"faceZone", faceZoneName_);
290 <<
indent <<
"Net mass flux [kg/s] : " << sumPhi_ <<
nl
291 <<
indent <<
"Total heat exchange [W] : " << Qt_ <<
nl
292 <<
indent <<
"Reference T [K] : " << Tref_
299 writeCurrentTime(
os);
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...
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.
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.
static autoPtr< heatExchangerModel > New(const fvMesh &mesh, const word &name, const dictionary &coeffs)
Return a reference to the selected heat exchanger model.
labelList faceSign_
List of +1/-1 representing face flip map (1 use as is, -1 negate).
A heat exchanger model where heat exchange is calculated via a temperature table.
virtual void initialise()
Initialise data members of the model.
virtual bool read(const dictionary &dict)
Read top-level dictionary.
referenceTemperature(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.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
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.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
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).
static constexpr const zero Zero
Global zero (0).
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.
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.