Galeri  Development
 All Classes Files Functions Variables Pages
Public Member Functions | List of all members
Galeri::FiniteElements::TriangleRectangleGrid Class Reference

Creates a grid composed by triangles, the domain is a rectangle. More...

#include <Galeri_TriangleRectangleGrid.h>

Inheritance diagram for Galeri::FiniteElements::TriangleRectangleGrid:
Inheritance graph
[legend]
Collaboration diagram for Galeri::FiniteElements::TriangleRectangleGrid:
Collaboration graph
[legend]

Public Member Functions

 TriangleRectangleGrid (const Epetra_Comm &Comm, const int nx, const int ny, const int mx, const int my, const double lx=1.0, const double ly=1.0)
 Constructor. More...
 
virtual int NumDimensions () const
 Returns the number of dimensions of the grid.
 
virtual int NumVerticesPerElement () const
 Returns the number of vertices contained in each element.
 
virtual int NumFacesPerElement () const
 Returns the number of faces contained in each element.
 
virtual int NumVerticesPerFace () const
 Returns the number of vertices contained in each face.
 
virtual std::string ElementType () const
 Returns GALERI_TRIANGLE.
 
virtual const Epetra_Comm & Comm () const
 Returns a reference to the communicator object.
 
virtual int NumMyElements () const
 Returns the number of finite elements on the calling process.
 
virtual int NumGlobalElements () const
 Returns the global number of finite elements.
 
virtual int NumMyVertices () const
 Returns the number of vertices on the calling process.
 
virtual int NumGlobalVertices () const
 Returns the global number of vertices.
 
virtual int NumMyBoundaryFaces () const
 Returns the number of boundary faces on the calling process.
 
virtual int NumGlobalBoundaryFaces () const
 Returns the global number of boundary faces.
 
virtual void VertexCoord (const int LocalID, double *coord) const
 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
 Returns the coordinates of specified local vertices. More...
 
virtual void ElementVertices (const int LocalID, int *elements) const
 Returns the local vertex IDs of the specified local finite element. More...
 
virtual double ElementMinLength (const int LocalElement) const
 Returns the volume of the specified local finite element.
 
virtual double ElementMaxLength (const int LocalElement) const
 Returns the volume of the specified local finite element.
 
virtual const RefCountPtr
< Epetra_Map > 
RCPVertexMap () const
 
virtual const RefCountPtr
< Epetra_Map > 
RCPElementMap () const
 
virtual const Epetra_Map & VertexMap () const
 Returns a reference to the map representing the vertex distribution.
 
virtual const Epetra_Map & ElementMap () const
 
virtual const Epetra_Map & RowMap () const
 Returns a reference to the map representing the distribution of rows.
 
virtual const Epetra_Import & Importer () const
 
virtual int ElementTag (const int ID) const
 
virtual int VertexTag (const int ID) const
 
virtual double ElementVolume () const
 
virtual void FaceVertices (const int LocalFace, int &tag, int *IDs) const
 Returns the local vertex IDs of vertices contained in the specified boundary face.
 
int FacePatch (const int LocalFace) const
 Returns the patch ID of the specified face. More...
 
int NumMyElementsX () const
 
int NumMyElementsY () const
 
int NumMyVerticesX () const
 
int NumMyVerticesY () const
 
int NumGlobalElementsX () const
 
int NumGlobalElementsY () const
 
int NumGlobalVerticesX () const
 
int NumGlobalVerticesY () const
 
double LengthX () const
 
double LengthY () const
 
double DeltaX () const
 
double DeltaY () const
 
virtual double ElementVolume (const int LocalElement) const
 Returns the volume of the specified local finite element. More...
 
virtual double FaceArea (const int LocalFace) const
 Returns the area of the specified local face. More...
 
virtual double MyVolume () const
 Returns the volume of all local elements.
 
virtual double GlobalVolume () const
 Returns the global volume of the grid.
 
int NumDomainsX () const
 
int NumDomainsY () const
 
void ExportToVertexMap (const Epetra_DistObject &RowObject, Epetra_DistObject &VertexObject) const
 Exports distributed object from RowMap() to VertexMap().
 
void ExportToRowMap (const Epetra_DistObject &VertexObject, Epetra_DistObject &RowObject) const
 Exports distributed object from VertexMap() to RowMap().
 
int NumNeighborsPerElement () const
 Returns the number of neighboring elements.
 
void ElementNeighbors (int, int *) const
 Returns the local IDs of neighboring elements.
 
- Public Member Functions inherited from Galeri::FiniteElements::AbstractGrid
virtual ~AbstractGrid ()
 Destructor.
 

Detailed Description

Creates a grid composed by triangles, the domain is a rectangle.

This class defined, on-the-fly, the triangulation of a 2D rectangular domain. The elements are all triangles. For parallel run, the rectangle is subdivided along the X- and Y-axis, as specified by the user.

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

Constructor & Destructor Documentation

Galeri::FiniteElements::TriangleRectangleGrid::TriangleRectangleGrid ( const Epetra_Comm &  Comm,
const int  nx,
const int  ny,
const int  mx,
const int  my,
const double  lx = 1.0,
const double  ly = 1.0 
)
inline

Constructor.

Parameters
Comm- (In) Communicator object.
nx- (In) number of elements along the X-axis.
ny- (In) number of elements along the Y-axis.
mx- (In) Number of subdomains along the X-axis.
my- (In) Number of subdomains along the Y-axis.
lx- (In) Length of the rectangle along the X-axis.
ly- (In) Length of the rectangle along the Y-axis.
Note
The total number of processors must equal mx * my.

Member Function Documentation

virtual void Galeri::FiniteElements::TriangleRectangleGrid::ElementVertices ( const int  LocalElement,
int *  elements 
) const
inlinevirtual

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.

Implements Galeri::FiniteElements::AbstractGrid.

virtual double Galeri::FiniteElements::TriangleRectangleGrid::ElementVolume ( const int  LocalElement) const
inlinevirtual

Returns the volume of the specified local finite element.

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

Implements Galeri::FiniteElements::AbstractGrid.

virtual double Galeri::FiniteElements::TriangleRectangleGrid::FaceArea ( const int  LocalFace) const
inlinevirtual

Returns the area of the specified local face.

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

Implements Galeri::FiniteElements::AbstractGrid.

int Galeri::FiniteElements::TriangleRectangleGrid::FacePatch ( const int  LocalFace) const
inlinevirtual

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.

Implements Galeri::FiniteElements::AbstractGrid.

virtual void Galeri::FiniteElements::TriangleRectangleGrid::VertexCoord ( const int  LocalVertex,
double *  coord 
) const
inlinevirtual

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.

Implements Galeri::FiniteElements::AbstractGrid.

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

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.

Implements Galeri::FiniteElements::AbstractGrid.


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