Galeri  Development
 All Classes Files Functions Variables Pages
List of all members
Galeri::FiniteElements::AbstractGrid Class Referenceabstract

Abstract interface to access finite element grids. More...

#include <Galeri_AbstractGrid.h>

Inheritance diagram for Galeri::FiniteElements::AbstractGrid:
Inheritance graph
[legend]

Public Member Functions

virtual ~AbstractGrid ()
 Destructor.
 
virtual int NumDimensions () const =0
 Returns the number of dimensions of the grid.
 
virtual int NumVerticesPerElement () const =0
 Returns the number of vertices contained in each element.
 
virtual int NumFacesPerElement () const =0
 Returns the number of faces contained in each element.
 
virtual int NumVerticesPerFace () const =0
 Returns the number of vertices contained in each face.
 
virtual std::string ElementType () const =0
 Returns a string containing the element type. More...
 
virtual int NumNeighborsPerElement () const =0
 Returns the number of neighboring elements.
 
virtual int NumMyElements () const =0
 Returns the number of finite elements on the calling process.
 
virtual int NumGlobalElements () const =0
 Returns the global number of finite elements.
 
virtual int NumMyVertices () const =0
 Returns the number of vertices on the calling process.
 
virtual int NumGlobalVertices () const =0
 Returns the global number of vertices.
 
virtual int NumMyBoundaryFaces () const =0
 Returns the number of boundary faces on the calling process.
 
virtual int NumGlobalBoundaryFaces () const =0
 Returns the global number of boundary faces.
 
virtual double MyVolume () const =0
 Returns the volume of all local elements.
 
virtual double GlobalVolume () const =0
 Returns the global volume of the grid.
 
virtual void VertexCoord (const int LocalVertex, double *coord) const =0
 Returns the coordinates of local vertex LocalVertex in vector coord. More...
 
virtual void VertexCoord (const int Length, const int *IDs, double *x, double *y, double *z) const =0
 Returns the coordinates of specified local vertices. More...
 
virtual void ElementVertices (const int LocalElement, int *elements) const =0
 Returns the local vertex IDs of the specified local finite element. More...
 
virtual void ElementNeighbors (const int LocalElement, int *elements) const =0
 Returns the local IDs of neighboring elements.
 
virtual void FaceVertices (const int LocalFace, int &tag, int *IDs) const =0
 Returns the local vertex IDs of vertices contained in the specified boundary face.
 
virtual int FacePatch (const int LocalFace) const =0
 Returns the patch ID of the specified face. More...
 
virtual double ElementMinLength (const int LocalElement) const =0
 Returns the volume of the specified local finite element.
 
virtual double ElementMaxLength (const int LocalElement) const =0
 Returns the volume of the specified local finite element.
 
virtual double ElementVolume (const int LocalElement) const =0
 Returns the volume of the specified local finite element. More...
 
virtual double FaceArea (const int LocalFace) const =0
 Returns the area of the specified local face. More...
 
virtual const Epetra_Map & VertexMap () const =0
 Returns a reference to the map representing the vertex distribution.
 
virtual const Epetra_Map & RowMap () const =0
 Returns a reference to the map representing the distribution of rows.
 
virtual void ExportToVertexMap (const Epetra_DistObject &RowObject, Epetra_DistObject &VertexObject) const =0
 Exports distributed object from RowMap() to VertexMap().
 
virtual void ExportToRowMap (const Epetra_DistObject &RowObject, Epetra_DistObject &VertexObject) const =0
 Exports distributed object from VertexMap() to RowMap().
 
virtual const Epetra_Comm & Comm () const =0
 Returns a reference to the communicator object.
 

Detailed Description

Abstract interface to access finite element grids.

AbstractGrid is a pure virtual function, that specifies the interface methods that a grid class must implement. Following an approach similar to getrow() for matrices, there is no grid format; instead, the user must implement a set of getelement() and similar methods. Therefore, it is possible to define grids that are built on-the-fly, a feature that is particularly convenient in testing phase.

This format is based on the following assumptions:

Two Epetra_Map's must be created:

We require two maps for the following reason: In parallel runs, vertices corresponding to internal boundary faces can are replicated over processors. This makes the contruction of the finite element matrix easier, since it is possible to work on local quantities only. However, such a distribution is not compatible with AztecOO and ML (and several other Trilinos packages), since these libraries require each row to belong to exactly one processor. Methods ExportToVertexMap() and ExportToRowMap() export objects from one map to the other. Typically, ExportToRowMap() is used in the assembly phase, whese local stiffness matrix and right-hand side are built using VertexMap(), then exported to RowMap(). ExportToVertexMap(), instead, is used after the solution of the linear system is computed, to update the values of the solution on the local vertices, so that norms can be computed, and the solution visualized. These methods are trivial in the serial case.

Author
Marzio Sala, SNL 9214
Date
Last updated on 12-Apr-05.

Member Function Documentation

virtual std::string Galeri::FiniteElements::AbstractGrid::ElementType ( ) const
pure virtual

Returns a string containing the element type.

Returns a string containing the type of element. This string is used in the quadrature class. Currently supported options are:

  • "GALERI_TRIANGLE"
  • "GALERI_QUAD"
  • "GALERI_HEX"
  • "GALERI_TET"

