compressibleTwoPhaseVoFMixture.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) 2013-2023 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 
27 
28 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
33 }
34 
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
39 (
40  const fvMesh& mesh
41 )
42 :
43  twoPhaseVoFMixture(mesh),
44 
45  totalInternalEnergy_
46  (
47  lookupOrDefault<Switch>("totalInternalEnergy", true)
48  ),
49 
50  p_
51  (
52  IOobject
53  (
54  "p",
55  mesh.time().name(),
56  mesh,
57  IOobject::MUST_READ,
58  IOobject::AUTO_WRITE
59  ),
60  mesh
61  ),
62  T_
63  (
64  IOobject
65  (
66  "T",
67  mesh.time().name(),
68  mesh,
69  IOobject::MUST_READ,
70  IOobject::AUTO_WRITE
71  ),
72  mesh
73  ),
74  thermo1_(nullptr),
75  thermo2_(nullptr),
76  rho_
77  (
78  IOobject
79  (
80  "rho",
81  mesh.time().name(),
82  mesh,
83  IOobject::NO_READ,
84  IOobject::AUTO_WRITE
85  ),
86  mesh,
87  dimensionedScalar("rho", dimDensity, 0)
88  ),
89  Alpha1_
90  (
91  IOobject
92  (
93  IOobject::groupName("Alpha", phase1Name()),
94  mesh.time().name(),
95  mesh
96  ),
97  alpha1(),
98  calculatedFvPatchScalarField::typeName
99  ),
100  Alpha2_
101  (
102  IOobject
103  (
104  IOobject::groupName("Alpha", phase2Name()),
105  mesh.time().name(),
106  mesh
107  ),
108  alpha2(),
109  calculatedFvPatchScalarField::typeName
110  )
111 {
112  {
113  volScalarField T1
114  (
115  IOobject
116  (
118  mesh.time().name(),
119  mesh
120  ),
121  T_,
122  calculatedFvPatchScalarField::typeName
123  );
124  T1.write();
125  }
126 
127  {
128  volScalarField T2
129  (
130  IOobject
131  (
133  mesh.time().name(),
134  mesh
135  ),
136  T_,
137  calculatedFvPatchScalarField::typeName
138  );
139  T2.write();
140  }
141 
142  // Note: we're writing files to be read in immediately afterwards.
143  // Avoid any thread-writing problems.
144  // fileHandler().flush();
145 
146  thermo1_ = rhoThermo::New(mesh, phase1Name());
147  thermo2_ = rhoThermo::New(mesh, phase2Name());
148 
149  // thermo1_->validate(phase1Name(), "e");
150  // thermo2_->validate(phase2Name(), "e");
151 
152  correct();
153 }
154 
155 
156 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
157 
159 {}
160 
161 
162 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
163 
165 {
166  thermo1_->T() = T_;
167  thermo1_->he() = thermo1_->he(p_, T_);
168  thermo1_->correct();
169 
170  thermo2_->T() = T_;
171  thermo2_->he() = thermo2_->he(p_, T_);
172  thermo2_->correct();
173 }
174 
175 
177 {
178  const volScalarField alphaRho1(alpha1()*thermo1_->rho());
179  const volScalarField alphaRho2(alpha2()*thermo2_->rho());
180 
181  rho_ = alphaRho1 + alphaRho2;
182  Alpha1_ = alphaRho1/rho_;
183  Alpha2_ = alphaRho2/rho_;
184 }
185 
186 
188 {
189  return (alpha1()*thermo1_->mu() + alpha2()*thermo2_->mu())/rho_;
190 }
191 
192 
194 (
195  const label patchi
196 ) const
197 {
198  return
199  (
200  alpha1().boundaryField()[patchi]*thermo1_->mu(patchi)
201  + alpha2().boundaryField()[patchi]*thermo2_->mu(patchi)
202  )/rho_.boundaryField()[patchi];
203 }
204 
205 
208 {
209  return
210  (
211  alpha1()*thermo1_->psi()/thermo1_->rho()
212  + alpha2()*thermo2_->psi()/thermo2_->rho()
213  );
214 }
215 
216 
218 {
220  {
221  totalInternalEnergy_ =
222  lookupOrDefault<Switch>("totalInternalEnergy", true);
223 
224  return true;
225  }
226  else
227  {
228  return false;
229  }
230 }
231 
232 
233 // ************************************************************************* //
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static word groupName(Name name, const word &group)
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
Class to represent a mixture of two rhoThermo-based phases.
virtual void correctThermo()
Correct the thermodynamics of each phase.
virtual void correct()
Update mixture properties.
compressibleTwoPhaseVoFMixture(const fvMesh &mesh)
Construct from a mesh.
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
virtual tmp< volScalarField > psiByRho() const
Return the mixture compressibility/density.
virtual bool read()
Read base phaseProperties dictionary.
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:406
virtual bool write(const bool write=true) const
Write using setting from DB.
static autoPtr< rhoThermo > New(const fvMesh &, const word &phaseName=word::null)
Standard selection based on fvMesh.
Definition: rhoThermo.C:64
A class for managing temporary objects.
Definition: tmp.H:55
const word & phase2Name() const
const word & phase1Name() const
Class to represent a VoF mixture.
virtual bool read()=0
Read base phaseProperties dictionary.
const fvMesh & mesh() const
Access the mesh.
Definition: twoPhases.H:71
label patchi
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const dimensionSet dimDensity
defineTypeNameAndDebug(combustionModel, 0)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47