PDRkEpsilon.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 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 "PDRkEpsilon.H"
27 #include "PDRDragModel.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace compressible
35 {
36 namespace RASModels
37 {
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 defineTypeNameAndDebug(PDRkEpsilon, 0);
42 addToRunTimeSelectionTable(RASModel, PDRkEpsilon, dictionary);
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 (
48  const geometricOneField& alpha,
49  const volScalarField& rho,
50  const volVectorField& U,
51  const surfaceScalarField& alphaRhoPhi,
52  const surfaceScalarField& phi,
53  const fluidThermo& thermophysicalModel,
54  const word& modelName
55 )
56 :
57  Foam::RASModels::kEpsilon<compressible::momentumTransportModel>
58  (
59  geometricOneField(),
60  rho,
61  U,
62  phi,
63  phi,
64  thermophysicalModel,
65  modelName
66  ),
67 
68  C4_
69  (
70  dimensioned<scalar>::lookupOrAddToDict
71  (
72  "C4",
73  coeffDict_,
74  0.1
75  )
76  )
77 {}
78 
79 
80 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
81 
82 PDRkEpsilon::~PDRkEpsilon()
83 {}
84 
85 
86 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
87 
88 bool PDRkEpsilon::read()
89 {
90  if (RASModel::read())
91  {
92  C4_.readIfPresent(coeffDict_);
93  return true;
94  }
95  else
96  {
97  return false;
98  }
99 }
100 
101 
103 {
104  if (!turbulence_)
105  {
106  // Re-calculate viscosity
107  nut_ = Cmu_*sqr(k_)/epsilon_;
108  nut_.correctBoundaryConditions();
109 
110  // Re-calculate thermal diffusivity
111  //***HGWalphat_ = mut_/Prt_;
112  // alphat_.correctBoundaryConditions();
113 
114  return;
115  }
116 
118 
120 
121  if (mesh_.moving())
122  {
123  divU += fvc::div(mesh_.phi());
124  }
125 
126  tmp<volTensorField> tgradU = fvc::grad(U_);
127  volScalarField G(GName(), rho_*nut_*(tgradU() && dev(twoSymm(tgradU()))));
128  tgradU.clear();
129 
130  // Update epsilon and G at the wall
131  epsilon_.boundaryFieldRef().updateCoeffs();
132 
133  // Add the blockage generation term so that it is included consistently
134  // in both the k and epsilon equations
135  const volScalarField& betav =
136  U_.db().lookupObject<volScalarField>("betav");
137 
138  const volScalarField& Lobs =
139  U_.db().lookupObject<volScalarField>("Lobs");
140 
141  const PDRDragModel& drag =
142  U_.db().lookupObject<PDRDragModel>("PDRDragModel");
143 
144  volScalarField GR(drag.Gk());
145 
146  volScalarField LI
147  (C4_*(Lobs + dimensionedScalar(dimLength, rootVSmall)));
148 
149  // Dissipation equation
150  tmp<fvScalarMatrix> epsEqn
151  (
152  betav*fvm::ddt(rho_, epsilon_)
153  + fvm::div(phi_, epsilon_)
154  - fvm::laplacian(rho_*DepsilonEff(), epsilon_)
155  ==
156  C1_*betav*G*epsilon_/k_
157  + 1.5*pow(Cmu_, 3.0/4.0)*GR*sqrt(k_)/LI
158  - fvm::SuSp(((2.0/3.0)*C1_)*betav*rho_*divU, epsilon_)
159  - fvm::Sp(C2_*betav*rho_*epsilon_/k_, epsilon_)
160  );
161 
162  epsEqn.ref().relax();
163 
164  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
165 
166  solve(epsEqn);
167  bound(epsilon_, epsilonMin_);
168 
169 
170  // Turbulent kinetic energy equation
171 
172  tmp<fvScalarMatrix> kEqn
173  (
174  betav*fvm::ddt(rho_, k_)
175  + fvm::div(phi_, k_)
176  - fvm::laplacian(rho_*DkEff(), k_)
177  ==
178  betav*G + GR
179  - fvm::SuSp((2.0/3.0)*betav*rho_*divU, k_)
180  - fvm::Sp(betav*rho_*epsilon_/k_, k_)
181  );
182 
183  kEqn.ref().relax();
184  solve(kEqn);
185  bound(k_, kMin_);
186 
187  // Re-calculate viscosity
188  nut_ = Cmu_*sqr(k_)/epsilon_;
189  nut_.correctBoundaryConditions();
190 
191  // Re-calculate thermal diffusivity
192  //***HGWalphat_ = mut_/Prt_;
193  // alphat_.correctBoundaryConditions();
194 }
195 
196 
197 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
198 
199 } // End namespace RASModels
200 } // End namespace compressible
201 } // End namespace Foam
202 
203 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
dimensionedScalar sqrt(const dimensionedScalar &ds)
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
virtual bool read()
Read model coefficients if they have changed.
Definition: RASModel.C:164
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
PDRkEpsilon(const geometricOneField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const fluidThermo &thermophysicalModel, const word &modelName=typeName)
Construct from components.
const compressible::RASModel & turbulence_
Definition: PDRDragModel.H:64
Macros for easy insertion into run-time selection tables.
const volScalarField & rho_
Definition: PDRDragModel.H:65
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
phi
Definition: pEqn.H:104
const volVectorField & U_
Definition: PDRDragModel.H:66
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
CompressibleMomentumTransportModel< fluidThermo > momentumTransportModel
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
const surfaceScalarField & phi_
Definition: PDRDragModel.H:67
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-8/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
defineTypeNameAndDebug(combustionModel, 0)
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.
PDRDragModel(const dictionary &PDRProperties, const compressible::RASModel &turbulence, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
U
Definition: pEqn.H:72
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
zeroField divU
Definition: alphaSuSp.H:3
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:322
const dimensionedScalar & G
Newtonian constant of gravitation.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RASModel.C:187
RASModel< momentumTransportModel > RASModel
Typedefs for turbulence, RAS and LES models for compressible flow based on the standard laminar trans...
Namespace for OpenFOAM.