interfaceCompression: New run-time selectable VoF interface compression scheme

A new run-time selectable interface compression scheme framework has been added
to the two-phase VoF solvers to provide greater flexibility, extensibility and
more consistent user-interface.  The previously built-in interface compression
is now in the standard run-time selectable surfaceInterpolationScheme
interfaceCompression:

Class
    Foam::interfaceCompression

Description
    Interface compression corrected scheme, based on counter-gradient
    transport, to maintain sharp interfaces during VoF simulations.

    The interface compression is applied to the face interpolated field from a
    suitable 2nd-order shape-preserving NVD or TVD scheme, e.g.  vanLeer or
    vanAlbada.  A coefficient is supplied to control the degree of compression,
    with a value of 1 suitable for most VoF cases to ensure interface integrity.
    A value larger than 1 can be used but the additional compression can bias
    the interface to follow the mesh more closely while a value smaller than 1
    can lead to interface smearing.

    Example:
    \verbatim
    divSchemes
    {
        .
        .
        div(phi,alpha)     Gauss interfaceCompression vanLeer 1;
        .
        .
    }
    \endverbatim

The separate scheme for the interface compression term "div(phirb,alpha)" is no
longer required or used nor is the compression coefficient cAlpha in fvSolution
as this is now part of the "div(phi,alpha)" scheme specification as shown above.

Backward-compatibility is provided by checking the specified "div(phi,alpha)"
scheme against the known interface compression schemes and if it is not one of
those the new interfaceCompression scheme is used with the cAlpha value
specified in fvSolution.

More details can be found here:
https://cfd.direct/openfoam/free-software/multiphase-interface-capturing

Henry G. Weller
CFD Direct Ltd.
This commit is contained in:
Henry Weller
2020-07-02 10:13:15 +01:00
parent 70cfa1a47c
commit fa79bab863
77 changed files with 325 additions and 290 deletions

View File

@ -1,5 +1,3 @@
#include "alphaControls.H"
tmp<surfaceScalarField> talphaPhi1(alphaPhi10); tmp<surfaceScalarField> talphaPhi1(alphaPhi10);
if (nAlphaSubCycles > 1) if (nAlphaSubCycles > 1)

View File

@ -105,6 +105,7 @@ int main(int argc, char *argv[])
// --- Pressure-velocity PIMPLE corrector loop // --- Pressure-velocity PIMPLE corrector loop
while (pimple.loop()) while (pimple.loop())
{ {
#include "alphaControls.H"
#include "compressibleAlphaEqnSubCycle.H" #include "compressibleAlphaEqnSubCycle.H"
turbulence.correctPhasePhi(); turbulence.correctPhasePhi();

View File

@ -138,6 +138,7 @@ int main(int argc, char *argv[])
divU = fvc::div(fvc::absolute(phi, U)); divU = fvc::div(fvc::absolute(phi, U));
#include "alphaControls.H"
#include "compressibleAlphaEqnSubCycle.H" #include "compressibleAlphaEqnSubCycle.H"
turbulence.correctPhasePhi(); turbulence.correctPhasePhi();

View File

@ -134,6 +134,7 @@ int main(int argc, char *argv[])
} }
} }
#include "alphaControls.H"
#include "alphaEqnSubCycle.H" #include "alphaEqnSubCycle.H"
mixture.correct(); mixture.correct();

View File

@ -1,21 +1,15 @@
const dictionary& alphaControls = mesh.solverDict(alpha1.name()); const dictionary& alphaControls = mesh.solverDict(alpha1.name());
label nAlphaCorr(alphaControls.lookup<label>("nAlphaCorr")); const label nAlphaCorr(alphaControls.lookup<label>("nAlphaCorr"));
label nAlphaSubCycles(alphaControls.lookup<label>("nAlphaSubCycles")); const label nAlphaSubCycles(alphaControls.lookup<label>("nAlphaSubCycles"));
bool MULESCorr(alphaControls.lookupOrDefault<Switch>("MULESCorr", false)); const bool MULESCorr(alphaControls.lookupOrDefault<Switch>("MULESCorr", false));
// Apply the compression correction from the previous iteration // Apply the compression correction from the previous iteration
// Improves efficiency for steady-simulations but can only be applied // Improves efficiency for steady-simulations but can only be applied
// once the alpha field is reasonably steady, i.e. fully developed // once the alpha field is reasonably steady, i.e. fully developed
//bool alphaApplyPrevCorr // const bool alphaApplyPrevCorr
//( // (
// alphaControls.lookupOrDefault<Switch>("alphaApplyPrevCorr", false) // alphaControls.lookupOrDefault<Switch>("alphaApplyPrevCorr", false)
//); // );
// compression coefficient
scalar cAlpha
(
alphaControls.lookupOrDefault<scalar>("cAlpha", 0)
);

View File

@ -1,8 +1,5 @@
{ {
word alphaScheme("div(phi,alpha)"); #include "alphaScheme.H"
word alpharScheme("div(phirb,alpha)");
surfaceScalarField phir("phir", phic*mixture.nHatf());
Pair<tmp<volScalarField>> vDotAlphal = mixture.vDotAlphal(); Pair<tmp<volScalarField>> vDotAlphal = mixture.vDotAlphal();
const volScalarField& vDotcAlphal = vDotAlphal[0](); const volScalarField& vDotcAlphal = vDotAlphal[0]();
@ -49,13 +46,7 @@
( (
phi, phi,
alpha1, alpha1,
alphaScheme compressionScheme.rewind()
)
+ fvc::flux
(
-fvc::flux(-phir, alpha2, alpharScheme),
alpha1,
alpharScheme
) )
); );

View File

@ -1,30 +1,33 @@
#include "alphaControls.H" if (nAlphaSubCycles > 1)
{ {
// Standard face-flux compression coefficient dimensionedScalar totalDeltaT = runTime.deltaT();
surfaceScalarField phic(cAlpha*mag(phi/mesh.magSf())); surfaceScalarField rhoPhiSum
(
if (nAlphaSubCycles > 1) IOobject
{
dimensionedScalar totalDeltaT = runTime.deltaT();
surfaceScalarField rhoPhiSum("rhoPhiSum", rhoPhi);
for
( (
subCycle<volScalarField> alphaSubCycle(alpha1, nAlphaSubCycles); "rhoPhiSum",
!(++alphaSubCycle).end(); runTime.timeName(),
) mesh
{ ),
#include "alphaEqn.H" mesh,
rhoPhiSum += (runTime.deltaT()/totalDeltaT)*rhoPhi; dimensionedScalar(rhoPhi.dimensions(), 0)
} );
rhoPhi = rhoPhiSum; for
} (
else subCycle<volScalarField> alphaSubCycle(alpha1, nAlphaSubCycles);
!(++alphaSubCycle).end();
)
{ {
#include "alphaEqn.H" #include "alphaEqn.H"
rhoPhiSum += (runTime.deltaT()/totalDeltaT)*rhoPhi;
} }
rho == alpha1*rho1 + alpha2*rho2; rhoPhi = rhoPhiSum;
} }
else
{
#include "alphaEqn.H"
}
rho == alpha1*rho1 + alpha2*rho2;

View File

@ -140,6 +140,7 @@ int main(int argc, char *argv[])
dimensionedScalar(dimMass/dimTime, 0) dimensionedScalar(dimMass/dimTime, 0)
); );
#include "alphaControls.H"
#include "alphaEqnSubCycle.H" #include "alphaEqnSubCycle.H"
mixture.correct(); mixture.correct();

View File

@ -1,5 +1,4 @@
interfaceProperties.C interfaceProperties.C
interfaceCompression/interfaceCompression.C
surfaceTensionModels/surfaceTensionModel/surfaceTensionModel.C surfaceTensionModels/surfaceTensionModel/surfaceTensionModel.C
surfaceTensionModels/surfaceTensionModel/newSurfaceTensionModel.C surfaceTensionModels/surfaceTensionModel/newSurfaceTensionModel.C

View File

@ -1,88 +0,0 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Class
Foam::interfaceCompressionLimiter
Description
Interface compression scheme currently based on the generic limited
scheme although it does not use the NVD/TVD functions.
SourceFiles
interfaceCompression.C
\*---------------------------------------------------------------------------*/
#ifndef interfaceCompression_H
#define interfaceCompression_H
#include "vector.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class interfaceCompressionWeight Declaration
\*---------------------------------------------------------------------------*/
class interfaceCompressionLimiter
{
public:
interfaceCompressionLimiter(Istream&)
{}
scalar limiter
(
const scalar cdWeight,
const scalar faceFlux,
const scalar phiP,
const scalar phiN,
const vector&,
const scalar
) const
{
// Quadratic compression scheme
// return min(max(4*min(phiP*(1 - phiP), phiN*(1 - phiN)), 0), 1);
// Quartic compression scheme
return
min(max(
1 - max(sqr(1 - 4*phiP*(1 - phiP)), sqr(1 - 4*phiN*(1 - phiN))),
0), 1);
}
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -1,3 +1,4 @@
twoPhaseMixture/twoPhaseMixture.C twoPhaseMixture/twoPhaseMixture.C
interfaceCompression/interfaceCompression.C
LIB = $(FOAM_LIBBIN)/libtwoPhaseMixture LIB = $(FOAM_LIBBIN)/libtwoPhaseMixture

View File

@ -1,27 +1,15 @@
const dictionary& alphaControls = mesh.solverDict(alpha1.name()); const dictionary& alphaControls = mesh.solverDict(alpha1.name());
label nAlphaCorr(alphaControls.lookup<label>("nAlphaCorr")); const label nAlphaCorr(alphaControls.lookup<label>("nAlphaCorr"));
label nAlphaSubCycles(alphaControls.lookup<label>("nAlphaSubCycles")); const label nAlphaSubCycles(alphaControls.lookup<label>("nAlphaSubCycles"));
bool MULESCorr(alphaControls.lookupOrDefault<Switch>("MULESCorr", false)); const bool MULESCorr(alphaControls.lookupOrDefault<Switch>("MULESCorr", false));
// Apply the compression correction from the previous iteration // Apply the compression correction from the previous iteration
// Improves efficiency for steady-simulations but can only be applied // Improves efficiency for steady-simulations but can only be applied
// once the alpha field is reasonably steady, i.e. fully developed // once the alpha field is reasonably steady, i.e. fully developed
bool alphaApplyPrevCorr const bool alphaApplyPrevCorr
( (
alphaControls.lookupOrDefault<Switch>("alphaApplyPrevCorr", false) alphaControls.lookupOrDefault<Switch>("alphaApplyPrevCorr", false)
); );
// compression coefficient
scalar cAlpha
(
alphaControls.lookupOrDefault<scalar>("cAlpha", 0)
);
// Shear compression coefficient
scalar scAlpha
(
alphaControls.lookupOrDefault<scalar>("scAlpha", 0)
);

View File

@ -1,6 +1,5 @@
{ {
word alphaScheme("div(phi,alpha)"); #include "alphaScheme.H"
word alpharScheme("div(phirb,alpha)");
// Set the off-centering coefficient according to ddt scheme // Set the off-centering coefficient according to ddt scheme
scalar ocCoeff = 0; scalar ocCoeff = 0;
@ -55,38 +54,16 @@
// Set the time blending factor, 1 for Euler // Set the time blending factor, 1 for Euler
scalar cnCoeff = 1.0/(1.0 + ocCoeff); scalar cnCoeff = 1.0/(1.0 + ocCoeff);
// Standard face-flux compression coefficient
surfaceScalarField phic(cAlpha*mag(phi/mesh.magSf()));
// Add the optional shear compression contribution
if (scAlpha > 0)
{
phic +=
scAlpha*mag(mesh.delta() & fvc::interpolate(symm(fvc::grad(U))));
}
surfaceScalarField::Boundary& phicBf =
phic.boundaryFieldRef();
// Do not compress interface at non-coupled boundary faces
// (inlets, outlets etc.)
forAll(phic.boundaryField(), patchi)
{
fvsPatchScalarField& phicp = phicBf[patchi];
if (!phicp.coupled())
{
phicp == 0;
}
}
tmp<surfaceScalarField> phiCN(phi); tmp<surfaceScalarField> phiCN(phi);
// Calculate the Crank-Nicolson off-centred volumetric flux // Calculate the Crank-Nicolson off-centred volumetric flux
if (ocCoeff > 0) if (ocCoeff > 0)
{ {
phiCN = cnCoeff*phi + (1.0 - cnCoeff)*phi.oldTime(); phiCN = surfaceScalarField::New
(
"phiCN",
cnCoeff*phi + (1.0 - cnCoeff)*phi.oldTime()
);
} }
if (MULESCorr) if (MULESCorr)
@ -147,26 +124,18 @@
mixture.correct(); mixture.correct();
} }
for (int aCorr=0; aCorr<nAlphaCorr; aCorr++) for (int aCorr=0; aCorr<nAlphaCorr; aCorr++)
{ {
#include "alphaSuSp.H" #include "alphaSuSp.H"
surfaceScalarField phir(phic*mixture.nHatf()); // Split operator
tmp<surfaceScalarField> talphaPhi1Un tmp<surfaceScalarField> talphaPhi1Un
( (
fvc::flux fvc::flux
( (
phiCN(), phiCN(),
cnCoeff*alpha1 + (1.0 - cnCoeff)*alpha1.oldTime(), (cnCoeff*alpha1 + (1.0 - cnCoeff)*alpha1.oldTime())(),
alphaScheme compressionScheme.rewind()
)
+ fvc::flux
(
-fvc::flux(-phir, alpha2, alpharScheme),
alpha1,
alpharScheme
) )
); );

View File

@ -1,5 +1,3 @@
#include "alphaControls.H"
if (nAlphaSubCycles > 1) if (nAlphaSubCycles > 1)
{ {
dimensionedScalar totalDeltaT = runTime.deltaT(); dimensionedScalar totalDeltaT = runTime.deltaT();

View File

@ -0,0 +1,29 @@
static const wordHashSet compressionSchemes
{
"interfaceCompression",
"PLIC",
"PLICU",
"MPLIC",
"MPLICU"
};
static const word divAlphaName("div(phi,alpha)");
const word alphaScheme(mesh.divScheme(divAlphaName)[1].wordToken());
ITstream compressionScheme
(
compressionSchemes.found(alphaScheme)
? mesh.divScheme(divAlphaName)
: ITstream
(
divAlphaName,
tokenList
{
word("Gauss"),
word("interfaceCompression"),
alphaScheme,
alphaControls.lookup<scalar>("cAlpha")
}
)
);

View File

@ -0,0 +1,3 @@
interfaceCompression.C
LIB = $(FOAM_USER_LIBBIN)/libinterfaceCompression

View File

@ -0,0 +1,7 @@
EXE_INC = \
-I$(LIB_SRC)/finiteVolume/lnInclude \
-I$(LIB_SRC)/meshTools/lnInclude
EXE_LIBS = \
-lfiniteVolume \
-lmeshTools

View File

@ -2,7 +2,7 @@
========= | ========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org \\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation \\ / A nd | Copyright (C) 2020 OpenFOAM Foundation
\\/ M anipulation | \\/ M anipulation |
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
@ -23,19 +23,62 @@ License
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/
#include "PhiScheme.H"
#include "interfaceCompression.H" #include "interfaceCompression.H"
#include "linear.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam namespace Foam
{ {
makePhiSurfaceInterpolationScheme defineTypeNameAndDebug(interfaceCompressionNew, 0);
(
interfaceCompression, surfaceInterpolationScheme<scalar>::
interfaceCompressionLimiter, addMeshFluxConstructorToTable<interfaceCompressionNew>
scalar addinterfaceCompressionScalarMeshFluxConstructorToTable_;
)
} }
Foam::tmp<Foam::surfaceScalarField>
Foam::interfaceCompressionNew::interpolate
(
const volScalarField& vf
) const
{
const surfaceScalarField& nHatf =
mesh().lookupObject<const surfaceScalarField>
(
"nHatf"
);
const surfaceScalarField vff
(
linear<scalar>(mesh()).interpolate(vf)
);
surfaceScalarField vfc
(
cAlpha_*sign(phi_)*vff*(1 - vff)*nHatf/mesh().magSf()
);
surfaceScalarField::Boundary& vfcBf = vfc.boundaryFieldRef();
// Do not compress interface at non-coupled boundary faces
// (inlets, outlets etc.)
forAll(vfc.boundaryField(), patchi)
{
fvsPatchScalarField& vfcp = vfcBf[patchi];
if (!vfcp.coupled())
{
vfcp == 0;
}
}
tmp<surfaceScalarField> tvff(tScheme_().interpolate(vf) + vfc);
tvff.ref().maxMin(0, 1);
return tvff;
}
// ************************************************************************* // // ************************************************************************* //

View File

@ -0,0 +1,154 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2020 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Class
Foam::interfaceCompression
Description
Interface compression corrected scheme, based on counter-gradient
transport, to maintain sharp interfaces during VoF simulations.
The interface compression is applied to the face interpolated field from a
suitable 2nd-order shape-preserving NVD or TVD scheme, e.g. vanLeer or
vanAlbada. A coefficient is supplied to control the degree of compression,
with a value of 1 suitable for most VoF cases to ensure interface integrity.
A value larger than 1 can be used but the additional compression can bias
the interface to follow the mesh more closely while a value smaller than 1
can lead to interface smearing.
Example:
\verbatim
divSchemes
{
.
.
div(phi,alpha) Gauss interfaceCompression vanLeer 1;
.
.
}
\endverbatim
See also
Foam::vanLeer
Foam::vanAlbada
Foam::PLIC
Foam::PLICU
Foam::MPLIC
Foam::MPLICU
SourceFiles
interfaceCompression.C
\*---------------------------------------------------------------------------*/
#ifndef interfaceCompression_H
#define interfaceCompression_H
#include "surfaceInterpolationScheme.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class interfaceCompression Declaration
\*---------------------------------------------------------------------------*/
class interfaceCompressionNew
:
public surfaceInterpolationScheme<scalar>
{
// Private member data
const surfaceScalarField& phi_;
//- Base scheme to which the compression is applied
tmp<surfaceInterpolationScheme<scalar>> tScheme_;
//- Compression factor
const scalar cAlpha_;
public:
//- Runtime type information
TypeName("interfaceCompression");
// Constructors
//- Construct from faceFlux and Istream
interfaceCompressionNew
(
const fvMesh& mesh,
const surfaceScalarField& faceFlux,
Istream& is
)
:
surfaceInterpolationScheme<scalar>(mesh),
phi_(faceFlux),
tScheme_
(
surfaceInterpolationScheme<scalar>::New(mesh, faceFlux, is)
),
cAlpha_(readScalar(is))
{}
// Member Functions
//- Return the interpolation weighting factors
virtual tmp<surfaceScalarField> weights
(
const volScalarField&
) const
{
NotImplemented;
return tmp<surfaceScalarField>(nullptr);
}
//- Return the face-interpolate of the given cell field
virtual tmp<surfaceScalarField> interpolate
(
const volScalarField& vf
) const;
// Member Operators
//- Disallow default bitwise assignment
void operator=(const interfaceCompressionNew&) = delete;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -30,8 +30,7 @@ divSchemes
div(rhoPhi,U) Gauss linearUpwindV grad(U); div(rhoPhi,U) Gauss linearUpwindV grad(U);
div(rhoPhi,T) Gauss linearUpwind grad(T); div(rhoPhi,T) Gauss linearUpwind grad(T);
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(phi,p) Gauss upwind; div(phi,p) Gauss upwind;
div(rhoPhi,K) Gauss upwind; div(rhoPhi,K) Gauss upwind;

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 2; nAlphaCorr 2;
nAlphaSubCycles 1; nAlphaSubCycles 1;
cAlpha 1;
MULESCorr yes; MULESCorr yes;
nLimiterIter 5; nLimiterIter 5;

View File

@ -27,8 +27,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(rhoPhi,U) Gauss upwind; div(rhoPhi,U) Gauss upwind;
div(rhoPhi,T) Gauss upwind; div(rhoPhi,T) Gauss upwind;

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 1; nAlphaCorr 1;
nAlphaSubCycles 2; nAlphaSubCycles 2;
cAlpha 1;
MULESCorr no; MULESCorr no;
nLimiterIter 5; nLimiterIter 5;

View File

@ -27,8 +27,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(rhoPhi,U) Gauss upwind; div(rhoPhi,U) Gauss upwind;
div(rhoPhi,T) Gauss upwind; div(rhoPhi,T) Gauss upwind;

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 1; nAlphaCorr 1;
nAlphaSubCycles 2; nAlphaSubCycles 2;
cAlpha 1;
MULESCorr no; MULESCorr no;
nLimiterIter 5; nLimiterIter 5;

View File

@ -27,8 +27,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(rhoPhi,U) Gauss vanLeerV; div(rhoPhi,U) Gauss vanLeerV;
div(rhoPhi,T) Gauss vanLeer; div(rhoPhi,T) Gauss vanLeer;

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 1; nAlphaCorr 1;
nAlphaSubCycles 3; nAlphaSubCycles 3;
cAlpha 1.5;
} }
".*(rho|rhoFinal)" ".*(rho|rhoFinal)"

View File

@ -28,8 +28,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss linear; div(rhoPhi,U) Gauss linear;
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(phi,k) Gauss limitedLinear 1; div(phi,k) Gauss limitedLinear 1;
div(phi,B) Gauss limitedLinear 1; div(phi,B) Gauss limitedLinear 1;
div(B) Gauss linear; div(B) Gauss linear;

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 1; nAlphaCorr 1;
nAlphaSubCycles 4; nAlphaSubCycles 4;
cAlpha 1;
} }
"pcorr.*" "pcorr.*"

View File

@ -29,8 +29,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss linearUpwind grad(U); div(rhoPhi,U) Gauss linearUpwind grad(U);
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(phi,k) Gauss linearUpwind limitedGrad; div(phi,k) Gauss linearUpwind limitedGrad;
div(phi,omega) Gauss linearUpwind limitedGrad; div(phi,omega) Gauss linearUpwind limitedGrad;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear; div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 2; nAlphaCorr 2;
nAlphaSubCycles 1; nAlphaSubCycles 1;
cAlpha 1;
MULESCorr yes; MULESCorr yes;
nLimiterIter 10; nLimiterIter 10;

View File

@ -29,8 +29,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss linearUpwind grad(U); div(rhoPhi,U) Gauss linearUpwind grad(U);
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(phi,k) Gauss linearUpwind limitedGrad; div(phi,k) Gauss linearUpwind limitedGrad;
div(phi,omega) Gauss linearUpwind limitedGrad; div(phi,omega) Gauss linearUpwind limitedGrad;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear; div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 3; nAlphaCorr 3;
nAlphaSubCycles 1; nAlphaSubCycles 1;
cAlpha 1;
MULESCorr yes; MULESCorr yes;
nLimiterIter 15; nLimiterIter 15;

View File

@ -29,8 +29,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss linearUpwind grad(U); div(rhoPhi,U) Gauss linearUpwind grad(U);
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(phi,k) Gauss linearUpwind limitedGrad; div(phi,k) Gauss linearUpwind limitedGrad;
div(phi,omega) Gauss linearUpwind limitedGrad; div(phi,omega) Gauss linearUpwind limitedGrad;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear; div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 3; nAlphaCorr 3;
nAlphaSubCycles 1; nAlphaSubCycles 1;
cAlpha 1;
MULESCorr yes; MULESCorr yes;
nLimiterIter 15; nLimiterIter 15;

View File

@ -28,8 +28,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss upwind; div(rhoPhi,U) Gauss upwind;
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(phi,k) Gauss upwind; div(phi,k) Gauss upwind;
div(phi,epsilon) Gauss upwind; div(phi,epsilon) Gauss upwind;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear; div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 2; nAlphaCorr 2;
nAlphaSubCycles 1; nAlphaSubCycles 1;
cAlpha 1;
MULESCorr yes; MULESCorr yes;
nLimiterIter 3; nLimiterIter 3;

View File

@ -28,8 +28,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss linearUpwind grad(U); div(rhoPhi,U) Gauss linearUpwind grad(U);
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(phi,k) Gauss upwind; div(phi,k) Gauss upwind;
div(phi,epsilon) Gauss upwind; div(phi,epsilon) Gauss upwind;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear; div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 2; nAlphaCorr 2;
nAlphaSubCycles 1; nAlphaSubCycles 1;
cAlpha 1;
MULESCorr yes; MULESCorr yes;
nLimiterIter 3; nLimiterIter 3;

View File

@ -28,8 +28,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss linearUpwind grad(U); div(rhoPhi,U) Gauss linearUpwind grad(U);
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(phi,k) Gauss upwind; div(phi,k) Gauss upwind;
div(phi,epsilon) Gauss upwind; div(phi,epsilon) Gauss upwind;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear; div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 2; nAlphaCorr 2;
nAlphaSubCycles 1; nAlphaSubCycles 1;
cAlpha 1;
MULESCorr yes; MULESCorr yes;
nLimiterIter 3; nLimiterIter 3;

View File

@ -28,8 +28,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss vanLeerV; div(rhoPhi,U) Gauss vanLeerV;
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(phi,k) Gauss upwind; div(phi,k) Gauss upwind;
div(phi,epsilon) Gauss upwind; div(phi,epsilon) Gauss upwind;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear; div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 2; nAlphaCorr 2;
nAlphaSubCycles 1; nAlphaSubCycles 1;
cAlpha 1;
MULESCorr yes; MULESCorr yes;
nLimiterIter 5; nLimiterIter 5;

View File

@ -29,8 +29,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss linearUpwind grad(U); div(rhoPhi,U) Gauss linearUpwind grad(U);
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(phid1,p_rgh) Gauss upwind; div(phid1,p_rgh) Gauss upwind;
div(phid2,p_rgh) Gauss upwind; div(phid2,p_rgh) Gauss upwind;
div(rhoPhi,T) Gauss linearUpwind unlimited; div(rhoPhi,T) Gauss linearUpwind unlimited;

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 1; nAlphaCorr 1;
nAlphaSubCycles 2; nAlphaSubCycles 2;
cAlpha 1;
} }
p_rgh p_rgh

