mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
gaussLaplacianScheme: Expose primitive uncorrected implementation
This commit is contained in:
@ -46,13 +46,10 @@ tmp<fvMatrix<Type> >
|
||||
gaussLaplacianScheme<Type, GType>::fvmLaplacianUncorrected
|
||||
(
|
||||
const surfaceScalarField& gammaMagSf,
|
||||
const surfaceScalarField& deltaCoeffs,
|
||||
const GeometricField<Type, fvPatchField, volMesh>& vf
|
||||
)
|
||||
{
|
||||
tmp<surfaceScalarField> tdeltaCoeffs =
|
||||
this->tsnGradScheme_().deltaCoeffs(vf);
|
||||
const surfaceScalarField& deltaCoeffs = tdeltaCoeffs();
|
||||
|
||||
tmp<fvMatrix<Type> > tfvm
|
||||
(
|
||||
new fvMatrix<Type>
|
||||
@ -66,14 +63,14 @@ gaussLaplacianScheme<Type, GType>::fvmLaplacianUncorrected
|
||||
fvm.upper() = deltaCoeffs.internalField()*gammaMagSf.internalField();
|
||||
fvm.negSumDiag();
|
||||
|
||||
forAll(vf.boundaryField(), patchI)
|
||||
forAll(vf.boundaryField(), patchi)
|
||||
{
|
||||
const fvPatchField<Type>& psf = vf.boundaryField()[patchI];
|
||||
const fvPatchField<Type>& psf = vf.boundaryField()[patchi];
|
||||
const fvsPatchScalarField& patchGamma =
|
||||
gammaMagSf.boundaryField()[patchI];
|
||||
gammaMagSf.boundaryField()[patchi];
|
||||
|
||||
fvm.internalCoeffs()[patchI] = patchGamma*psf.gradientInternalCoeffs();
|
||||
fvm.boundaryCoeffs()[patchI] = -patchGamma*psf.gradientBoundaryCoeffs();
|
||||
fvm.internalCoeffs()[patchi] = patchGamma*psf.gradientInternalCoeffs();
|
||||
fvm.boundaryCoeffs()[patchi] = -patchGamma*psf.gradientBoundaryCoeffs();
|
||||
}
|
||||
|
||||
return tfvm;
|
||||
@ -162,7 +159,12 @@ gaussLaplacianScheme<Type, GType>::fvmLaplacian
|
||||
);
|
||||
const surfaceVectorField SfGammaCorr(SfGamma - SfGammaSn*Sn);
|
||||
|
||||
tmp<fvMatrix<Type> > tfvm = fvmLaplacianUncorrected(SfGammaSn, vf);
|
||||
tmp<fvMatrix<Type> > tfvm = fvmLaplacianUncorrected
|
||||
(
|
||||
SfGammaSn,
|
||||
this->tsnGradScheme_().deltaCoeffs(vf),
|
||||
vf
|
||||
);
|
||||
fvMatrix<Type>& fvm = tfvm();
|
||||
|
||||
tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tfaceFluxCorrection
|
||||
|
||||
@ -58,12 +58,6 @@ class gaussLaplacianScheme
|
||||
{
|
||||
// Private Member Functions
|
||||
|
||||
tmp<fvMatrix<Type> > fvmLaplacianUncorrected
|
||||
(
|
||||
const surfaceScalarField& gammaMagSf,
|
||||
const GeometricField<Type, fvPatchField, volMesh>&
|
||||
);
|
||||
|
||||
tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > gammaSnGradCorr
|
||||
(
|
||||
const surfaceVectorField& SfGammaCorr,
|
||||
@ -116,6 +110,13 @@ public:
|
||||
|
||||
// Member Functions
|
||||
|
||||
static tmp<fvMatrix<Type> > fvmLaplacianUncorrected
|
||||
(
|
||||
const surfaceScalarField& gammaMagSf,
|
||||
const surfaceScalarField& deltaCoeffs,
|
||||
const GeometricField<Type, fvPatchField, volMesh>&
|
||||
);
|
||||
|
||||
tmp<GeometricField<Type, fvPatchField, volMesh> > fvcLaplacian
|
||||
(
|
||||
const GeometricField<Type, fvPatchField, volMesh>&
|
||||
|
||||
@ -53,7 +53,12 @@ Foam::fv::gaussLaplacianScheme<Foam::Type, Foam::scalar>::fvmLaplacian \
|
||||
gamma*mesh.magSf() \
|
||||
); \
|
||||
\
|
||||
tmp<fvMatrix<Type> > tfvm = fvmLaplacianUncorrected(gammaMagSf, vf); \
|
||||
tmp<fvMatrix<Type> > tfvm = fvmLaplacianUncorrected \
|
||||
( \
|
||||
gammaMagSf, \
|
||||
this->tsnGradScheme_().deltaCoeffs(vf), \
|
||||
vf \
|
||||
); \
|
||||
fvMatrix<Type>& fvm = tfvm(); \
|
||||
\
|
||||
if (this->tsnGradScheme_().corrected()) \
|
||||
|
||||
Reference in New Issue
Block a user