hydrostaticPressure.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) 2018-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 
28 #include "hydrostaticPressure.H"
29 #include "basicThermo.H"
31 #include "volFields.H"
32 #include "surfaceInterpolate.H"
33 #include "fvcDiv.H"
34 #include "fvmLaplacian.H"
35 #include "fvcSnGrad.H"
36 #include "constrainPressure.H"
37 
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 namespace functionObjects
45 {
46  defineTypeNameAndDebug(hydrostaticPressure, 0);
47 
49  (
50  functionObject,
51  hydrostaticPressure,
52  dictionary
53  );
54 }
55 }
56 
57 
60 {
61  if (pRefName_ == "none")
62  {
64  }
65  else if (pRefName_ == "pInf")
66  {
67  return dimensionedScalar("pRef", dimPressure, pRefValue_);
68  }
69  else
70  {
72  }
73 }
74 
75 
77 {
78  const auto& pRef = this->pRef();
79  const auto& U = mesh_.lookupObject<volVectorField>(UName_);
80  const auto& gh = mesh_.lookupObject<volScalarField>(ghName_);
81  const auto& ghf = mesh_.lookupObject<surfaceScalarField>(ghfName_);
82  auto& rho = mesh_.lookupObjectRef<volScalarField>(rhoName_);
83  auto& thermo = mesh_.lookupObjectRef<basicThermo>(basicThermo::dictName);
84  auto& p_rgh = mesh_.lookupObjectRef<volScalarField>(p_rghName_);
85  auto& ph_rgh = mesh_.lookupObjectRef<volScalarField>(ph_rghName_);
86 
87  auto& p = thermo.p();
88 
89  Info<< "Performing hydrostatic pressure initialisation";
90  if (mesh_.name() != polyMesh::defaultRegion)
91  {
92  Info<< "for region " << mesh_.name();
93  }
94 
95 
96  if (thermo.incompressible())
97  {
98  Info<< ": incompressible" << endl;
99 
100  // Constant density and temperature
101 
102  thermo.correct();
103  rho = thermo.rho();
104  p = ph_rgh + rho*gh + pRef;
105  p_rgh = ph_rgh;
106  }
107  else
108  {
109  Info<< ": compressible" << endl;
110 
111  p = ph_rgh + rho*gh + pRef;
112  thermo.correct();
113  rho = thermo.rho();
114 
115  for (label i=0; i<nCorrectors_; ++i)
116  {
118 
120  (
121  "phig",
122  -rhof*ghf*fvc::snGrad(rho)*mesh_.magSf()
123  );
124 
125  // Update the pressure BCs to ensure flux consistency
126  constrainPressure(ph_rgh, rho, U, phig, rhof);
127 
128  fvScalarMatrix ph_rghEqn
129  (
130  fvm::laplacian(rhof, ph_rgh) == fvc::div(phig)
131  );
132 
133  ph_rghEqn.relax();
134 
135  ph_rghEqn.solve();
136 
137  p = ph_rgh + rho*gh + pRef;
138  thermo.correct();
139  rho = thermo.rho();
140 
141  Info<< "Hydrostatic pressure variation "
142  << (max(ph_rgh) - min(ph_rgh)).value() << endl;
143  }
144 
145  p_rgh = ph_rgh;
146 
147  Log << " writing field " << ph_rgh.name() << nl;
148  ph_rgh.write();
149  }
150 
151  Log << " writing field " << rho.name() << nl;
152  rho.write();
153 
154  Log << " writing field " << p_rgh.name() << nl;
155  p_rgh.write();
156 
157  Log << " writing field " << p.name() << nl;
158  p.write();
159 }
160 
161 
162 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
163 
165 (
166  const word& name,
167  const Time& runTime,
168  const dictionary& dict
169 )
170 :
172  p_rghName_("p_rgh"),
173  ph_rghName_("ph_rgh"),
174  pRefName_("pRef"),
175  pRefValue_(0),
176  rhoName_("rho"),
177  UName_("U"),
178  ghName_("gh"),
179  ghfName_("ghf"),
180  nCorrectors_(5)
181 {
182  if (read(dict))
183  {
184  // Read and store the initial ph_rgh field
185  volScalarField* ph_rghPtr =
186  new volScalarField
187  (
188  IOobject
189  (
190  ph_rghName_,
191  runTime.timeName(),
192  mesh_,
194  IOobject::AUTO_WRITE // To enable restart
195  ),
196  mesh_
197  );
198 
199  mesh_.objectRegistry::store(ph_rghPtr);
200 
201  bool reInitialise = dict.getOrDefault("reInitialise", false);
202 
203  if (runTime.timeIndex() == 0 || reInitialise)
204  {
205  calculateAndWrite();
206  }
207  }
208 }
209 
210 
211 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
212 
214 {
216  {
217  dict.readIfPresent("p_rgh", p_rghName_);
218  dict.readIfPresent("ph_rgh", ph_rghName_);
219  dict.readIfPresent("pRef", pRefName_);
220  dict.readIfPresent("rho", rhoName_);
221  dict.readIfPresent("U", UName_);
222  dict.readIfPresent("gh", ghName_);
223  dict.readIfPresent("ghf", ghfName_);
224  dict.readIfPresent("nCorrectors", nCorrectors_);
225 
226  pRefValue_ = 0;
227  if (pRefName_ == "pInf")
228  {
229  pRefValue_ = dict.get<scalar>("pRefValue");
230  }
231 
232  return true;
233  }
234 
235  return false;
236 }
237 
238 
240 {
241  return true;
242 }
243 
244 
246 {
247  return true;
248 }
249 
250 
251 // ************************************************************************* //
Foam::constrainPressure
void constrainPressure(volScalarField &p, const RhoType &rho, const volVectorField &U, const surfaceScalarField &phiHbyA, const RAUType &rhorAU, const MRFType &MRF)
Definition: constrainPressure.C:39
volFields.H
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::dimPressure
const dimensionSet dimPressure
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:47
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
basicThermo.H
Log
#define Log
Definition: PDRblock.C:35
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
ghf
const surfaceScalarField & ghf
Definition: setRegionFluidFields.H:18
Foam::functionObjects::hydrostaticPressure::hydrostaticPressure
hydrostaticPressure(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: hydrostaticPressure.C:165
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::polyMesh::defaultRegion
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:318
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::functionObjects::hydrostaticPressure::calculateAndWrite
void calculateAndWrite()
Calculate the fields and write.
Definition: hydrostaticPressure.C:76
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
fvcSnGrad.H
Calculate the snGrad of the given volField.
constrainPressure.H
fvcDiv.H
Calculate the divergence of the given field.
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
Foam::basicThermo
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:63
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
gh
const volScalarField & gh
Definition: setRegionFluidFields.H:17
Foam::functionObjects::fvMeshFunctionObject
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Definition: fvMeshFunctionObject.H:64
rho
rho
Definition: readInitialConditions.H:88
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
hydrostaticPressure.H
Foam::functionObjects::hydrostaticPressure::read
virtual bool read(const dictionary &dict)
Read the hydrostaticPressure data.
Definition: hydrostaticPressure.C:213
Foam::UniformDimensionedField< scalar >
Foam::functionObjects::hydrostaticPressure::pRef
dimensionedScalar pRef() const
Helper function to return the reference pressure.
Definition: hydrostaticPressure.C:59
p_rgh
volScalarField & p_rgh
Definition: setRegionFluidFields.H:15
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
Foam::fvMatrix::relax
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1183
fvmLaplacian.H
Calculate the matrix for the laplacian of the field.
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::functionObjects::hydrostaticPressure::write
virtual bool write()
Write the p_rgh and derived fields.
Definition: hydrostaticPressure.C:245
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::functionObjects::regionFunctionObject::read
virtual bool read(const dictionary &dict)
Read optional controls.
Definition: regionFunctionObject.C:173
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::fvc::interpolate
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.
Foam::functionObjects::hydrostaticPressure::pRefName_
word pRefName_
Name of uniform pressure reference field, default is "pRef".
Definition: hydrostaticPressure.H:158
uniformDimensionedFields.H
U
U
Definition: pEqn.H:72
phig
surfaceScalarField phig("phig", -rhorAUf *ghf *fvc::snGrad(rho) *mesh.magSf())
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Foam::functionObjects::hydrostaticPressure::pRefValue_
scalar pRefValue_
Reference pressure if pRefName is set to "pInf".
Definition: hydrostaticPressure.H:161
Foam::functionObjects::fvMeshFunctionObject::mesh_
const fvMesh & mesh_
Reference to the fvMesh.
Definition: fvMeshFunctionObject.H:73
Foam::functionObjects::hydrostaticPressure::execute
virtual bool execute()
Calculate the p_rgh field.
Definition: hydrostaticPressure.C:239
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(ObukhovLength, 0)
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:319
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
rhof
surfaceScalarField rhof(fvc::interpolate(rho, "div(phi,rho)"))
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
Foam::dictionary::dictName
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:60