Implemented in Galeri::FiniteElements::FileGrid, Galeri::FiniteElements::TRIANGLEGrid, Galeri::FiniteElements::HexCubeGrid, Galeri::FiniteElements::TetCubeGrid, Galeri::FiniteElements::TriangleRectangleGrid, and Galeri::FiniteElements::QuadRectangleGrid.

Referenced by Galeri::FiniteElements::MEDITInterface::Write().

virtual void Galeri::FiniteElements::AbstractGrid::ElementVertices ( const int  LocalElement,
int *  elements 
) const
pure virtual

Returns the local vertex IDs of the specified local finite element.

Parameters
LocalElement- (In) ID of the required local element.
elements- (Out) array of length NumElementVertices(), in output will contain the local ID of the vertices of the specified element.

Implemented in Galeri::FiniteElements::FileGrid, Galeri::FiniteElements::TRIANGLEGrid, Galeri::FiniteElements::HexCubeGrid, Galeri::FiniteElements::TetCubeGrid, Galeri::FiniteElements::TriangleRectangleGrid, and Galeri::FiniteElements::QuadRectangleGrid.

Referenced by Galeri::FiniteElements::LinearProblem::Compute(), Galeri::FiniteElements::LinearProblem::ComputeNorms(), and Galeri::FiniteElements::MEDITInterface::Write().

virtual double Galeri::FiniteElements::AbstractGrid::ElementVolume ( const int  LocalElement) const
pure virtual

Returns the volume of the specified local finite element.

Returns the area (in 2D) or the volume (in 3D) of the specified local element

Implemented in Galeri::FiniteElements::FileGrid, Galeri::FiniteElements::HexCubeGrid, Galeri::FiniteElements::TetCubeGrid, Galeri::FiniteElements::QuadRectangleGrid, Galeri::FiniteElements::TriangleRectangleGrid, and Galeri::FiniteElements::TRIANGLEGrid.

virtual double Galeri::FiniteElements::AbstractGrid::FaceArea ( const int  LocalFace) const
pure virtual

Returns the area of the specified local face.

Returns the length (in 2D) or the area (in 3D) of the specified boundary face

Implemented in Galeri::FiniteElements::FileGrid, Galeri::FiniteElements::HexCubeGrid, Galeri::FiniteElements::TetCubeGrid, Galeri::FiniteElements::QuadRectangleGrid, Galeri::FiniteElements::TriangleRectangleGrid, and Galeri::FiniteElements::TRIANGLEGrid.

virtual int Galeri::FiniteElements::AbstractGrid::FacePatch ( const int  LocalFace) const
pure virtual

Returns the patch ID of the specified face.

Returns an integer ID that identifies the given boundary face as belonging to a given part of the domain. It can be used by the user to specify the value and the type of the boundary condition.

Implemented in Galeri::FiniteElements::FileGrid, Galeri::FiniteElements::QuadRectangleGrid, Galeri::FiniteElements::HexCubeGrid, Galeri::FiniteElements::TetCubeGrid, Galeri::FiniteElements::TriangleRectangleGrid, and Galeri::FiniteElements::TRIANGLEGrid.

Referenced by Galeri::FiniteElements::LinearProblem::Compute().

virtual void Galeri::FiniteElements::AbstractGrid::VertexCoord ( const int  LocalVertex,
double *  coord 
) const
pure virtual

Returns the coordinates of local vertex LocalVertex in vector coord.

Parameters
LocalVertex- (In) Local ID of the vertex for whic coordinates are required. Must be contained in the interval [0, NumMyVertices())
coord- (Out) double array of size 3. In output, contains the x-, y- and z-coordinate of the specified vertex.
Note
Parameter coord must be allocated of size 3 for both 2D and 3D problems.

Implemented in Galeri::FiniteElements::FileGrid, Galeri::FiniteElements::TRIANGLEGrid, Galeri::FiniteElements::HexCubeGrid, Galeri::FiniteElements::TetCubeGrid, Galeri::FiniteElements::TriangleRectangleGrid, and Galeri::FiniteElements::QuadRectangleGrid.

Referenced by Galeri::FiniteElements::LinearProblem::Compute(), Galeri::FiniteElements::LinearProblem::ComputeNorms(), and Galeri::FiniteElements::MEDITInterface::Write().

virtual void Galeri::FiniteElements::AbstractGrid::VertexCoord ( const int  Length,
const int *  IDs,
double *  x,
double *  y,
double *  z 
) const
pure virtual

Returns the coordinates of specified local vertices.

Parameters
Length- (In) Length of array IDs.
IDs- (In) Contains the list of vertices of which coordinates are required.
x- (Out) double array of size Length. In output, contains the x-coordinates of the specified vertices.
y- (Out) double array of size Length. In output, contains the y-coordinates of the specified vertices.
z- (Out) double array of size Length. In output, contains the z-coordinates of the specified vertices.
Note
The z array must be allocated for both 2D and 3D problems.

Implemented in Galeri::FiniteElements::FileGrid, Galeri::FiniteElements::TRIANGLEGrid, Galeri::FiniteElements::HexCubeGrid, Galeri::FiniteElements::TetCubeGrid, Galeri::FiniteElements::TriangleRectangleGrid, and Galeri::FiniteElements::QuadRectangleGrid.


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