functionObjects::scalarTransport: Added support for optional laminar and turbulent diffusion coefficients

Description
    Evolves a passive scalar transport equation.

    - To specify the field name set the \c field entry
    - To employ the same numerical schemes as another field set
      the \c schemesField entry,
    - A constant diffusivity may be specified with the \c D entry,

    - Alternatively if a turbulence model is available a turbulent diffusivity
      may be constructed from the laminar and turbulent viscosities using the
      optional diffusivity coefficients \c alphaD and \c alphaDt (which default
      to 1):
      \verbatim
          D = alphaD*nu + alphaDt*nut
      \endverbatim

Resolves feature request https://bugs.openfoam.org/view.php?id=2453
This commit is contained in:
Henry Weller
2017-02-12 17:19:27 +00:00
parent dd8b8bceac
commit ae9522f017
2 changed files with 23 additions and 12 deletions

View File

@ -89,7 +89,7 @@ Foam::tmp<Foam::volScalarField> Foam::functionObjects::scalarTransport::D
turbulenceModel::propertiesName
);
return model.nuEff();
return alphaD_*model.nu() + alphaDt_*model.nut();
}
else if (mesh_.foundObject<cmpModel>(turbulenceModel::propertiesName))
{
@ -98,7 +98,7 @@ Foam::tmp<Foam::volScalarField> Foam::functionObjects::scalarTransport::D
turbulenceModel::propertiesName
);
return model.muEff();
return alphaD_*model.mu() + alphaDt_*model.mut();
}
else
{
@ -169,11 +169,9 @@ bool Foam::functionObjects::scalarTransport::read(const dictionary& dict)
rhoName_ = dict.lookupOrDefault<word>("rho", "rho");
schemesField_ = dict.lookupOrDefault<word>("schemesField", fieldName_);
constantD_ = false;
if (dict.readIfPresent("D", D_))
{
constantD_ = true;
}
constantD_ = dict.readIfPresent("D", D_);
alphaD_ = dict.lookupOrDefault("alphat", 1.0);
alphaDt_ = dict.lookupOrDefault("alphaDt", 1.0);
dict.readIfPresent("nCorr", nCorr_);