plane Class Reference

A sampledSurface defined by a plane which 'cuts' the mesh using the cuttingPlane algorithm. The plane is triangulated by default. More...

Inheritance diagram for plane:
Collaboration diagram for plane:

Public Member Functions

 TypeName ("plane")
 Runtime type information. More...
 
 plane (const word &name, const polyMesh &mesh, const dictionary &dict)
 Construct from dictionary. More...
 
virtual ~plane ()
 Destructor. More...
 
virtual bool needsUpdate () const
 Does the surface need an update? More...
 
virtual bool expire ()
 Mark the surface as needing an update. More...
 
virtual bool update ()
 Update the surface as required. More...
 
virtual const pointFieldpoints () const
 Points of surface. More...
 
virtual const faceListfaces () const
 Faces of surface. More...
 
const labelListmeshCells () const
 For every face original cell in mesh. More...
 
virtual tmp< scalarFieldsample (const volScalarField &) const
 Sample field on surface. More...
 
virtual tmp< vectorFieldsample (const volVectorField &) const
 Sample field on surface. More...
 
virtual tmp< sphericalTensorFieldsample (const volSphericalTensorField &) const
 Sample field on surface. More...
 
virtual tmp< symmTensorFieldsample (const volSymmTensorField &) const
 Sample field on surface. More...
 
virtual tmp< tensorFieldsample (const volTensorField &) const
 Sample field on surface. More...
 
virtual tmp< scalarFieldinterpolate (const interpolation< scalar > &) const
 Interpolate field on surface. More...
 
virtual tmp< vectorFieldinterpolate (const interpolation< vector > &) const
 Interpolate field on surface. More...
 
virtual tmp< sphericalTensorFieldinterpolate (const interpolation< sphericalTensor > &) const
 Interpolate field on surface. More...
 
virtual tmp< symmTensorFieldinterpolate (const interpolation< symmTensor > &) const
 Interpolate field on surface. More...
 
virtual tmp< tensorFieldinterpolate (const interpolation< tensor > &) const
 Interpolate field on surface. More...
 
virtual void print (Ostream &) const
 Write. More...
 
template<class Type >
Foam::tmp< Foam::Field< Type > > sampleField (const GeometricField< Type, fvPatchField, volMesh > &vField) const
 
template<class Type >
Foam::tmp< Foam::Field< Type > > interpolateField (const interpolation< Type > &interpolator) const
 
- Public Member Functions inherited from sampledSurface
 TypeName ("sampledSurface")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, sampledSurface, word,(const word &name, const polyMesh &mesh, const dictionary &dict),(name, mesh, dict))
 Declare run-time constructor selection table. More...
 
 sampledSurface (const word &name, const polyMesh &, const bool interpolate=false)
 Construct from name, mesh. More...
 
 sampledSurface (const word &name, const polyMesh &, const dictionary &)
 Construct from dictionary. More...
 
autoPtr< sampledSurfaceclone () const
 Clone. More...
 
virtual ~sampledSurface ()
 Destructor. More...
 
const polyMeshmesh () const
 Access to the underlying mesh. More...
 
const wordname () const
 Name of surface. More...
 
bool interpolate () const
 Interpolation requested for surface. More...
 
virtual const vectorFieldSf () const
 Return face area vectors. More...
 
virtual const scalarFieldmagSf () const
 Return face area magnitudes. More...
 
virtual const vectorFieldCf () const
 Return face centres as vectorField. More...
 
scalar area () const
 The total surface area. More...
 
template<class Type >
Type integrate (const Field< Type > &) const
 Integration of a field across the surface. More...
 
template<class Type >
Type integrate (const tmp< Field< Type >> &) const
 Integration of a field across the surface. More...
 
template<class Type >
Type average (const Field< Type > &) const
 Area-averaged value of a field across the surface. More...
 
template<class Type >
Type average (const tmp< Field< Type >> &) const
 Area-averaged value of a field across the surface. More...
 
tmp< Field< scalar > > project (const Field< scalar > &) const
 Project field onto surface. More...
 
tmp< Field< scalar > > project (const Field< vector > &) const
 Project field onto surface. More...
 
tmp< Field< vector > > project (const Field< sphericalTensor > &) const
 Project field onto surface. More...
 
