phaseScalarTransport: Improved interface and documentation
This function now looks up an alphaPhi field by default and exits with an error if the field cannot be found. In order to solve for alphaPhi a new 'solveAlphaPhi' switch has to be set. The documentation has been updated to reflect the fact that the VoF solvers now store all the alphaPhi fluxes necessary for a in-phase equation to be constructed. The phaseScalarTransport function has been added to the damBreakLaminar tutorial.
This commit is contained in:
@ -90,7 +90,7 @@ Foam::volScalarField& Foam::functionObjects::phaseScalarTransport::Phi()
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"Phi" + s_.name(),
|
||||
typedName(IOobject::groupName("Phi", phaseName_)),
|
||||
time_.name(),
|
||||
mesh_,
|
||||
IOobject::READ_IF_PRESENT,
|
||||
@ -112,16 +112,11 @@ Foam::volScalarField& Foam::functionObjects::phaseScalarTransport::Phi()
|
||||
Foam::tmp<Foam::surfaceScalarField>
|
||||
Foam::functionObjects::phaseScalarTransport::alphaPhi()
|
||||
{
|
||||
// If alphaPhi exists then return it
|
||||
if (mesh_.foundObject<surfaceScalarField>(alphaPhiName_))
|
||||
if (!solveAlphaPhi_)
|
||||
{
|
||||
return mesh_.lookupObject<surfaceScalarField>(alphaPhiName_);
|
||||
}
|
||||
|
||||
// Otherwise generate it ...
|
||||
Info<< type() << ": " << surfaceScalarField::typeName << " "
|
||||
<< alphaPhiName_ << " was not found, so generating it" << endl;
|
||||
|
||||
const volScalarField& alpha =
|
||||
mesh_.lookupObject<volScalarField>(alphaName_);
|
||||
const surfaceScalarField& phi =
|
||||
@ -288,8 +283,6 @@ Foam::functionObjects::phaseScalarTransport::phaseScalarTransport
|
||||
fvMeshFunctionObject(name, runTime, dict),
|
||||
fieldName_(dict.lookup("field")),
|
||||
phaseName_(IOobject::group(fieldName_)),
|
||||
nCorr_(0),
|
||||
residualAlpha_(rootSmall),
|
||||
s_
|
||||
(
|
||||
IOobject
|
||||
@ -330,18 +323,20 @@ bool Foam::functionObjects::phaseScalarTransport::read(const dictionary& dict)
|
||||
{
|
||||
fvMeshFunctionObject::read(dict);
|
||||
|
||||
solveAlphaPhi_ = dict.lookupOrDefault<bool>("solveAlphaPhi", false);
|
||||
|
||||
alphaName_ =
|
||||
dict.lookupOrDefault<word>
|
||||
(
|
||||
"alpha",
|
||||
IOobject::groupName("alpha", phaseName_)
|
||||
);
|
||||
const word defaultAlphaPhiName =
|
||||
IOobject::groupName("alphaPhi", phaseName_);
|
||||
alphaPhiName_ =
|
||||
dict.lookupOrDefault<word>
|
||||
(
|
||||
"alphaPhi",
|
||||
IOobject::groupName("alphaPhi", phaseName_)
|
||||
);
|
||||
solveAlphaPhi_
|
||||
? typedName(defaultAlphaPhiName)
|
||||
: dict.lookupOrDefault<word>("alphaPhi", defaultAlphaPhiName);
|
||||
phiName_ = dict.lookupOrDefault<word>("phi", "phi");
|
||||
rhoName_ =
|
||||
dict.lookupOrDefault<word>
|
||||
@ -356,8 +351,8 @@ bool Foam::functionObjects::phaseScalarTransport::read(const dictionary& dict)
|
||||
alphaD_ = dict.lookupOrDefault<scalar>("alphaD", 1);
|
||||
alphaDt_ = dict.lookupOrDefault<scalar>("alphaDt", 1);
|
||||
|
||||
dict.readIfPresent("nCorr", nCorr_);
|
||||
dict.readIfPresent("residualAlpha", residualAlpha_);
|
||||
nCorr_ = dict.lookupOrDefault<label>("nCorr", 0);
|
||||
residualAlpha_ = dict.lookupOrDefault<scalar>("residualAlpha", rootSmall);
|
||||
writeAlphaField_ = dict.lookupOrDefault<bool>("writeAlphaField", true);
|
||||
|
||||
return true;
|
||||
|
||||
@ -37,27 +37,25 @@ Description
|
||||
detailed below. Note that the phase-name will be determined by stripping
|
||||
the extension from the supplied field name.
|
||||
|
||||
If \c alphaPhi is specified and found in the database then it will be used
|
||||
to transport the field. This is the likely mode of operation if this
|
||||
function is used with Euler-Euler solvers, in which \c alphaPhi -like
|
||||
fluxes are available. If \c alphaPhi is not found, then a pressure-like
|
||||
equation will be solved in order to construct it so that it exactly matches
|
||||
the time-derivative of the phase volume or mass. This is likely to be
|
||||
necessary in volume-of-fluid solvers where \c alphaPhi is not part of the
|
||||
solution procedure. The pressure field name will be required in this case.
|
||||
If the solver does not provide an \c alphaPhi flux, or that flux is for
|
||||
some reason unreliable, then the \c solveAlphaPhi switch can be used to
|
||||
make this function solve a pressure-like equation from which \c alphaPhi is
|
||||
recovered.
|
||||
|
||||
Usage
|
||||
\table
|
||||
Property | Description | Req'd? | Default
|
||||
alpha | Name of the volume-fraction field | no\\
|
||||
Property | Description | Req'd? | Default
|
||||
alpha | Name of the volume-fraction field | no\\
|
||||
| alpha.<phase-name>
|
||||
alphaPhi | Name of the phase-flux field | no\\
|
||||
alphaPhi | Name of the phase-flux field | no\\
|
||||
| alphaPhi.<phase-name>
|
||||
p | Name of the pressure field | no | p
|
||||
solveAlphaPhi | Solve for the alphaPhi flux, rather than looking it\\
|
||||
up | no | false
|
||||
p | Name of the pressure field | no | p
|
||||
residualAlpha | Small volume fraction used to stabilise the solution\\
|
||||
| no | rootSmall
|
||||
| no | rootSmall
|
||||
writeAlphaField | Also write out alpha multiplied by the field\\
|
||||
| no | true
|
||||
| no | true
|
||||
\endtable
|
||||
|
||||
Example specification for incompressibleVoF:
|
||||
@ -68,7 +66,6 @@ Usage
|
||||
libs ("libsolverFunctionObjects.so");
|
||||
|
||||
field s.water;
|
||||
p p_rgh;
|
||||
}
|
||||
\endverbatim
|
||||
|
||||
@ -123,6 +120,9 @@ class phaseScalarTransport
|
||||
//- Name of the phase in which to solve
|
||||
const word phaseName_;
|
||||
|
||||
//- Solve for the alphaPhi flux, rather than looking it up (optional)
|
||||
bool solveAlphaPhi_;
|
||||
|
||||
//- Name of phase volume-fraction field (optional)
|
||||
word alphaName_;
|
||||
|
||||
|
||||
Reference in New Issue
Block a user