momentumPredictor.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 
30 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
31 
33 {
35 
36  const surfaceVectorField phiUp
37  (
39  + (a_pos()*p_pos() + a_neg()*p_neg())*mesh.Sf()
40  );
41 
42  // Construct the divDevTau matrix first
43  // so that the maxwellSlipU BC can access the explicit part
44  tmp<fvVectorMatrix> divDevTau;
45  if (!inviscid)
46  {
47  divDevTau = momentumTransport->divDevTau(U);
48  }
49 
51  (
52  fvm::ddt(rho, U) + fvc::div(phiUp)
53  ==
54  fvModels().source(rho, U)
55  );
56 
57  if (!inviscid)
58  {
59  UEqn += divDevTau();
60  }
61 
62  UEqn.relax();
63 
65 
66  solve(UEqn);
67 
69  K = 0.5*magSqr(U);
70 
71  if (!inviscid)
72  {
73  devTau = divDevTau->flux();
74  }
75 }
76 
77 
78 // ************************************************************************* //
Generic GeometricField class.
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 surfaceVectorField & Sf() const
Return cell face area vectors.
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
volVectorField U_
Velocity field.
Definition: shockFluid.H:95
tmp< surfaceVectorField > devTau
Definition: shockFluid.H:149
tmp< surfaceScalarField > p_neg
Definition: shockFluid.H:139
tmp< surfaceScalarField > aphiv_pos
Definition: shockFluid.H:146
tmp< surfaceVectorField > rhoU_pos
Definition: shockFluid.H:132
volScalarField K
Kinetic energy field.
Definition: shockFluid.H:102
const volVectorField & U
Reference to the velocity field.
Definition: shockFluid.H:215
tmp< surfaceScalarField > p_pos
Definition: shockFluid.H:138
virtual void momentumPredictor()
Construct and optionally solve the momentum equation.
tmp< surfaceScalarField > a_pos
Definition: shockFluid.H:141
const volScalarField & rho
Reference to the continuity density field.
Definition: shockFluid.H:212
tmp< surfaceScalarField > a_neg
Definition: shockFluid.H:142
tmp< surfaceVectorField > rhoU_neg
Definition: shockFluid.H:133
tmp< surfaceScalarField > aphiv_neg
Definition: shockFluid.H:147
autoPtr< compressible::momentumTransportModel > momentumTransport
Pointer to the momentum transport model.
Definition: shockFluid.H:110
A class for managing temporary objects.
Definition: tmp.H:55
fvVectorMatrix & UEqn
Definition: UEqn.H:11
Calculate the divergence of the given field.
Calculate the matrix for the first temporal derivative.
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
dimensioned< scalar > magSqr(const dimensioned< Type > &)
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.