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-2016 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  <
48  RASModel<EddyDiffusivity<ThermalDiffusivity
49  <
50  PhaseCompressibleTurbulenceModel<phaseModel>
51  >>>
52  >
53  (
54  type,
55  alpha,
56  rho,
57  U,
58  alphaRhoPhi,
59  phi,
60  phase,
61  propertiesName
62  ),
63 
64  alphaMax_(coeffDict_.get<scalar>("alphaMax")),
65  preAlphaExp_(coeffDict_.get<scalar>("preAlphaExp")),
66  expMax_(coeffDict_.get<scalar>("expMax")),
67  g0_("g0", dimPressure, coeffDict_)
68 {
70 
71  if (type == typeName)
72  {
74  }
75 }
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
81 {
82  if
83  (
84  eddyViscosity
85  <
86  RASModel<EddyDiffusivity<ThermalDiffusivity
87  <
88  PhaseCompressibleTurbulenceModel<phaseModel>
89  >>>
90  >::read()
91  )
92  {
93  coeffDict().readEntry("alphaMax", alphaMax_);
94  coeffDict().readEntry("preAlphaExp", preAlphaExp_);
95  coeffDict().readEntry("expMax", expMax_);
96  g0_.readIfPresent(coeffDict());
97 
98  return true;
99  }
100 
101  return false;
102 }
103 
104 
107 {
109  return nut_;
110 }
111 
112 
115 {
117  return nut_;
118 }
119 
120 
123 {
125  return nullptr;
126 }
127 
128 
131 {
133  (
134  IOobject
135  (
136  IOobject::groupName("R", U_.group()),
137  runTime_.timeName(),
138  mesh_,
141  ),
142  mesh_,
143  dimensioned<symmTensor>(dimensionSet(0, 2, -2, 0, 0), Zero)
144  );
145 }
146 
147 
150 {
151  tmp<volScalarField> tpPrime
152  (
153  g0_
154  *min
155  (
156  exp(preAlphaExp_*(alpha_ - alphaMax_)),
157  expMax_
158  )
159  );
160 
161  volScalarField::Boundary& bpPrime =
162  tpPrime.ref().boundaryFieldRef();
163 
164  forAll(bpPrime, patchi)
165  {
166  if (!bpPrime[patchi].coupled())
167  {
168  bpPrime[patchi] == 0;
169  }
170  }
171 
172  return tpPrime;
173 }
174 
175 
178 {
179  tmp<surfaceScalarField> tpPrime
180  (
181  g0_
182  *min
183  (
184  exp(preAlphaExp_*(fvc::interpolate(alpha_) - alphaMax_)),
185  expMax_
186  )
187  );
188 
190  tpPrime.ref().boundaryFieldRef();
191 
192  forAll(bpPrime, patchi)
193  {
194  if (!bpPrime[patchi].coupled())
195  {
196  bpPrime[patchi] == 0;
197  }
198  }
199 
200  return tpPrime;
201 }
202 
203 
206 {
207  return devRhoReff(U_);
208 }
209 
210 
213 (
214  const volVectorField& U
215 ) const
216 {
218  (
219  IOobject
220  (
221  IOobject::groupName("devRhoReff", U.group()),
222  runTime_.timeName(),
223  mesh_,
226  ),
227  mesh_,
228  dimensioned<symmTensor>
229  (
230  rho_.dimensions()*dimensionSet(0, 2, -2, 0, 0), Zero
231  )
232  );
233 }
234 
235 
238 (
240 ) const
241 {
242  return tmp<fvVectorMatrix>
243  (
244  new fvVectorMatrix
245  (
246  U,
247  rho_.dimensions()*dimensionSet(0, 4, -2, 0, 0)
248  )
249  );
250 }
251 
252 
254 {}
255 
256 
257 // ************************************************************************* //
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
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
Ignore writing from objectRegistry::writeObject()
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
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
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
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.
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.
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.
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:40
Nothing to be read.
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:686
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.