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 "fluid.H"
27 #include "fvcDdt.H"
28 #include "fvmDiv.H"
29 
30 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
31 
33 {
35 
37  (
39  + fvc::ddt(rho, K) + fvc::div(phi, K)
40  + pressureWork
41  (
42  he.name() == "e"
43  ? fvc::div(phi, p/rho)()
44  : -dpdt
45  )
47  ==
48  (
49  buoyancy.valid()
50  ? fvModels().source(rho, he) + rho*(U & buoyancy->g)
51  : fvModels().source(rho, he)
52  )
53  );
54 
55  EEqn.relax();
56 
58 
59  EEqn.solve();
60 
62 
63  thermo_.correct();
64 }
65 
66 
67 // ************************************************************************* //
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.
virtual void correct()=0
Update properties.
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
Finite volume models.
Definition: fvModels.H:65
Foam::fvConstraints & fvConstraints() const
Return the fvConstraints that are created on demand.
Definition: solver.C:107
Buoyancy related data for the Foam::solvers::isothermalFluid solver module when solving buoyant cases...
Definition: buoyancy.H:70
uniformDimensionedVectorField g
Gravitational acceleration.
Definition: buoyancy.H:83
virtual void thermophysicalPredictor()
Construct and solve the energy equation,.
autoPtr< fluidThermoThermophysicalTransportModel > thermophysicalTransport
Definition: fluid.H:77
const surfaceScalarField & phi
Mass-flux field.
fluidThermo & thermo_
Reference to the fluid thermophysical properties.
volScalarField K
Kinetic energy field.
const volVectorField & U
Velocity field.
tmp< volScalarField::Internal > pressureWork(const tmp< volScalarField::Internal > &) const
Adds the mesh-motion work to the pressure work term provided.
const volScalarField & rho
Reference to the continuity density field.
volScalarField::Internal dpdt
Rate of change of the pressure.
const volScalarField & p
Reference to the pressure field.
Calculate the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
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 > > 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
thermo he()