49void Foam::fv::fanMomentumSource::writeProps
78void Foam::fv::fanMomentumSource::initUpstreamFaces()
81 const vectorField& cellCentres = mesh_.cellCentres();
82 const scalarField& cellVolumes = mesh_.cellVolumes();
85 scalar cellZoneVolume = 0;
86 for (
const label celli : cells_)
88 const scalar cellVolume = cellVolumes[celli];
89 centreGravityCellZone += cellCentres[celli]*cellVolume;
90 cellZoneVolume += cellVolume;
93 reduce(centreGravityCellZone, sumOp<vector>());
94 reduce(cellZoneVolume, sumOp<scalar>());
96 if (cellZoneVolume < SMALL)
99 <<
"Cannot initialize upstream faces because the total cell zone"
100 <<
" volume is zero or negative: " << cellZoneVolume <<
"." <<
nl
101 <<
"Check whether there are cells in fan momentum source cell zone."
105 centreGravityCellZone /= cellZoneVolume;
108 const faceZone& fZone = mesh_.faceZones()[surroundingFaceZoneID_];
109 const vectorField& faceCentreGravity = mesh_.faceCentres();
114 for (
const label facei : fZone)
118 (flowDir_ & (faceCentreGravity[facei] - centreGravityCellZone)) < 0
123 if (mesh_.isInternalFace(facei))
126 patchFaceInfo.second() = facei;
130 patchFaceInfo.first() = mesh_.boundaryMesh().whichPatch(facei);
131 const polyPatch&
pp = mesh_.boundaryMesh()[patchFaceInfo.first()];
136 patchFaceInfo.second() =
138 ?
pp.whichFace(facei)
143 patchFaceInfo.second() =
pp.whichFace(facei);
149 if (patchFaceInfo.second() >= 0)
151 upstreamPatchFaceInfo_[
count] = patchFaceInfo;
157 upstreamPatchFaceInfo_.setSize(count);
160 for (
const label celli : cells_)
162 cellsInZones_.insert(celli);
167 const labelUList& neighbours = mesh_.neighbour();
168 for (
const labelPair& patchFaceInfo : upstreamPatchFaceInfo_)
170 if (patchFaceInfo.first() == -1)
172 const label facei = patchFaceInfo.second();
173 const label own = owners[facei];
174 const label nei = neighbours[facei];
178 if (cellsInZones_.found(own) == cellsInZones_.found(nei))
181 <<
"It seems that the faceZone is not part of the cellZone "
191Foam::fv::fanMomentumSource::calcFlowRate
202 <<
"Incorrect usage of the function. " <<
nl
203 <<
"For incompressible flow, pass surfaceScalarField::null for rhof"
204 <<
" and volumetric flux for phi." <<
nl
205 <<
"For compressible flow, pass face-interpolated density as rhof"
206 <<
" and mass flux for phi."
211 scalarList phif(upstreamPatchFaceInfo_.size());
215 forAll(upstreamPatchFaceInfo_, i)
217 const labelPair& patchFaceInfo = upstreamPatchFaceInfo_[i];
223 patchFaceInfo.first() < 0
224 ?
phi[patchFaceInfo.second()]
225 :
phi.boundaryField()[patchFaceInfo];
231 patchFaceInfo.first() < 0
232 ?
phi[patchFaceInfo.second()]/
233 rhof.internalField()[patchFaceInfo.second()]
234 :
phi.boundaryField()[patchFaceInfo]/
235 rhof.boundaryField()[patchFaceInfo];
242 patchFaceInfo.first() < 0
243 && cellsInZones_.found(owners[patchFaceInfo.second()])
258 const word& sourceName,
259 const word& modelType,
265 upstreamPatchFaceInfo_(),
268 UName_(coeffs_.getOrDefault<
word>(
"U",
"U")),
269 flowDir_(coeffs_.get<
vector>(
"flowDir")),
270 thickness_(coeffs_.get<scalar>(
"thickness")),
272 surroundingFaceZoneID_(-1),
273 rho_(coeffs_.getOrDefault<scalar>(
"rho", SMALL)),
274 useRefRho_(coeffs_.
found(
"rho"))
285 if (surroundingFaceZoneID_ < 0)
288 <<
type() <<
" " << this->
name() <<
": "
289 <<
" Unknown face zone name: " << faceZoneName
295 if (
mag(flowDir_) < SMALL)
298 <<
"'flowDir' cannot be a zero-valued vector."
302 if (thickness_ < VSMALL)
305 <<
"'thickness' cannot be non-positive."
312 <<
"Found 'fields' entry, which will be ignored. This fvOption"
313 <<
" will operate only on field 'U' = " << UName_
323 mesh_.time().timePath()/
"uniform"/(
name_ +
"Properties")
326 if (propsFile.good())
328 Info<<
" Reading pressure gradient from file" <<
endl;
330 propsDict.readEntry(
"gradient", gradPFan_);
333 Info<<
" Initial pressure gradient = " << gradPFan_ <<
nl <<
endl;
338 <<
"'rho' cannot be zero or negative."
359 name_ + fieldNames_[fieldi] +
"Sup",
360 mesh_.time().timeName(),
375 <<
"You called incompressible variant of addSup for case with "
376 <<
"a mass flux and not volumetric flux. This is not allowed."
383 <<
"You called incompressible addSup without providing a "
384 <<
"reference density. Add 'rho' entry to the dictionary."
393 fanCurvePtr_->value(
max(flowRate, scalar(0)))/thickness_/rho_;
396 UIndirectList<vector>(
Su, cells_) = flowDir_*gradPFan_;
400 writeProps(gradPFan_, flowRate);
415 name_ + fieldNames_[fieldi] +
"Sup",
416 mesh_.time().timeName(),
431 <<
"You called compressible variant of addSup for case with "
432 <<
"a volumetric flux and not mass flux. This is not allowed."
437 const scalar flowRate = calcFlowRate(
phi,
rhof);
441 gradPFan_ = fanCurvePtr_->value(
max(flowRate, scalar(0)))/thickness_;
448 writeProps(gradPFan_, flowRate);
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
IOdictionary propsDict(dictIO)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
static const this_type & null() noexcept
DimensionedField< vector, volMesh > Internal
Input from file stream as an ISstream, normally using std::ifstream for the actual input.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
bool good() const noexcept
True if next operation might succeed.
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
bool writeTime() const noexcept
True if this is a write interval.
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
A List with indirect addressing. Like IndirectList but does not store addressing.
Vector< Cmpt > & normalise(const scalar tol=ROOTVSMALL)
Inplace normalise the vector by its magnitude.
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
wordList names() const
A list of the zone names.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
word getWord(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Same as get<word>(const word&, keyType::option).
A special matrix type and solver, designed for finite volume solutions of scalar equations....
const dimensionSet & dimensions() const noexcept
Mesh data needed to do the Finite Volume discretisation.
const Time & time() const
Return the top-level database.
Intermediate abstract class for handling cell-set options for the derived fvOptions.
labelList cells_
Set of cells to apply source to.
cellSetOption(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
This source term models the action of a fan on the flow. It calculates flow rate through a set of ups...
fanMomentumSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
virtual bool read(const dictionary &dict)
Read source dictionary.
virtual void addSup(fvMatrix< vector > &eqn, const label fieldi)
Add explicit contribution to momentum equation.
Base abstract class for handling finite volume options (i.e. fvOption).
const word & name() const noexcept
Return const access to the source name.
const fvMesh & mesh_
Reference to the mesh database.
static autoPtr< option > New(const word &name, const dictionary &dict, const fvMesh &mesh)
Return a reference to the selected fvOption model.
wordList fieldNames_
Field names to apply source to - populated by derived models.
dictionary coeffs_
Dictionary containing source coefficients.
virtual bool isActive()
Is the source active?
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
const word name_
Source name.
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
volVectorField gradP(fvc::grad(p))
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
surfaceScalarField rhof(fvc::interpolate(rho, "div(phi,rho)"))
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Namespace for finite-volume.
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.
Pair< label > labelPair
A pair of labels.
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.
Type gSum(const FieldField< Field, Type > &f)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
const dimensionSet dimVelocity
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.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
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).
errorManip< error > abort(error &err)
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
const dimensionSet dimVolume(pow3(dimLength))
errorManipArg< error, int > exit(error &err, const int errNo=1)
UList< label > labelUList
A UList of labels.
List< scalar > scalarList
List of scalar.
void Su(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
bool isNull(const T *ptr) noexcept
True if ptr is a pointer (of type T) to the nullObject.
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.