incompressibleVoF.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 "incompressibleVoF.H"
27 #include "localEulerDdtScheme.H"
28 #include "fvCorrectPhi.H"
29 #include "geometricZeroField.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace solvers
37 {
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 :
49  (
50  mesh,
52  ),
53 
54  mixture
55  (
57  ),
58 
59  p
60  (
61  IOobject
62  (
63  "p",
64  runTime.name(),
65  mesh,
66  IOobject::NO_READ,
67  IOobject::AUTO_WRITE
68  ),
69  p_rgh + rho*buoyancy.gh
70  ),
71 
72  pressureReference_
73  (
74  p,
75  p_rgh,
76  pimple.dict()
77  ),
78 
79  momentumTransport
80  (
81  U,
82  phi,
83  alphaPhi1,
84  alphaPhi2,
85  mixture
86  )
87 {
88  if (correctPhi || mesh.topoChanging())
89  {
90  rAU = new volScalarField
91  (
92  IOobject
93  (
94  "rAU",
95  runTime.name(),
96  mesh,
99  ),
100  mesh,
102  );
103  }
104 
105  if (!runTime.restart() || !divergent())
106  {
107  correctUphiBCs(U_, phi_, true);
108 
110  (
111  phi_,
112  U,
113  p_rgh,
114  rAU,
117  pimple
118  );
119  }
120 }
121 
122 
123 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
124 
126 {}
127 
128 
129 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
130 
132 {
134 
135  const dimensionedScalar& rho1 = mixture.rho1();
136  const dimensionedScalar& rho2 = mixture.rho2();
137 
138  // Calculate the mass-flux
139  rhoPhi = alphaPhi1*rho1 + alphaPhi2*rho2;
140 
141  if (pimple.predictTransport())
142  {
143  momentumTransport.predict();
144  }
145 }
146 
147 
149 {
150  incompressiblePressureCorrector(p);
151 }
152 
153 
155 {}
156 
157 
159 {
160  if (pimple.correctTransport())
161  {
162  momentumTransport.correct();
163  }
164 }
165 
166 
167 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
bool restart() const
Return true if the run is a restart, i.e. startTime != beginTime.
Definition: Time.H:432
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
bool topoChanging() const
Does the mesh topology change?
Definition: fvMesh.C:598
Class to represent a mixture of two constant density phases.
Abstract base class for run-time selectable region solvers.
Definition: solver.H:56
pimpleNoLoopControl pimple
PIMPLE inner-loop controls.
Definition: solver.H:107
const Time & runTime
Time.
Definition: solver.H:104
const fvMesh & mesh
Region mesh.
Definition: solver.H:101
volScalarField & p_rgh
Reference to the buoyant pressure for buoyant cases.
Definition: VoFSolver.H:107
volVectorField U_
Velocity field.
Definition: VoFSolver.H:94
autoPtr< volScalarField > rAU
Inverse momentum equation diagonal.
Definition: VoFSolver.H:129
const volVectorField & U
Reference to the velocity field.
Definition: VoFSolver.H:209
surfaceScalarField phi_
Volumetric flux field.
Definition: VoFSolver.H:97
Buoyancy related data for the Foam::solvers::isothermalFluid solver module when solving buoyant cases...
Definition: buoyancy.H:70
bool correctPhi
Switch to correct the flux after mesh change.
Definition: fluidSolver.H:101
Solver module for 2 incompressible, isothermal immiscible fluids using a VOF (volume of fluid) phase-...
virtual void thermophysicalPredictor()
Construct and solve the energy equation,.
virtual const Foam::pressureReference & pressureReference() const
Return the pressure reference.
virtual void prePredictor()
Called at the start of the PIMPLE loop.
virtual bool divergent() const
Is the flow divergent?
virtual ~incompressibleVoF()
Destructor.
virtual void pressureCorrector()
Construct and solve the pressure equation in the PISO loop.
virtual void postCorrector()
Correct the momentum and thermophysical transport modelling.
incompressibleVoF(fvMesh &mesh)
Construct from region mesh.
virtual void prePredictor()
Called at the start of the PIMPLE loop.
Solver module base-class for 2 immiscible fluids using a VOF (volume of fluid) phase-fraction based i...
Class to represent a VoF mixture.
pimpleControl pimple(mesh)
Flux correction functions to ensure continuity.
U
Definition: pEqn.H:72
void correctPhi(surfaceScalarField &phi, const volVectorField &U, const volScalarField &p, const autoPtr< volScalarField > &rAU, const autoPtr< volScalarField > &divU, const pressureReference &pressureReference, nonOrthogonalSolutionControl &pcorrControl)
Definition: fvCorrectPhi.C:32
addToRunTimeSelectionTable(solver, compressibleMultiphaseVoF, fvMesh)
defineTypeNameAndDebug(compressibleMultiphaseVoF, 0)
Namespace for OpenFOAM.
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:129
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimTime
const dimensionSet dimDensity
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
void correctUphiBCs(volVectorField &U, surfaceScalarField &phi, const bool evaluateUBCs)
If the mesh is moving correct the velocity BCs on the moving walls to.
dictionary dict
volScalarField & p