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-2020 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  faceId_.setSize(fZone.size());
80  facePatchId_.setSize(fZone.size());
81 
82  label count = 0;
83  forAll(fZone, i)
84  {
85  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  if (isA<coupledPolyPatch>(pp))
99  {
100  if (refCast<const coupledPolyPatch>(pp).owner())
101  {
102  faceId = pp.whichFace(faceI);
103  }
104  else
105  {
106  faceId = -1;
107  }
108  }
109  else if (!isA<emptyPolyPatch>(pp))
110  {
111  faceId = faceI - pp.start();
112  }
113  else
114  {
115  faceId = -1;
116  facePatchId = -1;
117  }
118  }
119 
120  if (faceId >= 0)
121  {
122  facePatchId_[count] = facePatchId;
123  faceId_[count] = faceId;
124  count++;
125  }
126  }
127  faceId_.setSize(count);
128  facePatchId_.setSize(count);
129 }
130 
131 
132 void Foam::fv::directionalPressureGradientExplicitSource::writeProps
133 (
134  const vectorField& gradP
135 ) const
136 {
137  // Only write on output time
138  if (mesh_.time().writeTime())
139  {
140  IOdictionary propsDict
141  (
142  IOobject
143  (
144  name_ + "Properties",
145  mesh_.time().timeName(),
146  "uniform",
147  mesh_,
150  )
151  );
152  propsDict.add("gradient", gradP);
153  propsDict.regIOobject::write();
154  }
155 }
156 
157 
158 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
159 
162 (
163  const word& sourceName,
164  const word& modelType,
165  const dictionary& dict,
166  const fvMesh& mesh
167 )
168 :
169  cellSetOption(sourceName, modelType, dict, mesh),
170  model_(pressureDropModelNames_.get("model", coeffs_)),
171  gradP0_(cells_.size(), Zero),
172  dGradP_(cells_.size(), Zero),
173  gradPporous_(cells_.size(), Zero),
174  flowDir_(coeffs_.get<vector>("flowDir")),
175  invAPtr_(nullptr),
176  D_(0),
177  I_(0),
178  length_(0),
179  pressureDrop_(0),
180  flowRate_(),
181  faceZoneName_(coeffs_.get<word>("faceZone")),
182  zoneID_(mesh_.faceZones().findZoneID(faceZoneName_)),
183  faceId_(),
184  facePatchId_(),
185  relaxationFactor_(coeffs_.getOrDefault<scalar>("relaxationFactor", 0.3)),
186  cellFaceMap_(cells_.size(), -1)
187 {
188  coeffs_.readEntry("fields", fieldNames_);
189 
190  flowDir_.normalise();
191 
192  if (fieldNames_.size() != 1)
193  {
195  << "Source can only be applied to a single field. Current "
196  << "settings are:" << fieldNames_ << exit(FatalError);
197  }
198 
199  if (zoneID_ < 0)
200  {
202  << type() << " " << this->name() << ": "
203  << " Unknown face zone name: " << faceZoneName_
204  << ". Valid face zones are: " << mesh_.faceZones().names()
205  << exit(FatalError);
206  }
207 
208  if (model_ == pVolumetricFlowRateTable)
209  {
210  flowRate_ = interpolationTable<scalar>(coeffs_);
211  }
212  else if (model_ == pConstant)
213  {
214  coeffs_.readEntry("pressureDrop", pressureDrop_);
215  }
216  else if (model_ == pDarcyForchheimer)
217  {
218  coeffs_.readEntry("D", D_);
219  coeffs_.readEntry("I", I_);
220  coeffs_.readEntry("length", length_);
221  }
222  else
223  {
225  << "Did not find mode " << model_ << nl
226  << "Please set 'model' to one of "
227  << pressureDropModelNames_
228  << exit(FatalError);
229  }
230 
231  applied_.setSize(fieldNames_.size(), false);
232 
233  // Read the initial pressure gradient from file if it exists
234  IFstream propsFile
235  (
236  mesh_.time().timePath()/"uniform"/(name_ + "Properties")
237  );
238 
239  if (propsFile.good())
240  {
241  Info<< " Reading pressure gradient from file" << endl;
243  propsDict.readEntry("gradient", gradP0_);
244  }
245 
246  Info<< " Initial pressure gradient = " << gradP0_ << nl << endl;
247 
248  initialise();
249 }
250 
251 
252 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
253 
255 (
257 )
258 {
259  const scalarField& rAU = invAPtr_().internalField();
260 
261  const scalarField magUn(mag(U), cells_);
262 
263  const auto& phi = mesh().lookupObject<surfaceScalarField>("phi");
264 
265  switch (model_)
266  {
267  case pDarcyForchheimer:
268  {
269  if (phi.dimensions() == dimVelocity*dimArea)
270  {
271  const auto& turbModel =
272  mesh().lookupObject<incompressible::turbulenceModel>
273  (
275  );
276 
277  const scalarField nu(turbModel.nu(), cells_);
278 
279  gradPporous_ = -flowDir_*(D_*nu + I_*0.5*magUn)*magUn*length_;
280  }
281  else
282  {
283  const auto& turbModel =
284  mesh().lookupObject<compressible::turbulenceModel>
285  (
287  );
288 
289  const scalarField mu(turbModel.mu(),cells_);
290 
291  const scalarField rho(turbModel.rho(),cells_);
292 
293  gradPporous_ =
294  - flowDir_*(D_*mu + I_*0.5*rho*magUn)*magUn*length_;
295  }
296  break;
297  }
298  case pConstant:
299  {
300  gradPporous_ = -flowDir_*pressureDrop_;
301  break;
302  }
303 
304  case pVolumetricFlowRateTable:
305  {
306  scalar volFlowRate = 0;
307  scalar totalphi = 0;
308 
309  forAll(faceId_, i)
310  {
311  label faceI = faceId_[i];
312  if (facePatchId_[i] != -1)
313  {
314  label patchI = facePatchId_[i];
315  totalphi += phi.boundaryField()[patchI][faceI];
316  }
317  else
318  {
319  totalphi += phi[faceI];
320  }
321  }
322  reduce(totalphi, sumOp<scalar>());
323 
324  if (phi.dimensions() == dimVelocity*dimArea)
325  {
326  volFlowRate = mag(totalphi);
327  }
328  else
329  {
330  const auto& turbModel =
331  mesh().lookupObject<compressible::turbulenceModel>
332  (
334  );
335  const scalarField rho(turbModel.rho(),cells_);
336  const scalarField cv(mesh_.V(), cells_);
337  scalar rhoAve = gSumProd(rho, cv)/gSum(cv);
338  volFlowRate = mag(totalphi)/rhoAve;
339  }
340 
341  gradPporous_ = -flowDir_*flowRate_(volFlowRate);
342  break;
343  }
344  }
345 
346  const faceZone& fZone = mesh_.faceZones()[zoneID_];
347 
348  labelList meshToLocal(mesh_.nCells(), -1);
349  forAll(cells_, i)
350  {
351  meshToLocal[cells_[i]] = i;
352  }
353 
354  labelList faceToCellIndex(fZone.size(), -1);
355  const labelList& mc = fZone.masterCells();
356  const labelList& sc = fZone.slaveCells();
357 
358  forAll(fZone, i)
359  {
360  label masterCellI = mc[i];
361 
362  if (meshToLocal[masterCellI] != -1 && masterCellI != -1)
363  {
364  faceToCellIndex[i] = meshToLocal[masterCellI];
365  }
366  else if (meshToLocal[masterCellI] == -1)
367  {
369  << "Did not find cell " << masterCellI
370  << "in cellZone :" << cellSetName()
371  << exit(FatalError);
372  }
373  }
374 
375  // Accumulate 'upstream' velocity into cells
376  vectorField UfCells(cells_.size(), Zero);
377  scalarField UfCellWeights(cells_.size(), Zero);
378 
379  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
380 
381  FieldField<Field, vector> upwindValues(pbm.size());
382 
383  forAll(U.boundaryField(), patchI)
384  {
385  const fvPatchVectorField& pf = U.boundaryField()[patchI];
386 
387  if (pf.coupled())
388  {
389  upwindValues.set(patchI, pf.patchNeighbourField());
390  }
391  else if (!isA<emptyFvPatchScalarField>(pf))
392  {
393  upwindValues.set(patchI, new vectorField(pf));
394  }
395  }
396 
397  forAll(fZone, i)
398  {
399  label faceI = fZone[i];
400  label cellId = faceToCellIndex[i];
401 
402  if (cellId != -1)
403  {
404  label sourceCellId = sc[i];
405  if (mesh_.isInternalFace(faceI))
406  {
407  scalar w = mesh_.magSf()[faceI];
408  UfCells[cellId] += U[sourceCellId]*w;
409  UfCellWeights[cellId] += w;
410  }
411  else if (fZone.flipMap()[i])
412  {
413  label patchI = pbm.patchID()[faceI-mesh_.nInternalFaces()];
414  label localFaceI = pbm[patchI].whichFace(faceI);
415 
416  scalar w = mesh_.magSf().boundaryField()[patchI][localFaceI];
417 
418  if (upwindValues.set(patchI))
419  {
420  UfCells[cellId] += upwindValues[patchI][localFaceI]*w;
421  UfCellWeights[cellId] += w;
422  }
423  }
424  }
425  }
426 
427  UfCells /= UfCellWeights;
428 
429  forAll(cells_, i)
430  {
431  label cellI = cells_[i];
432 
433  const vector Ufnorm(UfCells[i]/(mag(UfCells[i]) + SMALL));
434 
435  const tensor D(rotationTensor(Ufnorm, flowDir_));
436 
437  dGradP_[i] +=
438  relaxationFactor_*
439  (
440  (D & UfCells[i]) - U[cellI]
441  )/rAU[cellI];
442 
443 
444  if (debug)
445  {
446  Info<< "Difference mag(U) = "
447  << mag(UfCells[i]) - mag(U[cellI])
448  << endl;
449  Info<< "Pressure drop in flowDir direction : "
450  << gradPporous_[i] << endl;
451  Info<< "UfCell:= " << UfCells[i] << "U : " << U[cellI] << endl;
452  }
453  }
454  writeProps(gradP0_ + dGradP_);
455 }
456 
457 
459 (
460  fvMatrix<vector>& eqn,
461  const label fieldI
462 )
463 {
465  (
466  IOobject
467  (
468  name_ + fieldNames_[fieldI] + "Sup",
469  mesh_.time().timeName(),
470  mesh_,
473  ),
474  mesh_,
476  );
477 
478  UIndirectList<vector>(Su, cells_) = gradP0_ + dGradP_ + gradPporous_;
479 
480  eqn += Su;
481 }
482 
483 
485 (
486  const volScalarField& rho,
487  fvMatrix<vector>& eqn,
488  const label fieldI
489 )
490 {
491  this->addSup(eqn, fieldI);
492 }
493 
494 
496 (
497  fvMatrix<vector>& eqn,
498  const label
499 )
500 {
501  if (!invAPtr_)
502  {
503  invAPtr_.reset
504  (
505  new volScalarField
506  (
507  IOobject
508  (
509  name_ + ":invA",
510  mesh_.time().timeName(),
511  mesh_,
514  ),
515  1.0/eqn.A()
516  )
517  );
518  }
519  else
520  {
521  invAPtr_() = 1.0/eqn.A();
522  }
523 
524  gradP0_ += dGradP_;
525  dGradP_ = Zero;
526 }
527 
528 
530 (
531  Ostream& os
532 ) const
533 {
535 }
536 
537 
539 (
540  const dictionary& dict
541 )
542 {
543  const dictionary coeffs(dict.subDict(typeName + "Coeffs"));
544 
545  relaxationFactor_ =
546  coeffs.getOrDefault<scalar>("relaxationFactor", 0.3);
547 
548  coeffs.readEntry("flowDir", flowDir_);
549  flowDir_.normalise();
550 
551  if (model_ == pConstant)
552  {
553  coeffs.readEntry("pressureDrop", pressureDrop_);
554  }
555  else if (model_ == pDarcyForchheimer)
556  {
557  coeffs.readEntry("D", D_);
558  coeffs.readEntry("I", I_);
559  coeffs.readEntry("length", length_);
560  }
561 
562  return false;
563 }
564 
565 
566 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvPatchField< vector >
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::Tensor< scalar >
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
gradP
volVectorField gradP(fvc::grad(p))
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
rAU
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fv::cellSetOption
Intermediate abstract class for handling cell-set options for the derived fvOptions.
Definition: cellSetOption.H:163
Foam::FieldField
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:53
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFieldRefs.H:4
FieldField.H
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::IFstream
Input from file stream, using an ISstream.
Definition: IFstream.H:53
DimensionedField.H
Foam::fvMatrix::dimensions
const dimensionSet & dimensions() const
Definition: fvMatrix.H:292
Foam::fv::directionalPressureGradientExplicitSource::constrain
virtual void constrain(fvMatrix< vector > &eqn, const label fieldI)
Set 1/A coefficient.
Definition: directionalPressureGradientExplicitSource.C:496
turbulentTransportModel.H
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
rho
rho
Definition: readInitialConditions.H:88
Foam::fv::option::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:139
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:486
Foam::dictionary::null
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition: dictionary.H:385
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Su
zeroField Su
Definition: alphaSuSp.H:1
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:445
Foam::fvPatchField::coupled
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:331
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:60
Foam::Field< scalar >
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:65
vectorFieldIOField.H
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
directionalPressureGradientExplicitSource.H
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
faceId
label faceId(-1)
Foam::faceZone::slaveCells
const labelList & slaveCells() const
Return labels of slave cells.
Definition: faceZone.C:407
propsDict
IOdictionary propsDict(IOobject("particleTrackProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED))
Foam::polyBoundaryMesh::whichPatch
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Definition: polyBoundaryMesh.C:810
Foam::fv::directionalPressureGradientExplicitSource::directionalPressureGradientExplicitSource
directionalPressureGradientExplicitSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
Definition: directionalPressureGradientExplicitSource.C:162
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
emptyFvPatchFields.H
IFstream.H
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::gSumProd
scalarProduct< Type, Type >::type gSumProd(const UList< Type > &f1, const UList< Type > &f2, const label comm)
Definition: FieldFunctions.C:605
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMatrix::A
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:770
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::fv::directionalPressureGradientExplicitSource::pressureDropModel
pressureDropModel
Modes of pressure drop.
Definition: directionalPressureGradientExplicitSource.H:231
cellId
label cellId
Definition: interrogateWallPatches.H:67
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::fvPatchField::patchNeighbourField
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:434
fv
labelList fv(nPoints)
U
U
Definition: pEqn.H:72
Foam::faceZone::masterCells
const labelList & masterCells() const
Definition: faceZone.C:396
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::fv::directionalPressureGradientExplicitSource::addSup
virtual void addSup(fvMatrix< vector > &eqn, const label fieldI)
Add explicit contribution to momentum equation.
Definition: directionalPressureGradientExplicitSource.C:459
Foam::BitOps::count
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:77
Foam::Vector< scalar >
Foam::ThermalDiffusivity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Definition: phaseCompressibleTurbulenceModelFwd.H:47
Foam::List< label >
Foam::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
Definition: primitiveMeshI.H:102
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
D
const dimensionedScalar & D
Definition: solveBulkSurfactant.H:4
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:60
Foam::IncompressibleTurbulenceModel
Templated abstract base class for single-phase incompressible turbulence models.
Definition: IncompressibleTurbulenceModel.H:55
Foam::fv::defineTypeNameAndDebug
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
Foam::fv::directionalPressureGradientExplicitSource::writeData
virtual void writeData(Ostream &os) const
Write the source properties.
Definition: directionalPressureGradientExplicitSource.C:530
Foam::interpolationTable< scalar >
Foam::fv::addToRunTimeSelectionTable
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
transform.H
3D tensor transformation operations.
Foam::rotationTensor
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition: transform.H:51
Foam::IOstream::good
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:224
Foam::fv::directionalPressureGradientExplicitSource::correct
virtual void correct(volVectorField &U)
Correct the pressure gradient.
Definition: directionalPressureGradientExplicitSource.C:255
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:122
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::faceZone::flipMap
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:271
turbulenceModel.H
Foam::fv::directionalPressureGradientExplicitSource::read
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: directionalPressureGradientExplicitSource.C:539
turbulentFluidThermoModel.H
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54