00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include "gaussLaplacianScheme.H"
00027 #include <finiteVolume/surfaceInterpolate.H>
00028 #include <finiteVolume/fvcDiv.H>
00029 #include <finiteVolume/fvcGrad.H>
00030 #include <finiteVolume/fvMatrices.H>
00031
00032
00033
00034 namespace Foam
00035 {
00036
00037
00038
00039 namespace fv
00040 {
00041
00042
00043
00044 template<class Type, class GType>
00045 tmp<fvMatrix<Type> >
00046 gaussLaplacianScheme<Type, GType>::fvmLaplacianUncorrected
00047 (
00048 const surfaceScalarField& gammaMagSf,
00049 GeometricField<Type, fvPatchField, volMesh>& vf
00050 )
00051 {
00052 tmp<surfaceScalarField> tdeltaCoeffs =
00053 this->tsnGradScheme_().deltaCoeffs(vf);
00054 const surfaceScalarField& deltaCoeffs = tdeltaCoeffs();
00055
00056 tmp<fvMatrix<Type> > tfvm
00057 (
00058 new fvMatrix<Type>
00059 (
00060 vf,
00061 deltaCoeffs.dimensions()*gammaMagSf.dimensions()*vf.dimensions()
00062 )
00063 );
00064 fvMatrix<Type>& fvm = tfvm();
00065
00066 fvm.upper() = deltaCoeffs.internalField()*gammaMagSf.internalField();
00067 fvm.negSumDiag();
00068
00069 forAll(fvm.psi().boundaryField(), patchI)
00070 {
00071 const fvPatchField<Type>& psf = fvm.psi().boundaryField()[patchI];
00072 const fvsPatchScalarField& patchGamma =
00073 gammaMagSf.boundaryField()[patchI];
00074
00075 fvm.internalCoeffs()[patchI] = patchGamma*psf.gradientInternalCoeffs();
00076 fvm.boundaryCoeffs()[patchI] = -patchGamma*psf.gradientBoundaryCoeffs();
00077 }
00078
00079 return tfvm;
00080 }
00081
00082
00083 template<class Type, class GType>
00084 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
00085 gaussLaplacianScheme<Type, GType>::gammaSnGradCorr
00086 (
00087 const surfaceVectorField& SfGammaCorr,
00088 const GeometricField<Type, fvPatchField, volMesh>& vf
00089 )
00090 {
00091 const fvMesh& mesh = this->mesh();
00092
00093 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tgammaSnGradCorr
00094 (
00095 new GeometricField<Type, fvsPatchField, surfaceMesh>
00096 (
00097 IOobject
00098 (
00099 "gammaSnGradCorr("+vf.name()+')',
00100 vf.instance(),
00101 mesh,
00102 IOobject::NO_READ,
00103 IOobject::NO_WRITE
00104 ),
00105 mesh,
00106 SfGammaCorr.dimensions()
00107 *vf.dimensions()*mesh.deltaCoeffs().dimensions()
00108 )
00109 );
00110
00111 for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; cmpt++)
00112 {
00113 tgammaSnGradCorr().replace
00114 (
00115 cmpt,
00116 SfGammaCorr & fvc::interpolate(fvc::grad(vf.component(cmpt)))
00117 );
00118 }
00119
00120 return tgammaSnGradCorr;
00121 }
00122
00123
00124
00125
00126 template<class Type, class GType>
00127 tmp<GeometricField<Type, fvPatchField, volMesh> >
00128 gaussLaplacianScheme<Type, GType>::fvcLaplacian
00129 (
00130 const GeometricField<Type, fvPatchField, volMesh>& vf
00131 )
00132 {
00133 const fvMesh& mesh = this->mesh();
00134
00135 tmp<GeometricField<Type, fvPatchField, volMesh> > tLaplacian
00136 (
00137 fvc::div(this->tsnGradScheme_().snGrad(vf)*mesh.magSf())
00138 );
00139
00140 tLaplacian().rename("laplacian(" + vf.name() + ')');
00141
00142 return tLaplacian;
00143 }
00144
00145
00146 template<class Type, class GType>
00147 tmp<fvMatrix<Type> >
00148 gaussLaplacianScheme<Type, GType>::fvmLaplacian
00149 (
00150 const GeometricField<GType, fvsPatchField, surfaceMesh>& gamma,
00151 GeometricField<Type, fvPatchField, volMesh>& vf
00152 )
00153 {
00154 const fvMesh& mesh = this->mesh();
00155
00156 surfaceVectorField Sn = mesh.Sf()/mesh.magSf();
00157
00158 surfaceVectorField SfGamma = mesh.Sf() & gamma;
00159 GeometricField<scalar, fvsPatchField, surfaceMesh> SfGammaSn = SfGamma & Sn;
00160 surfaceVectorField SfGammaCorr = SfGamma - SfGammaSn*Sn;
00161
00162 tmp<fvMatrix<Type> > tfvm = fvmLaplacianUncorrected(SfGammaSn, vf);
00163 fvMatrix<Type>& fvm = tfvm();
00164
00165 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tfaceFluxCorrection
00166 = gammaSnGradCorr(SfGammaCorr, vf);
00167
00168 if (this->tsnGradScheme_().corrected())
00169 {
00170 tfaceFluxCorrection() +=
00171 SfGammaSn*this->tsnGradScheme_().correction(vf);
00172 }
00173
00174 fvm.source() -= mesh.V()*fvc::div(tfaceFluxCorrection())().internalField();
00175
00176 if (mesh.fluxRequired(vf.name()))
00177 {
00178 fvm.faceFluxCorrectionPtr() = tfaceFluxCorrection.ptr();
00179 }
00180
00181 return tfvm;
00182 }
00183
00184
00185 template<class Type, class GType>
00186 tmp<GeometricField<Type, fvPatchField, volMesh> >
00187 gaussLaplacianScheme<Type, GType>::fvcLaplacian
00188 (
00189 const GeometricField<GType, fvsPatchField, surfaceMesh>& gamma,
00190 const GeometricField<Type, fvPatchField, volMesh>& vf
00191 )
00192 {
00193 const fvMesh& mesh = this->mesh();
00194
00195 surfaceVectorField Sn = mesh.Sf()/mesh.magSf();
00196
00197 surfaceVectorField SfGamma = mesh.Sf() & gamma;
00198 GeometricField<scalar, fvsPatchField, surfaceMesh> SfGammaSn = SfGamma & Sn;
00199 surfaceVectorField SfGammaCorr = SfGamma - SfGammaSn*Sn;
00200
00201 tmp<GeometricField<Type, fvPatchField, volMesh> > tLaplacian
00202 (
00203 fvc::div
00204 (
00205 SfGammaSn*this->tsnGradScheme_().snGrad(vf)
00206 + gammaSnGradCorr(SfGammaCorr, vf)
00207 )
00208 );
00209
00210 tLaplacian().rename("laplacian(" + gamma.name() + ',' + vf.name() + ')');
00211
00212 return tLaplacian;
00213 }
00214
00215
00216
00217
00218 }
00219
00220
00221
00222 }
00223
00224