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 {
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  (
128  IOobject::groupName("R", U_.group()),
130  mesh_,
131  dimensioned<symmTensor>(dimensionSet(0, 2, -2, 0, 0), Zero)
132  );
133 }
134 
135 
138 {
139  tmp<volScalarField> tpPrime
140  (
141  g0_
142  *min
143  (
144  exp(preAlphaExp_*(alpha_ - alphaMax_)),
145  expMax_
146  )
147  );
148 
149  auto& bpPrime = tpPrime.ref().boundaryFieldRef();
150 
151  forAll(bpPrime, patchi)
152  {
153  if (!bpPrime[patchi].coupled())
154  {
155  bpPrime[patchi] == 0;
156  }
157  }
158 
159  return tpPrime;
160 }
161 
162 
165 {
167  (
168  g0_
169  *min
170  (
171  exp(preAlphaExp_*(fvc::interpolate(alpha_) - alphaMax_)),
172  expMax_
173  )
174  );
175 
176  auto& bpPrime = tpPrime.ref().boundaryFieldRef();
177 
178  forAll(bpPrime, patchi)
179  {
180  if (!bpPrime[patchi].coupled())
181  {
182  bpPrime[patchi] == 0;
183  }
184  }
185 
186  return tpPrime;
187 }
188 
189 
192 {
193  return devRhoReff(U_);
194 }
195 
196 
199 (
200  const volVectorField& U
201 ) const
202 {
204  (
205  IOobject::groupName("devRhoReff", U.group()),
207  mesh_,
209  (
210  rho_.dimensions()*dimensionSet(0, 2, -2, 0, 0),
211  Zero
212  )
213  );
214 }
215 
216 
219 (
221 ) const
222 {
224  (
225  U,
226  rho_.dimensions()*dimensionSet(0, 4, -2, 0, 0)
227  );
228 }
229 
230 
232 {}
233 
234 
235 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
type
Types of root.
Definition: Roots.H:52
ThermalDiffusivity< PhaseCompressibleTurbulenceModel< phaseModel > > phaseCompressibleTurbulenceModel
Typedef for phaseCompressibleTurbulenceModel.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:51
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:77
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
dimensionedScalar exp(const dimensionedScalar &ds)
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
virtual void correct()
Solve the kinetic theory equations and correct the viscosity.
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
const dimensionSet dimPressure
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure&#39;.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
static tmp< GeometricField< symmTensor, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< symmTensor >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
virtual tmp< volScalarField > omega() const
Return the specific dissipation rate.
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
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< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
bool coupled
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:696
Do not request registration (bool: false)
const dimensionSet & dimensions() const noexcept
Return dimensions.
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure&#39;.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
virtual bool read()
Re-read model coefficients if they have changed.