phasePressureModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2013-2015 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("zero", nut_.dimensions(), 0.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 {
136  return tmp<volSymmTensorField>
137  (
139  (
140  IOobject
141  (
142  IOobject::groupName("R", U_.group()),
143  runTime_.timeName(),
144  mesh_,
147  ),
148  mesh_,
149  dimensioned<symmTensor>
150  (
151  "R",
152  dimensionSet(0, 2, -2, 0, 0),
154  )
155  )
156  );
157 }
158 
159 
162 {
163  tmp<volScalarField> tpPrime
164  (
165  g0_
166  *min
167  (
168  exp(preAlphaExp_*(alpha_ - alphaMax_)),
169  expMax_
170  )
171  );
172 
173  volScalarField::GeometricBoundaryField& bpPrime = tpPrime().boundaryField();
174 
175  forAll(bpPrime, patchi)
176  {
177  if (!bpPrime[patchi].coupled())
178  {
179  bpPrime[patchi] == 0;
180  }
181  }
182 
183  return tpPrime;
184 }
185 
186 
189 {
190  tmp<surfaceScalarField> tpPrime
191  (
192  g0_
193  *min
194  (
195  exp(preAlphaExp_*(fvc::interpolate(alpha_) - alphaMax_)),
196  expMax_
197  )
198  );
199 
200  surfaceScalarField::GeometricBoundaryField& bpPrime =
201  tpPrime().boundaryField();
202 
203  forAll(bpPrime, patchi)
204  {
205  if (!bpPrime[patchi].coupled())
206  {
207  bpPrime[patchi] == 0;
208  }
209  }
210 
211  return tpPrime;
212 }
213 
214 
217 {
218  return tmp<volSymmTensorField>
219  (
221  (
222  IOobject
223  (
224  IOobject::groupName("devRhoReff", U_.group()),
225  runTime_.timeName(),
226  mesh_,
229  ),
230  mesh_,
231  dimensioned<symmTensor>
232  (
233  "R",
234  rho_.dimensions()*dimensionSet(0, 2, -2, 0, 0),
236  )
237  )
238  );
239 }
240 
241 
244 (
246 ) const
247 {
248  return tmp<fvVectorMatrix>
249  (
250  new fvVectorMatrix
251  (
252  U,
253  rho_.dimensions()*dimensionSet(0, 4, -2, 0, 0)
254  )
255  );
256 }
257 
258 
260 {}
261 
262 
263 // ************************************************************************* //
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
Foam::RASModels::phasePressureModel::divDevRhoReff
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::fvc::interpolate
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Definition: surfaceInterpolate.C:110
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::RASModels::phasePressureModel::phasePressureModel
phasePressureModel(const phasePressureModel &)
Disallow default bitwise copy construct.
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::RASModels::phasePressureModel::k
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Foam::RASModels::phasePressureModel::~phasePressureModel
virtual ~phasePressureModel()
Destructor.
U
U
Definition: pEqn.H:46
Foam::RASModels::phasePressureModel::pPrime
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure'.
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:365
Foam::fvVectorMatrix
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
Foam::RASModels::phasePressureModel::pPrimef
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure'.
Foam::volSymmTensorField
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::SymmTensor< scalar >::zero
static const SymmTensor zero
Definition: SymmTensor.H:77
Foam::RASModels::phasePressureModel::read
virtual bool read()
Re-read model coefficients if they have changed.
readScalar
#define readScalar
Definition: doubleScalar.C:38
rho
rho
Definition: pEqn.H:3
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Foam::RASModels::phasePressureModel::devRhoReff
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
Foam::RASModels::phasePressureModel::epsilon
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
RASModel
Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);singlePhaseTransportModel laminarTransport(U, phi);autoPtr< incompressible::RASModel > RASModel(incompressible::New< incompressible::RASModel >(U, phi, laminarTransport))
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::RASModels::phasePressureModel::correct
virtual void correct()
Solve the kinetic theory equations and correct the viscosity.
Foam::RASModels::phasePressureModel::R
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress