compressibleInterPhaseThermophysicalTransportModel.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) 2022-2024 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 
28 #include "fvcSnGrad.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
34 (
35  const compressibleInterPhaseTransportModel& momentumTransport
36 )
37 :
38  thermophysicalTransportModel(momentumTransport.mixture_.mesh(), word::null),
39  momentumTransport_(momentumTransport)
40 {}
41 
42 
43 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
44 
46 {
47  return true;
48 }
49 
50 
53 {
54  const compressibleTwoPhaseVoFMixture& mixture_ =
55  momentumTransport_.mixture_;
56 
57  if (momentumTransport_.twoPhaseTransport_)
58  {
59  return
60  mixture_.alpha1()
61  *(
62  mixture_.thermo1().kappa()
63  + mixture_.thermo1().rho()*mixture_.thermo1().Cp()
64  *momentumTransport_.momentumTransport1_->nut()
65  )
66  + mixture_.alpha2()
67  *(
68  mixture_.thermo2().kappa()
69  + mixture_.thermo2().rho()*mixture_.thermo2().Cp()
70  *momentumTransport_.momentumTransport2_->nut()
71  );
72  }
73  else
74  {
75  return
76  mixture_.alpha1()
77  *(
78  mixture_.thermo1().kappa()
79  + mixture_.thermo1().rho()*mixture_.thermo1().Cp()
80  *momentumTransport_.mixtureMomentumTransport_->nut()
81  )
82  + mixture_.alpha2()
83  *(
84  mixture_.thermo2().kappa()
85  + mixture_.thermo2().rho()*mixture_.thermo2().Cp()
86  *momentumTransport_.mixtureMomentumTransport_->nut()
87  );
88  }
89 }
90 
91 
94 (
95  const label patchi
96 ) const
97 {
98  const compressibleTwoPhaseVoFMixture& mixture_ =
99  momentumTransport_.mixture_;
100 
101  if (momentumTransport_.twoPhaseTransport_)
102  {
103  return
104  mixture_.alpha1().boundaryField()[patchi]
105  *(
106  mixture_.thermo1().kappa().boundaryField()[patchi]
107  + mixture_.thermo1().rho(patchi)
108  *mixture_.thermo1().Cp().boundaryField()[patchi]
109  *momentumTransport_.momentumTransport1_->nut(patchi)
110  )
111  + mixture_.alpha2().boundaryField()[patchi]
112  *(
113  mixture_.thermo2().kappa().boundaryField()[patchi]
114  + mixture_.thermo2().rho(patchi)
115  *mixture_.thermo2().Cp().boundaryField()[patchi]
116  *momentumTransport_.momentumTransport2_->nut(patchi)
117  );
118  }
119  else
120  {
121  return
122  mixture_.alpha1().boundaryField()[patchi]
123  *(
124  mixture_.thermo1().kappa().boundaryField()[patchi]
125  + mixture_.thermo1().rho(patchi)
126  *mixture_.thermo1().Cp().boundaryField()[patchi]
127  *momentumTransport_.mixtureMomentumTransport_->nut(patchi)
128  )
129  + mixture_.alpha2().boundaryField()[patchi]
130  *(
131  mixture_.thermo2().kappa().boundaryField()[patchi]
132  + mixture_.thermo2().rho(patchi)
133  *mixture_.thermo2().Cp().boundaryField()[patchi]
134  *momentumTransport_.mixtureMomentumTransport_->nut(patchi)
135  );
136  }
137 }
138 
139 
142 {
143  const compressibleTwoPhaseVoFMixture& mixture_ =
144  momentumTransport_.mixture_;
145 
146  if (momentumTransport_.twoPhaseTransport_)
147  {
148  return
149  mixture_.alpha1()
150  *(
151  mixture_.thermo1().kappa()
152  + mixture_.thermo1().rho()*mixture_.thermo1().Cp()
153  *momentumTransport_.momentumTransport1_->nut()
154  )/mixture_.thermo1().Cv()
155  + mixture_.alpha2()
156  *(
157  mixture_.thermo2().kappa()
158  + mixture_.thermo2().rho()*mixture_.thermo2().Cp()
159  *momentumTransport_.momentumTransport2_->nut()
160  )/mixture_.thermo2().Cv();
161  }
162  else
163  {
164  return
165  mixture_.alpha1()
166  *(
167  mixture_.thermo1().kappa()
168  + mixture_.thermo1().rho()*mixture_.thermo1().Cp()
169  *momentumTransport_.mixtureMomentumTransport_->nut()
170  )/mixture_.thermo1().Cv()
171  + mixture_.alpha2()
172  *(
173  mixture_.thermo2().kappa()
174  + mixture_.thermo2().rho()*mixture_.thermo2().Cp()
175  *momentumTransport_.mixtureMomentumTransport_->nut()
176  )/mixture_.thermo2().Cv();
177  }
178 }
179 
180 
183 {
185  (
186  "q",
187  -fvc::interpolate(kappaEff())
188  *fvc::snGrad(momentumTransport_.mixture_.T())
189  );
190 }
191 
192 
195 (
196  const label patchi
197 ) const
198 {
199  return
200  - kappaEff(patchi)
201  *momentumTransport_.mixture_.T().boundaryField()[patchi].snGrad();
202 }
203 
204 
207 (
208  const label patchi
209 ) const
210 {
211  return tmp<scalarField>(nullptr);
212 }
213 
214 
217 (
219 ) const
220 {
222 
223  return tmp<fvScalarMatrix>(nullptr);
224 }
225 
226 
228 {}
229 
230 
232 {}
233 
234 
235 // ************************************************************************* //
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
virtual const volScalarField & Cv() const =0
Heat capacity at constant volume [J/kg/K].
virtual const volScalarField & Cp() const =0
Heat capacity at constant pressure [J/kg/K].
virtual const volScalarField & kappa() const =0
Thermal conductivity of mixture [W/m/K].
virtual void correct()
Solve the thermophysical transport model equations.
virtual tmp< scalarField > qCorr(const label patchi) const
Return the patch heat flux correction [W/m^2].
virtual void predict()
Predict the transport coefficients if possible.
tmp< volScalarField > alphaEff() const
Return the effective temperature transport coefficient.
compressibleInterPhaseThermophysicalTransportModel(const compressibleInterPhaseTransportModel &momentumTransport)
Construct from components.
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
virtual tmp< volScalarField > kappaEff() const
Effective thermal turbulent conductivity.
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
Transport model selection class for the compressibleInterFoam family of solvers.
Class to represent a mixture of two rhoFluidThermo-based phases.
const rhoFluidThermo & thermo1() const
Return the thermo for phase 1.
const rhoFluidThermo & thermo2() const
Return the thermo for phase 2.
virtual tmp< volScalarField > rho() const =0
Density [kg/m^3].
Abstract base class for all fluid and solid thermophysical transport models.
A class for managing temporary objects.
Definition: tmp.H:55
const volScalarField & alpha1() const
Return the phase-fraction of phase 1.
const volScalarField & alpha2() const
Return the phase-fraction of phase 2.
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
Calculate the snGrad of the given volField.
label patchi
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
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
thermo he()