View File

@ -30,8 +30,7 @@ divSchemes
default none; default none;
div(rhoPhi,U) Gauss linearUpwind grad(U); div(rhoPhi,U) Gauss linearUpwind grad(U);
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
"div\(phi,(k|omega)\)" Gauss upwind; "div\(phi,(k|omega)\)" Gauss upwind;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear; div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 1; nAlphaCorr 1;
nAlphaSubCycles 1; nAlphaSubCycles 1;
cAlpha 1;
MULESCorr yes; MULESCorr yes;
nLimiterIter 3; nLimiterIter 3;

View File

@ -28,8 +28,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss linear; div(rhoPhi,U) Gauss linear;
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(phi,k) Gauss upwind; div(phi,k) Gauss upwind;
div(phi,epsilon) Gauss upwind; div(phi,epsilon) Gauss upwind;
div(phi,R) Gauss upwind; div(phi,R) Gauss upwind;

View File

@ -20,7 +20,6 @@ solvers
{ {
nAlphaCorr 1; nAlphaCorr 1;
nAlphaSubCycles 2; nAlphaSubCycles 2;
cAlpha 1;
} }
"pcorr.*" "pcorr.*"

View File

@ -28,8 +28,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss upwind; div(rhoPhi,U) Gauss upwind;
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear; div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;
} }

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 1; nAlphaCorr 1;
nAlphaSubCycles 2; nAlphaSubCycles 2;
cAlpha 1;
} }
"pcorr.*" "pcorr.*"

