Loading...
Searching...
No Matches
directionalPressureGradientExplicitSource.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) 2015-2025 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
29#include "fvMatrices.H"
30#include "DimensionedField.H"
31#include "IFstream.H"
33#include "transform.H"
34#include "surfaceInterpolate.H"
35#include "turbulenceModel.H"
38#include "vectorFieldIOField.H"
39#include "FieldField.H"
40#include "emptyFvPatchFields.H"
42// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
43
44namespace Foam
45{
46namespace fv
47{
50 (
51 option,
54 );
55}
56}
57
58
59// * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * * //
60
61const Foam::Enum
62<
64>
65Foam::fv::directionalPressureGradientExplicitSource::pressureDropModelNames_
66({
67 { pressureDropModel::pVolumetricFlowRateTable, "volumetricFlowRateTable" },
68 { pressureDropModel::pConstant, "constant" },
69 { pressureDropModel::pDarcyForchheimer, "DarcyForchheimer" },
70});
71
72
73// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
74
75void Foam::fv::directionalPressureGradientExplicitSource::initialise()
76{
77 const faceZone& fZone = mesh_.faceZones()[zoneID_];
78
79 // Total number of faces selected
80 label numFaces = fZone.size();
81
82 faceId_.resize_nocopy(numFaces);
83 facePatchId_.resize_nocopy(numFaces);
84
85 numFaces = 0;
86
87 // TDB: handle multiple zones
88 {
89 forAll(fZone, i)
90 {
91 const label meshFacei = fZone[i];
92
93 // Internal faces
94 label faceId = meshFacei;
95 label facePatchId = -1;
96
97 // Boundary faces
98 if (!mesh_.isInternalFace(meshFacei))
99 {
100 facePatchId = mesh_.boundaryMesh().whichPatch(meshFacei);
101 const polyPatch& pp = mesh_.boundaryMesh()[facePatchId];
102
104 {
105 continue; // Ignore empty patch
106 }
107
108 const auto* cpp = isA<coupledPolyPatch>(pp);
109
110 if (cpp && !cpp->owner())
111 {
112 continue; // Ignore neighbour side
113 }
114
115 faceId = pp.whichFace(meshFacei);
116 }
117
118 if (faceId >= 0)
119 {
120 faceId_[numFaces] = faceId;
121 facePatchId_[numFaces] = facePatchId;
122
123 ++numFaces;
124 }
125 }
126 }
127
128 // Shrink to size used
129 faceId_.resize(numFaces);
130 facePatchId_.resize(numFaces);
131}
132
133
134void Foam::fv::directionalPressureGradientExplicitSource::writeProps
135(
136 const vectorField& gradP
137) const
138{
139 // Only write on output time
140 if (mesh_.time().writeTime())
141 {
142 IOdictionary propsDict
143 (
144 IOobject
145 (
146 name_ + "Properties",
147 mesh_.time().timeName(),
148 "uniform",
149 mesh_,
153 )
154 );
155 propsDict.add("gradient", gradP);
156 propsDict.regIOobject::write();
157 }
158}
159
160
161// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
162
165(
166 const word& sourceName,
167 const word& modelType,
168 const dictionary& dict,
169 const fvMesh& mesh
170)
171:
172 fv::cellSetOption(sourceName, modelType, dict, mesh),
173 model_(pressureDropModelNames_.get("model", coeffs_)),
174 gradP0_(cells_.size(), Zero),
175 dGradP_(cells_.size(), Zero),
176 gradPporous_(cells_.size(), Zero),
177 flowDir_(coeffs_.get<vector>("flowDir")),
178 invAPtr_(nullptr),
179 D_(0),
180 I_(0),
181 length_(0),
182 pressureDrop_(0),
183 flowRate_(),
184 faceZoneName_(coeffs_.get<word>("faceZone")),
185 zoneID_(mesh_.faceZones().findZoneID(faceZoneName_)),
186 faceId_(),
187 facePatchId_(),
188 relaxationFactor_(coeffs_.getOrDefault<scalar>("relaxationFactor", 0.3)),
189 cellFaceMap_(cells_.size(), -1)
190{
191 coeffs_.readEntry("fields", fieldNames_);
192
193 flowDir_.normalise();
194
195 if (fieldNames_.size() != 1)
196 {
197 FatalErrorInFunction
198 << "Source can only be applied to a single field. Current "
199 << "settings are:" << fieldNames_ << exit(FatalError);
200 }
201
202 if (zoneID_ < 0)
203 {
205 << type() << " " << this->name() << ": "
206 << " Unknown face zone name: " << faceZoneName_
207 << ". Valid face zones are: " << mesh_.faceZones().names()
208 << exit(FatalError);
209 }
210
211 if (model_ == pVolumetricFlowRateTable)
212 {
213 flowRate_ = interpolationTable<scalar>(coeffs_);
214 }
215 else if (model_ == pConstant)
216 {
217 coeffs_.readEntry("pressureDrop", pressureDrop_);
218 }
219 else if (model_ == pDarcyForchheimer)
220 {
221 coeffs_.readEntry("D", D_);
222 coeffs_.readEntry("I", I_);
223 coeffs_.readEntry("length", length_);
224 }
225 else
226 {
228 << "Did not find mode " << model_ << nl
229 << "Please set 'model' to one of "
230 << pressureDropModelNames_
231 << exit(FatalError);
232 }
233
235
236 // Read the initial pressure gradient from file if it exists
237 IFstream propsFile
238 (
239 mesh_.time().timePath()/"uniform"/(name_ + "Properties")
240 );
241
242 if (propsFile.good())
243 {
244 Info<< " Reading pressure gradient from file" << endl;
245 dictionary propsDict(propsFile);
246 propsDict.readEntry("gradient", gradP0_);
247 }
248
249 Info<< " Initial pressure gradient = " << gradP0_ << nl << endl;
251 initialise();
252}
253
254
255// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
256
258(
260)
261{
262 const scalarField& rAU = invAPtr_().internalField();
263
264 const scalarField magUn(mag(U), cells_);
265
266 const auto& phi = mesh().lookupObject<surfaceScalarField>("phi");
267
268 switch (model_)
269 {
270 case pDarcyForchheimer:
271 {
272 if (phi.dimensions() == dimVolume/dimTime)
273 {
274 const auto& turbModel =
276 (
278 );
279
280 const scalarField nu(turbModel.nu(), cells_);
281
282 gradPporous_ = -flowDir_*(D_*nu + I_*0.5*magUn)*magUn*length_;
283 }
284 else
285 {
286 const auto& turbModel =
287 mesh().lookupObject<compressible::turbulenceModel>
288 (
290 );
291
292 const scalarField mu(turbModel.mu(),cells_);
293
294 const scalarField rho(turbModel.rho(),cells_);
295
296 gradPporous_ =
297 - flowDir_*(D_*mu + I_*0.5*rho*magUn)*magUn*length_;
298 }
299 break;
300 }
301 case pConstant:
302 {
303 gradPporous_ = -flowDir_*pressureDrop_;
304 break;
305 }
306
307 case pVolumetricFlowRateTable:
308 {
309 scalar volFlowRate = 0;
310 scalar totalphi = 0;
311
312 forAll(faceId_, i)
313 {
314 label faceI = faceId_[i];
315 if (facePatchId_[i] != -1)
316 {
317 label patchI = facePatchId_[i];
318 totalphi += phi.boundaryField()[patchI][faceI];
319 }
320 else
321 {
322 totalphi += phi[faceI];
323 }
324 }
325 reduce(totalphi, sumOp<scalar>());
326
327 if (phi.dimensions() == dimVolume/dimTime)
328 {
329 volFlowRate = mag(totalphi);
330 }
331 else
332 {
333 const auto& turbModel =
334 mesh().lookupObject<compressible::turbulenceModel>
335 (
337 );
338 const scalarField rho(turbModel.rho(),cells_);
339 const scalarField cv(mesh_.V(), cells_);
340 scalar rhoAve = gWeightedAverage(cv, rho);
341 volFlowRate = mag(totalphi)/rhoAve;
342 }
343
344 gradPporous_ = -flowDir_*flowRate_(volFlowRate);
345 break;
346 }
347 }
348
349 const faceZone& fZone = mesh_.faceZones()[zoneID_];
350
351 labelList meshToLocal(mesh_.nCells(), -1);
352 forAll(cells_, i)
353 {
354 meshToLocal[cells_[i]] = i;
355 }
356
357 labelList faceToCellIndex(fZone.size(), -1);
358 const labelList& mc = fZone.masterCells();
359 const labelList& sc = fZone.slaveCells();
360
361 forAll(fZone, i)
362 {
363 label masterCellI = mc[i];
364
365 if (meshToLocal[masterCellI] != -1 && masterCellI != -1)
366 {
367 faceToCellIndex[i] = meshToLocal[masterCellI];
368 }
369 else if (meshToLocal[masterCellI] == -1)
370 {
372 << "Did not find cell " << masterCellI
373 << "in cellZone :" << zoneName()
374 << exit(FatalError);
375 }
376 }
377
378 // Accumulate 'upstream' velocity into cells
379 vectorField UfCells(cells_.size(), Zero);
380 scalarField UfCellWeights(cells_.size(), Zero);
381
382 const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
383
384 FieldField<Field, vector> upwindValues(pbm.size());
385
386 forAll(U.boundaryField(), patchI)
387 {
388 const fvPatchVectorField& pf = U.boundaryField()[patchI];
389
390 if (pf.coupled())
391 {
392 upwindValues.set(patchI, pf.patchNeighbourField());
393 }
394 else if (!isA<emptyFvPatchScalarField>(pf))
395 {
396 upwindValues.set(patchI, new vectorField(pf));
397 }
398 }
399
400 forAll(fZone, i)
401 {
402 label faceI = fZone[i];
403 label cellId = faceToCellIndex[i];
404
405 if (cellId != -1)
406 {
407 label sourceCellId = sc[i];
408 if (mesh_.isInternalFace(faceI))
409 {
410 scalar w = mesh_.magSf()[faceI];
411 UfCells[cellId] += U[sourceCellId]*w;
412 UfCellWeights[cellId] += w;
413 }
414 else if (fZone.flipMap()[i])
415 {
416 const label patchI = pbm.patchID(faceI);
417 label localFaceI = pbm[patchI].whichFace(faceI);
418
419 scalar w = mesh_.magSf().boundaryField()[patchI][localFaceI];
420
421 if (upwindValues.set(patchI))
422 {
423 UfCells[cellId] += upwindValues[patchI][localFaceI]*w;
424 UfCellWeights[cellId] += w;
425 }
426 }
427 }
428 }
429
430 UfCells /= UfCellWeights;
431
432 forAll(cells_, i)
433 {
434 label cellI = cells_[i];
435
436 const vector Ufnorm(UfCells[i]/(mag(UfCells[i]) + SMALL));
437
438 const tensor D(rotationTensor(Ufnorm, flowDir_));
439
440 dGradP_[i] +=
441 relaxationFactor_*
442 (
443 (D & UfCells[i]) - U[cellI]
444 )/rAU[cellI];
445
446
447 if (debug)
448 {
449 Info<< "Difference mag(U) = "
450 << mag(UfCells[i]) - mag(U[cellI])
451 << endl;
452 Info<< "Pressure drop in flowDir direction : "
453 << gradPporous_[i] << endl;
454 Info<< "UfCell:= " << UfCells[i] << "U : " << U[cellI] << endl;
455 }
456 }
457 writeProps(gradP0_ + dGradP_);
458}
459
460
462(
463 fvMatrix<vector>& eqn,
464 const label fieldI
465)
466{
468 (
469 name_ + fieldNames_[fieldI] + "Sup",
471 mesh_,
472 dimensionedVector(eqn.dimensions()/dimVolume, Zero)
473 );
474 auto& Su = tSu.ref();
476 UIndirectList<vector>(Su, cells_) = gradP0_ + dGradP_ + gradPporous_;
477
478 eqn += tSu;
479}
480
481
483(
484 const volScalarField& rho,
485 fvMatrix<vector>& eqn,
486 const label fieldI
487)
488{
489 this->addSup(eqn, fieldI);
490}
491
492
494(
495 fvMatrix<vector>& eqn,
496 const label
497)
498{
499 if (!invAPtr_)
500 {
501 invAPtr_.reset
502 (
504 (
505 mesh_.newIOobject(IOobject::scopedName(name_, "invA")),
506 1.0/eqn.A()
507 )
508 );
509 }
510 else
511 {
512 invAPtr_() = 1.0/eqn.A();
514
515 gradP0_ += dGradP_;
516 dGradP_ = Zero;
517}
518
519
521(
523) const
524{
526}
527
528
530(
531 const dictionary& dict
532)
533{
534 const dictionary coeffs(dict.subDict(typeName + "Coeffs"));
535
536 relaxationFactor_ =
537 coeffs.getOrDefault<scalar>("relaxationFactor", 0.3);
538
539 coeffs.readEntry("flowDir", flowDir_);
540 flowDir_.normalise();
541
542 if (model_ == pConstant)
543 {
544 coeffs.readEntry("pressureDrop", pressureDrop_);
545 }
546 else if (model_ == pDarcyForchheimer)
547 {
548 coeffs.readEntry("D", D_);
549 coeffs.readEntry("I", I_);
550 coeffs.readEntry("length", length_);
551 }
552
553 return false;
554}
555
556
557// ************************************************************************* //
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)
if(patchID !=-1)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const Field< Type > &iField)
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions, copy of internal field....
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
A field of fields is a PtrList of fields with reference counting.
Definition FieldField.H:77
Input from file stream as an ISstream, normally using std::ifstream for the actual input.
Definition IFstream.H:55
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
bool good() const noexcept
True if next operation might succeed.
Definition IOstream.H:281
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Vector< Cmpt > & normalise(const scalar tol=ROOTVSMALL)
Inplace normalise the vector by its magnitude.
Definition VectorI.H:114
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
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
A subset of mesh faces organised as a primitive patch.
Definition faceZone.H:63
const labelList & slaveCells() const
Deprecated(2023-09) same as backCells.
Definition faceZone.H:559
const boolList & flipMap() const noexcept
Return face flip map.
Definition faceZone.H:389
const labelList & masterCells() const
Deprecated(2023-09) same as frontCells.
Definition faceZone.H:552
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
tmp< volScalarField > A() const
Return the central coefficient.
Definition fvMatrix.C:1306
const dimensionSet & dimensions() const noexcept
Definition fvMatrix.H:530
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
virtual bool coupled() const
True if the patch field is coupled.
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
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.
const wordRe & zoneName() const
Return const access to the first set/zone name.
Applies an explicit pressure gradient source in such a way to deflect the flow towards an specific di...
directionalPressureGradientExplicitSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
virtual void writeData(Ostream &os) const
Write the source properties.
virtual bool read(const dictionary &dict)
Read source dictionary.
virtual void addSup(fvMatrix< vector > &eqn, const label fieldI)
Add explicit contribution to momentum equation.
virtual void constrain(fvMatrix< vector > &eqn, const label fieldI)
Set 1/A coefficient.
virtual void correct(volVectorField &U)
Correct the pressure gradient.
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 dictionary & coeffs() const noexcept
Return dictionary.
Definition fvOptionI.H:36
const fvMesh & mesh_
Reference to the mesh database.
Definition fvOption.H:142
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
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
An interpolation/look-up table of scalar vs <Type> values. The reference scalar values must be monoto...
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition polyMesh.H:671
static const word propertiesName
Default name of the turbulence properties dictionary.
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
U
Definition pEqn.H:72
volVectorField gradP(fvc::grad(p))
dynamicFvMesh & mesh
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition error.H:688
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
A special matrix type and solver, designed for finite volume solutions of scalar equations.
tmp< volScalarField > rAU
zeroField Su
Definition alphaSuSp.H:1
label cellId
label faceId(-1)
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for finite-volume.
IncompressibleTurbulenceModel< transportModel > turbulenceModel
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Type gWeightedAverage(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted average of a field, using the mag() of the weights.
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition transform.H:47
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Tensor< scalar > tensor
Definition symmTensor.H:57
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.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
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).
Field< vector > vectorField
Specialisation of Field<T> for vector.
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
void Su(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
fvPatchField< vector > fvPatchVectorField
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & nu
dictionary dict
const dimensionedScalar & D
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
3D tensor transformation operations.