mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
1) Implementation of the compressibleIsoInterFOam solver 2) Implementation of a new PLIC interpolation scheme. 3) New tutorials associated with the solvers This implementation was carried out by Henning Scheufler (DLR) and Johan Roenby (DHI), following : \verbatim Henning Scheufler, Johan Roenby, Accurate and efficient surface reconstruction from volume fraction data on general meshes, Journal of Computational Physics, 2019, doi 10.1016/j.jcp.2019.01.009 \endverbatim The integration of the code was carried out by Andy Heather and Sergio Ferraris from OpenCFD Ltd.
111 lines
1.9 KiB
C
111 lines
1.9 KiB
C
#include "createRDeltaT.H"
|
|
|
|
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<< "Constructing twoPhaseMixtureThermo\n" << endl;
|
|
twoPhaseMixtureThermo mixture(U, phi);
|
|
|
|
|
|
volScalarField& alpha1(mixture.alpha1());
|
|
volScalarField& alpha2(mixture.alpha2());
|
|
|
|
Info<< "Reading thermophysical properties\n" << endl;
|
|
|
|
const volScalarField& rho1 = mixture.thermo1().rho();
|
|
const volScalarField& rho2 = mixture.thermo2().rho();
|
|
|
|
volScalarField rho
|
|
(
|
|
IOobject
|
|
(
|
|
"rho",
|
|
runTime.timeName(),
|
|
mesh,
|
|
IOobject::READ_IF_PRESENT,
|
|
IOobject::AUTO_WRITE
|
|
),
|
|
alpha1*rho1 + alpha2*rho2
|
|
);
|
|
|
|
|
|
dimensionedScalar pMin
|
|
(
|
|
"pMin",
|
|
dimPressure,
|
|
mixture
|
|
);
|
|
|
|
mesh.setFluxRequired(p_rgh.name());
|
|
mesh.setFluxRequired(alpha1.name());
|
|
|
|
|
|
#include "readGravitationalAcceleration.H"
|
|
#include "readhRef.H"
|
|
#include "gh.H"
|
|
|
|
|
|
// Mass flux
|
|
// Initialisation does not matter because rhoPhi is reset after the
|
|
// alpha1 solution before it is used in the U equation.
|
|
surfaceScalarField rhoPhi
|
|
(
|
|
IOobject
|
|
(
|
|
"rhoPhi",
|
|
runTime.timeName(),
|
|
mesh,
|
|
IOobject::NO_READ,
|
|
IOobject::NO_WRITE
|
|
),
|
|
fvc::interpolate(rho)*phi
|
|
);
|
|
|
|
volScalarField dgdt(alpha1*fvc::div(phi));
|
|
|
|
#include "createAlphaFluxes.H"
|
|
|
|
Foam::isoAdvection advector(alpha1,phi,U);
|
|
// Construct compressible turbulence model
|
|
compressibleInterPhaseTransportModel turbulence
|
|
(
|
|
rho,
|
|
U,
|
|
phi,
|
|
rhoPhi,
|
|
alphaPhi10,
|
|
mixture
|
|
);
|
|
|
|
#include "createK.H"
|
|
|
|
#include "createMRF.H"
|
|
#include "createFvOptions.H"
|