CoBlended: added support for mass-flux

but currently assuming that the density needed to calculate the volumetric flux is "rho"
This commit is contained in:
Henry
2014-01-23 15:18:54 +00:00
parent a7ad50ce7d
commit 9368035620

View File

@ -64,6 +64,7 @@ SourceFiles
#include "surfaceInterpolationScheme.H" #include "surfaceInterpolationScheme.H"
#include "blendedSchemeBase.H" #include "blendedSchemeBase.H"
#include "surfaceInterpolate.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -181,11 +182,6 @@ public:
} }
//- Destructor
virtual ~CoBlended()
{}
// Member Functions // Member Functions
//- Return the face-based blending factor //- Return the face-based blending factor
@ -195,39 +191,48 @@ public:
) const ) const
{ {
const fvMesh& mesh = this->mesh(); const fvMesh& mesh = this->mesh();
tmp<surfaceScalarField> tUflux = faceFlux_;
tmp<surfaceScalarField> tbf if (faceFlux_.dimensions() == dimDensity*dimVelocity*dimArea)
{
// Currently assume that the density field
// corresponding to the mass-flux is named "rho"
const volScalarField& rho =
mesh.objectRegistry::template lookupObject<volScalarField>
("rho");
tUflux = faceFlux_/fvc::interpolate(rho);
}
else if (faceFlux_.dimensions() != dimVelocity*dimArea)
{
FatalErrorIn
(
"CoBlended::blendingFactor()"
) << "dimensions of faceFlux are not correct"
<< exit(FatalError);
}
return tmp<surfaceScalarField>
( (
new surfaceScalarField new surfaceScalarField
( (
IOobject vf.name() + "BlendingFactor",
scalar(1)
- max
( (
vf.name() + "BlendingFactor", min
mesh.time().timeName(), (
mesh (
), mesh.time().deltaT()*mesh.deltaCoeffs()
mesh, *mag(tUflux)/mesh.magSf()
dimensionedScalar("blendingFactor", dimless, 0.0) - Co1_
)/(Co2_ - Co1_),
scalar(1)
),
scalar(0)
)
) )
); );
tbf() =
scalar(1)
- max
(
min
(
(
mesh.time().deltaT()*mesh.deltaCoeffs()
*mag(faceFlux_)/mesh.magSf()
- Co1_
)/(Co2_ - Co1_),
scalar(1)
),
scalar(0)
);
return tbf;
} }
@ -240,9 +245,6 @@ public:
{ {
surfaceScalarField bf(blendingFactor(vf)); surfaceScalarField bf(blendingFactor(vf));
Info<< "weights " << max(bf).value() << " " << min(bf).value()
<< endl;
return return
bf*tScheme1_().weights(vf) bf*tScheme1_().weights(vf)
+ (scalar(1.0) - bf)*tScheme2_().weights(vf); + (scalar(1.0) - bf)*tScheme2_().weights(vf);