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-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 
26 #include "multicomponentFluid.H"
27 #include "fvcDdt.H"
28 
29 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
30 
32 {
34  (
36  (
37  mesh,
38  fields,
39  phi,
40  mesh.schemes().div("div(phi,Yi_h)")
41  )
42  );
43 
44  reaction->correct();
45 
46  forAll(Y, i)
47  {
48  volScalarField& Yi = Y_[i];
49 
50  if (thermo_.solveSpecie(i))
51  {
52  fvScalarMatrix YiEqn
53  (
54  fvm::ddt(rho, Yi)
55  + mvConvection->fvmDiv(phi, Yi)
56  + thermophysicalTransport->divj(Yi)
57  ==
58  reaction->R(Yi)
59  + fvModels().source(rho, Yi)
60  );
61 
62  YiEqn.relax();
63 
64  fvConstraints().constrain(YiEqn);
65 
66  YiEqn.solve("Yi");
67 
69  }
70  else
71  {
73  }
74  }
75 
77 
78 
80 
82  (
83  fvm::ddt(rho, he) + mvConvection->fvmDiv(phi, he)
84  + fvc::ddt(rho, K) + fvc::div(phi, K)
85  + pressureWork
86  (
87  he.name() == "e"
88  ? mvConvection->fvcDiv(phi, p/rho)()
89  : -dpdt
90  )
92  ==
93  reaction->Qdot()
94  + (
95  buoyancy.valid()
96  ? fvModels().source(rho, he) + rho*(U & buoyancy->g)
97  : fvModels().source(rho, he)
98  )
99  );
100 
101  EEqn.relax();
102 
104 
105  EEqn.solve();
106 
108 
109  thermo_.correct();
110 }
111 
112 
113 // ************************************************************************* //
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))
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
tmp< fv::convectionScheme< scalar > > mvConvection(fv::convectionScheme< scalar >::New(mesh, fields, phi, mesh.schemes().div("div(phi,ft_b_ha_hau)")))
Generic GeometricField class.
void correctBoundaryConditions()
Correct boundary field.
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
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:603
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1751
Finite volume models.
Definition: fvModels.H:65
ITstream & div(const word &name) const
Definition: fvSchemes.C:479
Abstract base class for convection schemes.
bool solveSpecie(const label speciei) const
Should the given specie be solved for? I.e., is it active and.
void normaliseY()
Normalise the mass fractions by clipping positive and deriving.
Reaction base-class holding the specie names and coefficients.
Definition: reaction.H:57
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
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
const surfaceScalarField & phi
Mass-flux field.
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.
virtual void thermophysicalPredictor()
Construct and solve the energy equation,.
const PtrList< volScalarField > & Y
Reference to the composition.
fluidMulticomponentThermo & thermo_
PtrList< volScalarField > & Y_
autoPtr< fluidMulticomponentThermophysicalTransportModel > thermophysicalTransport
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
A class for managing temporary objects.
Definition: tmp.H:55
Calculate 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
thermo he()