|
[Sponsors] |
Error "Unknown null-constructable patchField type kqRWallFunction for patch cyclone " |
|
LinkBack | Thread Tools | Search this Thread | Display Modes |
December 8, 2023, 06:59 |
Error "Unknown null-constructable patchField type kqRWallFunction for patch cyclone "
|
#1 |
New Member
Ishak
Join Date: Dec 2023
Posts: 3
Rep Power: 3 |
Hello Everyone,
I'm trying to update a PANS turbulent model that I compiled back then in OpenFoam v7 but now with the new version of OpenFoam v11 I was obliged to change some new command.. Whatever, I correct and update the code with the new command and it compiled as well. but when I try to launch a simulation I got this error : Code:
--> FOAM FATAL ERROR: Unknown null-constructable patchField type kqRWallFunction for patch cyclone of type wall for field kU Valid null-constructable patchField types are : 24 ( calculated cyclic cyclicSlip empty extrapolatedCalculated fixedEnergy fixedUnburntEnthalpy fixedValue gradientEnergy gradientUnburntEnthalpy internal mixedEnergy mixedUnburntEnthalpy nonConformalCyclic nonConformalError nonConformalProcessorCyclic processor processorCyclic slip symmetry symmetryPlane wedge zeroGradient zeroInletOutlet ) From function static Foam::tmp<Foam::fvPatchField<Type> > Foam::fvPatchField<Type>::New(const Foam::word&, const Foam::word&, const Foam::fvPatch&, const Foam::DimensionedField<Type, Foam::volMesh>&) [with Type = double] in file lnInclude/fvPatchFieldNew.C at line 51. FOAM aborting The turbulence model called PANS and here's somme of the added constructor : Code:
fEpsilon_ ( dimensioned<scalar>::lookupOrAddToDict ( "fEpsilon", this->coeffDict_, 1.0 ) ), uLim_ ( dimensioned<scalar>::lookupOrAddToDict ( "fKupperLimit", this->coeffDict_, 1.0 ) ), loLim_ ( dimensioned<scalar>::lookupOrAddToDict ( "fKlowerLimit", this->coeffDict_, 0.1 ) ), fK_ ( IOobject ( this->groupName("fK"), this->runTime_.name(), this->mesh_, IOobject::NO_READ, IOobject::AUTO_WRITE ), this->mesh_, dimensionedScalar("zero",loLim_) ), C2U ( IOobject ( "C2U", this->runTime_.timeName(), this->mesh_ ), C1_ + (fK_/fEpsilon_)*(C2_ - C1_) ), delta_ ( LESdelta::New ( IOobject::groupName("delta", U.group()), *this, this->coeffDict_ ) ), k_ ( IOobject ( this->groupName("k"), this->runTime_.name(), this->mesh_, IOobject::MUST_READ, IOobject::AUTO_WRITE ), this->mesh_ ), kU_ ( IOobject ( this->groupName("kU"), this->runTime_.name(), this->mesh_, IOobject::NO_READ, IOobject::AUTO_WRITE ), //this->mesh_, k_*fK_, k_.boundaryField().types() ), epsilon_ ( IOobject ( this->groupName("epsilon"), this->runTime_.name(), this->mesh_, IOobject::MUST_READ, IOobject::AUTO_WRITE ), this->mesh_ ), epsilonU_ ( IOobject ( this->groupName("epsilonU"), this->runTime_.name(), this->mesh_, IOobject::NO_READ, IOobject::AUTO_WRITE ), //this->mesh_ epsilon_*fEpsilon_, epsilon_.boundaryField().types() ) |
|
July 2, 2024, 09:59 |
|
#3 |
New Member
Ishak
Join Date: Dec 2023
Posts: 3
Rep Power: 3 |
I fixed the problem.
I figured how to configure the PANS turbulence model in OF v 11, otherwise here's the code. The code worked well for me, but as I'm new in OF don't hesitate to correct me if I'm wrong please. C file: Code:
#include "kEpsilonPANS.H" #include "fvModels.H" #include "fvConstraints.H" #include "bound.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // namespace Foam { namespace RASModels { // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * // template<class BasicMomentumTransportModel> void kEpsilonPANS<BasicMomentumTransportModel>::correctNut() { this->nut_ = Cmu_*sqr(kU_)/epsilonU_; this->nut_.correctBoundaryConditions(); fvConstraints::New(this->mesh_).constrain(this->nut_); } template<class BasicMomentumTransportModel> tmp<fvScalarMatrix> kEpsilonPANS<BasicMomentumTransportModel>::kSource() const { return tmp<fvScalarMatrix> ( new fvScalarMatrix ( kU_, dimVolume*this->rho_.dimensions()*kU_.dimensions() /dimTime ) ); } template<class BasicMomentumTransportModel> tmp<fvScalarMatrix> kEpsilonPANS<BasicMomentumTransportModel>::epsilonSource() const { return tmp<fvScalarMatrix> ( new fvScalarMatrix ( epsilonU_, dimVolume*this->rho_.dimensions()*epsilonU_.dimensions() /dimTime ) ); } // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // template<class BasicMomentumTransportModel> kEpsilonPANS<BasicMomentumTransportModel>::kEpsilonPANS ( const alphaField& alpha, const rhoField& rho, const volVectorField& U, const surfaceScalarField& alphaRhoPhi, const surfaceScalarField& phi, const viscosity& viscosity, const word& type ) : eddyViscosity<RASModel<BasicMomentumTransportModel>> ( type, alpha, rho, U, alphaRhoPhi, phi, viscosity ), Cmu_ ( dimensioned<scalar>::lookupOrAddToDict ( "Cmu", this->coeffDict_, 0.09 ) ), C1_ ( dimensioned<scalar>::lookupOrAddToDict ( "C1", this->coeffDict_, 1.44 ) ), C2_ ( dimensioned<scalar>::lookupOrAddToDict ( "C2", this->coeffDict_, 1.92 ) ), C3_ ( dimensioned<scalar>::lookupOrAddToDict ( "C3", this->coeffDict_, -0.33 ) ), sigmak_ ( dimensioned<scalar>::lookupOrAddToDict ( "sigmak", this->coeffDict_, 1.0 ) ), sigmaEps_ ( dimensioned<scalar>::lookupOrAddToDict ( "sigmaEps", this->coeffDict_, 1.3 ) ), fEpsilon_ ( dimensioned<scalar>::lookupOrAddToDict ( "fEpsilon", this->coeffDict_, 1.0 ) ), uLim_ ( dimensioned<scalar>::lookupOrAddToDict ( "fKupperLimit", this->coeffDict_, 1.0 ) ), loLim_ ( dimensioned<scalar>::lookupOrAddToDict ( "fKlowerLimit", this->coeffDict_, 0.1 ) ), fK_ ( IOobject ( this->groupName("fK"), this->runTime_.name(), this->mesh_, IOobject::NO_READ, IOobject::AUTO_WRITE ), this->mesh_, dimensionedScalar("zero",loLim_) ), C2U ( IOobject ( "C2U", this->runTime_.name(), this->mesh_ ), C1_ + (fK_/fEpsilon_)*(C2_ - C1_) ), delta_ ( LESdelta::New ( IOobject::groupName("delta", U.group()), *this, this->coeffDict_ ) ), k_ ( IOobject ( this->groupName("k"), this->runTime_.name(), this->mesh_, IOobject::MUST_READ, IOobject::AUTO_WRITE ), this->mesh_ ), epsilon_ ( IOobject ( this->groupName("epsilon"), this->runTime_.name(), this->mesh_, IOobject::MUST_READ, IOobject::AUTO_WRITE ), this->mesh_ ), kU_ ( IOobject ( this->groupName("kU"), this->runTime_.name(), this->mesh_, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE ), //this->mesh_ this->k_*this->fK_ //this->k_.boundaryField().types() ), epsilonU_ ( IOobject ( this->groupName("epsilonU"), this->runTime_.name(), this->mesh_, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE ), //this->mesh_ this->epsilon_*this->fEpsilon_ //this->epsilon_.boundaryField().types() ) { bound(k_, this->kMin_); bound(epsilon_, this->epsilonMin_); bound(kU_, min(fK_)*this->kMin_); bound(epsilonU_, fEpsilon_*this->epsilonMin_); if (type == typeName) { this->printCoeffs(type); correctNut(); } } // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // template<class BasicMomentumTransportModel> bool kEpsilonPANS<BasicMomentumTransportModel>::read() { if (eddyViscosity<RASModel<BasicMomentumTransportModel>>::read()) { Cmu_.readIfPresent(this->coeffDict()); C1_.readIfPresent(this->coeffDict()); C2_.readIfPresent(this->coeffDict()); C3_.readIfPresent(this->coeffDict()); sigmak_.readIfPresent(this->coeffDict()); sigmaEps_.readIfPresent(this->coeffDict()); fEpsilon_.readIfPresent(this->coeffDict()); uLim_.readIfPresent(this->coeffDict()); loLim_.readIfPresent(this->coeffDict()); return true; } else { return false; } } template<class BasicMomentumTransportModel> void kEpsilonPANS<BasicMomentumTransportModel>::correct() { if (!this->turbulence_) { return; } // Local references const alphaField& alpha = this->alpha_; const rhoField& rho = this->rho_; const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_; const volVectorField& U = this->U_; volScalarField& nut = this->nut_; const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_)); const Foam::fvConstraints& fvConstraints ( Foam::fvConstraints::New(this->mesh_) ); eddyViscosity<RASModel<BasicMomentumTransportModel>>::correct(); volScalarField::Internal divU ( fvc::div(fvc::absolute(this->phi(), U))().v() ); tmp<volTensorField> tgradU = fvc::grad(U); volScalarField::Internal G ( this->GName(), nut.v()*(dev(twoSymm(tgradU().v())) && tgradU().v()) ); tgradU.clear(); // Update epsilon and G at the wall epsilonU_.boundaryFieldRef().updateCoeffs(); // Dissipation equation tmp<fvScalarMatrix> epsUEqn ( fvm::ddt(alpha, rho, epsilonU_) + fvm::div(alphaRhoPhi, epsilonU_) - fvm::laplacian(alpha*rho*DepsilonUEff(), epsilonU_) == C1_*alpha()*rho()*G*epsilonU_()/kU_() - fvm::SuSp(((2.0/3.0)*C1_ + C3_)*alpha()*rho()*divU, epsilonU_) - fvm::Sp(C2U*alpha()*rho()*epsilonU_()/kU_(), epsilonU_) + epsilonSource() + fvModels.source(alpha, rho, epsilonU_) ); epsUEqn.ref().relax(); fvConstraints.constrain(epsUEqn.ref()); epsUEqn.ref().boundaryManipulate(epsilonU_.boundaryFieldRef()); solve(epsUEqn); fvConstraints.constrain(epsilonU_); bound(epsilonU_, fEpsilon_*this->epsilonMin_); // Turbulent kinetic energy equation tmp<fvScalarMatrix> kUEqn ( fvm::ddt(alpha, rho, kU_) + fvm::div(alphaRhoPhi, kU_) - fvm::laplacian(alpha*rho*DkUEff(), kU_) == alpha()*rho()*G - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, kU_) - fvm::Sp(alpha()*rho()*epsilonU_()/kU_(), kU_) + kSource() + fvModels.source(alpha, rho, kU_) ); kUEqn.ref().relax(); fvConstraints.constrain(kUEqn.ref()); solve(kUEqn); fvConstraints.constrain(kU_); bound(kU_, min(fK_)*this->kMin_); // Calculation of Turbulent kinetic energy and Dissipation rate k_ = kU_/fK_; k_.correctBoundaryConditions(); epsilon_ = epsilonU_/fEpsilon_; epsilon_.correctBoundaryConditions(); bound(k_, this->kMin_); bound(epsilon_, this->epsilonMin_); correctNut(); // Recalculate fK and C2U with new kU and epsilonU // Calculate the turbulence integral length scale volScalarField::Internal Lambda ( pow(k_,1.5)/epsilon_ ); // update fK fK_.primitiveFieldRef() = min(max( sqrt(Cmu_.value())*pow(delta()/Lambda,2.0/3.0), loLim_), uLim_); // update C2U C2U = C1_ + (fK_/fEpsilon_)*(C2_ - C1_); } // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // } // End namespace RASModels } // End namespace Foam // ************************************************************************* // Code:
/* Description PANS implementation based on the Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distortion theory (RDT) based compression term. The fK parameter, responsible for switching from RAS to DNS is calculated based on [Girimaji,2005]: fK_ = min( max( sqrt(Cmu_.value())*(pow(delta()/ (pow(kU_,1.5)/epsilonU_),2.0/3.0)) , loLimVec ), uLimVec ); bounded between loLim (DNS -> 0.0) and uLim (RAS -> 1.0). The lower the uLim, defined in turbulenceProperties dict, more scales are resolved. Reference: \verbatim Standard model: Launder, B. E., & Spalding, D. B. (1972). Lectures in mathematical models of turbulence. Launder, B. E., & Spalding, D. B. (1974). The numerical computation of turbulent flows. Computer methods in applied mechanics and engineering, 3(2), 269-289. For the RDT-based compression term: El Tahry, S. H. (1983). k-epsilon equation for compressible reciprocating engine flows. Journal of Energy, 7(4), 345-353. PANS: Girimaji, S. S. & Abdol-Hamid K. S. (2005). Partially-averaged Navier-Stokes model for turbulence: Implementation and Validation. 43rd AIAA Aerospace Science Meeting and Exhibit. Girimaji, S. S. (2006). Partially-averaged Navier-Stokes method for turbulence: A Reynolds averaged Navier-Stokes to direct numerical simulation bridging method. Journal of Applied Mechanics, Vol 73, 413-421. Girimaji, S. & Jeong, E. & Srinivasan, R. (2006). Partially-averaged Navier-Stokes method for turbulence: Fixed point analysis and comparison with unsteady partially averaged Navier-Stokes. Journal of Applied Mechanics, Vol 73, 422-429. \endverbatim The default model coefficients are \verbatim kEpsilonPANSCoeffs { Cmu 0.09; C1 1.44; C2 1.92; C3 -0.33; sigmak 1.0; sigmaEps 1.3; fEpsilon 1.0; fKupperLimit 1.0; fKlowerLimit 0.1; // Delta must be specified for PANS e.g. delta cubeRootVol; cubeRootVolCoeffs {} } \endverbatim SourceFiles kEpsilonPANS.C */ \*---------------------------------------------------------------------------*/ #ifndef kEpsilonPANS_H #define kEpsilonPANS_H #include "RASModel.H" #include "eddyViscosity.H" #include "LESdelta.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // namespace Foam { namespace RASModels { /*---------------------------------------------------------------------------*\ Class kEpsilon Declaration \*---------------------------------------------------------------------------*/ template<class BasicMomentumTransportModel> class kEpsilonPANS : public eddyViscosity<RASModel<BasicMomentumTransportModel>> { // Private Member Functions /* // Disallow default bitwise copy construct and assignment kEpsilonPANS(const kEpsilonPANS&); kEpsilonPANS& operator=(const kEpsilonPANS&); */ protected: // Protected data // Model coefficients dimensionedScalar Cmu_; dimensionedScalar C1_; dimensionedScalar C2_; dimensionedScalar C3_; dimensionedScalar sigmak_; dimensionedScalar sigmaEps_; dimensionedScalar fEpsilon_; dimensionedScalar uLim_; dimensionedScalar loLim_; // Fields volScalarField fK_; volScalarField C2U; //- Run-time selectable delta model autoPtr<Foam::LESdelta> delta_; volScalarField k_; volScalarField epsilon_; volScalarField kU_; volScalarField epsilonU_; // Protected Member Functions virtual void correctNut(); virtual tmp<fvScalarMatrix> kSource() const; virtual tmp<fvScalarMatrix> epsilonSource() const; public: typedef typename BasicMomentumTransportModel::alphaField alphaField; typedef typename BasicMomentumTransportModel::rhoField rhoField; //- Runtime type information TypeName("kEpsilonPANS"); // Constructors //- Construct from components kEpsilonPANS ( const alphaField& alpha, const rhoField& rho, const volVectorField& U, const surfaceScalarField& alphaRhoPhi, const surfaceScalarField& phi, const viscosity& viscosity, const word& type = typeName ); //- Disallow default bitwise copy construction kEpsilonPANS(const kEpsilonPANS&) = delete; //- Destructor virtual ~kEpsilonPANS() {} // Member Functions //- Re-read model coefficients if they have changed virtual bool read(); //- Access function to filter width inline const volScalarField& delta() const { return delta_(); } //- Return the effective diffusivity for unresolved k tmp<volScalarField> DkUEff() const { return volScalarField::New ( "DkUEff", (this->nut_/(fK_*fK_*sigmak_/fEpsilon_) + this->nu()) ); } //- Return the effective diffusivity for unresolved epsilon tmp<volScalarField> DepsilonUEff() const { return tmp<volScalarField> ( new volScalarField ( "DepsilonUEff", (this->nut_/(fK_*fK_*sigmaEps_/fEpsilon_) + this->nu()) ) ); } //- Return the turbulence specific dissipation rate virtual tmp<volScalarField> omega() const { return volScalarField::New ( "omega", epsilon_/(Cmu_*k_) ); } //- Return the turbulence kinetic energy virtual tmp<volScalarField> k() const { return k_; } //- Return the turbulence kinetic energy dissipation rate virtual tmp<volScalarField> epsilon() const { return epsilon_; } //- Return the unresolved turbulence kinetic energy virtual tmp<volScalarField> kU() const { return kU_; } //- Return the unresolved turbulence kinetic energy dissipation rate virtual tmp<volScalarField> epsilonU() const { return epsilonU_; } //- Solve the turbulence equations and correct the turbulence viscosity virtual void correct(); // Member Operators //- Disallow default bitwise assignment void operator=(const kEpsilonPANS&) = delete; }; // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // } // End namespace RASModels } // End namespace Foam // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // #ifdef NoRepository #include "kEpsilonPANS.C" #endif // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // #endif // ************************************************************************* // |
|
Tags |
c++ code, openfoam 11, turbulence free |
|
|
Similar Threads | ||||
Thread | Thread Starter | Forum | Replies | Last Post |
convergence problem of steady 2D film cooling calculation using chtMultiRegionFoam | ruanyg968tf | OpenFOAM Running, Solving & CFD | 1 | April 10, 2024 03:23 |
[Commercial meshers] Problem with Mesh conversion from FLUENT Meshing to OpenFOAM | mn17jyf | OpenFOAM Meshing & Mesh Conversion | 3 | November 1, 2023 10:49 |
Instability in buoyantSimpleFoam | Avandri | OpenFOAM Running, Solving & CFD | 0 | August 7, 2020 17:13 |
Pressure instability with rhoSimpleFoam | daniel_mills | OpenFOAM Running, Solving & CFD | 44 | February 17, 2011 18:08 |
Problem with compile the setParabolicInlet | ivanyao | OpenFOAM Running, Solving & CFD | 6 | September 5, 2008 21:50 |