tmp< Field< vector > > project (const Field< symmTensor > &) const
 Project field onto surface. More...
 
tmp< Field< vector > > project (const Field< tensor > &) const
 Project field onto surface. More...
 
template<class Type >
tmp< GeometricField< Type, fvPatchField, volMesh > > pointAverage (const GeometricField< Type, pointPatchField, pointMesh > &pfld) const
 Interpolate from points to cell centre. More...
 
virtual tmp< scalarFieldsample (const surfaceScalarField &) const
 Surface sample field on surface. More...
 
virtual tmp< symmTensorFieldsample (const surfaceSymmTensorField &) const
 Surface sample field on surface. More...
 
virtual tmp< tensorFieldsample (const surfaceTensorField &) const
 Surface sample field on surface. More...
 
virtual void rename (const word &newName)
 Rename. More...
 
template<class ReturnType , class Type >
Foam::tmp< Foam::Field< ReturnType > > project (const tmp< Field< Type >> &field) const
 
template<class Type >
Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > pointAverage (const GeometricField< Type, pointPatchField, pointMesh > &pfld) const
 
- Public Member Functions inherited from cuttingPlane
 TypeName ("cuttingPlane")
 Runtime type information. More...
 
 cuttingPlane (const word &name, const polyMesh &mesh, const dictionary &dict)
 Construct from dictionary. More...
 
virtual ~cuttingPlane ()
 Destructor. More...
 
const isoSurfacesurface () const
 
template<class Type >
Foam::tmp< Foam::Field< Type > > sampleField (const GeometricField< Type, fvPatchField, volMesh > &vField) const
 
template<class Type >
Foam::tmp< Foam::Field< Type > > interpolateField (const interpolation< Type > &interpolator) const
 

Additional Inherited Members

- Static Public Member Functions inherited from sampledSurface
static autoPtr< sampledSurfaceNew (const word &name, const polyMesh &, const dictionary &)
 Return a reference to the selected surface. More...
 
- Protected Member Functions inherited from sampledSurface
virtual void clearGeom () const
 

Detailed Description

A sampledSurface defined by a plane which 'cuts' the mesh using the cuttingPlane algorithm. The plane is triangulated by default.

Note
Does not actually cut until update() called.
Source files

Definition at line 55 of file sampledPlane.H.

Constructor & Destructor Documentation

◆ plane()

plane ( const word name,
const polyMesh mesh,
const dictionary dict 
)

Construct from dictionary.

Definition at line 44 of file sampledPlane.C.

References polyMesh::cellZones(), Foam::endl(), ZoneMesh< ZoneType, MeshType >::findIndex(), dictionary::found(), Foam::Info, dictionary::readIfPresent(), and dictionary::subDict().

Here is the call graph for this function:

◆ ~plane()

~plane ( )
virtual

Destructor.

Definition at line 81 of file sampledPlane.C.

Member Function Documentation

◆ TypeName()

TypeName ( "plane"  )

Runtime type information.

◆ needsUpdate()

bool needsUpdate ( ) const
virtual

Does the surface need an update?

Reimplemented from cuttingPlane.

Definition at line 87 of file sampledPlane.C.

◆ expire()

bool expire ( )
virtual

Mark the surface as needing an update.

May also free up unneeded data. Return false if surface was already marked as expired.

Reimplemented from cuttingPlane.

Definition at line 93 of file sampledPlane.C.

References sampledSurface::clearGeom().

Here is the call graph for this function:

◆ update()

bool update ( )
virtual

Update the surface as required.

Do nothing (and return false) if no update was needed

Reimplemented from cuttingPlane.

Definition at line 108 of file sampledPlane.C.

References sampledSurface::clearGeom(), UList< T >::empty(), Foam::endl(), mesh, Foam::Pout, and plane::sample().

Here is the call graph for this function:

◆ points()

virtual const pointField& points ( ) const
inlinevirtual

Points of surface.

Reimplemented from cuttingPlane.

Definition at line 123 of file sampledPlane.H.

References cuttingPlane::points().

Here is the call graph for this function:

◆ faces()

virtual const faceList& faces ( ) const
inlinevirtual

Faces of surface.

Reimplemented from cuttingPlane.

Definition at line 129 of file sampledPlane.H.

References cuttingPlane::faces().

