Rationalise all correctPhi implementations

This commit is contained in:
Henry
2013-09-11 22:59:00 +01:00
parent 6f62380b65
commit e9b6c5f293
6 changed files with 74 additions and 68 deletions

View File

@ -1,30 +1,30 @@
if (mesh.changing())
{
if (mesh.changing())
forAll(U.boundaryField(), patchi)
{
forAll(U.boundaryField(), patchi)
if (U.boundaryField()[patchi].fixesValue())
{
if (U.boundaryField()[patchi].fixesValue())
{
U.boundaryField()[patchi].initEvaluate();
}
}
forAll(U.boundaryField(), patchi)
{
if (U.boundaryField()[patchi].fixesValue())
{
U.boundaryField()[patchi].evaluate();
phi.boundaryField()[patchi] =
rho.boundaryField()[patchi]
*(
U.boundaryField()[patchi]
& mesh.Sf().boundaryField()[patchi]
);
}
U.boundaryField()[patchi].initEvaluate();
}
}
forAll(U.boundaryField(), patchi)
{
if (U.boundaryField()[patchi].fixesValue())
{
U.boundaryField()[patchi].evaluate();
phi.boundaryField()[patchi] =
rho.boundaryField()[patchi]
*(
U.boundaryField()[patchi]
& mesh.Sf().boundaryField()[patchi]
);
}
}
}
{
volScalarField pcorr
(
IOobject
@ -40,13 +40,13 @@
pcorrTypes
);
dimensionedScalar Dp("Dp", dimTime, 1.0);
dimensionedScalar rAUf("rAUf", dimTime, 1.0);
while (pimple.correctNonOrthogonal())
{
fvScalarMatrix pcorrEqn
(
fvm::laplacian(Dp, pcorr) == fvc::div(phi) - divrhoU
fvm::laplacian(rAUf, pcorr) == fvc::div(phi) - divrhoU
);
pcorrEqn.solve();

View File

@ -37,13 +37,13 @@ if (mesh.changing())
pcorrTypes
);
dimensionedScalar Dp("Dp", dimTime, 1.0);
dimensionedScalar rAUf("rAUf", dimTime, 1.0);
while (pimple.correctNonOrthogonal())
{
fvScalarMatrix pcorrEqn
(
fvm::laplacian(Dp, pcorr) == fvc::div(phi)
fvm::laplacian(rAUf, pcorr) == fvc::div(phi)
);
pcorrEqn.setReference(pRefCell, pRefValue);
@ -54,6 +54,6 @@ if (mesh.changing())
phi -= pcorrEqn.flux();
}
}
}
#include "continuityErrs.H"
#include "continuityErrs.H"
}

View File

@ -1,3 +1,26 @@
if (mesh.changing())
{
forAll(U.boundaryField(), patchI)
{
if (U.boundaryField()[patchI].fixesValue())
{
U.boundaryField()[patchI].initEvaluate();
}
}
forAll(U.boundaryField(), patchI)
{
if (U.boundaryField()[patchI].fixesValue())
{
U.boundaryField()[patchI].evaluate();
phiv.boundaryField()[patchI] =
U.boundaryField()[patchI]
& mesh.Sf().boundaryField()[patchI];
}
}
}
{
volScalarField pcorr
(

View File

@ -1,28 +1,26 @@
if (mesh.changing())
{
if (mesh.changing())
forAll(U.boundaryField(), patchi)
{
forAll(U.boundaryField(), patchi)
if (U.boundaryField()[patchi].fixesValue())
{
if (U.boundaryField()[patchi].fixesValue())
{
U.boundaryField()[patchi].initEvaluate();
}
}
forAll(U.boundaryField(), patchi)
{
if (U.boundaryField()[patchi].fixesValue())
{
U.boundaryField()[patchi].evaluate();
phi.boundaryField()[patchi] =
U.boundaryField()[patchi] & mesh.Sf().boundaryField()[patchi];
}
U.boundaryField()[patchi].initEvaluate();
}
}
#include "continuityErrs.H"
forAll(U.boundaryField(), patchi)
{
if (U.boundaryField()[patchi].fixesValue())
{
U.boundaryField()[patchi].evaluate();
phi.boundaryField()[patchi] =
U.boundaryField()[patchi] & mesh.Sf().boundaryField()[patchi];
}
}
}
{
volScalarField pcorr
(
IOobject
@ -38,13 +36,13 @@
pcorrTypes
);
dimensionedScalar Dp("Dp", dimTime/rho.dimensions(), 1.0);
dimensionedScalar rAUf("rAUf", dimTime/rho.dimensions(), 1.0);
while (pimple.correctNonOrthogonal())
{
fvScalarMatrix pcorrEqn
(
fvm::laplacian(Dp, pcorr) == fvc::div(phi) - divU
fvm::laplacian(rAUf, pcorr) == fvc::div(phi) - divU
);
pcorrEqn.solve();
@ -54,8 +52,8 @@
phi -= pcorrEqn.flux();
}
}
phi.oldTime() = phi;
#include "continuityErrs.H"
}
phi.oldTime() = phi;
#include "continuityErrs.H"

View File

@ -1,20 +1,6 @@
{
#include "continuityErrs.H"
wordList pcorrTypes
(
p.boundaryField().size(),
zeroGradientFvPatchScalarField::typeName
);
forAll (p.boundaryField(), i)
{
if (p.boundaryField()[i].fixesValue())
{
pcorrTypes[i] = fixedValueFvPatchScalarField::typeName;
}
}
volScalarField pcorr
(
IOobject
@ -30,9 +16,7 @@
pcorrTypes
);
dimensionedScalar rAUf("(1|A(U))", dimTime/rho.dimensions(), 1.0);
//adjustPhi(phi, U, pcorr);
dimensionedScalar rAUf("rAUf", dimTime/rho.dimensions(), 1.0);
while (pimple.correctNonOrthogonal())
{

View File

@ -57,6 +57,7 @@ int main(int argc, char *argv[])
#include "createMRFZones.H"
#include "initContinuityErrs.H"
#include "readTimeControls.H"
#include "createPcorrTypes.H"
#include "correctPhi.H"
#include "CourantNo.H"
#include "setInitialDeltaT.H"