incompressibleDriftFlux.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-2025 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 
27 #include "fvCorrectPhi.H"
29 
30 #include "fvmDdt.H"
31 
32 #include "fvcDdt.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace solvers
39 {
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 void Foam::solvers::incompressibleDriftFlux::correctCoNum()
49 {
51 }
52 
53 
54 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
55 
57 (
58  volScalarField& rDeltaT
59 )
60 {}
61 
62 
64 {}
65 
66 
69 {
71  (
72  "surfaceTensionForce",
73  mesh,
75  );
76 }
77 
78 
79 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
80 
82 :
84  (
85  mesh,
87  ),
88 
89  mixture
90  (
92  .initialise(U)
93  ),
94 
95  p
96  (
97  IOobject
98  (
99  "p",
100  runTime.name(),
101  mesh,
102  IOobject::NO_READ,
103  IOobject::AUTO_WRITE
104  ),
105  p_rgh + rho*buoyancy.gh
106  ),
107 
108  pressureReference_
109  (
110  p,
111  p_rgh,
112  pimple.dict()
113  ),
114 
115  relativeVelocity
116  (
118  ),
119 
120  packingDispersion
121  (
122  packingDispersionModel::New(mixture, relativeVelocity)
123  ),
124 
125  momentumTransport
126  (
127  compressible::momentumTransportModel::New(rho, U, rhoPhi, mixture)
128  )
129 {
130  if (transient())
131  {
132  correctCoNum();
133  }
134 
135  if (correctPhi || mesh.topoChanging())
136  {
137  rAU = new volScalarField
138  (
139  IOobject
140  (
141  "rAU",
142  runTime.name(),
143  mesh,
146  ),
147  mesh,
149  );
150  }
151 
152  if (!runTime.restart() || !divergent())
153  {
154  correctUphiBCs(U_, phi_, true);
155 
157  (
158  phi_,
159  U,
160  p_rgh,
161  rAU,
164  pimple
165  );
166  }
167 }
168 
169 
170 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
171 
173 {}
174 
175 
176 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
177 
179 {
180  return fluidSolver::maxDeltaT();
181 }
182 
183 
185 {
186  alphaPredictor();
187 
188  // Apply the diffusion term separately to allow implicit solution
189  // and boundedness of the explicit advection
190  {
191  volScalarField nuEff(momentumTransport->nut());
192  nuEff += packingDispersion->Dd();
193 
194  fvScalarMatrix alpha1Eqn
195  (
196  fvm::ddt(alpha1) - fvc::ddt(alpha1)
197  - fvm::laplacian(nuEff, alpha1)
198  );
199 
200  alpha1Eqn.solve(alpha1.name() + "Diffusion");
201 
202  alphaPhi1 += alpha1Eqn.flux();
203 
204  alpha2 = scalar(1) - alpha1;
205  alphaPhi2 = phi - alphaPhi1;
206 
207  Info<< "Phase-1 volume fraction = "
208  << alpha1.weightedAverage(mesh.Vsc()).value()
209  << " Min(" << alpha1.name() << ") = " << min(alpha1).value()
210  << " Max(" << alpha1.name() << ") = " << max(alpha1).value()
211  << endl;
212  }
213 
214  mixture.correct();
215 
216  const dimensionedScalar& rho1 = mixture.rho1();
217  const dimensionedScalar& rho2 = mixture.rho2();
218 
219  // Calculate the mass-flux
220  rhoPhi = alphaPhi1*rho1 + alphaPhi2*rho2;
221 
222  relativeVelocity->correct();
223 }
224 
225 
227 {
228  momentumTransport->predict();
229 }
230 
231 
233 {}
234 
235 
237 {
238  incompressiblePressureCorrector(p);
239 }
240 
241 
243 {}
244 
245 
247 {
248  momentumTransport->correct();
249 }
250 
251 
253 {}
254 
255 
256 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
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:426
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.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
tmp< SurfaceField< Type > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:964
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
bool topoChanging() const
Does the mesh topology change?
Definition: fvMesh.C:705
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
Class to represent a mixture of two constant density phases.
Abstract base class for turbulence models (RAS, LES and laminar).
Packing dispersion model.
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
volVectorField U_
Velocity field.
Definition: VoFSolver.H:94
const volScalarField & p_rgh
Reference to the buoyant pressure for buoyant cases.
Definition: VoFSolver.H:210
autoPtr< volScalarField > rAU
Inverse momentum equation diagonal.
Definition: VoFSolver.H:129
const volVectorField & U
Reference to the velocity field.
Definition: VoFSolver.H:213
surfaceScalarField phi_
Volumetric flux field.
Definition: VoFSolver.H:97
virtual void correctCoNum()=0
Correct the cached Courant numbers.
Definition: VoFSolver.C:75
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
virtual scalar maxDeltaT() const
Return the current maximum time-step for stable solution.
Definition: fluidSolver.C:229
Solver module for 2 incompressible fluids using the mixture approach with the drift-flux approximatio...
virtual void thermophysicalPredictor()
Construct and solve the energy equation,.
virtual tmp< surfaceScalarField > surfaceTensionForce() const
Return the interface surface tension force for the momentum equation.
virtual const Foam::pressureReference & pressureReference() const
Return the pressure reference.
virtual void momentumTransportCorrector()
Correct the momentum transport.
virtual void prePredictor()
Called at the start of the PIMPLE loop.
virtual bool divergent() const
Is the flow divergent?
virtual void momentumTransportPredictor()
Predict the momentum transport.
incompressibleDriftFlux(fvMesh &mesh)
Construct from region mesh.
virtual scalar maxDeltaT() const
Return the current maximum time-step for stable solution.
virtual void pressureCorrector()
Construct and solve the pressure equation in the PISO loop.
virtual void setInterfaceRDeltaT(volScalarField &rDeltaT)
Adjust the rDeltaT in the vicinity of the interface.
virtual void thermophysicalTransportCorrector()
Correct the thermophysical transport.
virtual void correctInterface()
Correct the interface properties following mesh-change.
virtual void thermophysicalTransportPredictor()
Predict thermophysical transport.
Solver module base-class for 2 immiscible fluids, with optional mesh motion and mesh topology changes...
A class for managing temporary objects.
Definition: tmp.H:55
Class to represent a VoF mixture.
pimpleControl pimple(mesh)
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
const scalar nuEff
Flux correction functions to ensure continuity.
Calculate the first temporal derivative.
Calculate the matrix for the first temporal derivative.
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
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
addToRunTimeSelectionTable(solver, compressibleMultiphaseVoF, fvMesh)
defineTypeNameAndDebug(compressibleMultiphaseVoF, 0)
Namespace for OpenFOAM.
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:134
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
messageStream Info
const dimensionSet dimForce
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
const dimensionSet dimDensity
const dimensionSet dimVolume
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
void correctUphiBCs(volVectorField &U, surfaceScalarField &phi, const bool evaluateUBCs)
If the mesh is moving correct the velocity BCs on the moving walls to.
tmp< DimensionedField< TypeR, GeoMesh, Field > > New(const tmp< DimensionedField< TypeR, GeoMesh, Field >> &tdf1, const word &name, const dimensionSet &dimensions)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
volScalarField & p