hydrostaticInitialisation.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 
28 #include "fluidThermo.H"
29 #include "fvmLaplacian.H"
30 #include "fvcDiv.H"
31 #include "fvcSnGrad.H"
32 #include "surfaceInterpolate.H"
33 #include "constrainPressure.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
39 (
40  volScalarField& p_rgh,
43  const volVectorField& U,
44  const volScalarField& gh,
45  const surfaceScalarField& ghf,
48  const dictionary& dict
49 )
50 {
51  if (dict.lookupOrDefault<bool>("hydrostaticInitialisation", false))
52  {
53  const fvMesh& mesh = p_rgh.mesh();
54 
55  volScalarField& ph_rgh = regIOobject::store
56  (
57  new volScalarField
58  (
59  IOobject
60  (
61  "ph_rgh",
62  "0",
63  mesh,
64  IOobject::MUST_READ,
65  IOobject::NO_WRITE
66  ),
67  mesh
68  )
69  );
70 
71  if (!mesh.time().restart())
72  {
73  p = ph_rgh + rho*gh + pRef;
74  thermo.correct();
75  rho = thermo.rho();
76 
77  label nCorr
78  (
79  dict.lookupOrDefault<label>("nHydrostaticCorrectors", 5)
80  );
81 
82  for (label i=0; i<nCorr; i++)
83  {
85 
87  (
88  "phig",
89  -rhof*ghf*fvc::snGrad(rho)*mesh.magSf()
90  );
91 
92  // Update the pressure BCs to ensure flux consistency
93  constrainPressure(ph_rgh, rho, U, phig, rhof);
94 
95  fvScalarMatrix ph_rghEqn
96  (
97  fvm::laplacian(rhof, ph_rgh) == fvc::div(phig)
98  );
99 
100  ph_rghEqn.solve();
101 
102  p = ph_rgh + rho*gh + pRef;
103  thermo.correct();
104  rho = thermo.rho();
105 
106  Info<< "Hydrostatic pressure variation "
107  << (max(ph_rgh) - min(ph_rgh)).value() << endl;
108  }
109 
110  ph_rgh.write();
111 
112  p_rgh = ph_rgh;
113  }
114  }
115  else
116  {
117  // Force p_rgh to be consistent with p
118  p_rgh = p - rho*gh - pRef;
119  }
120 }
121 
122 
123 // ************************************************************************* //
const Mesh & mesh() const
Return mesh.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
bool restart() const
Return true if the run is a restart, i.e. startTime != beginTime.
Definition: Time.H:428
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:406
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
virtual bool write(const bool write=true) const
Write using setting from DB.
Calculate the divergence of the given field.
Calculate the snGrad of the given volField.
Calculate the matrix for the laplacian of the field.
U
Definition: pEqn.H:72
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvcLaplacian.C:45
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
void constrainPressure(volScalarField &p, const RhoType &rho, const volVectorField &U, const surfaceScalarField &phiHbyA, const RAUType &rhorAU, const MRFType &MRF)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
void hydrostaticInitialisation(volScalarField &p_rgh, volScalarField &p, volScalarField &rho, const volVectorField &U, const volScalarField &gh, const surfaceScalarField &ghf, const uniformDimensionedScalarField &pRef, fluidThermo &thermo, const dictionary &dict)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dictionary dict
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31