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) 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 "shockFluid.H"
27 #include "fvmDdt.H"
28 #include "fvcDiv.H"
29 #include "fvcDdt.H"
30 
31 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
32 
34 {
36 
37  const surfaceScalarField e_pos(interpolate(e, pos, thermo.T().name()));
38  const surfaceScalarField e_neg(interpolate(e, neg, thermo.T().name()));
39 
40  surfaceScalarField phiEp
41  (
42  "phiEp",
43  aphiv_pos()*(rho_pos()*(e_pos + 0.5*magSqr(U_pos())) + p_pos())
44  + aphiv_neg()*(rho_neg()*(e_neg + 0.5*magSqr(U_neg())) + p_neg())
45  + aSf()*(p_pos() - p_neg())
46  );
47 
48  // Make flux for pressure-work absolute
49  if (mesh.moving())
50  {
51  phiEp += mesh.phi()*(a_pos()*p_pos() + a_neg()*p_neg());
52  }
53 
55  (
56  fvm::ddt(rho, e) + fvc::div(phiEp)
57  + fvc::ddt(rho, K)
58  ==
59  fvModels().source(rho, e)
60  );
61 
62  if (!inviscid)
63  {
64  const surfaceScalarField devTauDotU
65  (
66  "devTauDotU",
67  devTau() & (a_pos()*U_pos() + a_neg()*U_neg())
68  );
69 
70  EEqn += thermophysicalTransport->divq(e) + fvc::div(devTauDotU);
71  }
72 
73  EEqn.relax();
74 
76 
77  EEqn.solve();
78 
80 
81  thermo_.correct();
82 }
83 
84 
85 // ************************************************************************* //
fvScalarMatrix EEqn(fvm::div(phi, he)+(he.name()=="e" ? fvc::div(phi, volScalarField("Ekp", 0.5 *magSqr(U)+p/rho)) :fvc::div(phi, volScalarField("K", 0.5 *magSqr(U))))+thermophysicalTransport->divq(he)==fvModels.source(rho, he))
Generic GeometricField class.
const word & name() const
Return name.
Definition: IOobject.H:310
virtual void correct()=0
Update properties.
virtual const volScalarField & T() const =0
Temperature [K].
virtual const volScalarField & he() const =0
Enthalpy/Internal energy [J/kg].
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
const surfaceScalarField & phi() const
Return cell face motion fluxes.
bool moving() const
Is mesh moving.
Definition: polyMesh.H:476
Foam::fvModels & fvModels() const
Return the fvModels that are created on demand.
Definition: solver.C:96
Foam::fvConstraints & fvConstraints() const
Return the fvConstraints that are created on demand.
Definition: solver.C:107
const fvMesh & mesh
Region mesh.
Definition: solver.H:101
virtual void thermophysicalPredictor()
Construct and solve the energy equation,.
tmp< surfaceVectorField > U_neg
Definition: shockFluid.H:136
tmp< surfaceVectorField > U_pos
Definition: shockFluid.H:135
tmp< surfaceVectorField > devTau
Definition: shockFluid.H:149
tmp< surfaceScalarField > p_neg
Definition: shockFluid.H:139
tmp< surfaceScalarField > aphiv_pos
Definition: shockFluid.H:146
tmp< surfaceScalarField > neg
Definition: shockFluid.H:127
const psiThermo & thermo
Reference to the fluid thermophysical properties.
Definition: shockFluid.H:206
volScalarField K
Kinetic energy field.
Definition: shockFluid.H:102
psiThermo & thermo_
Reference to the fluid thermophysical properties.
Definition: shockFluid.H:83
tmp< surfaceScalarField > p_pos
Definition: shockFluid.H:138
tmp< surfaceScalarField > rho_neg
Definition: shockFluid.H:130
tmp< surfaceScalarField > pos
Definition: shockFluid.H:126
tmp< surfaceScalarField > aSf
Definition: shockFluid.H:144
tmp< surfaceScalarField > a_pos
Definition: shockFluid.H:141
autoPtr< fluidThermoThermophysicalTransportModel > thermophysicalTransport
Definition: shockFluid.H:116
const volScalarField & rho
Reference to the continuity density field.
Definition: shockFluid.H:212
tmp< surfaceScalarField > a_neg
Definition: shockFluid.H:142
tmp< surfaceScalarField > rho_pos
Definition: shockFluid.H:129
tmp< surfaceScalarField > aphiv_neg
Definition: shockFluid.H:147
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the matrix for the first temporal derivative.
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< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
const doubleScalar e
Definition: doubleScalar.H:106
dimensioned< scalar > magSqr(const dimensioned< Type > &)