functionObjects/utilities/scalarTransport: Set UName_ during construction for boundaryTypes()

Resolves bug-report http://www.openfoam.org/mantisbt/view.php?id=1972
This commit is contained in:
Henry Weller
2016-01-15 14:47:56 +00:00
parent 36fc8595db
commit fa191d3456
2 changed files with 102 additions and 119 deletions

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2012-2015 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2012-2016 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -149,10 +149,9 @@ Foam::scalarTransport::scalarTransport
:
name_(name),
mesh_(refCast<const fvMesh>(obr)),
active_(true),
phiName_("phi"),
UName_("U"),
rhoName_("rho"),
phiName_(dict.lookupOrDefault<word>("phiName", "phi")),
UName_(dict.lookupOrDefault<word>("UName", "U")),
rhoName_(dict.lookupOrDefault<word>("rhoName", "rho")),
DT_(0.0),
userDT_(false),
resetOnStartUp_(false),
@ -193,137 +192,124 @@ Foam::scalarTransport::~scalarTransport()
void Foam::scalarTransport::read(const dictionary& dict)
{
if (active_)
Info<< type() << ":" << nl;
phiName_ = dict.lookupOrDefault<word>("phiName", "phi");
UName_ = dict.lookupOrDefault<word>("UName", "U");
rhoName_ = dict.lookupOrDefault<word>("rhoName", "rho");
userDT_ = false;
if (dict.readIfPresent("DT", DT_))
{
Info<< type() << ":" << nl;
phiName_ = dict.lookupOrDefault<word>("phiName", "phi");
UName_ = dict.lookupOrDefault<word>("UName", "U");
rhoName_ = dict.lookupOrDefault<word>("rhoName", "rho");
userDT_ = false;
if (dict.readIfPresent("DT", DT_))
{
userDT_ = true;
}
dict.lookup("resetOnStartUp") >> resetOnStartUp_;
dict.readIfPresent("nCorr", nCorr_);
dict.lookup("autoSchemes") >> autoSchemes_;
fvOptions_.reset(dict.subDict("fvOptions"));
userDT_ = true;
}
dict.lookup("resetOnStartUp") >> resetOnStartUp_;
dict.readIfPresent("nCorr", nCorr_);
dict.lookup("autoSchemes") >> autoSchemes_;
fvOptions_.reset(dict.subDict("fvOptions"));
}
void Foam::scalarTransport::execute()
{
if (active_)
Info<< type() << " output:" << endl;
const surfaceScalarField& phi =
mesh_.lookupObject<surfaceScalarField>(phiName_);
// calculate the diffusivity
volScalarField DT(this->DT(phi));
// set schemes
word schemeVar = T_.name();
if (autoSchemes_)
{
Info<< type() << " output:" << endl;
const surfaceScalarField& phi =
mesh_.lookupObject<surfaceScalarField>(phiName_);
// calculate the diffusivity
volScalarField DT(this->DT(phi));
// set schemes
word schemeVar = T_.name();
if (autoSchemes_)
{
schemeVar = UName_;
}
word divScheme("div(phi," + schemeVar + ")");
word laplacianScheme("laplacian(" + DT.name() + "," + schemeVar + ")");
// set under-relaxation coeff
scalar relaxCoeff = 0.0;
if (mesh_.relaxEquation(schemeVar))
{
relaxCoeff = mesh_.equationRelaxationFactor(schemeVar);
}
if (phi.dimensions() == dimMass/dimTime)
{
const volScalarField& rho =
mesh_.lookupObject<volScalarField>(rhoName_);
// solve
for (label i = 0; i <= nCorr_; i++)
{
fvScalarMatrix TEqn
(
fvm::ddt(rho, T_)
+ fvm::div(phi, T_, divScheme)
- fvm::laplacian(DT, T_, laplacianScheme)
==
fvOptions_(rho, T_)
);
TEqn.relax(relaxCoeff);
fvOptions_.constrain(TEqn);
TEqn.solve(mesh_.solverDict(schemeVar));
}
}
else if (phi.dimensions() == dimVolume/dimTime)
{
// solve
for (label i = 0; i <= nCorr_; i++)
{
fvScalarMatrix TEqn
(
fvm::ddt(T_)
+ fvm::div(phi, T_, divScheme)
- fvm::laplacian(DT, T_, laplacianScheme)
==
fvOptions_(T_)
);
TEqn.relax(relaxCoeff);
fvOptions_.constrain(TEqn);
TEqn.solve(mesh_.solverDict(schemeVar));
}
}
else
{
FatalErrorInFunction
<< "Incompatible dimensions for phi: " << phi.dimensions() << nl
<< "Dimensions should be " << dimMass/dimTime << " or "
<< dimVolume/dimTime << endl;
}
Info<< endl;
schemeVar = UName_;
}
word divScheme("div(phi," + schemeVar + ")");
word laplacianScheme("laplacian(" + DT.name() + "," + schemeVar + ")");
// set under-relaxation coeff
scalar relaxCoeff = 0.0;
if (mesh_.relaxEquation(schemeVar))
{
relaxCoeff = mesh_.equationRelaxationFactor(schemeVar);
}
if (phi.dimensions() == dimMass/dimTime)
{
const volScalarField& rho =
mesh_.lookupObject<volScalarField>(rhoName_);
// solve
for (label i = 0; i <= nCorr_; i++)
{
fvScalarMatrix TEqn
(
fvm::ddt(rho, T_)
+ fvm::div(phi, T_, divScheme)
- fvm::laplacian(DT, T_, laplacianScheme)
==
fvOptions_(rho, T_)
);
TEqn.relax(relaxCoeff);
fvOptions_.constrain(TEqn);
TEqn.solve(mesh_.solverDict(schemeVar));
}
}
else if (phi.dimensions() == dimVolume/dimTime)
{
// solve
for (label i = 0; i <= nCorr_; i++)
{
fvScalarMatrix TEqn
(
fvm::ddt(T_)
+ fvm::div(phi, T_, divScheme)
- fvm::laplacian(DT, T_, laplacianScheme)
==
fvOptions_(T_)
);
TEqn.relax(relaxCoeff);
fvOptions_.constrain(TEqn);
TEqn.solve(mesh_.solverDict(schemeVar));
}
}
else
{
FatalErrorInFunction
<< "Incompatible dimensions for phi: " << phi.dimensions() << nl
<< "Dimensions should be " << dimMass/dimTime << " or "
<< dimVolume/dimTime << endl;
}
Info<< endl;
}
void Foam::scalarTransport::end()
{
if (active_)
{
execute();
}
execute();
}
void Foam::scalarTransport::timeSet()
{
// Do nothing
}
{}
void Foam::scalarTransport::write()
{
// Do nothing
}
{}
// ************************************************************************* //

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2012-2013 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2012-2016 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -78,9 +78,6 @@ class scalarTransport
//- Reference to the mesh database
const fvMesh& mesh_;
//- On/off switch
bool active_;
//- Name of flux field (optional)
word phiName_;