View File

@ -28,8 +28,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss linearUpwind grad(U); div(rhoPhi,U) Gauss linearUpwind grad(U);
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear; div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;
} }

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 2; nAlphaCorr 2;
nAlphaSubCycles 1; nAlphaSubCycles 1;
cAlpha 1;
MULESCorr yes; MULESCorr yes;
nLimiterIter 5; nLimiterIter 5;

View File

@ -28,8 +28,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss upwind; div(rhoPhi,U) Gauss upwind;
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(phi,k) Gauss upwind; div(phi,k) Gauss upwind;
div(phi,omega) Gauss upwind; div(phi,omega) Gauss upwind;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear; div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 1; nAlphaCorr 1;
nAlphaSubCycles 3; nAlphaSubCycles 3;
cAlpha 1;
} }
p_rgh p_rgh

View File

@ -28,8 +28,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss linear; div(rhoPhi,U) Gauss linear;
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear; div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;
} }

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 1; nAlphaCorr 1;
nAlphaSubCycles 2; nAlphaSubCycles 2;
cAlpha 1;
} }
"pcorr.*" "pcorr.*"

View File

@ -29,8 +29,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss linearUpwindV grad(U); div(rhoPhi,U) Gauss linearUpwindV grad(U);
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear; div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;
} }

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 2; nAlphaCorr 2;
nAlphaSubCycles 1; nAlphaSubCycles 1;
cAlpha 1;
MULESCorr yes; MULESCorr yes;
nLimiterIter 8; nLimiterIter 8;

