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-2018 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 #include "twoPhaseSystem.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
32 (
33  const volScalarField& alpha,
34  const volScalarField& rho,
35  const volVectorField& U,
36  const surfaceScalarField& alphaRhoPhi,
37  const surfaceScalarField& phi,
38  const transportModel& phase,
39  const word& propertiesName,
40  const word& type
41 )
42 :
43  eddyViscosity
44  <
45  RASModel<EddyDiffusivity<ThermalDiffusivity
46  <
47  PhaseCompressibleTurbulenceModel<phaseModel>
48  >>>
49  >
50  (
51  type,
52  alpha,
53  rho,
54  U,
55  alphaRhoPhi,
56  phi,
57  phase,
58  propertiesName
59  ),
60 
61  phase_(phase),
62 
63  alphaMax_(readScalar(coeffDict_.lookup("alphaMax"))),
64  preAlphaExp_(readScalar(coeffDict_.lookup("preAlphaExp"))),
65  expMax_(readScalar(coeffDict_.lookup("expMax"))),
66  g0_
67  (
68  "g0",
69  dimensionSet(1, -1, -2, 0, 0),
70  coeffDict_.lookup("g0")
71  )
72 {
73  nut_ == dimensionedScalar(nut_.dimensions(), 0);
74 
75  if (type == typeName)
76  {
77  printCoeffs(type);
78  }
79 }
80 
81 
82 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
83 
85 {}
86 
87 
88 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
89 
91 {
92  if
93  (
94  eddyViscosity
95  <
96  RASModel<EddyDiffusivity<ThermalDiffusivity
97  <
98  PhaseCompressibleTurbulenceModel<phaseModel>
99  >>>
100  >::read()
101  )
102  {
103  coeffDict().lookup("alphaMax") >> alphaMax_;
104  coeffDict().lookup("preAlphaExp") >> preAlphaExp_;
105  coeffDict().lookup("expMax") >> expMax_;
106  g0_.readIfPresent(coeffDict());
107 
108  return true;
109  }
110  else
111  {
112  return false;
113  }
114 }
115 
116 
119 {
121  return nut_;
122 }
123 
124 
127 {
129  return nut_;
130 }
131 
132 
135 {
137  (
138  IOobject::groupName("R", U_.group()),
139  mesh_,
140  dimensioned<symmTensor>(dimensionSet(0, 2, -2, 0, 0), Zero)
141  );
142 }
143 
144 
147 {
148  tmp<volScalarField> tpPrime
149  (
150  g0_
151  *min
152  (
153  exp(preAlphaExp_*(alpha_ - alphaMax_)),
154  expMax_
155  )
156  );
157 
158  volScalarField::Boundary& bpPrime =
159  tpPrime.ref().boundaryFieldRef();
160 
161  forAll(bpPrime, patchi)
162  {
163  if (!bpPrime[patchi].coupled())
164  {
165  bpPrime[patchi] == 0;
166  }
167  }
168 
169  return tpPrime;
170 }
171 
172 
175 {
176  tmp<surfaceScalarField> tpPrime
177  (
178  g0_
179  *min
180  (
181  exp(preAlphaExp_*(fvc::interpolate(alpha_) - alphaMax_)),
182  expMax_
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 {
205  (
206  IOobject::groupName("devRhoReff", U_.group()),
207  mesh_,
208  dimensioned<symmTensor>
209  (
210  rho_.dimensions()*dimensionSet(0, 2, -2, 0, 0),
211  Zero
212  )
213  );
214 }
215 
216 
219 (
220  volVectorField& U
221 ) const
222 {
223  return tmp<fvVectorMatrix>
224  (
225  new fvVectorMatrix
226  (
227  U,
228  rho_.dimensions()*dimensionSet(0, 4, -2, 0, 0)
229  )
230  );
231 }
232 
233 
235 {}
236 
237 
238 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
static tmp< GeometricField< symmTensor, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< symmTensor >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure&#39;.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
stressControl lookup("compactNormalStress") >> compactNormalStress
dimensionedScalar exp(const dimensionedScalar &ds)
static word groupName(Name name, const word &group)
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
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.
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
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 tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
virtual ~phasePressureModel()
Destructor.
label patchi
phasePressureModel(const volScalarField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const phaseModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
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
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure&#39;.