fvModels: Improved interface for mass/volume sources
The interface for fvModels has been modified to improve its application
to "proxy" equations. That is, equations that are not straightforward
statements of conservation laws in OpenFOAM's usual convention.
A standard conservation law typically takes the following form:
fvMatrix<scalar> psiEqn
(
fvm::ddt(alpha, rho, psi)
+ <fluxes>
==
<sources>
);
A proxy equation, on the other hand, may be a derivation or
rearrangement of a law like this, and may be linearised in terms of a
different variable.
The pressure equation is the most common example of a proxy equation. It
represents a statement of the conservation of volume or mass, but it is
a rearrangement of the original continuity equation, and it has been
linearised in terms of a different variable; the pressure. Another
example is that in the pre-predictor of a VoF solver the
phase-continuity equation is constructed, but it is linearised in terms
of volume fraction rather than density.
In these situations, fvModels sources are now applied by calling:
fvModels().sourceProxy(<conserved-fields ...>, <equation-field>)
Where <conserved-fields ...> are (alpha, rho, psi), (rho, psi), just
(psi), or are omitted entirely (for volume continuity), and the
<equation-field> is the field associated with the proxy equation. This
produces a source term identical in value to the following call:
fvModels().source(<conserved-fields ...>)
It is only the linearisation in terms of <equation-field> that differs
between these two calls.
This change permits much greater flexibility in the handling of mass and
volume sources than the previous name-based system did. All the relevant
fields are available, dimensions can be used in the logic to determine
what sources are being constructed, and sources relating to a given
conservation law all share the same function.
This commit adds the functionality for injection-type sources in the
compressibleVoF solver. A following commit will add a volume source
model for use in incompressible solvers.
This commit is contained in:
@ -102,7 +102,7 @@ Foam::solvers::multiphaseEuler::compressibilityEqns
|
||||
// Option sources
|
||||
if (fvModels().addsSupToField(rho.name()))
|
||||
{
|
||||
pEqnComp -= (fvModels().source(alpha, rho) & rho)/rho;
|
||||
pEqnComp -= fvModels().sourceProxy(alpha, rho, p_rgh)/rho;
|
||||
}
|
||||
|
||||
// Mass transfer
|
||||
|
||||
@ -131,8 +131,8 @@ template<class RhoType>
|
||||
void Foam::fv::interfaceTurbulenceDamping::addRhoSup
|
||||
(
|
||||
const RhoType& rho,
|
||||
fvMatrix<scalar>& eqn,
|
||||
const word& fieldName
|
||||
const volScalarField& field,
|
||||
fvMatrix<scalar>& eqn
|
||||
) const
|
||||
{
|
||||
if (debug)
|
||||
@ -154,12 +154,12 @@ void Foam::fv::interfaceTurbulenceDamping::addRhoSup
|
||||
movingPhases[phasei]*sqr(movingPhases[phasei].thermo().nu()()());
|
||||
}
|
||||
|
||||
if (fieldName == "epsilon")
|
||||
if (field.name() == "epsilon")
|
||||
{
|
||||
eqn += rho*interfaceFraction(phase_)*C2_*aSqrnu*turbulence_.k()()
|
||||
/pow4(delta_);
|
||||
}
|
||||
else if (fieldName == "omega")
|
||||
else if (field.name() == "omega")
|
||||
{
|
||||
eqn += rho*interfaceFraction(phase_)*beta_*aSqrnu
|
||||
/(sqr(betaStar_)*pow4(delta_));
|
||||
@ -167,7 +167,7 @@ void Foam::fv::interfaceTurbulenceDamping::addRhoSup
|
||||
else
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "Support for field " << fieldName << " is not implemented"
|
||||
<< "Support for field " << field.name() << " is not implemented"
|
||||
<< exit(FatalError);
|
||||
}
|
||||
}
|
||||
@ -241,22 +241,22 @@ Foam::wordList Foam::fv::interfaceTurbulenceDamping::addSupFields() const
|
||||
|
||||
void Foam::fv::interfaceTurbulenceDamping::addSup
|
||||
(
|
||||
fvMatrix<scalar>& eqn,
|
||||
const word& fieldName
|
||||
const volScalarField& field,
|
||||
fvMatrix<scalar>& eqn
|
||||
) const
|
||||
{
|
||||
addRhoSup(one(), eqn, fieldName);
|
||||
addRhoSup(one(), field, eqn);
|
||||
}
|
||||
|
||||
|
||||
void Foam::fv::interfaceTurbulenceDamping::addSup
|
||||
(
|
||||
const volScalarField& rho,
|
||||
fvMatrix<scalar>& eqn,
|
||||
const word& fieldName
|
||||
const volScalarField& field,
|
||||
fvMatrix<scalar>& eqn
|
||||
) const
|
||||
{
|
||||
addRhoSup(rho(), eqn, fieldName);
|
||||
addRhoSup(rho(), field, eqn);
|
||||
}
|
||||
|
||||
|
||||
@ -264,8 +264,8 @@ void Foam::fv::interfaceTurbulenceDamping::addSup
|
||||
(
|
||||
const volScalarField& alpha,
|
||||
const volScalarField& rho,
|
||||
fvMatrix<scalar>& eqn,
|
||||
const word& fieldName
|
||||
const volScalarField& field,
|
||||
fvMatrix<scalar>& eqn
|
||||
) const
|
||||
{
|
||||
if (debug)
|
||||
@ -273,17 +273,17 @@ void Foam::fv::interfaceTurbulenceDamping::addSup
|
||||
Info<< type() << ": applying source to " << eqn.psi().name() << endl;
|
||||
}
|
||||
|
||||
volScalarField::Internal aSqrnu
|
||||
const volScalarField::Internal aSqrnu
|
||||
(
|
||||
alpha*sqr(phase_.thermo().nu()()())
|
||||
);
|
||||
|
||||
if (fieldName == IOobject::groupName("epsilon", phaseName_))
|
||||
if (field.name() == IOobject::groupName("epsilon", phaseName_))
|
||||
{
|
||||
eqn += rho()*interfaceFraction(alpha)
|
||||
*C2_*aSqrnu*turbulence_.k()()/pow4(delta_);
|
||||
}
|
||||
else if (fieldName == IOobject::groupName("omega", phaseName_))
|
||||
else if (field.name() == IOobject::groupName("omega", phaseName_))
|
||||
{
|
||||
eqn += rho()*interfaceFraction(alpha)
|
||||
*beta_*aSqrnu/(sqr(betaStar_)*pow4(delta_));
|
||||
@ -291,7 +291,7 @@ void Foam::fv::interfaceTurbulenceDamping::addSup
|
||||
else
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "Support for field " << fieldName << " is not implemented"
|
||||
<< "Support for field " << field.name() << " is not implemented"
|
||||
<< exit(FatalError);
|
||||
}
|
||||
}
|
||||
|
||||
@ -131,8 +131,8 @@ class interfaceTurbulenceDamping
|
||||
void addRhoSup
|
||||
(
|
||||
const RhoType& rho,
|
||||
fvMatrix<scalar>& eqn,
|
||||
const word& fieldName
|
||||
const volScalarField& field,
|
||||
fvMatrix<scalar>& eqn
|
||||
) const;
|
||||
|
||||
|
||||
@ -162,34 +162,38 @@ public:
|
||||
|
||||
// Member Functions
|
||||
|
||||
//- Return the list of fields for which the option adds source term
|
||||
// to the transport equation
|
||||
virtual wordList addSupFields() const;
|
||||
// Checks
|
||||
|
||||
//- Add explicit contribution to mixture epsilon or omega equation
|
||||
virtual void addSup
|
||||
(
|
||||
fvMatrix<scalar>& eqn,
|
||||
const word& fieldName
|
||||
) const;
|
||||
//- Return the list of fields for which the option adds source term
|
||||
// to the transport equation
|
||||
virtual wordList addSupFields() const;
|
||||
|
||||
//- Add explicit contribution to compressible
|
||||
// mixture epsilon or omega equation
|
||||
virtual void addSup
|
||||
(
|
||||
const volScalarField& rho,
|
||||
fvMatrix<scalar>& eqn,
|
||||
const word& fieldName
|
||||
) const;
|
||||
|
||||
//- Add explicit contribution to phase epsilon or omega equation
|
||||
virtual void addSup
|
||||
(
|
||||
const volScalarField& alpha,
|
||||
const volScalarField& rho,
|
||||
fvMatrix<scalar>& eqn,
|
||||
const word& fieldName
|
||||
) const;
|
||||
// Sources
|
||||
|
||||
//- Add source to mixture epsilon or omega equation
|
||||
virtual void addSup
|
||||
(
|
||||
const volScalarField& field,
|
||||
fvMatrix<scalar>& eqn
|
||||
) const;
|
||||
|
||||
//- Add source to compressible mixture epsilon or omega equation
|
||||
virtual void addSup
|
||||
(
|
||||
const volScalarField& rho,
|
||||
const volScalarField& field,
|
||||
fvMatrix<scalar>& eqn
|
||||
) const;
|
||||
|
||||
//- Add source to phase epsilon or omega equation
|
||||
virtual void addSup
|
||||
(
|
||||
const volScalarField& alpha,
|
||||
const volScalarField& rho,
|
||||
const volScalarField& field,
|
||||
fvMatrix<scalar>& eqn
|
||||
) const;
|
||||
|
||||
|
||||
// Mesh changes
|
||||
|
||||
@ -50,10 +50,11 @@ namespace Foam
|
||||
|
||||
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
|
||||
|
||||
void Foam::fv::phaseTurbulenceStabilisation::addSup
|
||||
void Foam::fv::phaseTurbulenceStabilisation::addAlphaRhoSup
|
||||
(
|
||||
const volScalarField& alpha,
|
||||
const volScalarField& rho,
|
||||
const volScalarField& field,
|
||||
fvMatrix<scalar>& eqn,
|
||||
tmp<volScalarField>
|
||||
(phaseCompressible::momentumTransportModel::*psi)() const
|
||||
@ -184,36 +185,39 @@ void Foam::fv::phaseTurbulenceStabilisation::addSup
|
||||
(
|
||||
const volScalarField& alpha,
|
||||
const volScalarField& rho,
|
||||
fvMatrix<scalar>& eqn,
|
||||
const word& fieldName
|
||||
const volScalarField& field,
|
||||
fvMatrix<scalar>& eqn
|
||||
) const
|
||||
{
|
||||
if (fieldName == IOobject::groupName("k", phaseName_))
|
||||
if (field.name() == IOobject::groupName("k", phaseName_))
|
||||
{
|
||||
addSup
|
||||
addAlphaRhoSup
|
||||
(
|
||||
alpha,
|
||||
rho,
|
||||
field,
|
||||
eqn,
|
||||
&phaseCompressible::momentumTransportModel::k
|
||||
);
|
||||
}
|
||||
else if (fieldName == IOobject::groupName("epsilon", phaseName_))
|
||||
else if (field.name() == IOobject::groupName("epsilon", phaseName_))
|
||||
{
|
||||
addSup
|
||||
addAlphaRhoSup
|
||||
(
|
||||
alpha,
|
||||
rho,
|
||||
field,
|
||||
eqn,
|
||||
&phaseCompressible::momentumTransportModel::epsilon
|
||||
);
|
||||
}
|
||||
else if (fieldName == IOobject::groupName("omega", phaseName_))
|
||||
else if (field.name() == IOobject::groupName("omega", phaseName_))
|
||||
{
|
||||
addSup
|
||||
addAlphaRhoSup
|
||||
(
|
||||
alpha,
|
||||
rho,
|
||||
field,
|
||||
eqn,
|
||||
&phaseCompressible::momentumTransportModel::omega
|
||||
);
|
||||
@ -221,7 +225,7 @@ void Foam::fv::phaseTurbulenceStabilisation::addSup
|
||||
else
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "Support for field " << fieldName << " is not implemented"
|
||||
<< "Support for field " << field.name() << " is not implemented"
|
||||
<< exit(FatalError);
|
||||
}
|
||||
}
|
||||
|
||||
@ -112,10 +112,11 @@ class phaseTurbulenceStabilisation
|
||||
// Private Member Functions
|
||||
|
||||
//- Add contribution to phase psi equation
|
||||
void addSup
|
||||
void addAlphaRhoSup
|
||||
(
|
||||
const volScalarField& alpha,
|
||||
const volScalarField& rho,
|
||||
const volScalarField& field,
|
||||
fvMatrix<scalar>& eqn,
|
||||
tmp<volScalarField>
|
||||
(phaseCompressible::momentumTransportModel::*psi)() const
|
||||
@ -148,20 +149,23 @@ public:
|
||||
|
||||
// Member Functions
|
||||
|
||||
//- Return the list of fields for which the option adds source term
|
||||
// to the transport equation
|
||||
virtual wordList addSupFields() const;
|
||||
// Checks
|
||||
|
||||
using fvModel::addSup;
|
||||
//- Return the list of fields for which the option adds source term
|
||||
// to the transport equation
|
||||
virtual wordList addSupFields() const;
|
||||
|
||||
//- Add contribution to phase k, epsilon or omega equation
|
||||
virtual void addSup
|
||||
(
|
||||
const volScalarField& alpha,
|
||||
const volScalarField& rho,
|
||||
fvMatrix<scalar>& eqn,
|
||||
const word& fieldName
|
||||
) const;
|
||||
|
||||
// Sources
|
||||
|
||||
//- Add contribution to phase k, epsilon or omega equation
|
||||
virtual void addSup
|
||||
(
|
||||
const volScalarField& alpha,
|
||||
const volScalarField& rho,
|
||||
const volScalarField& field,
|
||||
fvMatrix<scalar>& eqn
|
||||
) const;
|
||||
|
||||
|
||||
// Mesh changes
|
||||
|
||||
Reference in New Issue
Block a user