incompressibleInterPhaseTransportModel.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) 2020 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class Mixture>
35 (
36  const volScalarField& rho,
37  const volVectorField& U,
38  const surfaceScalarField& phi,
40  const Mixture& mixture
41 )
42 :
43  rhoType_(rhoType::UNIFORM),
44  phi_(phi),
45  rhoPhi_(rhoPhi)
46 {
47  {
48  IOdictionary turbulenceProperties
49  (
50  IOobject
51  (
52  turbulenceModel::propertiesName,
53  U.time().constant(),
54  U.db(),
55  IOobject::MUST_READ,
56  IOobject::NO_WRITE
57  )
58  );
59 
60  if (turbulenceProperties.found("density"))
61  {
62  const word densityMethod
63  (
64  turbulenceProperties.getWord("density")
65  );
66  if (densityMethod == "variable")
67  {
68  rhoType_ = rhoType::VARIABLE;
69  }
70  else if (densityMethod == "uniform")
71  {
72  rhoType_ = rhoType::UNIFORM;
73  }
74  else
75  {
77  << "The rho type provided is not correct " << nl
78  << " Available types are : " << nl
79  << " variable or uniform. " << nl
80  << nl << exit(FatalError);
81  }
82  }
83  }
84 
85  if (rhoType_ == rhoType::VARIABLE)
86  {
87  rhoIncTurbulence_ =
88  (
90  (
91  rho,
92  U,
93  rhoPhi,
94  phi,
95  mixture
96  )
97  );
98  }
99  else
100  {
101  incTurbulence_ = incompressible::turbulenceModel::New
102  (
103  U,
104  phi,
105  mixture
106  );
107 
108  incTurbulence_->validate();
109  }
110 }
111 
112 
113 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
114 
115 template<class Mixture>
118 (
119  const volScalarField& rho,
121 ) const
122 {
123  if (rhoType_ == rhoType::VARIABLE)
124  {
125  return rhoIncTurbulence_->divDevRhoReff(U);
126  }
127  else
128  {
129  return incTurbulence_->divDevRhoReff(rho, U);
130  }
131 }
132 
133 template<class Mixture>
135 {
136  if (rhoType_ == rhoType::VARIABLE)
137  {
138  rhoIncTurbulence_->correct();
139  }
140  else
141  {
142  incTurbulence_->correct();
143  }
144 }
145 
146 
147 // ************************************************************************* //
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
rhoPhi
rhoPhi
Definition: rhoEqn.H:10
incompressibleInterPhaseTransportModel.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
rho
rho
Definition: readInitialConditions.H:88
Foam::incompressibleInterPhaseTransportModel::correct
void correct()
Correct the phase or mixture transport models.
Definition: incompressibleInterPhaseTransportModel.C:134
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::FatalError
error FatalError
Foam::dictionary::getWord
word getWord(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Same as get< word >(const word&, keyType::option)
Definition: dictionary.H:1514
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
U
U
Definition: pEqn.H:72
Foam::incompressibleInterPhaseTransportModel::divDevRhoReff
tmp< fvVectorMatrix > divDevRhoReff(const volScalarField &rho, volVectorField &U) const
Return the effective momentum stress divergence.
Definition: incompressibleInterPhaseTransportModel.C:118
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::incompressibleInterPhaseTransportModel
Transport model selection class for the incompressibleInterFoam family of solvers.
Definition: incompressibleInterPhaseTransportModel.H:63
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::mixture
Definition: mixture.H:54