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) 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 "VoFSolver.H"
27 #include "fvmDiv.H"
28 #include "fvcSnGrad.H"
29 #include "fvcReconstruct.H"
30 
31 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
32 
34 {
35  volVectorField& U = U_;
36 
37  tUEqn =
38  (
40  + MRF.DDt(rho, U)
41  + divDevTau(U)
42  ==
43  fvModels().source(rho, U)
44  );
45  fvVectorMatrix& UEqn = tUEqn.ref();
46 
47  UEqn.relax();
48 
50 
52  {
53  solve
54  (
55  UEqn
56  ==
58  (
59  (
63  ) * mesh.magSf()
64  )
65  );
66 
68  }
69 }
70 
71 
72 // ************************************************************************* //
Generic GeometricField class.
tmp< volVectorField > DDt(const volVectorField &U) const
Return the Coriolis acceleration.
Definition: MRFZoneList.C:99
bool momentumPredictor() const
Flag to indicate to solve for momentum.
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 & magSf() const
Return cell face area magnitudes.
Foam::fvModels & fvModels() const
Return the fvModels that are created on demand.
Definition: solver.C:85
pimpleNoLoopControl pimple
PIMPLE inner-loop controls.
Definition: solver.H:100
Foam::fvConstraints & fvConstraints() const
Return the fvConstraints that are created on demand.
Definition: solver.C:96
const fvMesh & mesh
Region mesh.
Definition: solver.H:94
volScalarField & p_rgh
Reference to the buoyant pressure for buoyant cases.
Definition: VoFSolver.H:107
surfaceScalarField rhoPhi
Mass flux field.
Definition: VoFSolver.H:116
volVectorField U_
Velocity field.
Definition: VoFSolver.H:94
tmp< fvVectorMatrix > tUEqn
Cached momentum matrix.
Definition: VoFSolver.H:145
const volVectorField & U
Reference to the velocity field.
Definition: VoFSolver.H:209
virtual tmp< surfaceScalarField > surfaceTensionForce() const =0
Return the interface surface tension force for the momentum equation.
IOMRFZoneList MRF
MRF zone list.
Definition: VoFSolver.H:122
virtual void momentumPredictor()
Construct and optionally solve the momentum equation.
virtual tmp< fvVectorMatrix > divDevTau(volVectorField &U)=0
Return the momentum equation stress term.
const volScalarField & rho
Reference to the mixture continuity density field.
Definition: VoFSolver.H:110
Buoyancy related data for the Foam::solvers::isothermalFluid solver module when solving buoyant cases...
Definition: buoyancy.H:70
surfaceScalarField ghf
(g & hf) - ghRef
Definition: buoyancy.H:98
fvVectorMatrix & UEqn
Definition: UEqn.H:11
Reconstruct volField from a face flux field.
Calculate the snGrad of the given volField.
Calculate the matrix for the divergence of the given field and flux.
tmp< VolField< typename outerProduct< vector, Type >::type > > reconstruct(const SurfaceField< Type > &ssf)
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
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
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.