compressibleTwoPhaseMixture.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-2021 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 {
32  defineTypeNameAndDebug(compressibleTwoPhaseMixture, 0);
33 }
34 
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
39 (
40  const volVectorField& U,
41  const surfaceScalarField& phi
42 )
43 :
44  twoPhaseMixture(U.mesh()),
45  interfaceProperties(alpha1(), alpha2(), U, *this),
46 
47  totalInternalEnergy_
48  (
49  lookupOrDefault<Switch>("totalInternalEnergy", true)
50  ),
51 
52  p_
53  (
54  IOobject
55  (
56  "p",
57  U.mesh().time().timeName(),
58  U.mesh(),
59  IOobject::MUST_READ,
60  IOobject::AUTO_WRITE
61  ),
62  U.mesh()
63  ),
64  T_
65  (
66  IOobject
67  (
68  "T",
69  U.mesh().time().timeName(),
70  U.mesh(),
71  IOobject::MUST_READ,
72  IOobject::AUTO_WRITE
73  ),
74  U.mesh()
75  ),
76  thermo1_(nullptr),
77  thermo2_(nullptr),
78  rho_
79  (
80  IOobject
81  (
82  "thermo:rho",
83  U.mesh().time().timeName(),
84  U.mesh(),
85  IOobject::NO_READ,
86  IOobject::AUTO_WRITE
87  ),
88  U.mesh(),
89  dimensionedScalar("rho", dimDensity, 0)
90  ),
91  Alpha1_
92  (
93  IOobject
94  (
95  IOobject::groupName("Alpha", phase1Name()),
96  U.mesh().time().timeName(),
97  U.mesh()
98  ),
99  alpha1(),
100  calculatedFvPatchScalarField::typeName
101  ),
102  Alpha2_
103  (
104  IOobject
105  (
106  IOobject::groupName("Alpha", phase2Name()),
107  U.mesh().time().timeName(),
108  U.mesh()
109  ),
110  alpha2(),
111  calculatedFvPatchScalarField::typeName
112  )
113 {
114  {
115  volScalarField T1
116  (
117  IOobject
118  (
119  IOobject::groupName("T", phase1Name()),
120  U.mesh().time().timeName(),
121  U.mesh()
122  ),
123  T_,
124  calculatedFvPatchScalarField::typeName
125  );
126  T1.write();
127  }
128 
129  {
130  volScalarField T2
131  (
132  IOobject
133  (
134  IOobject::groupName("T", phase2Name()),
135  U.mesh().time().timeName(),
136  U.mesh()
137  ),
138  T_,
139  calculatedFvPatchScalarField::typeName
140  );
141  T2.write();
142  }
143 
144  // Note: we're writing files to be read in immediately afterwards.
145  // Avoid any thread-writing problems.
146  // fileHandler().flush();
147 
148  thermo1_ = rhoThermo::New(U.mesh(), phase1Name());
149  thermo2_ = rhoThermo::New(U.mesh(), phase2Name());
150 
151  // thermo1_->validate(phase1Name(), "e");
152  // thermo2_->validate(phase2Name(), "e");
153 
154  correct();
155 }
156 
157 
158 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
159 
161 {}
162 
163 
164 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
165 
167 {
168  thermo1_->T() = T_;
169  thermo1_->he() = thermo1_->he(p_, T_);
170  thermo1_->correct();
171 
172  thermo2_->T() = T_;
173  thermo2_->he() = thermo2_->he(p_, T_);
174  thermo2_->correct();
175 }
176 
177 
179 {
180  const volScalarField alphaRho1(alpha1()*thermo1_->rho());
181  const volScalarField alphaRho2(alpha2()*thermo2_->rho());
182 
183  rho_ = alphaRho1 + alphaRho2;
184  Alpha1_ = alphaRho1/rho_;
185  Alpha2_ = alphaRho2/rho_;
186 
188 }
189 
190 
192 {
193  return (alpha1()*thermo1_->mu() + alpha2()*thermo2_->mu())/rho_;
194 }
195 
196 
198 (
199  const label patchi
200 ) const
201 {
202  return
203  (
204  alpha1().boundaryField()[patchi]*thermo1_->mu(patchi)
205  + alpha2().boundaryField()[patchi]*thermo2_->mu(patchi)
206  )/rho_.boundaryField()[patchi];
207 }
208 
209 
211 {
212  if (twoPhaseMixture::read())
213  {
214  totalInternalEnergy_ =
215  lookupOrDefault<Switch>("totalInternalEnergy", true);
216 
217  return interfaceProperties::read();
218  }
219  else
220  {
221  return false;
222  }
223 }
224 
225 
226 // ************************************************************************* //
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
volScalarField & alpha1(mixture.alpha1())
static autoPtr< rhoThermo > New(const fvMesh &, const word &phaseName=word::null)
Standard selection based on fvMesh.
Definition: rhoThermo.C:92
bool read()
Read phaseProperties dictionary.
void correct()
Correct the phase or mixture transport models.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
fvMesh & mesh
virtual void correctThermo()
Correct the thermodynamics of each phase.
alpha2
Definition: alphaEqn.H:115
virtual void correct()
Update mixture properties.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
static word groupName(Name name, const word &group)
compressibleTwoPhaseMixture(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
const dimensionSet dimDensity
word timeName
Definition: getTimeIndex.H:3
virtual ~compressibleTwoPhaseMixture()
Destructor.
virtual bool read()
Read base phaseProperties dictionary.
defineTypeNameAndDebug(combustionModel, 0)
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual bool read()=0
Read base phaseProperties dictionary.
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Namespace for OpenFOAM.