OpenFOAM
v2512
The open source CFD toolbox
Loading...
Searching...
No Matches
alphaEqnSubCycle.H
Go to the documentation of this file.
1
if
(
pimple
.nCorrPIMPLE() > 1)
2
{
3
// If nOuterCorrectors > 1 then for all but the first loop the advection
4
// of alpha is done using an average, 0.5*phi+0.5*phiNew where phi is
5
// the flux at the beginning of the time step and phiNew is the flux
6
// estimate at the end of the time step from the previous outer
7
// iteration. Similarly we use 0.5*U + 0.5*UNew in later iterations.
8
if
(
pimple
.firstIter())
9
{
10
// To recalculate the alpha1 update in subsequent iterations, we
11
// must store its current value before overwriting with the new
12
// value
13
alpha1
.storePrevIter();
14
// Storing initial phi and U for use in later outer iterations.
15
phi
.storePrevIter();
16
U
.storePrevIter();
17
}
18
else
19
{
20
// Resetting alpha1 to value before advection in first PIMPLE
21
// iteration.
22
alpha1
=
alpha1
.prevIter();
23
24
// Temporarily setting U and phi with which to advect interface.
25
U
= 0.5*
U
.prevIter() + 0.5*
U
;
26
phi
= 0.5*
phi
.prevIter() + 0.5*
phi
;
27
}
28
}
29
30
if
(
nAlphaSubCycles
> 1)
31
{
32
dimensionedScalar totalDeltaT =
runTime
.deltaT();
33
surfaceScalarField rhoPhiSum
34
(
35
mesh
.newIOobject(
"rhoPhiSum"
),
36
mesh
,
37
dimensionedScalar(
rhoPhi
.dimensions(), Zero)
38
);
39
40
for
41
(
42
subCycle<volScalarField> alphaSubCycle(
alpha1
,
nAlphaSubCycles
);
43
!(++alphaSubCycle).end();
44
)
45
{
46
#include "
alphaEqn.H
"
47
rhoPhiSum += (
runTime
.deltaT()/totalDeltaT)*
rhoPhi
;
48
}
49
50
rhoPhi
= rhoPhiSum;
51
}
52
else
53
{
54
#include "
alphaEqn.H
"
55
}
56
57
if
(!
pimple
.firstIter())
58
{
59
// Resetting U and phi to value at latest iteration.
60
U
= 2.0*
U
-
U
.prevIter();
61
phi
= 2.0*
phi
-
phi
.prevIter();
62
}
63
64
rho
==
alpha1
*
rho1
+
alpha2
*
rho2
;
rhoPhi
rhoPhi
Definition
rhoEqn.H:10
pimple
pimpleControl & pimple
Definition
setRegionFluidFields.H:56
alpha1
const volScalarField & alpha1
Definition
setRegionFluidFields.H:6
rho2
volScalarField & rho2
Definition
setRegionFluidFields.H:30
alpha2
const volScalarField & alpha2
Definition
setRegionFluidFields.H:7
rho1
volScalarField & rho1
Definition
setRegionFluidFields.H:27
U
U
Definition
pEqn.H:72
phi
phi
Definition
correctPhiFaceMask.H:34
mesh
dynamicFvMesh & mesh
Definition
createDynamicFvMesh.H:6
runTime
engineTime & runTime
Definition
createEngineTime.H:13
alphaEqn.H
rho
rho
Definition
readInitialConditions.H:88
nAlphaSubCycles
label nAlphaSubCycles(alphaControls.get< label >("nAlphaSubCycles"))
applications
solvers
multiphase
interIsoFoam
alphaEqnSubCycle.H
Generated by
1.16.1