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 <finiteVolume/volFields.H>
00027 #include <finiteVolume/surfaceFields.H>
00028 #include <finiteVolume/fvMatrix.H>
00029 #include <finiteVolume/laplacianScheme.H>
00030
00031
00032
00033 namespace Foam
00034 {
00035
00036
00037
00038 namespace fvm
00039 {
00040
00041
00042
00043 template<class Type>
00044 tmp<fvMatrix<Type> >
00045 laplacian
00046 (
00047 GeometricField<Type, fvPatchField, volMesh>& vf,
00048 const word& name
00049 )
00050 {
00051 surfaceScalarField Gamma
00052 (
00053 IOobject
00054 (
00055 "1",
00056 vf.time().constant(),
00057 vf.mesh(),
00058 IOobject::NO_READ
00059 ),
00060 vf.mesh(),
00061 dimensionedScalar("1", dimless, 1.0)
00062 );
00063
00064 return fvm::laplacian(Gamma, vf, name);
00065 }
00066
00067
00068 template<class Type>
00069 tmp<fvMatrix<Type> >
00070 laplacian
00071 (
00072 GeometricField<Type, fvPatchField, volMesh>& vf
00073 )
00074 {
00075 surfaceScalarField Gamma
00076 (
00077 IOobject
00078 (
00079 "1",
00080 vf.time().constant(),
00081 vf.mesh(),
00082 IOobject::NO_READ
00083 ),
00084 vf.mesh(),
00085 dimensionedScalar("1", dimless, 1.0)
00086 );
00087
00088 return fvm::laplacian
00089 (
00090 Gamma,
00091 vf,
00092 "laplacian(" + vf.name() + ')'
00093 );
00094 }
00095
00096
00097 template<class Type>
00098 tmp<fvMatrix<Type> >
00099 laplacian
00100 (
00101 const zeroField&,
00102 GeometricField<Type, fvPatchField, volMesh>& vf,
00103 const word& name
00104 )
00105 {
00106 return tmp<fvMatrix<Type> >
00107 (
00108 new fvMatrix<Type>(vf, dimensionSet(0, 0, -2, 0, 0))
00109 );
00110 }
00111
00112
00113 template<class Type>
00114 tmp<fvMatrix<Type> >
00115 laplacian
00116 (
00117 const zeroField&,
00118 GeometricField<Type, fvPatchField, volMesh>& vf
00119 )
00120 {
00121 return tmp<fvMatrix<Type> >
00122 (
00123 new fvMatrix<Type>(vf, dimensionSet(0, 0, -2, 0, 0))
00124 );
00125 }
00126
00127
00128 template<class Type>
00129 tmp<fvMatrix<Type> >
00130 laplacian
00131 (
00132 const geometricOneField&,
00133 GeometricField<Type, fvPatchField, volMesh>& vf,
00134 const word& name
00135 )
00136 {
00137 return fvm::laplacian(vf, name);
00138 }
00139
00140
00141 template<class Type>
00142 tmp<fvMatrix<Type> >
00143 laplacian
00144 (
00145 const geometricOneField&,
00146 GeometricField<Type, fvPatchField, volMesh>& vf
00147 )
00148 {
00149 return fvm::laplacian(vf);
00150 }
00151
00152
00153 template<class Type, class GType>
00154 tmp<fvMatrix<Type> >
00155 laplacian
00156 (
00157 const dimensioned<GType>& gamma,
00158 GeometricField<Type, fvPatchField, volMesh>& vf,
00159 const word& name
00160 )
00161 {
00162 GeometricField<GType, fvsPatchField, surfaceMesh> Gamma
00163 (
00164 IOobject
00165 (
00166 gamma.name(),
00167 vf.instance(),
00168 vf.mesh(),
00169 IOobject::NO_READ
00170 ),
00171 vf.mesh(),
00172 gamma
00173 );
00174
00175 return fvm::laplacian(Gamma, vf, name);
00176 }
00177
00178
00179 template<class Type, class GType>
00180 tmp<fvMatrix<Type> >
00181 laplacian
00182 (
00183 const dimensioned<GType>& gamma,
00184 GeometricField<Type, fvPatchField, volMesh>& vf
00185 )
00186 {
00187 GeometricField<GType, fvsPatchField, surfaceMesh> Gamma
00188 (
00189 IOobject
00190 (
00191 gamma.name(),
00192 vf.instance(),
00193 vf.mesh(),
00194 IOobject::NO_READ
00195 ),
00196 vf.mesh(),
00197 gamma
00198 );
00199
00200 return fvm::laplacian(Gamma, vf);
00201 }
00202
00203
00204
00205
00206 template<class Type, class GType>
00207 tmp<fvMatrix<Type> >
00208 laplacian
00209 (
00210 const GeometricField<GType, fvPatchField, volMesh>& gamma,
00211 GeometricField<Type, fvPatchField, volMesh>& vf,
00212 const word& name
00213 )
00214 {
00215 return fv::laplacianScheme<Type, GType>::New
00216 (
00217 vf.mesh(),
00218 vf.mesh().laplacianScheme(name)
00219 )().fvmLaplacian(gamma, vf);
00220 }
00221
00222
00223 template<class Type, class GType>
00224 tmp<fvMatrix<Type> >
00225 laplacian
00226 (
00227 const tmp<GeometricField<GType, fvPatchField, volMesh> >& tgamma,
00228 GeometricField<Type, fvPatchField, volMesh>& vf,
00229 const word& name
00230 )
00231 {
00232 tmp<fvMatrix<Type> > Laplacian(fvm::laplacian(tgamma(), vf, name));
00233 tgamma.clear();
00234 return Laplacian;
00235 }
00236
00237
00238 template<class Type, class GType>
00239 tmp<fvMatrix<Type> >
00240 laplacian
00241 (
00242 const GeometricField<GType, fvPatchField, volMesh>& gamma,
00243 GeometricField<Type, fvPatchField, volMesh>& vf
00244 )
00245 {
00246 return fvm::laplacian
00247 (
00248 gamma,
00249 vf,
00250 "laplacian(" + gamma.name() + ',' + vf.name() + ')'
00251 );
00252 }
00253
00254
00255 template<class Type, class GType>
00256 tmp<fvMatrix<Type> >
00257 laplacian
00258 (
00259 const tmp<GeometricField<GType, fvPatchField, volMesh> >& tgamma,
00260 GeometricField<Type, fvPatchField, volMesh>& vf
00261 )
00262 {
00263 tmp<fvMatrix<Type> > Laplacian(fvm::laplacian(tgamma(), vf));
00264 tgamma.clear();
00265 return Laplacian;
00266 }
00267
00268
00269
00270
00271 template<class Type, class GType>
00272 tmp<fvMatrix<Type> >
00273 laplacian
00274 (
00275 const GeometricField<GType, fvsPatchField, surfaceMesh>& gamma,
00276 GeometricField<Type, fvPatchField, volMesh>& vf,
00277 const word& name
00278 )
00279 {
00280 return fv::laplacianScheme<Type, GType>::New
00281 (
00282 vf.mesh(),
00283 vf.mesh().laplacianScheme(name)
00284 )().fvmLaplacian(gamma, vf);
00285 }
00286
00287
00288 template<class Type, class GType>
00289 tmp<fvMatrix<Type> >
00290 laplacian
00291 (
00292 const tmp<GeometricField<GType, fvsPatchField, surfaceMesh> >& tgamma,
00293 GeometricField<Type, fvPatchField, volMesh>& vf,
00294 const word& name
00295 )
00296 {
00297 tmp<fvMatrix<Type> > tLaplacian = fvm::laplacian(tgamma(), vf, name);
00298 tgamma.clear();
00299 return tLaplacian;
00300 }
00301
00302
00303 template<class Type, class GType>
00304 tmp<fvMatrix<Type> >
00305 laplacian
00306 (
00307 const GeometricField<GType, fvsPatchField, surfaceMesh>& gamma,
00308 GeometricField<Type, fvPatchField, volMesh>& vf
00309 )
00310 {
00311 return fvm::laplacian
00312 (
00313 gamma,
00314 vf,
00315 "laplacian(" + gamma.name() + ',' + vf.name() + ')'
00316 );
00317 }
00318
00319
00320 template<class Type, class GType>
00321 tmp<fvMatrix<Type> >
00322 laplacian
00323 (
00324 const tmp<GeometricField<GType, fvsPatchField, surfaceMesh> >& tGamma,
00325 GeometricField<Type, fvPatchField, volMesh>& vf
00326 )
00327 {
00328 tmp<fvMatrix<Type> > tfvm(fvm::laplacian(tGamma(), vf));
00329 tGamma.clear();
00330 return tfvm;
00331 }
00332
00333
00334
00335
00336 }
00337
00338
00339
00340 }
00341
00342