mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
interDyMFoam, interPhaseChangeDyMFoam : Use 1/A from the previous time-step rather than 1 in the pcorr equation
This provides better convergence and consistency between p_rgh and pcorr
This commit is contained in:
@ -37,7 +37,8 @@ if (mesh.changing())
|
||||
pcorrTypes
|
||||
);
|
||||
|
||||
dimensionedScalar rAUf("rAUf", dimTime/rho.dimensions(), 1.0);
|
||||
// dimensionedScalar rAUf("rAUf", dimTime/rho.dimensions(), 1.0);
|
||||
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU));
|
||||
|
||||
while (pimple.correctNonOrthogonal())
|
||||
{
|
||||
|
||||
@ -57,6 +57,21 @@ int main(int argc, char *argv[])
|
||||
#include "createFields.H"
|
||||
#include "readTimeControls.H"
|
||||
#include "createPrghCorrTypes.H"
|
||||
|
||||
volScalarField rAU
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"rAU",
|
||||
runTime.timeName(),
|
||||
mesh,
|
||||
IOobject::READ_IF_PRESENT,
|
||||
IOobject::AUTO_WRITE
|
||||
),
|
||||
mesh,
|
||||
dimensionedScalar("rAUf", dimTime/rho.dimensions(), 1.0)
|
||||
);
|
||||
|
||||
#include "correctPhi.H"
|
||||
#include "createUf.H"
|
||||
#include "CourantNo.H"
|
||||
|
||||
@ -1,5 +1,5 @@
|
||||
{
|
||||
volScalarField rAU("rAU", 1.0/UEqn.A());
|
||||
rAU = 1.0/UEqn.A();
|
||||
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU));
|
||||
|
||||
volVectorField HbyA("HbyA", U);
|
||||
|
||||
@ -67,6 +67,21 @@ int main(int argc, char *argv[])
|
||||
#include "createFields.H"
|
||||
#include "readTimeControls.H"
|
||||
#include "createPcorrTypes.H"
|
||||
|
||||
volScalarField rAU
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"rAU",
|
||||
runTime.timeName(),
|
||||
mesh,
|
||||
IOobject::READ_IF_PRESENT,
|
||||
IOobject::AUTO_WRITE
|
||||
),
|
||||
mesh,
|
||||
dimensionedScalar("rAUf", dimTime/rho.dimensions(), 1.0)
|
||||
);
|
||||
|
||||
#include "../interFoam/interDyMFoam/correctPhi.H"
|
||||
#include "createUf.H"
|
||||
#include "CourantNo.H"
|
||||
|
||||
@ -1,5 +1,5 @@
|
||||
{
|
||||
volScalarField rAU("rAU", 1.0/UEqn.A());
|
||||
rAU = 1.0/UEqn.A();
|
||||
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU));
|
||||
|
||||
volVectorField HbyA("HbyA", U);
|
||||
|
||||
Reference in New Issue
Block a user