Skip to content
Snippets Groups Projects
Commit 30fb36f7 authored by andy's avatar andy
Browse files

ENH: Refactored blended convection schemes

parent 93a4ba56
Branches
Tags
No related merge requests found
...@@ -2,7 +2,7 @@ ...@@ -2,7 +2,7 @@
========= | ========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | \\ / O peration |
\\ / A nd | Copyright (C) 2011 OpenFOAM Foundation \\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation
\\/ M anipulation | \\/ M anipulation |
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
...@@ -36,6 +36,7 @@ SourceFiles ...@@ -36,6 +36,7 @@ SourceFiles
#define blended_H #define blended_H
#include "limitedSurfaceInterpolationScheme.H" #include "limitedSurfaceInterpolationScheme.H"
#include "blendedSchemeBase.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
...@@ -49,7 +50,8 @@ namespace Foam ...@@ -49,7 +50,8 @@ namespace Foam
template<class Type> template<class Type>
class blended class blended
: :
public limitedSurfaceInterpolationScheme<Type> public limitedSurfaceInterpolationScheme<Type>,
public blendedSchemeBase<Type>
{ {
// Private data // Private data
...@@ -111,8 +113,40 @@ public: ...@@ -111,8 +113,40 @@ public:
{} {}
//- Destructor
virtual ~blended()
{}
// Member Functions // 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 //- Return the interpolation limiter
virtual tmp<surfaceScalarField> limiter virtual tmp<surfaceScalarField> limiter
( (
......
...@@ -63,6 +63,7 @@ SourceFiles ...@@ -63,6 +63,7 @@ SourceFiles
#define CoBlended_H #define CoBlended_H
#include "surfaceInterpolationScheme.H" #include "surfaceInterpolationScheme.H"
#include "blendedSchemeBase.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
...@@ -76,7 +77,8 @@ namespace Foam ...@@ -76,7 +77,8 @@ namespace Foam
template<class Type> template<class Type>
class CoBlended class CoBlended
: :
public surfaceInterpolationScheme<Type> public surfaceInterpolationScheme<Type>,
public blendedSchemeBase<Type>
{ {
// Private data // Private data
...@@ -135,10 +137,7 @@ public: ...@@ -135,10 +137,7 @@ public:
), ),
faceFlux_ faceFlux_
( (
mesh.lookupObject<surfaceScalarField> mesh.lookupObject<surfaceScalarField>(word(is))
(
word(is)
)
) )
{ {
if (Co1_ < 0 || Co2_ < 0 || Co1_ >= Co2_) if (Co1_ < 0 || Co2_ < 0 || Co1_ >= Co2_)
...@@ -182,17 +181,39 @@ public: ...@@ -182,17 +181,39 @@ public:
} }
//- Destructor
virtual ~CoBlended()
{}
// Member Functions // Member Functions
//- Return the face-based Courant number blending factor //- Return the face-based blending factor
tmp<surfaceScalarField> blendingFactor() const 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) - new surfaceScalarField
max (
IOobject
(
vf.name() + "BlendingFactor",
mesh.time().timeName(),
mesh
),
mesh,
dimensionedScalar("blendingFactor", dimless, 0.0)
)
);
tbf() =
scalar(1)
- max
( (
min min
( (
...@@ -204,8 +225,9 @@ public: ...@@ -204,8 +225,9 @@ public:
scalar(1) scalar(1)
), ),
scalar(0) scalar(0)
) );
);
return tbf;
} }
...@@ -216,9 +238,10 @@ public: ...@@ -216,9 +238,10 @@ public:
const GeometricField<Type, fvPatchField, volMesh>& vf const GeometricField<Type, fvPatchField, volMesh>& vf
) const ) 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 return
bf*tScheme1_().weights(vf) bf*tScheme1_().weights(vf)
...@@ -234,7 +257,7 @@ public: ...@@ -234,7 +257,7 @@ public:
const GeometricField<Type, fvPatchField, volMesh>& vf const GeometricField<Type, fvPatchField, volMesh>& vf
) const ) const
{ {
surfaceScalarField bf(blendingFactor()); surfaceScalarField bf(blendingFactor(vf));
return return
bf*tScheme1_().interpolate(vf) bf*tScheme1_().interpolate(vf)
...@@ -257,7 +280,7 @@ public: ...@@ -257,7 +280,7 @@ public:
const GeometricField<Type, fvPatchField, volMesh>& vf const GeometricField<Type, fvPatchField, volMesh>& vf
) const ) const
{ {
surfaceScalarField bf(blendingFactor()); surfaceScalarField bf(blendingFactor(vf));
if (tScheme1_().corrected()) if (tScheme1_().corrected())
{ {
......
...@@ -2,7 +2,7 @@ ...@@ -2,7 +2,7 @@
========= | ========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | \\ / O peration |
\\ / A nd | Copyright (C) 2011 OpenFOAM Foundation \\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation
\\/ M anipulation | \\/ M anipulation |
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
...@@ -36,6 +36,7 @@ SourceFiles ...@@ -36,6 +36,7 @@ SourceFiles
#define localBlended_H #define localBlended_H
#include "surfaceInterpolationScheme.H" #include "surfaceInterpolationScheme.H"
#include "blendedSchemeBase.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
...@@ -49,7 +50,8 @@ namespace Foam ...@@ -49,7 +50,8 @@ namespace Foam
template<class Type> template<class Type>
class localBlended class localBlended
: :
public surfaceInterpolationScheme<Type> public surfaceInterpolationScheme<Type>,
public blendedSchemeBase<Type>
{ {
// Private Member Functions // Private Member Functions
...@@ -115,8 +117,27 @@ public: ...@@ -115,8 +117,27 @@ public:
{} {}
//- Destructor
virtual ~localBlended()
{}
// Member Functions // 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 //- Return the interpolation weighting factors
tmp<surfaceScalarField> weights tmp<surfaceScalarField> weights
( (
...@@ -125,10 +146,10 @@ public: ...@@ -125,10 +146,10 @@ public:
{ {
const surfaceScalarField& blendingFactor = const surfaceScalarField& blendingFactor =
this->mesh().objectRegistry::template this->mesh().objectRegistry::template
lookupObject<const surfaceScalarField> lookupObject<const surfaceScalarField>
( (
word(vf.name() + "BlendingFactor") word(vf.name() + "BlendingFactor")
); );
return return
blendingFactor*tScheme1_().weights(vf) blendingFactor*tScheme1_().weights(vf)
......
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment