Go to the documentation of this file.
44 effectivenessHeatExchangerSource,
64 label faceI = fZone[i];
66 label facePatchId = -1;
76 if (isA<coupledPolyPatch>(pp))
78 if (refCast<const coupledPolyPatch>(pp).owner())
87 else if (!isA<emptyPolyPatch>(pp))
129 label faceI = faceId_[i];
130 if (facePatchId_[i] != -1)
132 label patchI = facePatchId_[i];
133 area += mesh_.magSf().boundaryField()[patchI][faceI];
137 area += mesh_.magSf()[faceI];
149 const word& modelType,
155 secondaryMassFlowRate_(
readScalar(coeffs_.lookup(
"secondaryMassFlowRate"))),
156 secondaryInletT_(
readScalar(coeffs_.lookup(
"secondaryInletT"))),
157 primaryInletT_(
readScalar(coeffs_.lookup(
"primaryInletT"))),
159 UName_(coeffs_.lookupOrDefault<
word>(
"UName",
"U")),
160 TName_(coeffs_.lookupOrDefault<
word>(
"TName",
"T")),
161 phiName_(coeffs_.lookupOrDefault<
word>(
"phiName",
"phi")),
162 faceZoneName_(coeffs_.lookup(
"faceZone")),
163 zoneID_(mesh_.faceZones().findZoneID(faceZoneName_)),
172 <<
type() <<
" " << this->
name() <<
": "
173 <<
" Unknown face zone name: " << faceZoneName_
174 <<
". Valid face zones are: " << mesh_.faceZones().names()
184 fieldNames_.setSize(1,
thermo.he().name());
186 applied_.setSize(1,
false);
215 label faceI = faceId_[i];
216 if (facePatchId_[i] != -1)
218 label patchI = facePatchId_[i];
219 totalphi +=
phi.boundaryField()[patchI][faceI]*faceSign_[i];
222 Cpf.boundaryField()[patchI][faceI]
223 *mesh_.magSf().boundaryField()[patchI][faceI];
227 totalphi +=
phi[faceI]*faceSign_[i];
228 CpfMean += Cpf[faceI]*mesh_.magSf()[faceI];
235 eTable_()(
mag(totalphi), secondaryMassFlowRate_)
236 *(secondaryInletT_ - primaryInletT_)
237 *(CpfMean/faceZoneArea_)*
mag(totalphi);
258 deltaTCells[i] =
max(Tref - TCells[i], 0.0);
262 deltaTCells[i] =
max(TCells[i] - Tref, 0.0);
268 scalar sumWeight = 0;
271 sumWeight += V[cells_[i]]*
mag(
U[cells_[i]])*deltaTCells[i];
275 if (this->V() > VSMALL &&
mag(Qt) > VSMALL)
281 heSource[cells_[i]] -=
282 Qt*V[cells_[i]]*
mag(
U[cells_[i]])*deltaTCells[i]/sumWeight;
288 Info<<
indent <<
"Net mass flux [Kg/s] = " << totalphi <<
nl;
289 Info<<
indent <<
"Total energy exchange [W] = " << Qt <<
nl;
292 << eTable_()(
mag(totalphi), secondaryMassFlowRate_) <<
endl;
301 coeffs_.lookup(
"secondaryMassFlowRate") >> secondaryMassFlowRate_;
302 coeffs_.lookup(
"secondaryInletT") >> secondaryInletT_;
303 coeffs_.lookup(
"primaryInletT") >> primaryInletT_;
A class for handling words, derived from string.
Cell-set options abtract base class. Provides a base set of controls, e.g.
#define forAll(list, i)
Loop across all elements in list.
virtual bool read(const dictionary &dict)
Read dictionary.
static const word dictName
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Abstract base-class for fluid and solid thermodynamic properties.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Ostream & endl(Ostream &os)
Add newline and flush stream.
void initialise()
Initialise heat exchanger source model.
dimensioned< scalar > mag(const dimensioned< Type > &)
const fvMesh & mesh_
Reference to the mesh database.
const faceZoneMesh & faceZones() const
Return face zone mesh.
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldI)
Scalar.
2D table interpolation. The data must be in ascending order in both dimensions x and y.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
labelList faceId_
Local list of face IDs.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
A subset of mesh faces organised as a primitive patch.
A patch is a list of labels that address the faces in the global face list.
virtual bool read(const dictionary &dict)
Read source dictionary.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
labelList faceSign_
List of +1/-1 representing face flip map (1 use as is, -1 negate)
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Macros for easy insertion into run-time selection tables.
Mesh data needed to do the Finite Volume discretisation.
Ostream & indent(Ostream &os)
Indent stream.
label start() const
Return start label of this patch in the polyMesh face list.
scalar faceZoneArea_
Area of the face zone.
errorManipArg< error, int > exit(error &err, const int errNo=1)
static bool master(const label communicator=0)
Am I the master process.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
defineTypeNameAndDebug(option, 0)
void setSize(const label)
Reset size of List.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
void calculateTotalArea(scalar &area)
Calculate total area of faceZone accross processors.
label whichFace(const label l) const
Return label of face in patch from global face label.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
effectivenessHeatExchangerSource(const effectivenessHeatExchangerSource &)
Disallow default bitwise copy construct.
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
labelList facePatchId_
Local list of patch ID per face.
label zoneID_
Id for the face zone.
void size(const label)
Override size to be inconsistent with allocated storage.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Generic GeometricField class.
const boolList & flipMap() const
Return face flip map.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.