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  alphaMax_(coeffDict_.get<scalar>("alphaMax")),
62  preAlphaExp_(coeffDict_.get<scalar>("preAlphaExp")),
63  expMax_(coeffDict_.get<scalar>("expMax")),
64  g0_("g0", dimPressure, coeffDict_)
65 {
66  nut_ == dimensionedScalar(nut_.dimensions());
67 
68  if (type == typeName)
69  {
71  }
72 }
73 
74 
75 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
76 
78 {
79  if
80  (
82  <
84  >::read()
85  )
86  {
87  coeffDict().readEntry("alphaMax", alphaMax_);
88  coeffDict().readEntry("preAlphaExp", preAlphaExp_);
89  coeffDict().readEntry("expMax", expMax_);
90  g0_.readIfPresent(coeffDict());
91 
92  return true;
93  }
94 
95  return false;
96 }
97 
98 
101 {
103  return nut_;
104 }
105 
106 
109 {
111  return nut_;
112 }
113 
114 
117 {
119  return nullptr;
120 }
121 
122 
125 {
127  (
129  (
130  IOobject
131  (
132  IOobject::groupName("R", U_.group()),
133  runTime_.timeName(),
134  mesh_,
137  ),
138  mesh_,
139  dimensioned<symmTensor>(dimensionSet(0, 2, -2, 0, 0), Zero)
140  )
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 {
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 {
204  return devRhoReff(U_);
205 }
206 
207 
210 (
211  const volVectorField& U
212 ) const
213 {
215  (
217  (
218  IOobject
219  (
220  IOobject::groupName("devRhoReff", U.group()),
221  runTime_.timeName(),
222  mesh_,
225  ),
226  mesh_,
228  (
229  "R",
230  rho_.dimensions()*dimensionSet(0, 2, -2, 0, 0),
231  Zero
232  )
233  )
234  );
235 }
236 
237 
240 (
242 ) const
243 {
244  return tmp<fvVectorMatrix>
245  (
246  new fvVectorMatrix
247  (
248  U,
249  rho_.dimensions()*dimensionSet(0, 4, -2, 0, 0)
250  )
251  );
252 }
253 
254 
256 {}
257 
258 
259 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::dimPressure
const dimensionSet dimPressure
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
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:52
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:108
Foam::RASModels::phasePressureModel::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: phasePressureModel.C:77
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:108
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::dimensionSet
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:108
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:174
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
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:146
Foam::RASModels::phasePressureModel::divDevRhoReff
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
Definition: phasePressureModel.C:240
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
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:62
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:124
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:68
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::PtrListOps::get
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
Foam::RASModels::phasePressureModel::correct
virtual void correct()
Solve the kinetic theory equations and correct the viscosity.
Definition: phasePressureModel.C:255
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:202
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
coupled
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
Foam::GeometricField< symmTensor, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::RASModels::phasePressureModel::k
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: phasePressureModel.C:100
Foam::RASModels::phasePressureModel::omega
virtual tmp< volScalarField > omega() const
Return the specific dissipation rate.
Definition: phasePressureModel.C:116