Loading...
Searching...
No Matches
fanMomentumSource.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2022 Louis Vittoz, SimScale GmbH
9 Copyright (C) 2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "fanMomentumSource.H"
30#include "IFstream.H"
35// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace fv
40{
43}
44}
45
46
47// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48
49void Foam::fv::fanMomentumSource::writeProps
50(
51 const scalar gradP,
52 const scalar flowRate
53) const
54{
55 // Only write on output time
56 if (mesh_.time().writeTime())
57 {
59 (
61 (
62 name_ + "Properties",
64 "uniform",
65 mesh_,
69 )
70 );
71 propsDict.add("gradient", gradP);
72 propsDict.add("flow_rate", flowRate);
73 propsDict.regIOobject::write();
74 }
75}
76
77
78void Foam::fv::fanMomentumSource::initUpstreamFaces()
79{
80 // First calculate the centre of gravity of the cell zone
81 const vectorField& cellCentres = mesh_.cellCentres();
82 const scalarField& cellVolumes = mesh_.cellVolumes();
83
84 vector centreGravityCellZone(Zero);
85 scalar cellZoneVolume = 0;
86 for (const label celli : cells_)
87 {
88 const scalar cellVolume = cellVolumes[celli];
89 centreGravityCellZone += cellCentres[celli]*cellVolume;
90 cellZoneVolume += cellVolume;
91 }
92
93 reduce(centreGravityCellZone, sumOp<vector>());
94 reduce(cellZoneVolume, sumOp<scalar>());
95
96 if (cellZoneVolume < SMALL)
97 {
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."
102 << abort(FatalError);
103 }
104
105 centreGravityCellZone /= cellZoneVolume;
106
107 // Collect faces upstream of the centre of gavity
108 const faceZone& fZone = mesh_.faceZones()[surroundingFaceZoneID_];
109 const vectorField& faceCentreGravity = mesh_.faceCentres();
110
111 upstreamPatchFaceInfo_.resize_nocopy(fZone.size());
112
113 label count = 0;
114 for (const label facei : fZone)
115 {
116 if
117 (
118 (flowDir_ & (faceCentreGravity[facei] - centreGravityCellZone)) < 0
119 )
120 {
121 labelPair patchFaceInfo(-1, -1);
122
123 if (mesh_.isInternalFace(facei))
124 {
125 // Patch ID already set to -1, set only the face ID
126 patchFaceInfo.second() = facei;
127 }
128 else
129 {
130 patchFaceInfo.first() = mesh_.boundaryMesh().whichPatch(facei);
131 const polyPatch& pp = mesh_.boundaryMesh()[patchFaceInfo.first()];
132 const auto* cpp = isA<coupledPolyPatch>(pp);
133
134 if (cpp)
135 {
136 patchFaceInfo.second() =
137 cpp->owner()
138 ? pp.whichFace(facei)
139 : -1;
140 }
141 else if (!isA<emptyPolyPatch>(pp))
142 {
143 patchFaceInfo.second() = pp.whichFace(facei);
144 }
145 // else both face ID and patch ID remain at -1
146 }
147
148 // If this is an upstream face, set it in the list
149 if (patchFaceInfo.second() >= 0)
150 {
151 upstreamPatchFaceInfo_[count] = patchFaceInfo;
152 count++;
153 }
154 }
155 }
156
157 upstreamPatchFaceInfo_.setSize(count);
158
159 // Fill cellsInZones_ with all cell IDs
160 for (const label celli : cells_)
161 {
162 cellsInZones_.insert(celli);
163 }
164
165 // Sanity check
166 const labelUList& owners = mesh_.owner();
167 const labelUList& neighbours = mesh_.neighbour();
168 for (const labelPair& patchFaceInfo : upstreamPatchFaceInfo_)
169 {
170 if (patchFaceInfo.first() == -1)
171 {
172 const label facei = patchFaceInfo.second();
173 const label own = owners[facei];
174 const label nei = neighbours[facei];
175
176 // To be valid: one cell has to be inside the cellZone and the other
177 // one, outside
178 if (cellsInZones_.found(own) == cellsInZones_.found(nei))
179 {
181 << "It seems that the faceZone is not part of the cellZone "
182 << "boundaries."
183 << abort(FatalError);
184 }
185 }
186 }
187}
188
189
190Foam::scalar
191Foam::fv::fanMomentumSource::calcFlowRate
192(
193 const surfaceScalarField& phi,
195) const
196{
197 // Sanity check to make sure we're not left in an inconsistent state because
198 // here we're operating on primitive (non-dimensional) fields
199 if (isNull(rhof) != (phi.dimensions() == dimVolume/dimTime))
200 {
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."
207 << exit(FatalError);
208 }
209
210 // Calculate the flow rate through the upstream faces
211 scalarList phif(upstreamPatchFaceInfo_.size());
212
213 const labelUList& owners = mesh_.owner();
214
215 forAll(upstreamPatchFaceInfo_, i)
216 {
217 const labelPair& patchFaceInfo = upstreamPatchFaceInfo_[i];
218
219 if (isNull(rhof))
220 {
221 // Incompressible variant (phi = volumetric flux)
222 phif[i] =
223 patchFaceInfo.first() < 0
224 ? phi[patchFaceInfo.second()]
225 : phi.boundaryField()[patchFaceInfo];
226 }
227 else
228 {
229 // Compressible variant (phi = mass flux)
230 phif[i] =
231 patchFaceInfo.first() < 0
232 ? phi[patchFaceInfo.second()]/
233 rhof.internalField()[patchFaceInfo.second()]
234 : phi.boundaryField()[patchFaceInfo]/
235 rhof.boundaryField()[patchFaceInfo];
236 }
237
238 // Sign of the flux needs to be flipped if this is an internal face
239 // whose owner is found in the cell zone
240 if
241 (
242 patchFaceInfo.first() < 0
243 && cellsInZones_.found(owners[patchFaceInfo.second()])
244 )
245 {
246 phif[i] *= -1;
247 }
248 }
250 return gSum(phif);
251}
252
253
254// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
255
257(
258 const word& sourceName,
259 const word& modelType,
260 const dictionary& dict,
261 const fvMesh& mesh
262)
263:
264 fv::cellSetOption(sourceName, modelType, dict, mesh),
265 upstreamPatchFaceInfo_(),
266 cellsInZones_(),
267 fanCurvePtr_(Function1<scalar>::New("fanCurve", coeffs_, &mesh)),
268 UName_(coeffs_.getOrDefault<word>("U", "U")),
269 flowDir_(coeffs_.get<vector>("flowDir")),
270 thickness_(coeffs_.get<scalar>("thickness")),
271 gradPFan_(0.0),
272 surroundingFaceZoneID_(-1),
273 rho_(coeffs_.getOrDefault<scalar>("rho", SMALL)),
274 useRefRho_(coeffs_.found("rho"))
275{
276 // Skip all the checks if the source term has been deactivated
277 // because there are no selected cells.
278 // Such a situation typically occurs for multiple fluid regions
280 {
281 const word faceZoneName = coeffs_.getWord("faceZone");
282
283 surroundingFaceZoneID_ = mesh_.faceZones().findZoneID(faceZoneName);
284
285 if (surroundingFaceZoneID_ < 0)
286 {
288 << type() << " " << this->name() << ": "
289 << " Unknown face zone name: " << faceZoneName
290 << ". Valid face zones are: " << mesh_.faceZones().names()
291 << exit(FatalError);
292 }
293
294 flowDir_.normalise();
295 if (mag(flowDir_) < SMALL)
296 {
298 << "'flowDir' cannot be a zero-valued vector."
299 << exit(FatalIOError);
300 }
301
302 if (thickness_ < VSMALL)
303 {
305 << "'thickness' cannot be non-positive."
306 << exit(FatalIOError);
307 }
308
309 if (coeffs_.found("fields"))
310 {
312 << "Found 'fields' entry, which will be ignored. This fvOption"
313 << " will operate only on field 'U' = " << UName_
314 << endl;
315 }
316 fieldNames_.resize(1, UName_);
317
319
320 // Read the initial pressure gradient from file if it exists
321 IFstream propsFile
322 (
323 mesh_.time().timePath()/"uniform"/(name_ + "Properties")
324 );
325
326 if (propsFile.good())
327 {
328 Info<< " Reading pressure gradient from file" << endl;
329 dictionary propsDict(propsFile);
330 propsDict.readEntry("gradient", gradPFan_);
331 }
332
333 Info<< " Initial pressure gradient = " << gradPFan_ << nl << endl;
334
335 if (rho_ < SMALL)
336 {
338 << "'rho' cannot be zero or negative."
339 << exit(FatalIOError);
340 }
341
342 initUpstreamFaces();
343 }
344}
345
346
347// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
348
350(
351 fvMatrix<vector>& eqn,
352 const label fieldi
353)
354{
356 (
358 (
359 name_ + fieldNames_[fieldi] + "Sup",
360 mesh_.time().timeName(),
361 mesh_,
365 ),
366 mesh_,
368 );
369
370 const auto& phi = mesh().lookupObject<surfaceScalarField>("phi");
371
372 if (phi.dimensions() != dimVelocity*dimArea)
373 {
375 << "You called incompressible variant of addSup for case with "
376 << "a mass flux and not volumetric flux. This is not allowed."
377 << abort(FatalError);
378 }
379
380 if (!useRefRho_)
381 {
383 << "You called incompressible addSup without providing a "
384 << "reference density. Add 'rho' entry to the dictionary."
385 << abort(FatalError);
386 }
387
388 const scalar flowRate = calcFlowRate(phi, surfaceScalarField::null());
389
390 // Pressure drop for this flow rate
391 // if flow rate is negative, pressure is clipped at the static pressure
392 gradPFan_ =
393 fanCurvePtr_->value(max(flowRate, scalar(0)))/thickness_/rho_;
394
395 // Create the source term
396 UIndirectList<vector>(Su, cells_) = flowDir_*gradPFan_;
398 eqn += Su;
399
400 writeProps(gradPFan_, flowRate);
401}
402
403
405(
406 const volScalarField& rho,
407 fvMatrix<vector>& eqn,
408 const label fieldi
409)
410{
412 (
414 (
415 name_ + fieldNames_[fieldi] + "Sup",
416 mesh_.time().timeName(),
417 mesh_,
421 ),
422 mesh_,
424 );
425
426 const auto& phi = mesh().lookupObject<surfaceScalarField>("phi");
427
428 if (phi.dimensions() != dimMass/dimTime)
429 {
431 << "You called compressible variant of addSup for case with "
432 << "a volumetric flux and not mass flux. This is not allowed."
433 << abort(FatalError);
434 }
435
437 const scalar flowRate = calcFlowRate(phi, rhof);
438
439 // Pressure drop for this flow rate
440 // if flow rate is negative, pressure is clipped at the static pressure
441 gradPFan_ = fanCurvePtr_->value(max(flowRate, scalar(0)))/thickness_;
442
443 // Create the source term
444 UIndirectList<vector>(Su, cells_) = flowDir_*gradPFan_;
446 eqn += Su;
447
448 writeProps(gradPFan_, flowRate);
449}
450
451
453{
455
456 return false;
457}
458
459
460// ************************************************************************* //
bool found
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...
Definition Function1.H:92
DimensionedField< vector, volMesh > Internal
Input from file stream as an ISstream, normally using std::ifstream for the actual input.
Definition IFstream.H:55
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,...
Definition IOobject.H:191
bool good() const noexcept
True if next operation might succeed.
Definition IOstream.H:281
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition ListI.H:171
bool writeTime() const noexcept
True if this is a write interval.
Definition TimeStateI.H:73
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.
Definition Time.C:714
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.
Definition VectorI.H:114
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition ZoneMesh.C:757
wordList names() const
A list of the zone names.
Definition ZoneMesh.C:463
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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....
Definition fvMatrix.H:118
const dimensionSet & dimensions() const noexcept
Definition fvMatrix.H:530
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
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).
Definition fvOption.H:124
const word & name() const noexcept
Return const access to the source name.
Definition fvOptionI.H:24
const fvMesh & mesh_
Reference to the mesh database.
Definition fvOption.H:142
static autoPtr< option > New(const word &name, const dictionary &dict, const fvMesh &mesh)
Return a reference to the selected fvOption model.
Definition fvOption.C:78
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition fvOption.H:157
dictionary coeffs_
Dictionary containing source coefficients.
Definition fvOption.H:152
virtual bool isActive()
Is the source active?
Definition fvOption.C:115
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition fvOption.C:41
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
Definition fvOptionI.H:30
const word name_
Source name.
Definition fvOption.H:132
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition polyMesh.H:671
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volVectorField gradP(fvc::grad(p))
dynamicFvMesh & mesh
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition error.H:688
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
zeroField Su
Definition alphaSuSp.H:1
#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.
Definition BitOps.H:73
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.
Namespace for OpenFOAM.
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
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.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
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)
Definition errorManip.H:139
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).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
const dimensionSet dimVolume(pow3(dimLength))
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
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.
Definition nullObject.H:248
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299