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-2024 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 
66 
67 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
68 
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
76 {
77  if
78  (
80  read()
81  )
82  {
83  coeffDict().lookup("preAlphaExp") >> preAlphaExp_;
84  coeffDict().lookup("expMax") >> expMax_;
85  g0_.readIfPresent(coeffDict());
86 
87  return true;
88  }
89  else
90  {
91  return false;
92  }
93 }
94 
95 
98 {
100  return nut_;
101 }
102 
103 
106 {
108  return nut_;
109 }
110 
111 
114 {
116  return nut_;
117 }
118 
119 
122 {
124  (
126  (
127  IOobject::groupName("R", U_.group()),
128  mesh_,
130  (
131  "R",
132  dimensionSet(0, 2, -2, 0, 0),
133  Zero
134  )
135  )
136  );
137 }
138 
139 
141 Foam::RASModels::phasePressureModel::pPrime() const
142 {
143  tmp<volScalarField> tpPrime
144  (
146  (
147  IOobject::groupName("pPrime", U_.group()),
148  g0_
149  *min
150  (
151  exp(preAlphaExp_*(alpha_ - phase_.alphaMax())),
152  expMax_
153  )
154  )
155  );
156 
157  volScalarField::Boundary& bpPrime =
158  tpPrime.ref().boundaryFieldRef();
159 
160  forAll(bpPrime, patchi)
161  {
162  if (!bpPrime[patchi].coupled())
163  {
164  bpPrime[patchi] == 0;
165  }
166  }
167 
168  return tpPrime;
169 }
170 
171 
172 // Foam::tmp<Foam::surfaceScalarField>
173 // Foam::RASModels::phasePressureModel::pPrimef() const
174 // {
175 // tmp<surfaceScalarField> tpPrime
176 // (
177 // surfaceScalarField::New
178 // (
179 // IOobject::groupName("pPrimef", U_.group()),
180 // g0_
181 // *min
182 // (
183 // exp(preAlphaExp_
184 // *(fvc::interpolate(alpha_, "hmm") - phase_.alphaMax())),
185 // expMax_
186 // )
187 // )
188 // );
189 
190 // surfaceScalarField::Boundary& bpPrime =
191 // tpPrime.ref().boundaryFieldRef();
192 
193 // forAll(bpPrime, patchi)
194 // {
195 // if (!bpPrime[patchi].coupled())
196 // {
197 // bpPrime[patchi] == 0;
198 // }
199 // }
200 
201 // return tpPrime;
202 // }
203 
204 
207 {
209  (
210  IOobject::groupName("pPrimef", U_.group()),
211  fvc::interpolate(pPrime())
212  );
213 }
214 
215 
218 {
220  (
222  (
223  IOobject::groupName("devTau", U_.group()),
224  mesh_,
226  (
227  "devTau",
228  rho_.dimensions()*dimensionSet(0, 2, -2, 0, 0),
229  Zero
230  )
231  )
232  );
233 }
234 
235 
238 (
240 ) const
241 {
242  return tmp<fvVectorMatrix>
243  (
244  new fvVectorMatrix
245  (
246  U,
247  rho_.dimensions()*dimensionSet(0, 4, -2, 0, 0)
248  )
249  );
250 }
251 
252 
254 {}
255 
256 
257 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricBoundaryField class.
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
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 tmp< volSymmTensorField > R() const
Return the stress tensor (0) [m^2/s^2].
virtual tmp< surfaceVectorField > devTau() const
Return the effective stress.
virtual void correct()
Solve the kinetic theory equations and correct the viscosity.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
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< 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< 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:125
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
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:381
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.
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:134
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488