Skip to content
Snippets Groups Projects
Commit 6b2fb466 authored by Henry's avatar Henry
Browse files

twoPhaseEulerFoam: Update only the fixed-value phi patch fields before...

twoPhaseEulerFoam: Update only the fixed-value phi patch fields before constructing the pressure eqn
Avoids small continuity error in parallel
parent 3ebeea2b
Branches
Tags
No related merge requests found
...@@ -22,8 +22,6 @@ volScalarField rAU2 ...@@ -22,8 +22,6 @@ volScalarField rAU2
) )
); );
//surfaceScalarField alpharAUf1(fvc::interpolate(alpha1*rAU1));
//surfaceScalarField alpharAUf2(fvc::interpolate(alpha2*rAU2));
surfaceScalarField alpharAUf1 surfaceScalarField alpharAUf1
( (
fvc::interpolate(max(alpha1, fluid.residualAlpha(phase1))*rAU1) fvc::interpolate(max(alpha1, fluid.residualAlpha(phase1))*rAU1)
...@@ -94,10 +92,26 @@ while (pimple.correct()) ...@@ -94,10 +92,26 @@ while (pimple.correct())
#include "correctContErrs.H" #include "correctContErrs.H"
// Correct flux BCs to be consistent with the velocity BCs // Correct flux BCs to be consistent with the velocity BCs
phi1.boundaryField() == forAll(mesh.boundary(), patchi)
mrfZones.relative(mesh.Sf().boundaryField() & U1.boundaryField()); {
phi2.boundaryField() == if (isA<fixedValueFvsPatchScalarField>(phi1.boundaryField()[patchi]))
mrfZones.relative(mesh.Sf().boundaryField() & U2.boundaryField()); {
phi1.boundaryField()[patchi] ==
mrfZones.relative
(
mesh.Sf().boundaryField() & U1.boundaryField()
)()[patchi];
}
if (isA<fixedValueFvsPatchScalarField>(phi2.boundaryField()[patchi]))
{
phi2.boundaryField()[patchi] ==
mrfZones.relative
(
mesh.Sf().boundaryField() & U2.boundaryField()
)()[patchi];
}
}
volVectorField HbyA1 volVectorField HbyA1
( (
...@@ -227,12 +241,9 @@ while (pimple.correct()) ...@@ -227,12 +241,9 @@ while (pimple.correct())
p_rgh.boundaryField(), p_rgh.boundaryField(),
( (
phiHbyA.boundaryField() phiHbyA.boundaryField()
- mrfZones.relative - (
( alphaf1.boundaryField()*phi1.boundaryField()
alphaf1.boundaryField() + alphaf2.boundaryField()*phi2.boundaryField()
*(mesh.Sf().boundaryField() & U1.boundaryField())
+ alphaf2.boundaryField()
*(mesh.Sf().boundaryField() & U2.boundaryField())
) )
)/(mesh.magSf().boundaryField()*rAUf.boundaryField()) )/(mesh.magSf().boundaryField()*rAUf.boundaryField())
); );
......
...@@ -98,11 +98,27 @@ while (pimple.correct()) ...@@ -98,11 +98,27 @@ while (pimple.correct())
surfaceScalarField rhof1(fvc::interpolate(rho1)); surfaceScalarField rhof1(fvc::interpolate(rho1));
surfaceScalarField rhof2(fvc::interpolate(rho2)); surfaceScalarField rhof2(fvc::interpolate(rho2));
// Correct flux BCs to be consistent with the velocity BCs // Correct fixed-flux BCs to be consistent with the velocity BCs
phi1.boundaryField() == forAll(mesh.boundary(), patchi)
mrfZones.relative(mesh.Sf().boundaryField() & U1.boundaryField()); {
phi2.boundaryField() == if (isA<fixedValueFvsPatchScalarField>(phi1.boundaryField()[patchi]))
mrfZones.relative(mesh.Sf().boundaryField() & U2.boundaryField()); {
phi1.boundaryField()[patchi] ==
mrfZones.relative
(
mesh.Sf().boundaryField() & U1.boundaryField()
)()[patchi];
}
if (isA<fixedValueFvsPatchScalarField>(phi2.boundaryField()[patchi]))
{
phi2.boundaryField()[patchi] ==
mrfZones.relative
(
mesh.Sf().boundaryField() & U2.boundaryField()
)()[patchi];
}
}
surfaceScalarField alpharAUf1 surfaceScalarField alpharAUf1
( (
......
...@@ -30,6 +30,7 @@ License ...@@ -30,6 +30,7 @@ License
#include "PhaseCompressibleTurbulenceModel.H" #include "PhaseCompressibleTurbulenceModel.H"
#include "dragModel.H" #include "dragModel.H"
#include "heatTransferModel.H" #include "heatTransferModel.H"
#include "fixedValueFvsPatchFields.H"
#include "fixedValueFvPatchFields.H" #include "fixedValueFvPatchFields.H"
#include "slipFvPatchFields.H" #include "slipFvPatchFields.H"
#include "partialSlipFvPatchFields.H" #include "partialSlipFvPatchFields.H"
...@@ -151,7 +152,7 @@ Foam::phaseModel::phaseModel ...@@ -151,7 +152,7 @@ Foam::phaseModel::phaseModel
|| isA<partialSlipFvPatchVectorField>(U_.boundaryField()[i]) || isA<partialSlipFvPatchVectorField>(U_.boundaryField()[i])
) )
{ {
phiTypes[i] = fixedValueFvPatchScalarField::typeName; phiTypes[i] = fixedValueFvsPatchScalarField::typeName;
} }
} }
......
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment