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-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 
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 transportModel& transport,
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  transport
50  ),
51 
52  alphaMax_(coeffDict_.lookup<scalar>("alphaMax")),
53  preAlphaExp_(coeffDict_.lookup<scalar>("preAlphaExp")),
54  expMax_(coeffDict_.lookup<scalar>("expMax")),
55  g0_
56  (
57  "g0",
58  dimensionSet(1, -1, -2, 0, 0),
59  coeffDict_.lookup("g0")
60  )
61 {
62  nut_ == dimensionedScalar(nut_.dimensions(), 0);
63 
64  if (type == typeName)
65  {
66  printCoeffs(type);
67  }
68 }
69 
70 
71 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
72 
74 {}
75 
76 
77 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78 
80 {
81  if
82  (
83  eddyViscosity<RASModel<phaseCompressible::momentumTransportModel>>::
84  read()
85  )
86  {
87  coeffDict().lookup("alphaMax") >> alphaMax_;
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 {
120  return tmp<volSymmTensorField>
121  (
123  (
124  IOobject::groupName("R", U_.group()),
125  mesh_,
126  dimensioned<symmTensor>
127  (
128  "R",
129  dimensionSet(0, 2, -2, 0, 0),
130  Zero
131  )
132  )
133  );
134 }
135 
136 
139 {
140  tmp<volScalarField> tpPrime
141  (
143  (
144  IOobject::groupName("pPrime", U_.group()),
145  g0_
146  *min
147  (
148  exp(preAlphaExp_*(alpha_ - alphaMax_)),
149  expMax_
150  )
151  )
152  );
153 
154  volScalarField::Boundary& bpPrime =
155  tpPrime.ref().boundaryFieldRef();
156 
157  forAll(bpPrime, patchi)
158  {
159  if (!bpPrime[patchi].coupled())
160  {
161  bpPrime[patchi] == 0;
162  }
163  }
164 
165  return tpPrime;
166 }
167 
168 
171 {
172  tmp<surfaceScalarField> tpPrime
173  (
175  (
176  IOobject::groupName("pPrimef", U_.group()),
177  g0_
178  *min
179  (
180  exp(preAlphaExp_*(fvc::interpolate(alpha_) - alphaMax_)),
181  expMax_
182  )
183  )
184  );
185 
186  surfaceScalarField::Boundary& bpPrime =
187  tpPrime.ref().boundaryFieldRef();
188 
189  forAll(bpPrime, patchi)
190  {
191  if (!bpPrime[patchi].coupled())
192  {
193  bpPrime[patchi] == 0;
194  }
195  }
196 
197  return tpPrime;
198 }
199 
200 
203 {
204  return tmp<volSymmTensorField>
205  (
207  (
208  IOobject::groupName("devTau", U_.group()),
209  mesh_,
210  dimensioned<symmTensor>
211  (
212  "R",
213  rho_.dimensions()*dimensionSet(0, 2, -2, 0, 0),
214  Zero
215  )
216  )
217  );
218 }
219 
220 
223 (
224  volVectorField& U
225 ) const
226 {
227  return tmp<fvVectorMatrix>
228  (
229  new fvVectorMatrix
230  (
231  U,
232  rho_.dimensions()*dimensionSet(0, 4, -2, 0, 0)
233  )
234  );
235 }
236 
237 
239 {}
240 
241 
242 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual tmp< volSymmTensorField > devTau() const
Return the effective stress tensor.
virtual tmp< volSymmTensorField > sigma() const
Return the stress tensor [m^2/s^2].
virtual tmp< fvVectorMatrix > divDevTau(volVectorField &U) const
Return the source term for the momentum equation.
static tmp< GeometricField< symmTensor, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< symmTensor >> &)
Return a temporary field constructed from name,.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure&#39;.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
CompressibleMomentumTransportModel< dynamicTransportModel > momentumTransportModel
stressControl lookup("compactNormalStress") >> compactNormalStress
dimensionedScalar exp(const dimensionedScalar &ds)
phasePressureModel(const volScalarField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
static word groupName(Name name, const word &group)
static const zero Zero
Definition: zero.H:97
virtual void correct()
Solve the kinetic theory equations and correct the viscosity.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual bool read()
Re-read model coefficients if they have changed.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
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.
virtual ~phasePressureModel()
Destructor.
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:370
RASModel< momentumTransportModel > RASModel
Typedefs for turbulence, RAS and LES models for compressible flow based on the standard laminar trans...
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure&#39;.