phasePressureModel.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) 2013-2023 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 
26 #include "phasePressureModel.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
31 (
32  const volScalarField& alpha,
33  const volScalarField& rho,
34  const volVectorField& U,
35  const surfaceScalarField& alphaRhoPhi,
36  const surfaceScalarField& phi,
37  const viscosity& viscosity,
38  const word& type
39 )
40 :
41  eddyViscosity<RASModel<phaseCompressible::momentumTransportModel>>
42  (
43  type,
44  alpha,
45  rho,
46  U,
47  alphaRhoPhi,
48  phi,
49  viscosity
50  ),
51 
52  phase_(refCast<const phaseModel>(viscosity)),
53 
54  preAlphaExp_(coeffDict_.lookup<scalar>("preAlphaExp")),
55  expMax_(coeffDict_.lookup<scalar>("expMax")),
56  g0_
57  (
58  "g0",
59  dimensionSet(1, -1, -2, 0, 0),
60  coeffDict_.lookup("g0")
61  )
62 {
64 
65  if (type == typeName)
66  {
67  printCoeffs(type);
68  }
69 }
70 
71 
72 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
73 
75 {}
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
81 {
82  if
83  (
85  read()
86  )
87  {
88  coeffDict().lookup("preAlphaExp") >> preAlphaExp_;
89  coeffDict().lookup("expMax") >> expMax_;
90  g0_.readIfPresent(coeffDict());
91 
92  return true;
93  }
94  else
95  {
96  return false;
97  }
98 }
99 
100 
103 {
105  return nut_;
106 }
107 
108 
111 {
113  return nut_;
114 }
115 
116 
119 {
121  return nut_;
122 }
123 
124 
127 {
129  (
131  (
132  IOobject::groupName("R", U_.group()),
133  mesh_,
135  (
136  "R",
137  dimensionSet(0, 2, -2, 0, 0),
138  Zero
139  )
140  )
141  );
142 }
143 
144 
147 {
148  tmp<volScalarField> tpPrime
149  (
151  (
152  IOobject::groupName("pPrime", U_.group()),
153  g0_
154  *min
155  (
156  exp(preAlphaExp_*(alpha_ - phase_.alphaMax())),
157  expMax_
158  )
159  )
160  );
161 
162  volScalarField::Boundary& bpPrime =
163  tpPrime.ref().boundaryFieldRef();
164 
165  forAll(bpPrime, patchi)
166  {
167  if (!bpPrime[patchi].coupled())
168  {
169  bpPrime[patchi] == 0;
170  }
171  }
172 
173  return tpPrime;
174 }
175 
176 
179 {
181  (
183  (
184  IOobject::groupName("pPrimef", U_.group()),
185  g0_
186  *min
187  (
188  exp(preAlphaExp_
189  *(fvc::interpolate(alpha_) - phase_.alphaMax())),
190  expMax_
191  )
192  )
193  );
194 
196  tpPrime.ref().boundaryFieldRef();
197 
198  forAll(bpPrime, patchi)
199  {
200  if (!bpPrime[patchi].coupled())
201  {
202  bpPrime[patchi] == 0;
203  }
204  }
205 
206  return tpPrime;
207 }
208 
209 
212 {
214  (
216  (
217  IOobject::groupName("devTau", U_.group()),
218  mesh_,
220  (
221  "R",
222  rho_.dimensions()*dimensionSet(0, 2, -2, 0, 0),
223  Zero
224  )
225  )
226  );
227 }
228 
229 
232 (
234 ) const
235 {
236  return tmp<fvVectorMatrix>
237  (
238  new fvVectorMatrix
239  (
240  U,
241  rho_.dimensions()*dimensionSet(0, 4, -2, 0, 0)
242  )
243  );
244 }
245 
246 
248 {}
249 
250 
251 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricBoundaryField class.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
static word groupName(Name name, const word &group)
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:56
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
virtual void correct()
Solve the kinetic theory equations and correct the viscosity.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual tmp< volSymmTensorField > devTau() const
Return the effective stress tensor.
phasePressureModel(const volScalarField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure'.
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure'.
virtual tmp< fvVectorMatrix > divDevTau(volVectorField &U) const
Return the source term for the momentum equation.
virtual tmp< volSymmTensorField > sigma() const
Return the stress tensor [m^2/s^2].
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
virtual bool read()
Re-read model coefficients if they have changed.
Dimension set for the base types.
Definition: dimensionSet.H:122
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Abstract base class for turbulence models (RAS, LES and laminar).
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:353
label patchi
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
static const zero Zero
Definition: zero.H:97
dimensionedScalar exp(const dimensionedScalar &ds)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:111
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488