VoF Solvers: Relocate the correction, sub-cycling and compressions controls from PIMPLE to the alpha1 sub-dict

MULES: Add support for explicit correction
interPhaseChangeFoam: Add option for explicit MULES or as correction to an upwind solution
Deprecate implicit form of MULES
This commit is contained in:
Henry
2013-05-03 15:46:29 +01:00
parent 24e36ddf3f
commit d13c9810e8
24 changed files with 635 additions and 157 deletions

View File

@ -48,6 +48,25 @@ void Foam::MULES::explicitSolve
}
void Foam::MULES::correct
(
volScalarField& psi,
surfaceScalarField& phiPsiCorr,
const scalar psiMax,
const scalar psiMin
)
{
correct
(
geometricOneField(),
psi,
phiPsiCorr,
zeroField(), zeroField(),
psiMax, psiMin
);
}
void Foam::MULES::limitSum(UPtrList<scalarField>& phiPsiCorrs)
{
forAll(phiPsiCorrs[0], facei)

View File

@ -92,6 +92,38 @@ void explicitSolve
const scalar psiMin
);
template<class RhoType, class SpType, class SuType>
void correct
(
const RhoType& rho,
volScalarField& psi,
const surfaceScalarField& phiCorr,
const SpType& Sp,
const SuType& Su
);
template<class RhoType, class SpType, class SuType>
void correct
(
const RhoType& rho,
volScalarField& psi,
surfaceScalarField& phiCorr,
const SpType& Sp,
const SuType& Su,
const scalar psiMax,
const scalar psiMin
);
void correct
(
volScalarField& psi,
surfaceScalarField& phiCorr,
const scalar psiMax,
const scalar psiMin
);
template<class RhoType, class SpType, class SuType>
void limiter
(
@ -122,6 +154,35 @@ void limit
const bool returnCorr
);
template<class RhoType, class SpType, class SuType>
void limiterCorr
(
scalarField& allLambda,
const RhoType& rho,
const volScalarField& psi,
const surfaceScalarField& phiCorr,
const SpType& Sp,
const SuType& Su,
const scalar psiMax,
const scalar psiMin,
const label nLimiterIter
);
template<class RhoType, class SpType, class SuType>
void limitCorr
(
const RhoType& rho,
const volScalarField& psi,
surfaceScalarField& phiCorr,
const SpType& Sp,
const SuType& Su,
const scalar psiMax,
const scalar psiMin,
const label nLimiterIter
);
void limitSum(UPtrList<scalarField>& phiPsiCorrs);
template<class SurfaceScalarFieldList>

View File

@ -96,6 +96,74 @@ void Foam::MULES::explicitSolve
}
template<class RhoType, class SpType, class SuType>
void Foam::MULES::correct
(
const RhoType& rho,
volScalarField& psi,
const surfaceScalarField& phiCorr,
const SpType& Sp,
const SuType& Su
)
{
Info<< "MULES: Correcting " << psi.name() << endl;
const fvMesh& mesh = psi.mesh();
const scalar deltaT = mesh.time().deltaTValue();
scalarField psiIf(psi.size(), 0);
fvc::surfaceIntegrate(psiIf, phiCorr);
if (mesh.moving())
{
psi.internalField() =
(
rho.field()*psi.internalField()/deltaT
+ Su.field()
- psiIf
)/(rho.field()/deltaT - Sp.field());
}
else
{
psi.internalField() =
(
rho.field()*psi.internalField()/deltaT
+ Su.field()
- psiIf
)/(rho.field()/deltaT - Sp.field());
}
psi.correctBoundaryConditions();
}
template<class RhoType, class SpType, class SuType>
void Foam::MULES::correct
(
const RhoType& rho,
volScalarField& psi,
surfaceScalarField& phiCorr,
const SpType& Sp,
const SuType& Su,
const scalar psiMax,
const scalar psiMin
)
{
const fvMesh& mesh = psi.mesh();
const dictionary& MULEScontrols = mesh.solverDict(psi.name());
label nLimiterIter
(
readLabel(MULEScontrols.lookup("nLimiterIter"))
);
limitCorr(rho, psi, phiCorr, Sp, Su, psiMax, psiMin, nLimiterIter);
correct(rho, psi, phiCorr, Sp, Su);
}
template<class RhoType, class SpType, class SuType>
void Foam::MULES::limiter
(
@ -511,6 +579,364 @@ void Foam::MULES::limit
}
template<class RhoType, class SpType, class SuType>
void Foam::MULES::limiterCorr
(
scalarField& allLambda,
const RhoType& rho,
const volScalarField& psi,
const surfaceScalarField& phiCorr,
const SpType& Sp,
const SuType& Su,
const scalar psiMax,
const scalar psiMin,
const label nLimiterIter
)
{
const scalarField& psiIf = psi;
const volScalarField::GeometricBoundaryField& psiBf = psi.boundaryField();
const fvMesh& mesh = psi.mesh();
const labelUList& owner = mesh.owner();
const labelUList& neighb = mesh.neighbour();
tmp<volScalarField::DimensionedInternalField> tVsc = mesh.Vsc();
const scalarField& V = tVsc();
const scalar deltaT = mesh.time().deltaTValue();
const scalarField& phiCorrIf = phiCorr;
const surfaceScalarField::GeometricBoundaryField& phiCorrBf =
phiCorr.boundaryField();
slicedSurfaceScalarField lambda
(
IOobject
(
"lambda",
mesh.time().timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE,
false
),
mesh,
dimless,
allLambda,
false // Use slices for the couples
);
scalarField& lambdaIf = lambda;
surfaceScalarField::GeometricBoundaryField& lambdaBf =
lambda.boundaryField();
scalarField psiMaxn(psiIf.size(), psiMin);
scalarField psiMinn(psiIf.size(), psiMax);
scalarField sumPhip(psiIf.size(), VSMALL);
scalarField mSumPhim(psiIf.size(), VSMALL);
forAll(phiCorrIf, facei)
{
label own = owner[facei];
label nei = neighb[facei];
psiMaxn[own] = max(psiMaxn[own], psiIf[nei]);
psiMinn[own] = min(psiMinn[own], psiIf[nei]);
psiMaxn[nei] = max(psiMaxn[nei], psiIf[own]);
psiMinn[nei] = min(psiMinn[nei], psiIf[own]);
scalar phiCorrf = phiCorrIf[facei];
if (phiCorrf > 0.0)
{
sumPhip[own] += phiCorrf;
mSumPhim[nei] += phiCorrf;
}
else
{
mSumPhim[own] -= phiCorrf;
sumPhip[nei] -= phiCorrf;
}
}
forAll(phiCorrBf, patchi)
{
const fvPatchScalarField& psiPf = psiBf[patchi];
const scalarField& phiCorrPf = phiCorrBf[patchi];
const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
if (psiPf.coupled())
{
const scalarField psiPNf(psiPf.patchNeighbourField());
forAll(phiCorrPf, pFacei)
{
label pfCelli = pFaceCells[pFacei];
psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPNf[pFacei]);
psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPNf[pFacei]);
}
}
else
{
forAll(phiCorrPf, pFacei)
{
label pfCelli = pFaceCells[pFacei];
psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPf[pFacei]);
psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPf[pFacei]);
}
}
forAll(phiCorrPf, pFacei)
{
label pfCelli = pFaceCells[pFacei];
scalar phiCorrf = phiCorrPf[pFacei];
if (phiCorrf > 0.0)
{
sumPhip[pfCelli] += phiCorrf;
}
else
{
mSumPhim[pfCelli] -= phiCorrf;
}
}
}
psiMaxn = min(psiMaxn, psiMax);
psiMinn = max(psiMinn, psiMin);
// scalar smooth = 0.5;
// psiMaxn = min((1.0 - smooth)*psiIf + smooth*psiMaxn, psiMax);
// psiMinn = max((1.0 - smooth)*psiIf + smooth*psiMinn, psiMin);
psiMaxn =
V
*(
(rho.field()/deltaT - Sp.field())*psiMaxn
- Su.field()
- rho.field()*psi.internalField()/deltaT
);
psiMinn =
V
*(
Su.field()
- (rho.field()/deltaT - Sp.field())*psiMinn
+ rho.field()*psi.internalField()/deltaT
);
scalarField sumlPhip(psiIf.size());
scalarField mSumlPhim(psiIf.size());
for (int j=0; j<nLimiterIter; j++)
{
sumlPhip = 0.0;
mSumlPhim = 0.0;
forAll(lambdaIf, facei)
{
label own = owner[facei];
label nei = neighb[facei];
scalar lambdaPhiCorrf = lambdaIf[facei]*phiCorrIf[facei];
if (lambdaPhiCorrf > 0.0)
{
sumlPhip[own] += lambdaPhiCorrf;
mSumlPhim[nei] += lambdaPhiCorrf;
}
else
{
mSumlPhim[own] -= lambdaPhiCorrf;
sumlPhip[nei] -= lambdaPhiCorrf;
}
}
forAll(lambdaBf, patchi)
{
scalarField& lambdaPf = lambdaBf[patchi];
const scalarField& phiCorrfPf = phiCorrBf[patchi];
const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
forAll(lambdaPf, pFacei)
{
label pfCelli = pFaceCells[pFacei];
scalar lambdaPhiCorrf = lambdaPf[pFacei]*phiCorrfPf[pFacei];
if (lambdaPhiCorrf > 0.0)
{
sumlPhip[pfCelli] += lambdaPhiCorrf;
}
else
{
mSumlPhim[pfCelli] -= lambdaPhiCorrf;
}
}
}
forAll(sumlPhip, celli)
{
sumlPhip[celli] =
max(min
(
(sumlPhip[celli] + psiMaxn[celli])/mSumPhim[celli],
1.0), 0.0
);
mSumlPhim[celli] =
max(min
(
(mSumlPhim[celli] + psiMinn[celli])/sumPhip[celli],
1.0), 0.0
);
}
const scalarField& lambdam = sumlPhip;
const scalarField& lambdap = mSumlPhim;
forAll(lambdaIf, facei)
{
if (phiCorrIf[facei] > 0.0)
{
lambdaIf[facei] = min
(
lambdaIf[facei],
min(lambdap[owner[facei]], lambdam[neighb[facei]])
);
}
else
{
lambdaIf[facei] = min
(
lambdaIf[facei],
min(lambdam[owner[facei]], lambdap[neighb[facei]])
);
}
}
forAll(lambdaBf, patchi)
{
fvsPatchScalarField& lambdaPf = lambdaBf[patchi];
const scalarField& phiCorrfPf = phiCorrBf[patchi];
const fvPatchScalarField& psiPf = psiBf[patchi];
if (isA<wedgeFvPatch>(mesh.boundary()[patchi]))
{
lambdaPf = 0;
}
else if (psiPf.coupled())
{
const labelList& pFaceCells =
mesh.boundary()[patchi].faceCells();
forAll(lambdaPf, pFacei)
{
label pfCelli = pFaceCells[pFacei];
if (phiCorrfPf[pFacei] > 0.0)
{
lambdaPf[pFacei] =
min(lambdaPf[pFacei], lambdap[pfCelli]);
}
else
{
lambdaPf[pFacei] =
min(lambdaPf[pFacei], lambdam[pfCelli]);
}
}
}
else
{
const labelList& pFaceCells =
mesh.boundary()[patchi].faceCells();
const scalarField& phiCorrPf = phiCorrBf[patchi];
forAll(lambdaPf, pFacei)
{
// Limit outlet faces only
if (phiCorrPf[pFacei] > 0)
{
label pfCelli = pFaceCells[pFacei];
if (phiCorrfPf[pFacei] > 0.0)
{
lambdaPf[pFacei] =
min(lambdaPf[pFacei], lambdap[pfCelli]);
}
else
{
lambdaPf[pFacei] =
min(lambdaPf[pFacei], lambdam[pfCelli]);
}
}
}
}
}
syncTools::syncFaceList(mesh, allLambda, minEqOp<scalar>());
}
}
template<class RhoType, class SpType, class SuType>
void Foam::MULES::limitCorr
(
const RhoType& rho,
const volScalarField& psi,
surfaceScalarField& phiCorr,
const SpType& Sp,
const SuType& Su,
const scalar psiMax,
const scalar psiMin,
const label nLimiterIter
)
{
const fvMesh& mesh = psi.mesh();
scalarField allLambda(mesh.nFaces(), 1.0);
slicedSurfaceScalarField lambda
(
IOobject
(
"lambda",
mesh.time().timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE,
false
),
mesh,
dimless,
allLambda,
false // Use slices for the couples
);
limiterCorr
(
allLambda,
rho,
psi,
phiCorr,
Sp,
Su,
psiMax,
psiMin,
nLimiterIter
);
phiCorr *= lambda;
}
template<class SurfaceScalarFieldList>
void Foam::MULES::limitSum(SurfaceScalarFieldList& phiPsiCorrs)
{

View File

@ -159,7 +159,7 @@ Foam::interfaceProperties::interfaceProperties
(
readScalar
(
alpha1.mesh().solutionDict().subDict("PIMPLE").lookup("cAlpha")
alpha1.mesh().solverDict(alpha1.name()).lookup("cAlpha")
)
),
sigma_(dict.lookup("sigma")),