Here is the call graph for this function:

◆ meshCells()

const labelList& meshCells ( ) const
inline

For every face original cell in mesh.

Definition at line 135 of file sampledPlane.H.

References cuttingPlane::cutCells(), sampledSurface::interpolate(), plane::print(), and plane::sample().

Here is the call graph for this function:

◆ sample() [1/5]

Foam::tmp< Foam::scalarField > sample ( const volScalarField vField) const
virtual

Sample field on surface.

Reimplemented from cuttingPlane.

Definition at line 140 of file sampledPlane.C.

Referenced by plane::meshCells(), plane::sample(), and plane::update().

Here is the caller graph for this function:

◆ sample() [2/5]

Foam::tmp< Foam::vectorField > sample ( const volVectorField vField) const
virtual

Sample field on surface.

Reimplemented from cuttingPlane.

Definition at line 149 of file sampledPlane.C.

References plane::sample().

Here is the call graph for this function:

◆ sample() [3/5]

Foam::tmp< Foam::sphericalTensorField > sample ( const volSphericalTensorField vField) const
virtual

Sample field on surface.

Reimplemented from cuttingPlane.

Definition at line 158 of file sampledPlane.C.

References plane::sample().

Here is the call graph for this function:

◆ sample() [4/5]

Foam::tmp< Foam::symmTensorField > sample ( const volSymmTensorField vField) const
virtual

Sample field on surface.

Reimplemented from cuttingPlane.

Definition at line 167 of file sampledPlane.C.

References plane::sample().

Here is the call graph for this function:

◆ sample() [5/5]

Foam::tmp< Foam::tensorField > sample ( const volTensorField vField) const
virtual

Sample field on surface.

Reimplemented from cuttingPlane.

Definition at line 176 of file sampledPlane.C.

References sampledSurface::interpolate().

Here is the call graph for this function:

◆ interpolate() [1/5]

Foam::tmp< Foam::scalarField > interpolate ( const interpolation< scalar > &  interpolator) const
virtual

Interpolate field on surface.

Reimplemented from cuttingPlane.

Definition at line 185 of file sampledPlane.C.

References sampledSurface::interpolate().

Here is the call graph for this function:

◆ interpolate() [2/5]

Foam::tmp< Foam::vectorField > interpolate ( const interpolation< vector > &  interpolator) const
virtual

Interpolate field on surface.

Reimplemented from cuttingPlane.

Definition at line 194 of file sampledPlane.C.

References sampledSurface::interpolate().

Here is the call graph for this function:

◆ interpolate() [3/5]

Foam::tmp< Foam::sphericalTensorField > interpolate ( const interpolation< sphericalTensor > &  interpolator) const
virtual

Interpolate field on surface.

Reimplemented from cuttingPlane.

Definition at line 202 of file sampledPlane.C.

References sampledSurface::interpolate().

Here is the call graph for this function:

◆ interpolate() [4/5]

Foam::tmp< Foam::symmTensorField > interpolate ( const interpolation< symmTensor > &  interpolator) const
virtual

Interpolate field on surface.

Reimplemented from cuttingPlane.

Definition at line 211 of file sampledPlane.C.

References sampledSurface::interpolate().

Here is the call graph for this function:

◆ interpolate() [5/5]

Foam::tmp< Foam::tensorField > interpolate ( const interpolation< tensor > &  interpolator) const
virtual

Interpolate field on surface.

Reimplemented from cuttingPlane.

Definition at line 220 of file sampledPlane.C.

◆ print()

void print ( Ostream os) const
virtual

Write.

Reimplemented from cuttingPlane.

Definition at line 228 of file sampledPlane.C.

References functionObject::name(), and points.

Referenced by plane::meshCells().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ sampleField()

Foam::tmp<Foam::Field<Type> > sampleField ( const GeometricField< Type, fvPatchField, volMesh > &  vField) const

Definition at line 33 of file sampledPlaneTemplates.C.

◆ interpolateField()

Foam::tmp<Foam::Field<Type> > interpolateField ( const interpolation< Type > &  interpolator) const

Definition at line 44 of file sampledPlaneTemplates.C.

References f(), forAll, interpolation< Type >::interpolate(), points, and tmp< T >::ref().

Here is the call graph for this function:

The documentation for this class was generated from the following files: