thermophysicalPredictor.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-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 
26 #include "compressibleVoF.H"
27 #include "fvcMeshPhi.H"
28 #include "fvcDdt.H"
29 #include "fvmDiv.H"
30 #include "fvmSup.H"
31 #include "fvmLaplacian.H"
32 
33 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
34 
36 {
37  const volScalarField& rho1(mixture.rho1());
38  const volScalarField& rho2(mixture.rho2());
39 
40  const volScalarField& e1(mixture.thermo1().he());
41  const volScalarField& e2(mixture.thermo2().he());
42 
43  const fvScalarMatrix e1Source(fvModels().source(alpha1, rho1, e1));
44  const fvScalarMatrix e2Source(fvModels().source(alpha2, rho2, e2));
45 
47 
48  fvScalarMatrix TEqn
49  (
51  (
52  mixture.thermo1().Cv()()
53  *(
55  - (
56  e1Source.hasDiag()
57  ? fvm::Sp(contErr1(), T) + fvm::Sp(e1Source.A(), T)
58  : fvm::Sp(contErr1(), T)
59  )
60  )
61  + mixture.thermo2().Cv()()
62  *(
64  - (
65  e2Source.hasDiag()
66  ? fvm::Sp(contErr2(), T) + fvm::Sp(e2Source.A(), T)
67  : fvm::Sp(contErr2(), T)
68  )
69  )
70  )
71 
72  + fvc::ddt(alpha1, rho1, e1) + fvc::div(alphaRhoPhi1, e1)
73  - contErr1()*e1
74  + fvc::ddt(alpha2, rho2, e2) + fvc::div(alphaRhoPhi2, e2)
75  - contErr2()*e2
76 
78 
79  + (
80  mixture.totalInternalEnergy()
81  ?
82  fvc::div(fvc::absolute(phi, U), p)()()
83  + (fvc::ddt(rho, K) + fvc::div(rhoPhi, K))()()
84  - (U()&(fvModels().source(rho, U)&U)()) - (contErr1() + contErr2())*K
85  :
87  )
88  ==
89  (e1Source&e1)
90  + (e2Source&e2)
91  );
92 
93  TEqn.relax();
94 
95  fvConstraints().constrain(TEqn);
96 
97  TEqn.solve();
98 
100 
102  mixture_.correct();
103 }
104 
105 
106 // ************************************************************************* //
Generic GeometricField class.
virtual tmp< volScalarField > kappaEff() const
Effective thermal turbulent conductivity.
virtual void correctThermo()
Correct the thermodynamics of each phase.
virtual void correct()
Update mixture properties.
volScalarField & T()
Return mixture temperature [K].
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:604
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:809
bool hasDiag() const
Definition: lduMatrix.H:580
Foam::fvModels & fvModels() const
Return the fvModels that are created on demand.
Definition: solver.C:85
Foam::fvConstraints & fvConstraints() const
Return the fvConstraints that are created on demand.
Definition: solver.C:96
const surfaceScalarField & phi
Reference to the mass-flux field.
Definition: VoFSolver.H:212
surfaceScalarField rhoPhi
Mass flux field.
Definition: VoFSolver.H:116
const volVectorField & U
Reference to the velocity field.
Definition: VoFSolver.H:209
const volScalarField & rho
Reference to the mixture continuity density field.
Definition: VoFSolver.H:110
virtual void thermophysicalPredictor()
Construct and solve the energy equation,.
volScalarField & p
Reference to the mixture static pressure field.
tmp< volScalarField::Internal > contErr2
Phase-2 continuity error.
volScalarField K
Kinetic energy field.
compressibleTwoPhaseVoFMixture & mixture_
The compressible two-phase mixture.
tmp< volScalarField::Internal > contErr1
Phase-1 continuity error.
compressibleInterPhaseThermophysicalTransportModel thermophysicalTransport
Thermophysical transport model.
volScalarField & alpha1
Reference to the phase1-fraction.
volScalarField & alpha2
Reference to the phase2-fraction.
Calculate the first temporal derivative.
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for the laplacian of the field.
Calculate the matrix for implicit and explicit sources.
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)