From 30fb36f7e3d74f7e7646c4871a256432b775c5ad Mon Sep 17 00:00:00 2001 From: andy <andy> Date: Thu, 12 Dec 2013 12:46:42 +0000 Subject: [PATCH] ENH: Refactored blended convection schemes --- .../limitedSchemes/blended/blended.H | 38 ++++++++++++- .../schemes/CoBlended/CoBlended.H | 57 +++++++++++++------ .../schemes/localBlended/localBlended.H | 33 +++++++++-- 3 files changed, 103 insertions(+), 25 deletions(-) diff --git a/src/finiteVolume/interpolation/surfaceInterpolation/limitedSchemes/blended/blended.H b/src/finiteVolume/interpolation/surfaceInterpolation/limitedSchemes/blended/blended.H index 6d207fcdcde..f58b8c2393b 100644 --- a/src/finiteVolume/interpolation/surfaceInterpolation/limitedSchemes/blended/blended.H +++ b/src/finiteVolume/interpolation/surfaceInterpolation/limitedSchemes/blended/blended.H @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation + \\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation \\/ M anipulation | ------------------------------------------------------------------------------- License @@ -36,6 +36,7 @@ SourceFiles #define blended_H #include "limitedSurfaceInterpolationScheme.H" +#include "blendedSchemeBase.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // @@ -49,7 +50,8 @@ namespace Foam template<class Type> class blended : - public limitedSurfaceInterpolationScheme<Type> + public limitedSurfaceInterpolationScheme<Type>, + public blendedSchemeBase<Type> { // Private data @@ -111,8 +113,40 @@ public: {} + //- Destructor + virtual ~blended() + {} + + // Member Functions + //- Return the face-based blending factor + virtual tmp<surfaceScalarField> blendingFactor + ( + const GeometricField<Type, fvPatchField, volMesh>& vf + ) const + { + return tmp<surfaceScalarField> + ( + new surfaceScalarField + ( + IOobject + ( + vf.name() + "BlendingFactor", + this->mesh().time().timeName(), + this->mesh() + ), + this->mesh(), + dimensionedScalar + ( + "blendingFactor", + dimless, + blendingFactor_ + ) + ) + ); + } + //- Return the interpolation limiter virtual tmp<surfaceScalarField> limiter ( diff --git a/src/finiteVolume/interpolation/surfaceInterpolation/schemes/CoBlended/CoBlended.H b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/CoBlended/CoBlended.H index 8d420950b0a..3f15520061c 100644 --- a/src/finiteVolume/interpolation/surfaceInterpolation/schemes/CoBlended/CoBlended.H +++ b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/CoBlended/CoBlended.H @@ -63,6 +63,7 @@ SourceFiles #define CoBlended_H #include "surfaceInterpolationScheme.H" +#include "blendedSchemeBase.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // @@ -76,7 +77,8 @@ namespace Foam template<class Type> class CoBlended : - public surfaceInterpolationScheme<Type> + public surfaceInterpolationScheme<Type>, + public blendedSchemeBase<Type> { // Private data @@ -135,10 +137,7 @@ public: ), faceFlux_ ( - mesh.lookupObject<surfaceScalarField> - ( - word(is) - ) + mesh.lookupObject<surfaceScalarField>(word(is)) ) { if (Co1_ < 0 || Co2_ < 0 || Co1_ >= Co2_) @@ -182,17 +181,39 @@ public: } + //- Destructor + virtual ~CoBlended() + {} + + // Member Functions - //- Return the face-based Courant number blending factor - tmp<surfaceScalarField> blendingFactor() const + //- Return the face-based blending factor + virtual tmp<surfaceScalarField> blendingFactor + ( + const GeometricField<Type, fvPatchField, volMesh>& vf + ) const { - const fvMesh& mesh = faceFlux_.mesh(); + const fvMesh& mesh = this->mesh(); - return + tmp<surfaceScalarField> tbf ( - scalar(1) - - max + new surfaceScalarField + ( + IOobject + ( + vf.name() + "BlendingFactor", + mesh.time().timeName(), + mesh + ), + mesh, + dimensionedScalar("blendingFactor", dimless, 0.0) + ) + ); + + tbf() = + scalar(1) + - max ( min ( @@ -204,8 +225,9 @@ public: scalar(1) ), scalar(0) - ) - ); + ); + + return tbf; } @@ -216,9 +238,10 @@ public: const GeometricField<Type, fvPatchField, volMesh>& vf ) const { - surfaceScalarField bf(blendingFactor()); + surfaceScalarField bf(blendingFactor(vf)); - Info<< "weights " << max(bf) << " " << min(bf) << endl; + Info<< "weights " << max(bf).value() << " " << min(bf).value() + << endl; return bf*tScheme1_().weights(vf) @@ -234,7 +257,7 @@ public: const GeometricField<Type, fvPatchField, volMesh>& vf ) const { - surfaceScalarField bf(blendingFactor()); + surfaceScalarField bf(blendingFactor(vf)); return bf*tScheme1_().interpolate(vf) @@ -257,7 +280,7 @@ public: const GeometricField<Type, fvPatchField, volMesh>& vf ) const { - surfaceScalarField bf(blendingFactor()); + surfaceScalarField bf(blendingFactor(vf)); if (tScheme1_().corrected()) { diff --git a/src/finiteVolume/interpolation/surfaceInterpolation/schemes/localBlended/localBlended.H b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/localBlended/localBlended.H index 5f49be98fd2..0e85f0ea762 100644 --- a/src/finiteVolume/interpolation/surfaceInterpolation/schemes/localBlended/localBlended.H +++ b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/localBlended/localBlended.H @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation + \\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation \\/ M anipulation | ------------------------------------------------------------------------------- License @@ -36,6 +36,7 @@ SourceFiles #define localBlended_H #include "surfaceInterpolationScheme.H" +#include "blendedSchemeBase.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // @@ -49,7 +50,8 @@ namespace Foam template<class Type> class localBlended : - public surfaceInterpolationScheme<Type> + public surfaceInterpolationScheme<Type>, + public blendedSchemeBase<Type> { // Private Member Functions @@ -115,8 +117,27 @@ public: {} + //- Destructor + virtual ~localBlended() + {} + + // Member Functions + //- Return the face-based blending factor + virtual tmp<surfaceScalarField> blendingFactor + ( + const GeometricField<Type, fvPatchField, volMesh>& vf + ) const + { + return + this->mesh().objectRegistry::template + lookupObject<const surfaceScalarField> + ( + word(vf.name() + "BlendingFactor") + ); + } + //- Return the interpolation weighting factors tmp<surfaceScalarField> weights ( @@ -125,10 +146,10 @@ public: { const surfaceScalarField& blendingFactor = this->mesh().objectRegistry::template - lookupObject<const surfaceScalarField> - ( - word(vf.name() + "BlendingFactor") - ); + lookupObject<const surfaceScalarField> + ( + word(vf.name() + "BlendingFactor") + ); return blendingFactor*tScheme1_().weights(vf) -- GitLab