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-2021 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"
41
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 faceId_.setSize(fZone.size());
80 facePatchId_.setSize(fZone.size());
81
82 label count = 0;
83 forAll(fZone, i)
84 {
85 const label faceI = fZone[i];
86
87 label faceId = -1;
88 label facePatchId = -1;
89 if (mesh_.isInternalFace(faceI))
90 {
91 faceId = faceI;
92 facePatchId = -1;
93 }
94 else
95 {
96 facePatchId = mesh_.boundaryMesh().whichPatch(faceI);
97 const polyPatch& pp = mesh_.boundaryMesh()[facePatchId];
98 const auto* cpp = isA<coupledPolyPatch>(pp);
99
100 if (cpp)
101 {
102 faceId = (cpp->owner() ? pp.whichFace(faceI) : -1);
103 }
104 else if (!isA<emptyPolyPatch>(pp))
105 {
106 faceId = pp.whichFace(faceI);
107 }
108 else
109 {
110 faceId = -1;
111 facePatchId = -1;
112 }
113 }
114
115 if (faceId >= 0)
116 {
117 facePatchId_[count] = facePatchId;
118 faceId_[count] = faceId;
119 count++;
120 }
121 }
122 faceId_.setSize(count);
123 facePatchId_.setSize(count);
124}
125
126
127void Foam::fv::directionalPressureGradientExplicitSource::writeProps
128(
129 const vectorField& gradP
130) const
131{
132 // Only write on output time
133 if (mesh_.time().writeTime())
134 {
135 IOdictionary propsDict
136 (
137 IOobject
138 (
139 name_ + "Properties",
140 mesh_.time().timeName(),
141 "uniform",
142 mesh_,
145 )
146 );
147 propsDict.add("gradient", gradP);
148 propsDict.regIOobject::write();
149 }
150}
151
152
153// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
154
155Foam::fv::directionalPressureGradientExplicitSource::
156directionalPressureGradientExplicitSource
157(
158 const word& sourceName,
159 const word& modelType,
160 const dictionary& dict,
161 const fvMesh& mesh
162)
163:
164 fv::cellSetOption(sourceName, modelType, dict, mesh),
165 model_(pressureDropModelNames_.get("model", coeffs_)),
166 gradP0_(cells_.size(), Zero),
167 dGradP_(cells_.size(), Zero),
168 gradPporous_(cells_.size(), Zero),
169 flowDir_(coeffs_.get<vector>("flowDir")),
170 invAPtr_(nullptr),
171 D_(0),
172 I_(0),
173 length_(0),
174 pressureDrop_(0),
175 flowRate_(),
176 faceZoneName_(coeffs_.get<word>("faceZone")),
177 zoneID_(mesh_.faceZones().findZoneID(faceZoneName_)),
178 faceId_(),
179 facePatchId_(),
180 relaxationFactor_(coeffs_.getOrDefault<scalar>("relaxationFactor", 0.3)),
181 cellFaceMap_(cells_.size(), -1)
182{
183 coeffs_.readEntry("fields", fieldNames_);
184
185 flowDir_.normalise();
186
187 if (fieldNames_.size() != 1)
188 {
190 << "Source can only be applied to a single field. Current "
191 << "settings are:" << fieldNames_ << exit(FatalError);
192 }
193
194 if (zoneID_ < 0)
195 {
197 << type() << " " << this->name() << ": "
198 << " Unknown face zone name: " << faceZoneName_
199 << ". Valid face zones are: " << mesh_.faceZones().names()
200 << exit(FatalError);
201 }
202
203 if (model_ == pVolumetricFlowRateTable)
204 {
206 }
207 else if (model_ == pConstant)
208 {
209 coeffs_.readEntry("pressureDrop", pressureDrop_);
210 }
211 else if (model_ == pDarcyForchheimer)
212 {
213 coeffs_.readEntry("D", D_);
214 coeffs_.readEntry("I", I_);
215 coeffs_.readEntry("length", length_);
216 }
217 else
218 {
220 << "Did not find mode " << model_ << nl
221 << "Please set 'model' to one of "
222 << pressureDropModelNames_
223 << exit(FatalError);
224 }
225
227
228 // Read the initial pressure gradient from file if it exists
229 IFstream propsFile
230 (
231 mesh_.time().timePath()/"uniform"/(name_ + "Properties")
232 );
233
234 if (propsFile.good())
235 {
236 Info<< " Reading pressure gradient from file" << endl;
237 dictionary propsDict(propsFile);
238 propsDict.readEntry("gradient", gradP0_);
239 }
240
241 Info<< " Initial pressure gradient = " << gradP0_ << nl << endl;
242
243 initialise();
244}
245
246
247// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
248
250(
252)
253{
254 const scalarField& rAU = invAPtr_().internalField();
255
256 const scalarField magUn(mag(U), cells_);
257
258 const auto& phi = mesh().lookupObject<surfaceScalarField>("phi");
259
260 switch (model_)
261 {
262 case pDarcyForchheimer:
263 {
264 if (phi.dimensions() == dimVelocity*dimArea)
265 {
266 const auto& turbModel =
268 (
270 );
271
272 const scalarField nu(turbModel.nu(), cells_);
273
274 gradPporous_ = -flowDir_*(D_*nu + I_*0.5*magUn)*magUn*length_;
275 }
276 else
277 {
278 const auto& turbModel =
279 mesh().lookupObject<compressible::turbulenceModel>
280 (
282 );
283
284 const scalarField mu(turbModel.mu(),cells_);
285
286 const scalarField rho(turbModel.rho(),cells_);
287
288 gradPporous_ =
289 - flowDir_*(D_*mu + I_*0.5*rho*magUn)*magUn*length_;
290 }
291 break;
292 }
293 case pConstant:
294 {
295 gradPporous_ = -flowDir_*pressureDrop_;
296 break;
297 }
298
299 case pVolumetricFlowRateTable:
300 {
301 scalar volFlowRate = 0;
302 scalar totalphi = 0;
303
304 forAll(faceId_, i)
305 {
306 label faceI = faceId_[i];
307 if (facePatchId_[i] != -1)
308 {
309 label patchI = facePatchId_[i];
310 totalphi += phi.boundaryField()[patchI][faceI];
311 }
312 else
313 {
314 totalphi += phi[faceI];
315 }
316 }
317 reduce(totalphi, sumOp<scalar>());
318
319 if (phi.dimensions() == dimVelocity*dimArea)
320 {
321 volFlowRate = mag(totalphi);
322 }
323 else
324 {
325 const auto& turbModel =
326 mesh().lookupObject<compressible::turbulenceModel>
327 (
329 );
330 const scalarField rho(turbModel.rho(),cells_);
331 const scalarField cv(mesh_.V(), cells_);
332 scalar rhoAve = gSumProd(rho, cv)/gSum(cv);
333 volFlowRate = mag(totalphi)/rhoAve;
334 }
335
336 gradPporous_ = -flowDir_*flowRate_(volFlowRate);
337 break;
338 }
339 }
340
341 const faceZone& fZone = mesh_.faceZones()[zoneID_];
342
343 labelList meshToLocal(mesh_.nCells(), -1);
344 forAll(cells_, i)
345 {
346 meshToLocal[cells_[i]] = i;
347 }
348
349 labelList faceToCellIndex(fZone.size(), -1);
350 const labelList& mc = fZone.masterCells();
351 const labelList& sc = fZone.slaveCells();
352
353 forAll(fZone, i)
354 {
355 label masterCellI = mc[i];
356
357 if (meshToLocal[masterCellI] != -1 && masterCellI != -1)
358 {
359 faceToCellIndex[i] = meshToLocal[masterCellI];
360 }
361 else if (meshToLocal[masterCellI] == -1)
362 {
364 << "Did not find cell " << masterCellI
365 << "in cellZone :" << zoneName()
366 << exit(FatalError);
367 }
368 }
369
370 // Accumulate 'upstream' velocity into cells
371 vectorField UfCells(cells_.size(), Zero);
372 scalarField UfCellWeights(cells_.size(), Zero);
373
374 const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
375
376 FieldField<Field, vector> upwindValues(pbm.size());
377
378 forAll(U.boundaryField(), patchI)
379 {
380 const fvPatchVectorField& pf = U.boundaryField()[patchI];
381
382 if (pf.coupled())
383 {
384 upwindValues.set(patchI, pf.patchNeighbourField());
385 }
386 else if (!isA<emptyFvPatchScalarField>(pf))
387 {
388 upwindValues.set(patchI, new vectorField(pf));
389 }
390 }
391
392 forAll(fZone, i)
393 {
394 label faceI = fZone[i];
395 label cellId = faceToCellIndex[i];
396
397 if (cellId != -1)
398 {
399 label sourceCellId = sc[i];
400 if (mesh_.isInternalFace(faceI))
401 {
402 scalar w = mesh_.magSf()[faceI];
403 UfCells[cellId] += U[sourceCellId]*w;
404 UfCellWeights[cellId] += w;
405 }
406 else if (fZone.flipMap()[i])
407 {
408 label patchI = pbm.patchID()[faceI-mesh_.nInternalFaces()];
409 label localFaceI = pbm[patchI].whichFace(faceI);
410
411 scalar w = mesh_.magSf().boundaryField()[patchI][localFaceI];
412
413 if (upwindValues.set(patchI))
414 {
415 UfCells[cellId] += upwindValues[patchI][localFaceI]*w;
416 UfCellWeights[cellId] += w;
417 }
418 }
419 }
420 }
421
422 UfCells /= UfCellWeights;
423
424 forAll(cells_, i)
425 {
426 label cellI = cells_[i];
427
428 const vector Ufnorm(UfCells[i]/(mag(UfCells[i]) + SMALL));
429
430 const tensor D(rotationTensor(Ufnorm, flowDir_));
431
432 dGradP_[i] +=
433 relaxationFactor_*
434 (
435 (D & UfCells[i]) - U[cellI]
436 )/rAU[cellI];
437
438
439 if (debug)
440 {
441 Info<< "Difference mag(U) = "
442 << mag(UfCells[i]) - mag(U[cellI])
443 << endl;
444 Info<< "Pressure drop in flowDir direction : "
445 << gradPporous_[i] << endl;
446 Info<< "UfCell:= " << UfCells[i] << "U : " << U[cellI] << endl;
447 }
448 }
449 writeProps(gradP0_ + dGradP_);
450}
451
452
454(
455 fvMatrix<vector>& eqn,
456 const label fieldI
457)
458{
460 (
462 (
463 name_ + fieldNames_[fieldI] + "Sup",
464 mesh_.time().timeName(),
465 mesh_,
468 ),
469 mesh_,
471 );
472
473 UIndirectList<vector>(Su, cells_) = gradP0_ + dGradP_ + gradPporous_;
474
475 eqn += Su;
476}
477
478
480(
481 const volScalarField& rho,
482 fvMatrix<vector>& eqn,
483 const label fieldI
484)
485{
486 this->addSup(eqn, fieldI);
487}
488
489
491(
492 fvMatrix<vector>& eqn,
493 const label
494)
495{
496 if (!invAPtr_)
497 {
498 invAPtr_.reset
499 (
501 (
503 (
504 name_ + ":invA",
505 mesh_.time().timeName(),
506 mesh_,
509 ),
510 1.0/eqn.A()
511 )
512 );
513 }
514 else
515 {
516 invAPtr_() = 1.0/eqn.A();
517 }
518
519 gradP0_ += dGradP_;
520 dGradP_ = Zero;
521}
522
523
525(
526 Ostream& os
527) const
528{
530}
531
532
534(
535 const dictionary& dict
536)
537{
538 const dictionary coeffs(dict.subDict(typeName + "Coeffs"));
539
540 relaxationFactor_ =
541 coeffs.getOrDefault<scalar>("relaxationFactor", 0.3);
542
543 coeffs.readEntry("flowDir", flowDir_);
544 flowDir_.normalise();
545
546 if (model_ == pConstant)
547 {
548 coeffs.readEntry("pressureDrop", pressureDrop_);
549 }
550 else if (model_ == pDarcyForchheimer)
551 {
552 coeffs.readEntry("D", D_);
553 coeffs.readEntry("I", I_);
554 coeffs.readEntry("length", length_);
555 }
556
557 return false;
558}
559
560
561// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
surfaceScalarField & phi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: Enum.H:61
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:80
Input from file stream, using an ISstream.
Definition: IFstream.H:57
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
bool good() const noexcept
True if next operation might succeed.
Definition: IOstream.H:233
Templated abstract base class for single-phase incompressible turbulence models.
void setSize(const label n)
Alias for resize()
Definition: List.H:218
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
const T * set(const label i) const
Definition: PtrList.H:138
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual bool read()
Re-read model coefficients if they have changed.
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
fileName timePath() const
Return current time path.
Definition: Time.H:375
A List with indirect addressing. Like IndirectList but does not store addressing.
Definition: IndirectList.H:79
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
wordList names() const
A list of the zone names.
Definition: ZoneMesh.C:304
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:67
const labelList & masterCells() const
Definition: faceZone.C:389
const boolList & flipMap() const noexcept
Return face flip map.
Definition: faceZone.H:272
const labelList & slaveCells() const
Return labels of slave cells.
Definition: faceZone.C:400
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:1303
const dimensionSet & dimensions() const noexcept
Definition: fvMatrix.H:453
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:453
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:350
Intermediate abstract class for handling cell-set options for the derived fvOptions.
Applies an explicit pressure gradient source in such a way to deflect the flow towards an specific di...
virtual void writeData(Ostream &os) const
Write the source properties.
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.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition: fvOption.H:127
const word & name() const noexcept
Return const access to the source name.
Definition: fvOptionI.H:31
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:139
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition: fvOption.H:148
dictionary coeffs_
Dictionary containing source coefficients.
Definition: fvOption.H:145
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:48
const word name_
Source name.
Definition: fvOption.H:133
An interpolation/look-up table of scalar vs <Type> values. The reference scalar values must be monoto...
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
const labelList & patchID() const
Per boundary face label the patch index.
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:498
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
const volScalarField & mu
volVectorField gradP(fvc::grad(p))
dynamicFvMesh & mesh
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
A special matrix type and solver, designed for finite volume solutions of scalar equations.
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
zeroField Su
Definition: alphaSuSp.H:1
label cellId
label faceId(-1)
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:78
Namespace for OpenFOAM.
scalarProduct< Type, Type >::type gSumProd(const UList< Type > &f1, const UList< Type > &f2, const label comm)
Type gSum(const FieldField< Field, Type > &f)
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition: transform.H:51
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList fv(nPoints)
volScalarField & nu
dictionary dict
const dimensionedScalar & D
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
3D tensor transformation operations.
IOdictionary propsDict(dictIO)