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-2020 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& phase,
38  const word& type
39 )
40 :
42  (
43  type,
44  alpha,
45  rho,
46  U,
47  alphaRhoPhi,
48  phi,
49  phase
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<phaseCompressibleMomentumTransportModel>>
84  ::read())
85  {
86  coeffDict().lookup("alphaMax") >> alphaMax_;
87  coeffDict().lookup("preAlphaExp") >> preAlphaExp_;
88  coeffDict().lookup("expMax") >> expMax_;
89  g0_.readIfPresent(coeffDict());
90 
91  return true;
92  }
93  else
94  {
95  return false;
96  }
97 }
98 
99 
102 {
104  return nut_;
105 }
106 
107 
110 {
112  return nut_;
113 }
114 
115 
118 {
119  return tmp<volSymmTensorField>
120  (
122  (
123  IOobject
124  (
125  IOobject::groupName("R", U_.group()),
126  runTime_.timeName(),
127  mesh_,
130  ),
131  mesh_,
132  dimensioned<symmTensor>
133  (
134  "R",
135  dimensionSet(0, 2, -2, 0, 0),
136  Zero
137  )
138  )
139  );
140 }
141 
142 
145 {
146  tmp<volScalarField> tpPrime
147  (
148  g0_
149  *min
150  (
151  exp(preAlphaExp_*(alpha_ - alphaMax_)),
152  expMax_
153  )
154  );
155 
156  volScalarField::Boundary& bpPrime =
157  tpPrime.ref().boundaryFieldRef();
158 
159  forAll(bpPrime, patchi)
160  {
161  if (!bpPrime[patchi].coupled())
162  {
163  bpPrime[patchi] == 0;
164  }
165  }
166 
167  return tpPrime;
168 }
169 
170 
173 {
174  tmp<surfaceScalarField> tpPrime
175  (
176  g0_
177  *min
178  (
179  exp(preAlphaExp_*(fvc::interpolate(alpha_) - alphaMax_)),
180  expMax_
181  )
182  );
183 
184  surfaceScalarField::Boundary& bpPrime =
185  tpPrime.ref().boundaryFieldRef();
186 
187  forAll(bpPrime, patchi)
188  {
189  if (!bpPrime[patchi].coupled())
190  {
191  bpPrime[patchi] == 0;
192  }
193  }
194 
195  return tpPrime;
196 }
197 
198 
201 {
202  return tmp<volSymmTensorField>
203  (
205  (
206  IOobject
207  (
208  IOobject::groupName("devTau", U_.group()),
209  runTime_.timeName(),
210  mesh_,
213  ),
214  mesh_,
215  dimensioned<symmTensor>
216  (
217  "R",
218  rho_.dimensions()*dimensionSet(0, 2, -2, 0, 0),
219  Zero
220  )
221  )
222  );
223 }
224 
225 
228 (
229  volVectorField& U
230 ) const
231 {
232  return tmp<fvVectorMatrix>
233  (
234  new fvVectorMatrix
235  (
236  U,
237  rho_.dimensions()*dimensionSet(0, 4, -2, 0, 0)
238  )
239  );
240 }
241 
242 
244 {}
245 
246 
247 // ************************************************************************* //
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:61
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual tmp< volSymmTensorField > devTau() const
Return the effective stress tensor.
PhaseCompressibleMomentumTransportModel< phaseModel > phaseCompressibleMomentumTransportModel
Typedef for phaseCompressibleMomentumTransportModel.
phasePressureModel(const volScalarField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const phaseModel &transport, const word &type=typeName)
Construct from components.
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.
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
stressControl lookup("compactNormalStress") >> compactNormalStress
dimensionedScalar exp(const dimensionedScalar &ds)
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:366
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;.