openfoam/applications/solvers/multiphase/interPhaseChangeFoam/createFields.H
sergio 499933dbab ENH: Adding features for phase change solvers
1) Adding interfaceHeight FO
2) Adding interfaceHeatResistance mass transfer model to
   interCondensatingEvaporatingFoam with spread source approach
3) Reworking framework for icoReactingMultiphaseInterFoam
2020-02-19 14:18:25 +00:00

128 lines
2.2 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<< "Creating phaseChangeTwoPhaseMixture\n" << endl;
autoPtr<phaseChangeTwoPhaseMixture> mixture =
phaseChangeTwoPhaseMixture::New(U, phi);
volScalarField& alpha1(mixture->alpha1());
volScalarField& alpha2(mixture->alpha2());
const dimensionedScalar& rho1 = mixture->rho1();
const dimensionedScalar& rho2 = mixture->rho2();
// Need to store rho for ddt(rho, U)
volScalarField rho
(
IOobject
(
"rho",
runTime.timeName(),
mesh,
IOobject::READ_IF_PRESENT
),
alpha1*rho1 + alpha2*rho2
);
rho.oldTime();
// Construct interface from alpha1 distribution
interfaceProperties interface(alpha1, U, mixture());
// Construct incompressible turbulence model
autoPtr<incompressible::turbulenceModel> turbulence
(
incompressible::turbulenceModel::New(U, phi, 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
);
label pRefCell = 0;
scalar pRefValue = 0.0;
setRefCell
(
p,
p_rgh,
pimple.dict(),
pRefCell,
pRefValue
);
if (p_rgh.needReference())
{
p += dimensionedScalar
(
"p",
p.dimensions(),
pRefValue - getRefCellValue(p, pRefCell)
);
p_rgh = p - rho*gh;
}
mesh.setFluxRequired(p_rgh.name());
mesh.setFluxRequired(alpha1.name());
#include "createFvOptions.H"
IOobject alphaPhi10Header
(
IOobject::groupName("alphaPhi0", alpha1.group()),
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
);
// MULES flux from previous time-step
surfaceScalarField alphaPhi10
(
alphaPhi10Header,
phi*fvc::interpolate(alpha1)
);