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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2018 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "phasePressureModel.H"
30 #include "twoPhaseSystem.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 Foam::RASModels::phasePressureModel::phasePressureModel
35 (
36  const volScalarField& alpha,
37  const volScalarField& rho,
38  const volVectorField& U,
39  const surfaceScalarField& alphaRhoPhi,
40  const surfaceScalarField& phi,
41  const transportModel& phase,
42  const word& propertiesName,
43  const word& type
44 )
45 :
46  eddyViscosity
47  <
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_(coeffDict_.get<scalar>("alphaMax")),
64  preAlphaExp_(coeffDict_.get<scalar>("preAlphaExp")),
65  expMax_(coeffDict_.get<scalar>("expMax")),
66  g0_("g0", dimPressure, coeffDict_)
67 {
68  nut_ == dimensionedScalar(nut_.dimensions());
69 
70  if (type == typeName)
71  {
73  }
74 }
75 
76 
77 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78 
80 {
81  if
82  (
84  <
86  >::read()
87  )
88  {
89  coeffDict().readEntry("alphaMax", alphaMax_);
90  coeffDict().readEntry("preAlphaExp", preAlphaExp_);
91  coeffDict().readEntry("expMax", expMax_);
92  g0_.readIfPresent(coeffDict());
93 
94  return true;
95  }
96 
97  return false;
98 }
99 
100 
103 {
105  return nut_;
106 }
107 
108 
111 {
113  return nut_;
114 }
115 
116 
119 {
121  return nullptr;
122 }
123 
124 
127 {
129  (
131  (
132  IOobject
133  (
134  IOobject::groupName("R", U_.group()),
135  runTime_.timeName(),
136  mesh_,
139  ),
140  mesh_,
141  dimensioned<symmTensor>(dimensionSet(0, 2, -2, 0, 0), Zero)
142  )
143  );
144 }
145 
146 
149 {
150  tmp<volScalarField> tpPrime
151  (
152  g0_
153  *min
154  (
155  exp(preAlphaExp_*(alpha_ - alphaMax_)),
156  expMax_
157  )
158  );
159 
160  volScalarField::Boundary& bpPrime =
161  tpPrime.ref().boundaryFieldRef();
162 
163  forAll(bpPrime, patchi)
164  {
165  if (!bpPrime[patchi].coupled())
166  {
167  bpPrime[patchi] == 0;
168  }
169  }
170 
171  return tpPrime;
172 }
173 
174 
177 {
179  (
180  g0_
181  *min
182  (
183  exp(preAlphaExp_*(fvc::interpolate(alpha_) - alphaMax_)),
184  expMax_
185  )
186  );
187 
188  surfaceScalarField::Boundary& bpPrime =
189  tpPrime.ref().boundaryFieldRef();
190 
191  forAll(bpPrime, patchi)
192  {
193  if (!bpPrime[patchi].coupled())
194  {
195  bpPrime[patchi] == 0;
196  }
197  }
198 
199  return tpPrime;
200 }
201 
202 
205 {
207  (
209  (
210  IOobject
211  (
212  IOobject::groupName("devRhoReff", U_.group()),
213  runTime_.timeName(),
214  mesh_,
217  ),
218  mesh_,
220  (
221  "R",
222  rho_.dimensions()*dimensionSet(0, 2, -2, 0, 0),
223  Zero
224  )
225  )
226  );
227 }
228 
229 
232 (
234 ) const
235 {
236  return tmp<fvVectorMatrix>
237  (
238  new fvVectorMatrix
239  (
240  U,
241  rho_.dimensions()*dimensionSet(0, 4, -2, 0, 0)
242  )
243  );
244 }
245 
246 
248 {}
249 
250 
251 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::dimPressure
const dimensionSet dimPressure
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::RASModel
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:51
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:107
Foam::RASModels::phasePressureModel::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: phasePressureModel.C:79
Foam::phaseCompressibleTurbulenceModel
ThermalDiffusivity< PhaseCompressibleTurbulenceModel< phaseModel > > phaseCompressibleTurbulenceModel
Typedef for phaseCompressibleTurbulenceModel.
Definition: phaseCompressibleTurbulenceModel.H:47
Foam::RASModels::phasePressureModel::epsilon
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: phasePressureModel.C:110
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:65
rho
rho
Definition: readInitialConditions.H:88
phasePressureModel.H
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::EddyDiffusivity< phaseCompressibleTurbulenceModel >
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::RASModels::phasePressureModel::pPrimef
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure'.
Definition: phasePressureModel.C:176
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:436
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
Foam::dimensioned::readIfPresent
bool readIfPresent(const dictionary &dict)
Definition: dimensionedType.C:483
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:314
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::RASModels::phasePressureModel::pPrime
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure'.
Definition: phasePressureModel.C:148
Foam::RASModels::phasePressureModel::divDevRhoReff
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
Definition: phasePressureModel.C:232
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:43
Foam::fvc::interpolate
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.
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:60
Foam::RASModel< EddyDiffusivity< phaseCompressibleTurbulenceModel > >::coeffDict
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: RASModel.H:205
U
U
Definition: pEqn.H:72
Foam::compressible::RASModel
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
Definition: turbulentFluidThermoModel.H:66
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::RASModels::phasePressureModel::R
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
Definition: phasePressureModel.C:126
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::RASModel< EddyDiffusivity< phaseCompressibleTurbulenceModel > >::printCoeffs
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:34
Foam::eddyViscosity
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:55
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
Foam::RASModels::phasePressureModel::correct
virtual void correct()
Solve the kinetic theory equations and correct the viscosity.
Definition: phasePressureModel.C:247
Foam::eddyViscosity< RASModel< EddyDiffusivity< phaseCompressibleTurbulenceModel > > >::nut_
volScalarField nut_
Definition: eddyViscosity.H:66
Foam::RASModels::phasePressureModel::devRhoReff
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
Definition: phasePressureModel.C:204
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::GeometricField< symmTensor, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::RASModels::phasePressureModel::k
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: phasePressureModel.C:102
Foam::RASModels::phasePressureModel::omega
virtual tmp< volScalarField > omega() const
Return the specific dissipation rate.
Definition: phasePressureModel.C:118