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-2023 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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 
44 namespace Foam
45 {
46 namespace fv
47 {
48  defineTypeNameAndDebug(directionalPressureGradientExplicitSource, 0);
50  (
51  option,
52  directionalPressureGradientExplicitSource,
53  dictionary
54  );
55 }
56 }
57 
58 
59 // * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * * //
60 
61 const Foam::Enum
62 <
64 >
65 Foam::fv::directionalPressureGradientExplicitSource::pressureDropModelNames_
66 ({
67  { pressureDropModel::pVolumetricFlowRateTable, "volumetricFlowRateTable" },
68  { pressureDropModel::pConstant, "constant" },
69  { pressureDropModel::pDarcyForchheimer, "DarcyForchheimer" },
70 });
71 
72 
73 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
74 
75 void 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 
103  if (isA<emptyPolyPatch>(pp))
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 
134 void 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_,
152  )
153  );
154  propsDict.add("gradient", gradP);
155  propsDict.regIOobject::write();
156  }
157 }
158 
159 
160 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
161 
164 (
165  const word& sourceName,
166  const word& modelType,
167  const dictionary& dict,
168  const fvMesh& mesh
169 )
170 :
171  fv::cellSetOption(sourceName, modelType, dict, mesh),
172  model_(pressureDropModelNames_.get("model", coeffs_)),
173  gradP0_(cells_.size(), Zero),
174  dGradP_(cells_.size(), Zero),
175  gradPporous_(cells_.size(), Zero),
176  flowDir_(coeffs_.get<vector>("flowDir")),
177  invAPtr_(nullptr),
178  D_(0),
179  I_(0),
180  length_(0),
181  pressureDrop_(0),
182  flowRate_(),
183  faceZoneName_(coeffs_.get<word>("faceZone")),
184  zoneID_(mesh_.faceZones().findZoneID(faceZoneName_)),
185  faceId_(),
186  facePatchId_(),
187  relaxationFactor_(coeffs_.getOrDefault<scalar>("relaxationFactor", 0.3)),
188  cellFaceMap_(cells_.size(), -1)
189 {
190  coeffs_.readEntry("fields", fieldNames_);
191 
192  flowDir_.normalise();
193 
194  if (fieldNames_.size() != 1)
195  {
197  << "Source can only be applied to a single field. Current "
198  << "settings are:" << fieldNames_ << exit(FatalError);
199  }
200 
201  if (zoneID_ < 0)
202  {
204  << type() << " " << this->name() << ": "
205  << " Unknown face zone name: " << faceZoneName_
206  << ". Valid face zones are: " << mesh_.faceZones().names()
207  << exit(FatalError);
208  }
209 
210  if (model_ == pVolumetricFlowRateTable)
211  {
212  flowRate_ = interpolationTable<scalar>(coeffs_);
213  }
214  else if (model_ == pConstant)
215  {
216  coeffs_.readEntry("pressureDrop", pressureDrop_);
217  }
218  else if (model_ == pDarcyForchheimer)
219  {
220  coeffs_.readEntry("D", D_);
221  coeffs_.readEntry("I", I_);
222  coeffs_.readEntry("length", length_);
223  }
224  else
225  {
227  << "Did not find mode " << model_ << nl
228  << "Please set 'model' to one of "
229  << pressureDropModelNames_
230  << exit(FatalError);
231  }
232 
234 
235  // Read the initial pressure gradient from file if it exists
236  IFstream propsFile
237  (
238  mesh_.time().timePath()/"uniform"/(name_ + "Properties")
239  );
240 
241  if (propsFile.good())
242  {
243  Info<< " Reading pressure gradient from file" << endl;
244  dictionary propsDict(propsFile);
245  propsDict.readEntry("gradient", gradP0_);
246  }
247 
248  Info<< " Initial pressure gradient = " << gradP0_ << nl << endl;
249 
250  initialise();
251 }
252 
253 
254 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
255 
257 (
259 )
260 {
261  const scalarField& rAU = invAPtr_().internalField();
262 
263  const scalarField magUn(mag(U), cells_);
264 
265  const auto& phi = mesh().lookupObject<surfaceScalarField>("phi");
266 
267  switch (model_)
268  {
269  case pDarcyForchheimer:
270  {
271  if (phi.dimensions() == dimVolume/dimTime)
272  {
273  const auto& turbModel =
274  mesh().lookupObject<incompressible::turbulenceModel>
275  (
277  );
278 
279  const scalarField nu(turbModel.nu(), cells_);
280 
281  gradPporous_ = -flowDir_*(D_*nu + I_*0.5*magUn)*magUn*length_;
282  }
283  else
284  {
285  const auto& turbModel =
286  mesh().lookupObject<compressible::turbulenceModel>
287  (
289  );
290 
291  const scalarField mu(turbModel.mu(),cells_);
292 
293  const scalarField rho(turbModel.rho(),cells_);
294 
295  gradPporous_ =
296  - flowDir_*(D_*mu + I_*0.5*rho*magUn)*magUn*length_;
297  }
298  break;
299  }
300  case pConstant:
301  {
302  gradPporous_ = -flowDir_*pressureDrop_;
303  break;
304  }
305 
306  case pVolumetricFlowRateTable:
307  {
308  scalar volFlowRate = 0;
309  scalar totalphi = 0;
310 
311  forAll(faceId_, i)
312  {
313  label faceI = faceId_[i];
314  if (facePatchId_[i] != -1)
315  {
316  label patchI = facePatchId_[i];
317  totalphi += phi.boundaryField()[patchI][faceI];
318  }
319  else
320  {
321  totalphi += phi[faceI];
322  }
323  }
324  reduce(totalphi, sumOp<scalar>());
325 
326  if (phi.dimensions() == dimVolume/dimTime)
327  {
328  volFlowRate = mag(totalphi);
329  }
330  else
331  {
332  const auto& turbModel =
333  mesh().lookupObject<compressible::turbulenceModel>
334  (
336  );
337  const scalarField rho(turbModel.rho(),cells_);
338  const scalarField cv(mesh_.V(), cells_);
339  scalar rhoAve = gSumProd(rho, cv)/gSum(cv);
340  volFlowRate = mag(totalphi)/rhoAve;
341  }
342 
343  gradPporous_ = -flowDir_*flowRate_(volFlowRate);
344  break;
345  }
346  }
347 
348  const faceZone& fZone = mesh_.faceZones()[zoneID_];
349 
350  labelList meshToLocal(mesh_.nCells(), -1);
351  forAll(cells_, i)
352  {
353  meshToLocal[cells_[i]] = i;
354  }
355 
356  labelList faceToCellIndex(fZone.size(), -1);
357  const labelList& mc = fZone.masterCells();
358  const labelList& sc = fZone.slaveCells();
359 
360  forAll(fZone, i)
361  {
362  label masterCellI = mc[i];
363 
364  if (meshToLocal[masterCellI] != -1 && masterCellI != -1)
365  {
366  faceToCellIndex[i] = meshToLocal[masterCellI];
367  }
368  else if (meshToLocal[masterCellI] == -1)
369  {
371  << "Did not find cell " << masterCellI
372  << "in cellZone :" << zoneName()
373  << exit(FatalError);
374  }
375  }
376 
377  // Accumulate 'upstream' velocity into cells
378  vectorField UfCells(cells_.size(), Zero);
379  scalarField UfCellWeights(cells_.size(), Zero);
380 
381  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
382 
383  FieldField<Field, vector> upwindValues(pbm.size());
384 
385  forAll(U.boundaryField(), patchI)
386  {
387  const fvPatchVectorField& pf = U.boundaryField()[patchI];
388 
389  if (pf.coupled())
390  {
391  upwindValues.set(patchI, pf.patchNeighbourField());
392  }
393  else if (!isA<emptyFvPatchScalarField>(pf))
394  {
395  upwindValues.set(patchI, new vectorField(pf));
396  }
397  }
398 
399  forAll(fZone, i)
400  {
401  label faceI = fZone[i];
402  label cellId = faceToCellIndex[i];
403 
404  if (cellId != -1)
405  {
406  label sourceCellId = sc[i];
407  if (mesh_.isInternalFace(faceI))
408  {
409  scalar w = mesh_.magSf()[faceI];
410  UfCells[cellId] += U[sourceCellId]*w;
411  UfCellWeights[cellId] += w;
412  }
413  else if (fZone.flipMap()[i])
414  {
415  const label patchI = pbm.patchID(faceI);
416  label localFaceI = pbm[patchI].whichFace(faceI);
417 
418  scalar w = mesh_.magSf().boundaryField()[patchI][localFaceI];
419 
420  if (upwindValues.set(patchI))
421  {
422  UfCells[cellId] += upwindValues[patchI][localFaceI]*w;
423  UfCellWeights[cellId] += w;
424  }
425  }
426  }
427  }
428 
429  UfCells /= UfCellWeights;
430 
431  forAll(cells_, i)
432  {
433  label cellI = cells_[i];
434 
435  const vector Ufnorm(UfCells[i]/(mag(UfCells[i]) + SMALL));
436 
437  const tensor D(rotationTensor(Ufnorm, flowDir_));
438 
439  dGradP_[i] +=
440  relaxationFactor_*
441  (
442  (D & UfCells[i]) - U[cellI]
443  )/rAU[cellI];
444 
445 
446  if (debug)
447  {
448  Info<< "Difference mag(U) = "
449  << mag(UfCells[i]) - mag(U[cellI])
450  << endl;
451  Info<< "Pressure drop in flowDir direction : "
452  << gradPporous_[i] << endl;
453  Info<< "UfCell:= " << UfCells[i] << "U : " << U[cellI] << endl;
454  }
455  }
456  writeProps(gradP0_ + dGradP_);
457 }
458 
459 
461 (
462  fvMatrix<vector>& eqn,
463  const label fieldI
464 )
465 {
466  DimensionedField<vector, volMesh> Su
467  (
468  IOobject
469  (
470  name_ + fieldNames_[fieldI] + "Sup",
471  mesh_.time().timeName(),
472  mesh_,
475  ),
476  mesh_,
477  dimensionedVector(eqn.dimensions()/dimVolume, Zero)
478  );
479 
480  UIndirectList<vector>(Su, cells_) = gradP0_ + dGradP_ + gradPporous_;
481 
482  eqn += Su;
483 }
484 
485 
487 (
488  const volScalarField& rho,
489  fvMatrix<vector>& eqn,
490  const label fieldI
491 )
492 {
493  this->addSup(eqn, fieldI);
494 }
495 
496 
498 (
499  fvMatrix<vector>& eqn,
500  const label
501 )
502 {
503  if (!invAPtr_)
504  {
505  invAPtr_.reset
506  (
507  new volScalarField
508  (
509  IOobject
510  (
511  name_ + ":invA",
512  mesh_.time().timeName(),
513  mesh_,
516  ),
517  1.0/eqn.A()
518  )
519  );
520  }
521  else
522  {
523  invAPtr_() = 1.0/eqn.A();
524  }
526  gradP0_ += dGradP_;
527  dGradP_ = Zero;
528 }
529 
530 
532 (
533  Ostream& os
534 ) const
535 {
537 }
538 
539 
541 (
542  const dictionary& dict
543 )
544 {
545  const dictionary coeffs(dict.subDict(typeName + "Coeffs"));
546 
547  relaxationFactor_ =
548  coeffs.getOrDefault<scalar>("relaxationFactor", 0.3);
549 
550  coeffs.readEntry("flowDir", flowDir_);
551  flowDir_.normalise();
552 
553  if (model_ == pConstant)
554  {
555  coeffs.readEntry("pressureDrop", pressureDrop_);
556  }
557  else if (model_ == pDarcyForchheimer)
558  {
559  coeffs.readEntry("D", D_);
560  coeffs.readEntry("I", I_);
561  coeffs.readEntry("length", length_);
562  }
563 
564  return false;
565 }
566 
567 
568 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
const polyBoundaryMesh & pbm
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
label faceId(-1)
fvPatchField< vector > fvPatchVectorField
zeroField Su
Definition: alphaSuSp.H:1
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition: fvOption.H:157
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:160
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
fileName timePath() const
Return current time path = path/timeName.
Definition: Time.H:520
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:142
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
scalarProduct< Type, Type >::type gSumProd(const UList< Type > &f1, const UList< Type > &f2, const label comm)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Tensor< scalar > tensor
Definition: symmTensor.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition: ListI.H:175
Ignore writing from objectRegistry::writeObject()
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition: transform.H:47
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
Macros for easy insertion into run-time selection tables.
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
IOdictionary propsDict(dictIO)
virtual void correct(volVectorField &U)
Correct the pressure gradient.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
volVectorField gradP(fvc::grad(p))
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...
Templated abstract base class for single-phase incompressible turbulence models.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
const word name_
Source name.
Definition: fvOption.H:132
static const word propertiesName
Default name of the turbulence properties dictionary.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
3D tensor transformation operations.
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
const word & name() const noexcept
Return const access to the source name.
Definition: fvOptionI.H:24
Vector< scalar > vector
Definition: vector.H:57
label whichPatch(const label meshFacei) const
Return patch index for a given mesh face index. Uses binary search.
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
virtual void writeData(Ostream &os) const
Write the source properties.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:670
virtual bool read(const dictionary &dict)
Read source dictionary.
const dimensionedScalar mu
Atomic mass unit.
virtual void constrain(fvMatrix< vector > &eqn, const label fieldI)
Set 1/A coefficient.
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:1313
U
Definition: pEqn.H:72
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
directionalPressureGradientExplicitSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:56
label cellId
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:41
Nothing to be read.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
A special matrix type and solver, designed for finite volume solutions of scalar equations.
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
wordList names() const
A list of the zone names.
Definition: ZoneMesh.C:362
messageStream Info
Information stream (stdout output on master, null elsewhere)
Field< vector > vectorField
Specialisation of Field<T> for vector.
Intermediate abstract class for handling cell-set options for the derived fvOptions.
const dimensionedScalar & D
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
List< label > labelList
A List of labels.
Definition: List.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:686
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
dictionary coeffs_
Dictionary containing source coefficients.
Definition: fvOption.H:152
volScalarField & nu
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
virtual void addSup(fvMatrix< vector > &eqn, const label fieldI)
Add explicit contribution to momentum equation.