Basic second-order laplacian using face-gradients and Gauss' theorem. More...


Public Member Functions | |
| TypeName ("Gauss") | |
| Runtime type information. More... | |
| gaussLaplacianScheme (const fvMesh &mesh) | |
| Construct null. More... | |
| gaussLaplacianScheme (const fvMesh &mesh, Istream &is) | |
| Construct from Istream. More... | |
| gaussLaplacianScheme (const fvMesh &mesh, const tmp< surfaceInterpolationScheme< GType >> &igs, const tmp< snGradScheme< Type >> &sngs) | |
| Construct from mesh, interpolation and snGradScheme schemes. More... | |
| virtual | ~gaussLaplacianScheme () |
| Destructor. More... | |
| tmp< VolField< Type > > | fvcLaplacian (const VolField< Type > &) |
| tmp< fvMatrix< Type > > | fvmLaplacian (const SurfaceField< GType > &, const VolField< Type > &) |
| tmp< VolField< Type > > | fvcLaplacian (const SurfaceField< GType > &, const VolField< Type > &) |
| tmp< fvMatrix< scalar > > | fvmLaplacian (const SurfaceField< scalar > &, const VolField< scalar > &) |
| tmp< VolField< scalar > > | fvcLaplacian (const SurfaceField< scalar > &, const VolField< scalar > &) |
| tmp< fvMatrix< vector > > | fvmLaplacian (const SurfaceField< scalar > &, const VolField< vector > &) |
| tmp< VolField< vector > > | fvcLaplacian (const SurfaceField< scalar > &, const VolField< vector > &) |
| tmp< fvMatrix< sphericalTensor > > | fvmLaplacian (const SurfaceField< scalar > &, const VolField< sphericalTensor > &) |
| tmp< VolField< sphericalTensor > > | fvcLaplacian (const SurfaceField< scalar > &, const VolField< sphericalTensor > &) |
| tmp< fvMatrix< symmTensor > > | fvmLaplacian (const SurfaceField< scalar > &, const VolField< symmTensor > &) |
| tmp< VolField< symmTensor > > | fvcLaplacian (const SurfaceField< scalar > &, const VolField< symmTensor > &) |
| tmp< fvMatrix< tensor > > | fvmLaplacian (const SurfaceField< scalar > &, const VolField< tensor > &) |
| tmp< VolField< tensor > > | fvcLaplacian (const SurfaceField< scalar > &, const VolField< tensor > &) |
| Foam::tmp< Foam::VolField< Foam::scalar > > | fvcLaplacian (const SurfaceField< scalar > &gamma, const VolField< scalar > &vf) |
| Foam::tmp< Foam::fvMatrix< Foam::vector > > | fvmLaplacian (const SurfaceField< scalar > &gamma, const VolField< vector > &vf) |
| Foam::tmp< Foam::VolField< Foam::vector > > | fvcLaplacian (const SurfaceField< scalar > &gamma, const VolField< vector > &vf) |
| Foam::tmp< Foam::fvMatrix< Foam::sphericalTensor > > | fvmLaplacian (const SurfaceField< scalar > &gamma, const VolField< sphericalTensor > &vf) |
| Foam::tmp< Foam::VolField< Foam::sphericalTensor > > | fvcLaplacian (const SurfaceField< scalar > &gamma, const VolField< sphericalTensor > &vf) |
| Foam::tmp< Foam::fvMatrix< Foam::symmTensor > > | fvmLaplacian (const SurfaceField< scalar > &gamma, const VolField< symmTensor > &vf) |
| Foam::tmp< Foam::VolField< Foam::symmTensor > > | fvcLaplacian (const SurfaceField< scalar > &gamma, const VolField< symmTensor > &vf) |
| Foam::tmp< Foam::fvMatrix< Foam::tensor > > | fvmLaplacian (const SurfaceField< scalar > &gamma, const VolField< tensor > &vf) |
| Foam::tmp< Foam::VolField< Foam::tensor > > | fvcLaplacian (const SurfaceField< scalar > &gamma, const VolField< tensor > &vf) |
Public Member Functions inherited from laplacianScheme< Type, GType > | |
| virtual const word & | type () const =0 |
| Runtime type information. More... | |
| declareRunTimeSelectionTable (tmp, laplacianScheme, Istream,(const fvMesh &mesh, Istream &schemeData),(mesh, schemeData)) | |
| laplacianScheme (const fvMesh &mesh) | |
| Construct from mesh. More... | |
| laplacianScheme (const fvMesh &mesh, Istream &is) | |
| Construct from mesh and Istream. More... | |
| laplacianScheme (const fvMesh &mesh, const tmp< surfaceInterpolationScheme< GType >> &igs, const tmp< snGradScheme< Type >> &sngs) | |
| Construct from mesh, interpolation and snGradScheme schemes. More... | |
| laplacianScheme (const laplacianScheme &)=delete | |
| Disallow default bitwise copy construction. More... | |
| virtual | ~laplacianScheme () |
| Destructor. More... | |
| const fvMesh & | mesh () const |
| Return mesh reference. More... | |
| virtual tmp< fvMatrix< Type > > | fvmLaplacian (const VolField< GType > &, const VolField< Type > &) |
| virtual tmp< VolField< Type > > | fvcLaplacian (const VolField< GType > &, const VolField< Type > &) |
| void | operator= (const laplacianScheme &)=delete |
| Disallow default bitwise assignment. More... | |
Public Member Functions inherited from refCount | |
| int | count () const |
| Return the current reference count. More... | |
| bool | unique () const |
| Return true if the reference count is zero. More... | |
| void | operator++ () |
| Increment the reference count. More... | |
| void | operator++ (int) |
| Increment the reference count. More... | |
| void | operator-- () |
| Decrement the reference count. More... | |
| void | operator-- (int) |
| Decrement the reference count. More... | |
| void | operator= (const refCount &)=delete |
| Disallow bitwise assignment. More... | |
Static Public Member Functions | |
| static tmp< fvMatrix< Type > > | fvmLaplacianUncorrected (const surfaceScalarField &gammaMagSf, const surfaceScalarField &deltaCoeffs, const VolField< Type > &) |
Static Public Member Functions inherited from laplacianScheme< Type, GType > | |
| static tmp< laplacianScheme< Type, GType > > | New (const fvMesh &mesh, Istream &schemeData) |
| Return a pointer to a new laplacianScheme created on freestore. More... | |
Additional Inherited Members | |
Protected Member Functions inherited from refCount | |
| refCount () | |
| Construct null initialising count to 0. More... | |
| refCount (const refCount &)=delete | |
| Disallow copy. More... | |
Protected Attributes inherited from laplacianScheme< Type, GType > | |
| const fvMesh & | mesh_ |
| tmp< surfaceInterpolationScheme< GType > > | tinterpGammaScheme_ |
| tmp< snGradScheme< Type > > | tsnGradScheme_ |
Basic second-order laplacian using face-gradients and Gauss' theorem.
Definition at line 54 of file gaussLaplacianScheme.H.
|
inline |
Construct null.
Definition at line 82 of file gaussLaplacianScheme.H.
|
inline |
Construct from Istream.
Definition at line 88 of file gaussLaplacianScheme.H.
|
inline |
Construct from mesh, interpolation and snGradScheme schemes.
Definition at line 94 of file gaussLaplacianScheme.H.
|
inlinevirtual |
Destructor.
Definition at line 106 of file gaussLaplacianScheme.H.
| TypeName | ( | "Gauss" | ) |
Runtime type information.
|
static |
Definition at line 46 of file gaussLaplacianScheme.C.
References fvMatrix< Type >::boundaryCoeffs(), GeometricField< Type, GeoMesh, PrimitiveField >::boundaryField(), fvPatchField< Type >::coupled(), DimensionedField< Type, GeoMesh, PrimitiveField >::dimensions(), forAll, fvPatchField< Type >::gradientBoundaryCoeffs(), fvPatchField< Type >::gradientInternalCoeffs(), fvMatrix< Type >::internalCoeffs(), lduMatrix::negSumDiag(), patchi, GeometricField< Type, GeoMesh, PrimitiveField >::primitiveField(), tmp< T >::ref(), and lduMatrix::upper().

Implements laplacianScheme< Type, GType >.
Definition at line 129 of file gaussLaplacianScheme.C.
References Foam::fvc::div(), fvMesh::magSf(), mesh, IOobject::name(), tmp< T >::ref(), and Foam::fvc::snGrad().

|
virtual |
Implements laplacianScheme< Type, GType >.
Definition at line 149 of file gaussLaplacianScheme.C.
References Foam::fvc::div(), fvMatrix< Type >::faceFluxCorrectionPtr(), fvSchemes::fluxRequired(), fvMesh::magSf(), mesh, IOobject::name(), tmp< T >::ptr(), tmp< T >::ref(), fvMesh::schemes(), fvMesh::Sf(), Foam::Sn(), fvMatrix< Type >::source(), and fvMesh::V().

|
virtual |
Implements laplacianScheme< Type, GType >.
Definition at line 196 of file gaussLaplacianScheme.C.
References Foam::fvc::div(), fvMesh::magSf(), mesh, IOobject::name(), tmp< T >::ref(), fvMesh::Sf(), Foam::Sn(), and Foam::fvc::snGrad().

| tmp< fvMatrix< scalar > > fvmLaplacian | ( | const SurfaceField< scalar > & | , |
| const VolField< scalar > & | |||
| ) |
| tmp< VolField< scalar > > fvcLaplacian | ( | const SurfaceField< scalar > & | , |
| const VolField< scalar > & | |||
| ) |
| tmp< fvMatrix< vector > > fvmLaplacian | ( | const SurfaceField< scalar > & | , |
| const VolField< vector > & | |||
| ) |
| tmp< VolField< vector > > fvcLaplacian | ( | const SurfaceField< scalar > & | , |
| const VolField< vector > & | |||
| ) |
| tmp< fvMatrix< sphericalTensor > > fvmLaplacian | ( | const SurfaceField< scalar > & | , |
| const VolField< sphericalTensor > & | |||
| ) |
| tmp< VolField< sphericalTensor > > fvcLaplacian | ( | const SurfaceField< scalar > & | , |
| const VolField< sphericalTensor > & | |||
| ) |
| tmp< fvMatrix< symmTensor > > fvmLaplacian | ( | const SurfaceField< scalar > & | , |
| const VolField< symmTensor > & | |||
| ) |
| tmp< VolField< symmTensor > > fvcLaplacian | ( | const SurfaceField< scalar > & | , |
| const VolField< symmTensor > & | |||
| ) |
| tmp< fvMatrix< tensor > > fvmLaplacian | ( | const SurfaceField< scalar > & | , |
| const VolField< tensor > & | |||
| ) |
| tmp< VolField< tensor > > fvcLaplacian | ( | const SurfaceField< scalar > & | , |
| const VolField< tensor > & | |||
| ) |
| Foam::tmp< Foam::VolField< Foam::scalar > > fvcLaplacian | ( | const SurfaceField< scalar > & | gamma, |
| const VolField< scalar > & | vf | ||
| ) |
Definition at line 114 of file gaussLaplacianSchemes.C.
| Foam::tmp< Foam::fvMatrix< Foam::vector > > fvmLaplacian | ( | const SurfaceField< scalar > & | gamma, |
| const VolField< vector > & | vf | ||
| ) |
Definition at line 115 of file gaussLaplacianSchemes.C.
| Foam::tmp< Foam::VolField< Foam::vector > > fvcLaplacian | ( | const SurfaceField< scalar > & | gamma, |
| const VolField< vector > & | vf | ||
| ) |
Definition at line 115 of file gaussLaplacianSchemes.C.
| Foam::tmp< Foam::fvMatrix< Foam::sphericalTensor > > fvmLaplacian | ( | const SurfaceField< scalar > & | gamma, |
| const VolField< sphericalTensor > & | vf | ||
| ) |
Definition at line 116 of file gaussLaplacianSchemes.C.
| Foam::tmp< Foam::VolField< Foam::sphericalTensor > > fvcLaplacian | ( | const SurfaceField< scalar > & | gamma, |
| const VolField< sphericalTensor > & | vf | ||
| ) |
Definition at line 116 of file gaussLaplacianSchemes.C.
| Foam::tmp< Foam::fvMatrix< Foam::symmTensor > > fvmLaplacian | ( | const SurfaceField< scalar > & | gamma, |
| const VolField< symmTensor > & | vf | ||
| ) |
Definition at line 117 of file gaussLaplacianSchemes.C.
| Foam::tmp< Foam::VolField< Foam::symmTensor > > fvcLaplacian | ( | const SurfaceField< scalar > & | gamma, |
| const VolField< symmTensor > & | vf | ||
| ) |
Definition at line 117 of file gaussLaplacianSchemes.C.
| Foam::tmp< Foam::fvMatrix< Foam::tensor > > fvmLaplacian | ( | const SurfaceField< scalar > & | gamma, |
| const VolField< tensor > & | vf | ||
| ) |
Definition at line 118 of file gaussLaplacianSchemes.C.
| Foam::tmp< Foam::VolField< Foam::tensor > > fvcLaplacian | ( | const SurfaceField< scalar > & | gamma, |
| const VolField< tensor > & | vf | ||
| ) |
Definition at line 118 of file gaussLaplacianSchemes.C.