mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
turbulenceModels/LES/dynamicLagrangian: Added support for fvOptions
This commit is contained in:
@ -24,6 +24,7 @@ License
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "dynamicLagrangian.H"
|
||||
#include "fvOptions.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
@ -42,6 +43,9 @@ void dynamicLagrangian<BasicTurbulenceModel>::correctNut
|
||||
{
|
||||
this->nut_ = (flm_/fmm_)*sqr(this->delta())*mag(dev(symm(gradU)));
|
||||
this->nut_.correctBoundaryConditions();
|
||||
fv::options::New(this->mesh_).correct(this->nut_);
|
||||
|
||||
BasicTurbulenceModel::correctNut();
|
||||
}
|
||||
|
||||
|
||||
@ -157,6 +161,7 @@ void dynamicLagrangian<BasicTurbulenceModel>::correct()
|
||||
// Local references
|
||||
const surfaceScalarField& phi = this->phi_;
|
||||
const volVectorField& U = this->U_;
|
||||
fv::options& fvOptions(fv::options::New(this->mesh_));
|
||||
|
||||
LESeddyViscosity<BasicTurbulenceModel>::correct();
|
||||
|
||||
@ -190,11 +195,13 @@ void dynamicLagrangian<BasicTurbulenceModel>::correct()
|
||||
==
|
||||
invT*LM
|
||||
- fvm::Sp(invT, flm_)
|
||||
+ fvOptions(flm_)
|
||||
);
|
||||
|
||||
flmEqn.relax();
|
||||
fvOptions.constrain(flmEqn);
|
||||
flmEqn.solve();
|
||||
|
||||
fvOptions.correct(flm_);
|
||||
bound(flm_, flm0_);
|
||||
|
||||
volScalarField MM(M && M);
|
||||
@ -206,11 +213,13 @@ void dynamicLagrangian<BasicTurbulenceModel>::correct()
|
||||
==
|
||||
invT*MM
|
||||
- fvm::Sp(invT, fmm_)
|
||||
+ fvOptions(fmm_)
|
||||
);
|
||||
|
||||
fmmEqn.relax();
|
||||
fvOptions.constrain(fmmEqn);
|
||||
fmmEqn.solve();
|
||||
|
||||
fvOptions.correct(fmm_);
|
||||
bound(fmm_, fmm0_);
|
||||
|
||||
correctNut(gradU);
|
||||
|
||||
Reference in New Issue
Block a user