continuousGasKEqn.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) 2019-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 "continuousGasKEqn.H"
30 #include "twoPhaseSystem.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace LESModels
37 {
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
41 template<class BasicTurbulenceModel>
43 (
44  const alphaField& alpha,
45  const rhoField& rho,
46  const volVectorField& U,
47  const surfaceScalarField& alphaRhoPhi,
48  const surfaceScalarField& phi,
49  const transportModel& transport,
50  const word& propertiesName,
51  const word& type
52 )
53 :
55  (
56  alpha,
57  rho,
58  U,
59  alphaRhoPhi,
60  phi,
61  transport,
62  propertiesName,
63  type
64  ),
65 
66  liquidTurbulencePtr_(nullptr),
67 
68  alphaInversion_
69  (
71  (
72  "alphaInversion",
73  this->coeffDict_,
74  0.7
75  )
76  )
77 {
78  if (type == typeName)
79  {
80  this->printCoeffs(type);
81  }
82 }
83 
84 
85 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86 
87 template<class BasicTurbulenceModel>
89 {
91  {
92  alphaInversion_.readIfPresent(this->coeffDict());
93 
94  return true;
95  }
96 
97  return false;
98 }
99 
100 
101 template<class BasicTurbulenceModel>
102 const turbulenceModel&
104 {
105  if (!liquidTurbulencePtr_)
106  {
107  const volVectorField& U = this->U_;
108 
109  const transportModel& gas = this->transport();
110  const twoPhaseSystem& fluid =
111  refCast<const twoPhaseSystem>(gas.fluid());
112  const transportModel& liquid = fluid.otherPhase(gas);
113 
114  liquidTurbulencePtr_ =
115  &U.db().lookupObject<turbulenceModel>
116  (
118  (
120  liquid.name()
121  )
122  );
123  }
125  return *liquidTurbulencePtr_;
126 }
127 
128 
129 template<class BasicTurbulenceModel>
132 {
133  const volVectorField& U = this->U_;
134  const alphaField& alpha = this->alpha_;
135  const rhoField& rho = this->rho_;
136 
137  const turbulenceModel& liquidTurbulence = this->liquidTurbulence();
138 
139  return
140  (
141  max(alphaInversion_ - alpha, scalar(0))
142  *rho
143  *min
144  (
145  this->Ce_*sqrt(liquidTurbulence.k())/this->delta(),
146  1.0/U.time().deltaT()
147  )
148  );
149 }
150 
151 
152 template<class BasicTurbulenceModel>
155 {
156  const turbulenceModel& liquidTurbulence = this->liquidTurbulence();
157  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
158 
159  return
160  phaseTransferCoeff*liquidTurbulence.k()
161  - fvm::Sp(phaseTransferCoeff, this->k_);
162 }
163 
164 
165 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166 
167 } // End namespace LESModels
168 } // End namespace Foam
169 
170 // ************************************************************************* //
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:50
twoPhaseSystem & fluid
scalar delta
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedScalar sqrt(const dimensionedScalar &ds)
Generic dimensioned Type class.
Abstract base class for turbulence models (RAS, LES and laminar).
One equation eddy-viscosity model.
Definition: kEqn.H:73
const turbulenceModel & liquidTurbulence() const
Return the turbulence model for the liquid phase.
Class which solves the volume fraction equations for two phases.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
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.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
BasicTurbulenceModel::alphaField alphaField
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition: IOobject.C:456
virtual tmp< fvScalarMatrix > kSource() const
virtual bool read()
Read model coefficients if they have changed.
Definition: kEqn.C:122
tmp< volScalarField > phaseTransferCoeff() const
virtual bool read()
Read model coefficients if they have changed.
One-equation SGS model for the gas-phase in a two-phase system supporting phase-inversion.
U
Definition: pEqn.H:72
BasicTurbulenceModel::transportModel transportModel
Base-class for all transport models used by the incompressible turbulence models. ...
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
BasicTurbulenceModel::rhoField rhoField
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:61
Namespace for OpenFOAM.