Files
OpenFOAM-12/applications/solvers/multiphase/driftFluxFoam/createFields.H
Henry Weller e8ff92cd67 CorrectPhi: Added pressureReference argument to set the reference cell
so that the same reference cell is used for pcorr and p or p_rgh to improve
consistency between flux and flux correction.
2021-04-30 21:07:15 +01:00

125 lines
2.1 KiB
C++

Info<< "Reading field p_rgh\n" << endl;
volScalarField p_rgh
(
IOobject
(
"p_rgh",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
);
Info<< "Reading field U\n" << endl;
volVectorField U
(
IOobject
(
"U",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
);
#include "createPhi.H"
Info<< "Reading incompressibleTwoPhaseInteractingMixture\n" << endl;
incompressibleTwoPhaseInteractingMixture mixture(U, phi);
volScalarField& alpha1(mixture.alpha1());
// Mixture density
volScalarField rho
(
IOobject
(
"rho",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
mixture.rho()
);
// Mass flux
surfaceScalarField rhoPhi
(
IOobject
(
"rhoPhi",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
fvc::interpolate(rho)*phi
);
// Relative Velocity
autoPtr<relativeVelocityModel> UdmModelPtr
(
relativeVelocityModel::New
(
mixture,
mixture
)
);
// Construct compressible turbulence model
autoPtr<compressible::momentumTransportModel> turbulence
(
compressible::momentumTransportModel::New(rho, U, rhoPhi, mixture)
);
#include "readGravitationalAcceleration.H"
#include "readhRef.H"
#include "gh.H"
volScalarField p
(
IOobject
(
"p",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
p_rgh + rho*gh
);
pressureReference pressureReference(p, p_rgh, pimple.dict());
if (p_rgh.needReference())
{
p += dimensionedScalar
(
"p",
p.dimensions(),
pressureReference.refValue()
- getRefCellValue(p, pressureReference.refCell())
);
p_rgh = p - rho*gh;
}
mesh.setFluxRequired(p_rgh.name());
mesh.setFluxRequired(alpha1.name());
// MULES Correction
tmp<surfaceScalarField> talphaPhiCorr0;
#include "createMRF.H"
#include "createFvModels.H"
#include "createFvConstraints.H"