View File

@ -28,8 +28,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss vanLeerV; div(rhoPhi,U) Gauss vanLeerV;
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear; div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;
} }

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 1; nAlphaCorr 1;
nAlphaSubCycles 3; nAlphaSubCycles 3;
cAlpha 1.5;
} }
"pcorr.*" "pcorr.*"

View File

@ -28,8 +28,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss vanLeerV; div(rhoPhi,U) Gauss vanLeerV;
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear; div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;
} }

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 1; nAlphaCorr 1;
nAlphaSubCycles 3; nAlphaSubCycles 3;
cAlpha 1.5;
} }
"pcorr.*" "pcorr.*"

View File

@ -28,8 +28,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss vanLeerV; div(rhoPhi,U) Gauss vanLeerV;
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear; div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;
} }

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 1; nAlphaCorr 1;
nAlphaSubCycles 3; nAlphaSubCycles 3;
cAlpha 1.5;
} }
"pcorr.*" "pcorr.*"

View File

@ -28,8 +28,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss vanLeerV; div(rhoPhi,U) Gauss vanLeerV;
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear; div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;
} }

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 1; nAlphaCorr 1;
nAlphaSubCycles 3; nAlphaSubCycles 3;
cAlpha 1.5;
} }
"pcorr.*" "pcorr.*"

View File

@ -28,8 +28,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss vanLeerV; div(rhoPhi,U) Gauss vanLeerV;
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear; div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;
} }

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 1; nAlphaCorr 1;
nAlphaSubCycles 3; nAlphaSubCycles 3;
cAlpha 1.5;
} }
"pcorr.*" "pcorr.*"

View File

@ -28,8 +28,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss vanLeerV; div(rhoPhi,U) Gauss vanLeerV;
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear; div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;
} }

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 1; nAlphaCorr 1;
nAlphaSubCycles 3; nAlphaSubCycles 3;
cAlpha 1;
} }
"pcorr.*" "pcorr.*"

View File

@ -34,8 +34,7 @@ gradSchemes
divSchemes divSchemes
{ {
div(rhoPhi,U) Gauss linearUpwindV grad(U); div(rhoPhi,U) Gauss linearUpwindV grad(U);
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear; div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;
} }

View File

@ -21,7 +21,6 @@ solvers
{ {
nAlphaCorr 2; nAlphaCorr 2;
nAlphaSubCycles 1; nAlphaSubCycles 1;
cAlpha 1;
MULESCorr yes; MULESCorr yes;
nLimiterIter 3; nLimiterIter 3;

View File

@ -30,8 +30,7 @@ divSchemes
div(rhoPhi,U) Gauss linearUpwind grad(U); div(rhoPhi,U) Gauss linearUpwind grad(U);
div(phi,omega) Gauss linearUpwind grad(omega); div(phi,omega) Gauss linearUpwind grad(omega);
div(phi,k) Gauss linearUpwind grad(k); div(phi,k) Gauss linearUpwind grad(k);
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear; div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;
} }

View File

@ -18,7 +18,6 @@ solvers
{ {
"alpha.water.*" "alpha.water.*"
{ {
cAlpha 0;
nAlphaCorr 2; nAlphaCorr 2;
nAlphaSubCycles 1; nAlphaSubCycles 1;

View File

@ -30,8 +30,7 @@ divSchemes
{ {
default none; default none;
div(phi,alpha) Gauss vanLeer; div(phi,alpha) Gauss interfaceCompression vanLeer 1;
div(phirb,alpha) Gauss linear;
div(rhoPhi,U) Gauss linearUpwind grad(U); div(rhoPhi,U) Gauss linearUpwind grad(U);
div(phi,k) Gauss upwind; div(phi,k) Gauss upwind;
div(phi,epsilon) Gauss upwind; div(phi,epsilon) Gauss upwind;

View File

@ -18,7 +18,6 @@ solvers
{ {
"alpha.water.*" "alpha.water.*"
{ {
cAlpha 0;
nAlphaCorr 2; nAlphaCorr 2;
nAlphaSubCycles 1; nAlphaSubCycles 1;