This repository has been archived on 2023-11-20. You can view files and clone it, but cannot push or open issues or pull requests.
foam-extend4.1-coherent-io/applications/solvers/multiphase/interDyMFoam/createFields.H

150 lines
3.1 KiB
C
Raw Normal View History

Info<< "Reading field pd\n" << endl;
volScalarField pd
(
IOobject
(
"pd",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
);
2010-09-22 18:13:13 +00:00
Info<< "Reading field alpha1\n" << endl;
volScalarField alpha1
(
IOobject
(
2010-09-22 18:13:13 +00:00
"alpha1",
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 transportProperties\n" << endl;
2010-09-22 18:13:13 +00:00
twoPhaseMixture twoPhaseProperties(U, phi, "alpha1");
const dimensionedScalar& rho1 = twoPhaseProperties.rho1();
const dimensionedScalar& rho2 = twoPhaseProperties.rho2();
// Need to store rho for ddt(rho, U)
volScalarField rho
(
IOobject
(
"rho",
runTime.timeName(),
mesh,
IOobject::READ_IF_PRESENT
),
2010-09-22 18:13:13 +00:00
alpha1*rho1 + (scalar(1) - alpha1)*rho2,
alpha1.boundaryField().types()
);
rho.oldTime();
// Mass flux
// Initialisation does not matter because rhoPhi is reset after the
2010-09-22 18:13:13 +00:00
// alpha1 solution before it is used in the U equation.
surfaceScalarField rhoPhi
(
IOobject
(
"rho*phi",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
rho1*phi
);
2010-09-22 18:13:13 +00:00
// Construct interface from alpha1 distribution
interfaceProperties interface(alpha1, U, twoPhaseProperties);
2010-09-22 18:13:13 +00:00
// Construct incompressible turbulence model
autoPtr<incompressible::turbulenceModel> turbulence
(
2010-09-22 18:13:13 +00:00
incompressible::turbulenceModel::New(U, phi, twoPhaseProperties)
);
2010-09-22 18:13:13 +00:00
wordList pcorrTypes
(
pd.boundaryField().size(),
zeroGradientFvPatchScalarField::typeName
);
2010-09-22 18:13:13 +00:00
for (label i = 0; i < pd.boundaryField().size(); i++)
{
if (pd.boundaryField()[i].fixesValue())
{
pcorrTypes[i] = fixedValueFvPatchScalarField::typeName;
}
}
volScalarField p
(
IOobject
(
"p",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
pd + rho*(g & mesh.C())
);
label pdRefCell = 0;
scalar pdRefValue = 0.0;
setRefCell
(
pd,
mesh.solutionDict().subDict("PIMPLE"),
pdRefCell,
pdRefValue
);
scalar pRefValue = 0.0;
if (pd.needReference())
{
pRefValue = readScalar
(
mesh.solutionDict().subDict("PIMPLE").lookup("pRefValue")
);
p += dimensionedScalar
(
"p",
p.dimensions(),
pRefValue - getRefCellValue(p, pdRefCell)
);
}