mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
BUG: blendingFactor function object updated to handle bounded schemes. Fixes #205
This commit is contained in:
@ -3,7 +3,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) 2013-2015 OpenFOAM Foundation
|
\\ / A nd | Copyright (C) 2013-2015 OpenFOAM Foundation
|
||||||
\\/ M anipulation |
|
\\/ M anipulation | Copyright (C) 2016 OpenCFD Ltd.
|
||||||
-------------------------------------------------------------------------------
|
-------------------------------------------------------------------------------
|
||||||
License
|
License
|
||||||
This file is part of OpenFOAM.
|
This file is part of OpenFOAM.
|
||||||
@ -88,6 +88,8 @@ SourceFiles
|
|||||||
#include "OFstream.H"
|
#include "OFstream.H"
|
||||||
#include "Switch.H"
|
#include "Switch.H"
|
||||||
|
|
||||||
|
#include "convectionScheme.H"
|
||||||
|
|
||||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
namespace Foam
|
namespace Foam
|
||||||
@ -140,6 +142,14 @@ class blendingFactor
|
|||||||
//- Disallow default bitwise assignment
|
//- Disallow default bitwise assignment
|
||||||
void operator=(const blendingFactor&);
|
void operator=(const blendingFactor&);
|
||||||
|
|
||||||
|
//- Helper function to calculate the blending factor for the scheme
|
||||||
|
template<class Type>
|
||||||
|
void calcScheme
|
||||||
|
(
|
||||||
|
const GeometricField<Type, fvPatchField, volMesh>& field,
|
||||||
|
const typename fv::convectionScheme<Type>& cs
|
||||||
|
);
|
||||||
|
|
||||||
//- Calculate the blending factor
|
//- Calculate the blending factor
|
||||||
template<class Type>
|
template<class Type>
|
||||||
void calc();
|
void calc();
|
||||||
|
|||||||
@ -3,7 +3,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) 2013-2016 OpenFOAM Foundation
|
\\ / A nd | Copyright (C) 2013-2016 OpenFOAM Foundation
|
||||||
\\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
|
\\/ M anipulation | Copyright (C) 2015-2016 OpenCFD Ltd.
|
||||||
-------------------------------------------------------------------------------
|
-------------------------------------------------------------------------------
|
||||||
License
|
License
|
||||||
This file is part of OpenFOAM.
|
This file is part of OpenFOAM.
|
||||||
@ -24,6 +24,7 @@ License
|
|||||||
\*---------------------------------------------------------------------------*/
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
#include "gaussConvectionScheme.H"
|
#include "gaussConvectionScheme.H"
|
||||||
|
#include "boundedConvectionScheme.H"
|
||||||
#include "blendedSchemeBase.H"
|
#include "blendedSchemeBase.H"
|
||||||
#include "fvcCellReduce.H"
|
#include "fvcCellReduce.H"
|
||||||
#include "zeroGradientFvPatchFields.H"
|
#include "zeroGradientFvPatchFields.H"
|
||||||
@ -31,39 +32,35 @@ License
|
|||||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
template<class Type>
|
template<class Type>
|
||||||
void Foam::blendingFactor::calc()
|
void Foam::blendingFactor::calcScheme
|
||||||
|
(
|
||||||
|
const GeometricField<Type, fvPatchField, volMesh>& field,
|
||||||
|
const typename fv::convectionScheme<Type>& cs
|
||||||
|
)
|
||||||
{
|
{
|
||||||
typedef GeometricField<Type, fvPatchField, volMesh> fieldType;
|
if (!isA<fv::gaussConvectionScheme<Type>>(cs))
|
||||||
|
{
|
||||||
|
WarningInFunction
|
||||||
|
<< "Scheme for field " << field.name() << " is not a "
|
||||||
|
<< fv::gaussConvectionScheme<Type>::typeName
|
||||||
|
<< " scheme. Not calculating " << resultName_ << endl;
|
||||||
|
|
||||||
if (!obr_.foundObject<fieldType>(fieldName_))
|
|
||||||
{
|
|
||||||
return;
|
return;
|
||||||
}
|
}
|
||||||
|
|
||||||
const fvMesh& mesh = refCast<const fvMesh>(obr_);
|
|
||||||
|
|
||||||
const fieldType& field = mesh.lookupObject<fieldType>(fieldName_);
|
|
||||||
|
|
||||||
const word divScheme("div(" + phiName_ + ',' + fieldName_ + ')');
|
|
||||||
ITstream& its = mesh.divScheme(divScheme);
|
|
||||||
|
|
||||||
const surfaceScalarField& phi =
|
|
||||||
mesh.lookupObject<surfaceScalarField>(phiName_);
|
|
||||||
|
|
||||||
tmp<fv::convectionScheme<Type>> cs =
|
|
||||||
fv::convectionScheme<Type>::New(mesh, phi, its);
|
|
||||||
|
|
||||||
const fv::gaussConvectionScheme<Type>& gcs =
|
const fv::gaussConvectionScheme<Type>& gcs =
|
||||||
refCast<const fv::gaussConvectionScheme<Type>>(cs());
|
refCast<const fv::gaussConvectionScheme<Type>>(cs);
|
||||||
|
|
||||||
const surfaceInterpolationScheme<Type>& interpScheme =
|
|
||||||
gcs.interpScheme();
|
const surfaceInterpolationScheme<Type>& interpScheme = gcs.interpScheme();
|
||||||
|
|
||||||
if (!isA<blendedSchemeBase<Type>>(interpScheme))
|
if (!isA<blendedSchemeBase<Type>>(interpScheme))
|
||||||
{
|
{
|
||||||
FatalErrorInFunction
|
WarningInFunction
|
||||||
<< interpScheme.typeName << " is not a blended scheme"
|
<< interpScheme.type() << " is not a blended scheme"
|
||||||
<< exit(FatalError);
|
<< ". Not calculating " << resultName_ << endl;
|
||||||
|
|
||||||
|
return;
|
||||||
}
|
}
|
||||||
|
|
||||||
// Retrieve the face-based blending factor
|
// Retrieve the face-based blending factor
|
||||||
@ -123,4 +120,46 @@ void Foam::blendingFactor::calc()
|
|||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
template<class Type>
|
||||||
|
void Foam::blendingFactor::calc()
|
||||||
|
{
|
||||||
|
typedef GeometricField<Type, fvPatchField, volMesh> fieldType;
|
||||||
|
|
||||||
|
if (!obr_.foundObject<fieldType>(fieldName_))
|
||||||
|
{
|
||||||
|
return;
|
||||||
|
}
|
||||||
|
|
||||||
|
const fvMesh& mesh = refCast<const fvMesh>(obr_);
|
||||||
|
|
||||||
|
const fieldType& field = mesh.lookupObject<fieldType>(fieldName_);
|
||||||
|
|
||||||
|
const word divScheme("div(" + phiName_ + ',' + fieldName_ + ')');
|
||||||
|
ITstream& its = mesh.divScheme(divScheme);
|
||||||
|
|
||||||
|
const surfaceScalarField& phi =
|
||||||
|
mesh.lookupObject<surfaceScalarField>(phiName_);
|
||||||
|
|
||||||
|
tmp<fv::convectionScheme<Type>> tcs =
|
||||||
|
fv::convectionScheme<Type>::New(mesh, phi, its);
|
||||||
|
|
||||||
|
const fv::convectionScheme<Type>& cs = tcs();
|
||||||
|
|
||||||
|
if (isA<fv::boundedConvectionScheme<Type>>(cs))
|
||||||
|
{
|
||||||
|
const fv::boundedConvectionScheme<Type>& bcs =
|
||||||
|
refCast<const fv::boundedConvectionScheme<Type>>(cs);
|
||||||
|
|
||||||
|
calcScheme(field, bcs.scheme());
|
||||||
|
}
|
||||||
|
else
|
||||||
|
{
|
||||||
|
const fv::gaussConvectionScheme<Type>& gcs =
|
||||||
|
refCast<const fv::gaussConvectionScheme<Type>>(cs);
|
||||||
|
|
||||||
|
calcScheme(field, gcs);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
// ************************************************************************* //
|
// ************************************************************************* //
|
||||||
|
|||||||
Reference in New Issue
Block a user