View Issue Details
ID | Project | Category | View Status | Date Submitted | Last Update |
---|---|---|---|---|---|
0001710 | OpenFOAM | Bug | public | 2015-05-26 23:03 | 2015-06-15 17:45 |
Reporter | dhora | Assigned To | henry | ||
Priority | low | Severity | feature | Reproducibility | always |
Status | closed | Resolution | suspended | ||
Summary | 0001710: Output of SST blending functions | ||||
Description | I know that it's probably not of importance for most users but it were nice if the SST models would write out the values of the blending functions. Just for comparisons with other CFD codes. I added a modified incompressible version. | ||||
Tags | No tags attached. | ||||
|
kOmegaSST.H (8,851 bytes)
/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2011-2014 OpenFOAM Foundation \\/ M anipulation | ------------------------------------------------------------------------------- License This file is part of OpenFOAM. OpenFOAM is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. OpenFOAM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>. Class Foam::incompressible::RASModels::kOmegaSST Description Implementation of the k-omega-SST turbulence model for incompressible flows. Turbulence model described in \verbatim Menter, F., Esch, T., "Elements of Industrial Heat Transfer Prediction", 16th Brazilian Congress of Mechanical Engineering (COBEM), Nov. 2001. \endverbatim with updated coefficients from \verbatim Menter, F. R., Kuntz, M., and Langtry, R., "Ten Years of Industrial Experience with the SST Turbulence Model", Turbulence, Heat and Mass Transfer 4, 2003, pp. 625 - 632. \endverbatim but with the consistent production terms from the 2001 paper as form in the 2003 paper is a typo, see \verbatim http://turbmodels.larc.nasa.gov/sst.html \endverbatim and the addition of the optional F3 term for rough walls from \verbatim Hellsten, A. "Some Improvements in Menter’s k-omega-SST turbulence model" 29th AIAA Fluid Dynamics Conference, AIAA-98-2554, June 1998. \endverbatim Note that this implementation is written in terms of alpha diffusion coefficients rather than the more traditional sigma (alpha = 1/sigma) so that the blending can be applied to all coefficuients in a consistent manner. The paper suggests that sigma is blended but this would not be consistent with the blending of the k-epsilon and k-omega models. Also note that the error in the last term of equation (2) relating to sigma has been corrected. The default model coefficients correspond to the following: \verbatim kOmegaSSTCoeffs { alphaK1 0.85; alphaK2 1.0; alphaOmega1 0.5; alphaOmega2 0.856; beta1 0.075; beta2 0.0828; betaStar 0.09; gamma1 5/9; gamma2 0.44; a1 0.31; b1 1.0; c1 10.0; F3 no; } \endverbatim SourceFiles kOmegaSST.C \*---------------------------------------------------------------------------*/ #ifndef kOmegaSST_H #define kOmegaSST_H #include "RASModel.H" #include "wallDist.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // namespace Foam { namespace incompressible { namespace RASModels { /*---------------------------------------------------------------------------*\ Class kOmegaSST Declaration \*---------------------------------------------------------------------------*/ class kOmegaSST : public RASModel { protected: // Protected data // Model coefficients dimensionedScalar alphaK1_; dimensionedScalar alphaK2_; dimensionedScalar alphaOmega1_; dimensionedScalar alphaOmega2_; dimensionedScalar gamma1_; dimensionedScalar gamma2_; dimensionedScalar beta1_; dimensionedScalar beta2_; dimensionedScalar betaStar_; dimensionedScalar a1_; dimensionedScalar b1_; dimensionedScalar c1_; Switch F3_; //- Wall distance field // Note: different to wall distance in parent RASModel wallDist y_; // Fields volScalarField k_; volScalarField omega_; volScalarField nut_; volScalarField f1_; volScalarField f23_; // Protected Member Functions tmp<volScalarField> F1(const volScalarField& CDkOmega) const; tmp<volScalarField> F2() const; tmp<volScalarField> F3() const; tmp<volScalarField> F23() const; tmp<volScalarField> blend ( const volScalarField& F1, const dimensionedScalar& psi1, const dimensionedScalar& psi2 ) const { return F1*(psi1 - psi2) + psi2; } tmp<volScalarField> alphaK(const volScalarField& F1) const { return blend(F1, alphaK1_, alphaK2_); } tmp<volScalarField> alphaOmega(const volScalarField& F1) const { return blend(F1, alphaOmega1_, alphaOmega2_); } tmp<volScalarField> beta(const volScalarField& F1) const { return blend(F1, beta1_, beta2_); } tmp<volScalarField> gamma(const volScalarField& F1) const { return blend(F1, gamma1_, gamma2_); } public: //- Runtime type information TypeName("kOmegaSST"); // Constructors //- Construct from components kOmegaSST ( const volVectorField& U, const surfaceScalarField& phi, transportModel& transport, const word& turbulenceModelName = turbulenceModel::typeName, const word& modelName = typeName ); //- Destructor virtual ~kOmegaSST() {} // Member Functions //- Return the turbulence viscosity virtual tmp<volScalarField> nut() const { return nut_; } //- Return the effective diffusivity for k tmp<volScalarField> DkEff(const volScalarField& F1) const { return tmp<volScalarField> ( new volScalarField("DkEff", alphaK(F1)*nut_ + nu()) ); } //- Return the effective diffusivity for omega tmp<volScalarField> DomegaEff(const volScalarField& F1) const { return tmp<volScalarField> ( new volScalarField("DomegaEff", alphaOmega(F1)*nut_ + nu()) ); } //- Return the turbulence kinetic energy virtual tmp<volScalarField> k() const { return k_; } //- Return the turbulence specific dissipation rate virtual tmp<volScalarField> omega() const { return omega_; } //- Return the turbulence kinetic energy dissipation rate virtual tmp<volScalarField> epsilon() const { return tmp<volScalarField> ( new volScalarField ( IOobject ( "epsilon", mesh_.time().timeName(), mesh_ ), betaStar_*k_*omega_, omega_.boundaryField().types() ) ); } //- Return the Reynolds stress tensor virtual tmp<volSymmTensorField> R() const; //- Return the effective stress tensor including the laminar stress virtual tmp<volSymmTensorField> devReff() const; //- Return the source term for the momentum equation virtual tmp<fvVectorMatrix> divDevReff(volVectorField& U) const; //- Return the source term for the momentum equation virtual tmp<fvVectorMatrix> divDevRhoReff ( const volScalarField& rho, volVectorField& U ) const; //- Solve the turbulence equations and correct the turbulence viscosity virtual void correct(); //- Read RASProperties dictionary virtual bool read(); }; // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // } // End namespace RASModels } // namespace incompressible } // End namespace Foam // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // #endif // ************************************************************************* // |
|
kOmegaSST.C (10,571 bytes)
/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2011-2014 OpenFOAM Foundation \\/ M anipulation | ------------------------------------------------------------------------------- License This file is part of OpenFOAM. OpenFOAM is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. OpenFOAM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>. \*---------------------------------------------------------------------------*/ #include "kOmegaSST.H" #include "addToRunTimeSelectionTable.H" #include "backwardsCompatibilityWallFunctions.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // namespace Foam { namespace incompressible { namespace RASModels { // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // defineTypeNameAndDebug(kOmegaSST, 0); addToRunTimeSelectionTable(RASModel, kOmegaSST, dictionary); // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * // tmp<volScalarField> kOmegaSST::F1(const volScalarField& CDkOmega) const { tmp<volScalarField> CDkOmegaPlus = max ( CDkOmega, dimensionedScalar("1.0e-10", dimless/sqr(dimTime), 1.0e-10) ); tmp<volScalarField> arg1 = min ( min ( max ( (scalar(1)/betaStar_)*sqrt(k_)/(omega_*y_), scalar(500)*nu()/(sqr(y_)*omega_) ), (4*alphaOmega2_)*k_/(CDkOmegaPlus*sqr(y_)) ), scalar(10) ); return tanh(pow4(arg1)); } tmp<volScalarField> kOmegaSST::F2() const { tmp<volScalarField> arg2 = min ( max ( (scalar(2)/betaStar_)*sqrt(k_)/(omega_*y_), scalar(500)*nu()/(sqr(y_)*omega_) ), scalar(100) ); return tanh(sqr(arg2)); } tmp<volScalarField> kOmegaSST::F3() const { tmp<volScalarField> arg3 = min ( 150*nu()/(omega_*sqr(y_)), scalar(10) ); return 1 - tanh(pow4(arg3)); } tmp<volScalarField> kOmegaSST::F23() const { tmp<volScalarField> f23(F2()); if (F3_) { f23() *= F3(); } return f23; } // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // kOmegaSST::kOmegaSST ( const volVectorField& U, const surfaceScalarField& phi, transportModel& transport, const word& turbulenceModelName, const word& modelName ) : RASModel(modelName, U, phi, transport, turbulenceModelName), alphaK1_ ( dimensioned<scalar>::lookupOrAddToDict ( "alphaK1", coeffDict_, 0.85 ) ), alphaK2_ ( dimensioned<scalar>::lookupOrAddToDict ( "alphaK2", coeffDict_, 1.0 ) ), alphaOmega1_ ( dimensioned<scalar>::lookupOrAddToDict ( "alphaOmega1", coeffDict_, 0.5 ) ), alphaOmega2_ ( dimensioned<scalar>::lookupOrAddToDict ( "alphaOmega2", coeffDict_, 0.856 ) ), gamma1_ ( dimensioned<scalar>::lookupOrAddToDict ( "gamma1", coeffDict_, 5.0/9.0 ) ), gamma2_ ( dimensioned<scalar>::lookupOrAddToDict ( "gamma2", coeffDict_, 0.44 ) ), beta1_ ( dimensioned<scalar>::lookupOrAddToDict ( "beta1", coeffDict_, 0.075 ) ), beta2_ ( dimensioned<scalar>::lookupOrAddToDict ( "beta2", coeffDict_, 0.0828 ) ), betaStar_ ( dimensioned<scalar>::lookupOrAddToDict ( "betaStar", coeffDict_, 0.09 ) ), a1_ ( dimensioned<scalar>::lookupOrAddToDict ( "a1", coeffDict_, 0.31 ) ), b1_ ( dimensioned<scalar>::lookupOrAddToDict ( "b1", coeffDict_, 1.0 ) ), c1_ ( dimensioned<scalar>::lookupOrAddToDict ( "c1", coeffDict_, 10.0 ) ), F3_ ( Switch::lookupOrAddToDict ( "F3", coeffDict_, false ) ), y_(mesh_), k_ ( IOobject ( "k", runTime_.timeName(), mesh_, IOobject::NO_READ, IOobject::AUTO_WRITE ), autoCreateK("k", mesh_) ), omega_ ( IOobject ( "omega", runTime_.timeName(), mesh_, IOobject::NO_READ, IOobject::AUTO_WRITE ), autoCreateOmega("omega", mesh_) ), nut_ ( IOobject ( "nut", runTime_.timeName(), mesh_, IOobject::NO_READ, IOobject::AUTO_WRITE ), autoCreateNut("nut", mesh_) ), f1_ ( IOobject ( "F1", runTime_.timeName(), mesh_, IOobject::NO_READ, IOobject::AUTO_WRITE ), F1( (2*alphaOmega2_)*(fvc::grad(k_) & fvc::grad(omega_))/omega_ ) ), f23_ ( IOobject ( "F23", runTime_.timeName(), mesh_, IOobject::NO_READ, IOobject::AUTO_WRITE ), F23() ) { bound(k_, kMin_); bound(omega_, omegaMin_); nut_ = ( a1_*k_ / max ( a1_*omega_, b1_*F23()*sqrt(2.0)*mag(symm(fvc::grad(U_))) ) ); nut_.correctBoundaryConditions(); printCoeffs(); } // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // tmp<volSymmTensorField> kOmegaSST::R() const { return tmp<volSymmTensorField> ( new volSymmTensorField ( IOobject ( "R", runTime_.timeName(), mesh_, IOobject::NO_READ, IOobject::NO_WRITE ), ((2.0/3.0)*I)*k_ - nut_*twoSymm(fvc::grad(U_)), k_.boundaryField().types() ) ); } tmp<volSymmTensorField> kOmegaSST::devReff() const { return tmp<volSymmTensorField> ( new volSymmTensorField ( IOobject ( "devRhoReff", runTime_.timeName(), mesh_, IOobject::NO_READ, IOobject::NO_WRITE ), -nuEff()*dev(twoSymm(fvc::grad(U_))) ) ); } tmp<fvVectorMatrix> kOmegaSST::divDevReff(volVectorField& U) const { return ( - fvm::laplacian(nuEff(), U) - fvc::div(nuEff()*dev(T(fvc::grad(U)))) ); } tmp<fvVectorMatrix> kOmegaSST::divDevRhoReff ( const volScalarField& rho, volVectorField& U ) const { volScalarField muEff("muEff", rho*nuEff()); return ( - fvm::laplacian(muEff, U) - fvc::div(muEff*dev(T(fvc::grad(U)))) ); } bool kOmegaSST::read() { if (RASModel::read()) { alphaK1_.readIfPresent(coeffDict()); alphaK2_.readIfPresent(coeffDict()); alphaOmega1_.readIfPresent(coeffDict()); alphaOmega2_.readIfPresent(coeffDict()); gamma1_.readIfPresent(coeffDict()); gamma2_.readIfPresent(coeffDict()); beta1_.readIfPresent(coeffDict()); beta2_.readIfPresent(coeffDict()); betaStar_.readIfPresent(coeffDict()); a1_.readIfPresent(coeffDict()); b1_.readIfPresent(coeffDict()); c1_.readIfPresent(coeffDict()); F3_.readIfPresent("F3", coeffDict()); return true; } else { return false; } } void kOmegaSST::correct() { RASModel::correct(); if (!turbulence_) { return; } if (mesh_.changing()) { y_.correct(); } const volScalarField S2(2*magSqr(symm(fvc::grad(U_)))); volScalarField G(GName(), nut_*S2); // Update omega and G at the wall omega_.boundaryField().updateCoeffs(); const volScalarField CDkOmega ( (2*alphaOmega2_)*(fvc::grad(k_) & fvc::grad(omega_))/omega_ ); //const volScalarField F1(this->F1(CDkOmega)); f1_ = this->F1(CDkOmega); // Turbulent frequency equation tmp<fvScalarMatrix> omegaEqn ( fvm::ddt(omega_) + fvm::div(phi_, omega_) - fvm::laplacian(DomegaEff(f1_), omega_) == gamma(f1_) *min(S2, (c1_/a1_)*betaStar_*omega_*max(a1_*omega_, b1_*F23()*sqrt(S2))) - fvm::Sp(beta(f1_)*omega_, omega_) - fvm::SuSp ( (f1_ - scalar(1))*CDkOmega/omega_, omega_ ) ); omegaEqn().relax(); omegaEqn().boundaryManipulate(omega_.boundaryField()); solve(omegaEqn); bound(omega_, omegaMin_); // Turbulent kinetic energy equation tmp<fvScalarMatrix> kEqn ( fvm::ddt(k_) + fvm::div(phi_, k_) - fvm::laplacian(DkEff(f1_), k_) == min(G, c1_*betaStar_*k_*omega_) - fvm::Sp(betaStar_*omega_, k_) ); kEqn().relax(); solve(kEqn); bound(k_, kMin_); // Re-calculate viscosity f23_ = F23(); nut_ = a1_*k_/max(a1_*omega_, b1_*f23_*sqrt(S2)); nut_.correctBoundaryConditions(); } // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // } // End namespace RASModels } // End namespace incompressible } // End namespace Foam // ************************************************************************* // |
|
The problem with your proposal is that two additional fields are stored permanently and always written irrespective of need; generally this information is not needed so the overhead would be unwelcome. It would be better if the information could be generated and written on demand, via a functionObject for example, for which I am considering what kind of API would be needed and how many kinds of fields might be requested from which models and what the commonality between them is. To do this well will take quite a bit of thought and time so funding will be required. |
|
if (runTime_.outputTime() && somethingelse) { volScalarField f1 ( IOobject ( "F1", runTime_.timeName(), mesh_, IOobject::NO_READ, IOobject::AUTO_WRITE ), F1 ); f1.write(); volScalarField f23 ( IOobject ( "F23", runTime_.timeName(), mesh_, IOobject::NO_READ, IOobject::AUTO_WRITE ), F23() ); f23.write(); } |
|
Yes this is the code which would selectively write these fields for this particular model it does not provide the API or the functionality to write similar fields for all the other models. Rather than hacking a specific solution to your particular need I would like to see a solution to the more general need for optional diagnostics for all the turbulence models, including blending functions, switching functions, distance-to-wall functions etc. |
Date Modified | Username | Field | Change |
---|---|---|---|
2015-05-26 23:03 | dhora | New Issue | |
2015-05-26 23:03 | dhora | File Added: kOmegaSST.H | |
2015-05-26 23:03 | dhora | File Added: kOmegaSST.C | |
2015-05-26 23:18 | henry | Note Added: 0004816 | |
2015-05-26 23:52 | dhora | Note Added: 0004817 | |
2015-05-27 07:00 | henry | Note Added: 0004819 | |
2015-06-15 17:45 | henry | Status | new => closed |
2015-06-15 17:45 | henry | Assigned To | => henry |
2015-06-15 17:45 | henry | Resolution | open => suspended |