MFEM  v3.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Classes | Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Protected Member Functions | Static Protected Member Functions | Protected Attributes | Static Protected Attributes | Friends | List of all members
mfem::Mesh Class Reference

#include <mesh.hpp>

Inheritance diagram for mfem::Mesh:
[legend]
Collaboration diagram for mfem::Mesh:
[legend]

Classes

class  FaceInfo
 

Public Types

enum  { NORMAL, TWO_LEVEL_COARSE, TWO_LEVEL_FINE }
 

Public Member Functions

 Mesh ()
 
 Mesh (int _Dim, int NVert, int NElem, int NBdrElem=0, int _spaceDim=-1)
 
ElementNewElement (int geom)
 
void AddVertex (const double *)
 
void AddTri (const int *vi, int attr=1)
 
void AddTriangle (const int *vi, int attr=1)
 
void AddQuad (const int *vi, int attr=1)
 
void AddTet (const int *vi, int attr=1)
 
void AddHex (const int *vi, int attr=1)
 
void AddHexAsTets (const int *vi, int attr=1)
 
void AddElement (Element *elem)
 
void AddBdrSegment (const int *vi, int attr=1)
 
void AddBdrTriangle (const int *vi, int attr=1)
 
void AddBdrQuad (const int *vi, int attr=1)
 
void AddBdrQuadAsTriangles (const int *vi, int attr=1)
 
void GenerateBoundaryElements ()
 
void FinalizeTriMesh (int generate_edges=0, int refine=0, bool fix_orientation=true)
 
void FinalizeQuadMesh (int generate_edges=0, int refine=0, bool fix_orientation=true)
 
void FinalizeTetMesh (int generate_edges=0, int refine=0, bool fix_orientation=true)
 
void FinalizeHexMesh (int generate_edges=0, int refine=0, bool fix_orientation=true)
 
void SetAttributes ()
 
 Mesh (int nx, int ny, int nz, Element::Type type, int generate_edges=0, double sx=1.0, double sy=1.0, double sz=1.0)
 
 Mesh (int nx, int ny, Element::Type type, int generate_edges=0, double sx=1.0, double sy=1.0)
 
 Mesh (int n, double sx=1.0)
 
 Mesh (std::istream &input, int generate_edges=0, int refine=1, bool fix_orientation=true)
 
 Mesh (Mesh *mesh_array[], int num_pieces)
 Create a disjoint mesh from the given mesh array. More...
 
void Load (std::istream &input, int generate_edges=0, int refine=1, bool fix_orientation=true)
 
int MeshGenerator ()
 Truegrid or NetGen? More...
 
int GetNV () const
 
int GetNE () const
 Returns number of elements. More...
 
int GetNBE () const
 Returns number of boundary elements. More...
 
int GetNEdges () const
 Return the number of edges. More...
 
int GetNFaces () const
 Return the number of faces in a 3D mesh. More...
 
int GetNumFaces () const
 Return the number of faces (3D), edges (2D) or vertices (1D). More...
 
int EulerNumber () const
 Equals 1 + num_holes - num_loops. More...
 
int EulerNumber2D () const
 Equals 1 - num_holes. More...
 
int Dimension () const
 
int SpaceDimension () const
 
const double * GetVertex (int i) const
 Return pointer to vertex i's coordinates. More...
 
double * GetVertex (int i)
 
const ElementGetElement (int i) const
 
ElementGetElement (int i)
 
const ElementGetBdrElement (int i) const
 
ElementGetBdrElement (int i)
 
const ElementGetFace (int i) const
 
int GetFaceBaseGeometry (int i) const
 
int GetElementBaseGeometry (int i) const
 
int GetBdrElementBaseGeometry (int i) const
 
void GetElementVertices (int i, Array< int > &dofs) const
 Returns the indices of the dofs of element i. More...
 
void GetBdrElementVertices (int i, Array< int > &dofs) const
 Returns the indices of the dofs of boundary element i. More...
 
void GetElementEdges (int i, Array< int > &, Array< int > &) const
 Return the indices and the orientations of all edges of element i. More...
 
void GetBdrElementEdges (int i, Array< int > &, Array< int > &) const
 Return the indices and the orientations of all edges of bdr element i. More...
 
void GetFaceEdges (int i, Array< int > &, Array< int > &) const
 
void GetFaceVertices (int i, Array< int > &vert) const
 Returns the indices of the vertices of face i. More...
 
void GetEdgeVertices (int i, Array< int > &vert) const
 Returns the indices of the vertices of edge i. More...
 
TableGetFaceEdgeTable () const
 Returns the face-to-edge Table (3D) More...
 
TableGetEdgeVertexTable () const
 Returns the edge-to-vertex Table (3D) More...
 
void GetElementFaces (int i, Array< int > &, Array< int > &) const
 Return the indices and the orientations of all faces of element i. More...
 
void GetBdrElementFace (int i, int *, int *) const
 Return the index and the orientation of the face of bdr element i. (3D) More...
 
int GetBdrElementEdgeIndex (int i) const
 
int GetElementType (int i) const
 Returns the type of element i. More...
 
int GetBdrElementType (int i) const
 Returns the type of boundary element i. More...
 
void GetPointMatrix (int i, DenseMatrix &pointmat) const
 
void GetBdrPointMatrix (int i, DenseMatrix &pointmat) const
 
void GetElementTransformation (int i, IsoparametricTransformation *ElTr)
 
ElementTransformationGetElementTransformation (int i)
 Returns the transformation defining the i-th element. More...
 
void GetElementTransformation (int i, const Vector &nodes, IsoparametricTransformation *ElTr)
 
ElementTransformationGetBdrElementTransformation (int i)
 Returns the transformation defining the i-th boundary element. More...
 
void GetBdrElementTransformation (int i, IsoparametricTransformation *ElTr)
 
void GetFaceTransformation (int i, IsoparametricTransformation *FTr)
 
ElementTransformationGetFaceTransformation (int FaceNo)
 Returns the transformation defining the given face element. More...
 
void GetEdgeTransformation (int i, IsoparametricTransformation *EdTr)
 
ElementTransformationGetEdgeTransformation (int EdgeNo)
 Returns the transformation defining the given face element. More...
 
FaceElementTransformationsGetFaceElementTransformations (int FaceNo, int mask=31)
 
FaceElementTransformationsGetInteriorFaceTransformations (int FaceNo)
 
FaceElementTransformationsGetBdrFaceTransformations (int BdrElemNo)
 
bool FaceIsInterior (int FaceNo) const
 Return true if the given face is interior. More...
 
void GetFaceElements (int Face, int *Elem1, int *Elem2)
 
void GetFaceInfos (int Face, int *Inf1, int *Inf2)
 
void CheckElementOrientation (bool fix_it=true)
 Check the orientation of the elements. More...
 
void CheckBdrElementOrientation (bool fix_it=true)
 Check the orientation of the boundary elements. More...
 
int GetAttribute (int i) const
 Return the attribute of element i. More...
 
int GetBdrAttribute (int i) const
 Return the attribute of boundary element i. More...
 
const TableElementToElementTable ()
 
const TableElementToFaceTable () const
 
const TableElementToEdgeTable () const
 
TableGetVertexToElementTable ()
 The returned Table must be destroyed by the caller. More...
 
TableGetFaceToElementTable () const
 
virtual void ReorientTetMesh ()
 
int * CartesianPartitioning (int nxyz[])
 
int * GeneratePartitioning (int nparts, int part_method=1)
 
void CheckPartitioning (int *partitioning)
 
void CheckDisplacements (const Vector &displacements, double &tmax)
 
void MoveVertices (const Vector &displacements)
 
void GetVertices (Vector &vert_coord) const
 
void SetVertices (const Vector &vert_coord)
 
void GetNode (int i, double *coord)
 
void SetNode (int i, const double *coord)
 
void MoveNodes (const Vector &displacements)
 
void GetNodes (Vector &node_coord) const
 
void SetNodes (const Vector &node_coord)
 
GridFunctionGetNodes ()
 Return a pointer to the internal node GridFunction (may be NULL). More...
 
void NewNodes (GridFunction &nodes, bool make_owner=false)
 Replace the internal node GridFunction with the given GridFunction. More...
 
void SwapNodes (GridFunction *&nodes, int &own_nodes_)
 
void GetNodes (GridFunction &nodes) const
 Return the mesh nodes/vertices projected on the given GridFunction. More...
 
void SetNodalFESpace (FiniteElementSpace *nfes)
 
void SetNodalGridFunction (GridFunction *nodes, bool make_owner=false)
 
const FiniteElementSpaceGetNodalFESpace ()
 
void UniformRefinement ()
 
void GeneralRefinement (Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
 
void GeneralRefinement (Array< int > &el_to_refine, int nonconforming=-1, int nc_limit=0)
 
void KnotInsert (Array< KnotVector * > &kv)
 
void DegreeElevate (int t)
 
void UseTwoLevelState (int use)
 
void SetState (int s)
 Change the mesh state to NORMAL, TWO_LEVEL_COARSE, TWO_LEVEL_FINE. More...
 
int GetState () const
 
int GetNumFineElems (int i)
 
int GetRefinementType (int i)
 For a given coarse element i returns its refinement type. More...
 
int GetFineElem (int i, int j)
 
ElementTransformationGetFineElemTrans (int i, int j)
 
virtual void PrintXG (std::ostream &out=std::cout) const
 Print the mesh to the given stream using Netgen/Truegrid format. More...
 
virtual void Print (std::ostream &out=std::cout) const
 Print the mesh to the given stream using the default MFEM mesh format. More...
 
void PrintVTK (std::ostream &out)
 Print the mesh in VTK format (linear and quadratic meshes only). More...
 
void PrintVTK (std::ostream &out, int ref, int field_data=0)
 
void GetElementColoring (Array< int > &colors, int el0=0)
 
void PrintWithPartitioning (int *partitioning, std::ostream &out, int elem_attr=0) const
 
void PrintElementsWithPartitioning (int *partitioning, std::ostream &out, int interior_faces=0)
 
void PrintSurfaces (const Table &Aface_face, std::ostream &out) const
 Print set of disjoint surfaces: More...
 
void ScaleSubdomains (double sf)
 
void ScaleElements (double sf)
 
void Transform (void(*f)(const Vector &, Vector &))
 
double GetElementSize (int i, int type=0)
 
double GetElementSize (int i, const Vector &dir)
 
double GetElementVolume (int i)
 
void PrintCharacteristics (Vector *Vh=NULL, Vector *Vk=NULL)
 
void MesquiteSmooth (const int mesquite_option=0)
 
virtual ~Mesh ()
 Destroys mesh. More...
 

Static Public Member Functions

static FiniteElementGetTransformationFEforElementType (int)
 

Public Attributes

Array< int > attributes
 
Array< int > bdr_attributes
 
NURBSExtensionNURBSext
 
NCMeshncmesh
 

Protected Member Functions

void Init ()
 
void InitTables ()
 
void DeleteTables ()
 
void DeleteCoarseTables ()
 
ElementReadElementWithoutAttr (std::istream &)
 
ElementReadElement (std::istream &)
 
void SetMeshGen ()
 
double GetLength (int i, int j) const
 Return the length of the segment from node i to node j. More...
 
void GetElementJacobian (int i, DenseMatrix &J)
 
void MarkForRefinement ()
 
void MarkTriMeshForRefinement ()
 
void GetEdgeOrdering (DSTable &v_to_v, Array< int > &order)
 
void MarkTetMeshForRefinement ()
 
void PrepareNodeReorder (DSTable **old_v_to_v, Table **old_elem_vert)
 
void DoNodeReorder (DSTable *old_v_to_v, Table *old_elem_vert)
 
STable3DGetFacesTable ()
 
STable3DGetElementToFaceTable (int ret_ftbl=0)
 
void RedRefinement (int i, const DSTable &v_to_v, int *edge1, int *edge2, int *middle)
 
void GreenRefinement (int i, const DSTable &v_to_v, int *edge1, int *edge2, int *middle)
 
void Bisection (int i, const DSTable &, int *, int *, int *)
 
void Bisection (int i, const DSTable &, int *)
 
void UniformRefinement (int i, const DSTable &, int *, int *, int *)
 
void AverageVertices (int *indexes, int n, int result)
 
void UpdateNodes ()
 Update the nodes of a curved mesh after refinement. More...
 
virtual void QuadUniformRefinement ()
 Refine quadrilateral mesh. More...
 
virtual void HexUniformRefinement ()
 Refine hexahedral mesh. More...
 
virtual void NURBSUniformRefinement ()
 Refine NURBS mesh. More...
 
virtual void LocalRefinement (const Array< int > &marked_el, int type=3)
 This function is not public anymore. Use GeneralRefinement instead. More...
 
void NonconformingRefinement (const Array< Refinement > &refinements, int nc_limit=0)
 This function is not public anymore. Use GeneralRefinement instead. More...
 
void LoadPatchTopo (std::istream &input, Array< int > &edge_to_knot)
 Read NURBS patch/macro-element mesh. More...
 
void UpdateNURBS ()
 
void PrintTopo (std::ostream &out, const Array< int > &e_to_k) const
 
void BisectTriTrans (DenseMatrix &pointmat, Triangle *tri, int child)
 
void BisectTetTrans (DenseMatrix &pointmat, Tetrahedron *tet, int child)
 
int GetFineElemPath (int i, int j)
 
int GetBisectionHierarchy (Element *E)
 
void GetLocalPtToSegTransformation (IsoparametricTransformation &, int)
 Used in GetFaceElementTransformations (...) More...
 
void GetLocalSegToTriTransformation (IsoparametricTransformation &loc, int i)
 
void GetLocalSegToQuadTransformation (IsoparametricTransformation &loc, int i)
 
void GetLocalTriToTetTransformation (IsoparametricTransformation &loc, int i)
 Used in GetFaceElementTransformations (...) More...
 
void GetLocalQuadToHexTransformation (IsoparametricTransformation &loc, int i)
 Used in GetFaceElementTransformations (...) More...
 
void GetVertexToVertexTable (DSTable &) const
 
int GetElementToEdgeTable (Table &, Array< int > &)
 
void AddPointFaceElement (int lf, int gf, int el)
 Used in GenerateFaces() More...
 
void AddSegmentFaceElement (int lf, int gf, int el, int v0, int v1)
 
void AddTriangleFaceElement (int lf, int gf, int el, int v0, int v1, int v2)
 
void AddQuadFaceElement (int lf, int gf, int el, int v0, int v1, int v2, int v3)
 
bool FaceIsTrueInterior (int FaceNo) const
 
void FreeElement (Element *E)
 
void GenerateFaces ()
 
void InitMesh (int _Dim, int _spaceDim, int NVert, int NElem, int NBdrElem)
 Begin construction of a mesh. More...
 
void Make3D (int nx, int ny, int nz, Element::Type type, int generate_edges, double sx, double sy, double sz)
 
void Make2D (int nx, int ny, Element::Type type, int generate_edges, double sx, double sy)
 
void Make1D (int n, double sx=1.0)
 Creates a 1D mesh for the interval [0,sx] divided into n equal intervals. More...
 
 Mesh (NCMesh &ncmesh)
 Create from a nonconforming mesh. More...
 
void Swap (Mesh &other, bool non_geometry=false)
 

Static Protected Member Functions

static void PrintElementWithoutAttr (const Element *, std::ostream &)
 
static void PrintElement (const Element *, std::ostream &)
 
static int GetTriOrientation (const int *base, const int *test)
 Returns the orientation of "test" relative to "base". More...
 
static int GetQuadOrientation (const int *base, const int *test)
 Returns the orientation of "test" relative to "base". More...
 
static void GetElementArrayEdgeTable (const Array< Element * > &elem_array, const DSTable &v_to_v, Table &el_to_edge)
 
static void ShiftL2R (int &, int &, int &)
 
static void Rotate3 (int &, int &, int &)
 

Protected Attributes

int Dim
 
int spaceDim
 
int NumOfVertices
 
int NumOfElements
 
int NumOfBdrElements
 
int NumOfEdges
 
int NumOfFaces
 
int State
 
int WantTwoLevelState
 
int meshgen
 
int c_NumOfVertices
 
int c_NumOfElements
 
int c_NumOfBdrElements
 
int f_NumOfVertices
 
int f_NumOfElements
 
int f_NumOfBdrElements
 
int c_NumOfEdges
 
int c_NumOfFaces
 
int f_NumOfEdges
 
int f_NumOfFaces
 
Array< Element * > elements
 
Array< Vertexvertices
 
Array< Element * > boundary
 
Array< Element * > faces
 
Array< FaceInfofaces_info
 
Tableel_to_edge
 
Tableel_to_face
 
Tableel_to_el
 
Array< int > be_to_edge
 
Tablebel_to_edge
 
Array< int > be_to_face
 
Tableface_edge
 
Tableedge_vertex
 
Tablec_el_to_edge
 
Tablef_el_to_edge
 
Tablec_bel_to_edge
 
Tablef_bel_to_edge
 
Array< int > fc_be_to_edge
 
Tablec_el_to_face
 
Tablef_el_to_face
 
Array< FaceInfofc_faces_info
 
IsoparametricTransformation Transformation
 
IsoparametricTransformation Transformation2
 
IsoparametricTransformation FaceTransformation
 
IsoparametricTransformation EdgeTransformation
 
FaceElementTransformations FaceElemTr
 
GridFunctionNodes
 
int own_nodes
 
Meshnc_coarse_level
 
MemAlloc< Tetrahedron, 1024 > TetMemory
 
MemAlloc< BisectedElement, 1024 > BEMemory
 

Static Protected Attributes

static const int tet_faces [4][3]
 
static const int hex_faces [6][4]
 
static const int tri_orientations [6][3]
 
static const int quad_orientations [8][4]
 

Friends

class ParMesh
 
class NURBSExtension
 
class Tetrahedron
 

Detailed Description

Definition at line 40 of file mesh.hpp.

Member Enumeration Documentation

anonymous enum
Enumerator
NORMAL 
TWO_LEVEL_COARSE 
TWO_LEVEL_FINE 

Definition at line 302 of file mesh.hpp.

Constructor & Destructor Documentation

mfem::Mesh::Mesh ( NCMesh ncmesh)
protected

Create from a nonconforming mesh.

Definition at line 5798 of file mesh.cpp.

mfem::Mesh::Mesh ( )
inline

Definition at line 310 of file mesh.hpp.

mfem::Mesh::Mesh ( int  _Dim,
int  NVert,
int  NElem,
int  NBdrElem = 0,
int  _spaceDim = -1 
)
inline

Definition at line 312 of file mesh.hpp.

mfem::Mesh::Mesh ( int  nx,
int  ny,
int  nz,
Element::Type  type,
int  generate_edges = 0,
double  sx = 1.0,
double  sy = 1.0,
double  sz = 1.0 
)
inline

Creates mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into nx*ny*nz hexahedrals if type=HEXAHEDRON or into 6*nx*ny*nz tetrahedrons if type=TETRAHEDRON. If generate_edges = 0 (default) edges are not generated, if 1 edges are generated.

Definition at line 350 of file mesh.hpp.

mfem::Mesh::Mesh ( int  nx,
int  ny,
Element::Type  type,
int  generate_edges = 0,
double  sx = 1.0,
double  sy = 1.0 
)
inline

Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny quadrilaterals if type = QUADRILATERAL or into 2*nx*ny triangles if type = TRIANGLE. If generate_edges = 0 (default) edges are not generated, if 1 edges are generated.

Definition at line 360 of file mesh.hpp.

mfem::Mesh::Mesh ( int  n,
double  sx = 1.0 
)
inlineexplicit

Creates 1D mesh , divided into n equal intervals.

Definition at line 367 of file mesh.hpp.

mfem::Mesh::Mesh ( std::istream &  input,
int  generate_edges = 0,
int  refine = 1,
bool  fix_orientation = true 
)

Creates mesh by reading data stream in MFEM, netgen, or VTK format. If generate_edges = 0 (default) edges are not generated, if 1 edges are generated.

Definition at line 1663 of file mesh.cpp.

mfem::Mesh::Mesh ( Mesh mesh_array[],
int  num_pieces 
)

Create a disjoint mesh from the given mesh array.

Definition at line 2701 of file mesh.cpp.

mfem::Mesh::~Mesh ( )
virtual

Destroys mesh.

Definition at line 8221 of file mesh.cpp.

Member Function Documentation

void mfem::Mesh::AddBdrQuad ( const int *  vi,
int  attr = 1 
)

Definition at line 827 of file mesh.cpp.

void mfem::Mesh::AddBdrQuadAsTriangles ( const int *  vi,
int  attr = 1 
)

Definition at line 832 of file mesh.cpp.

void mfem::Mesh::AddBdrSegment ( const int *  vi,
int  attr = 1 
)

Definition at line 817 of file mesh.cpp.

void mfem::Mesh::AddBdrTriangle ( const int *  vi,
int  attr = 1 
)

Definition at line 822 of file mesh.cpp.

void mfem::Mesh::AddElement ( Element elem)
inline

Definition at line 329 of file mesh.hpp.

void mfem::Mesh::AddHex ( const int *  vi,
int  attr = 1 
)

Definition at line 797 of file mesh.cpp.

void mfem::Mesh::AddHexAsTets ( const int *  vi,
int  attr = 1 
)

Definition at line 802 of file mesh.cpp.

void mfem::Mesh::AddPointFaceElement ( int  lf,
int  gf,
int  el 
)
protected

Used in GenerateFaces()

Definition at line 3862 of file mesh.cpp.

void mfem::Mesh::AddQuad ( const int *  vi,
int  attr = 1 
)

Definition at line 779 of file mesh.cpp.

void mfem::Mesh::AddQuadFaceElement ( int  lf,
int  gf,
int  el,
int  v0,
int  v1,
int  v2,
int  v3 
)
protected

Definition at line 3930 of file mesh.cpp.

void mfem::Mesh::AddSegmentFaceElement ( int  lf,
int  gf,
int  el,
int  v0,
int  v1 
)
protected

Definition at line 3884 of file mesh.cpp.

void mfem::Mesh::AddTet ( const int *  vi,
int  attr = 1 
)

Definition at line 784 of file mesh.cpp.

void mfem::Mesh::AddTri ( const int *  vi,
int  attr = 1 
)

Definition at line 769 of file mesh.cpp.

void mfem::Mesh::AddTriangle ( const int *  vi,
int  attr = 1 
)

Definition at line 774 of file mesh.cpp.

void mfem::Mesh::AddTriangleFaceElement ( int  lf,
int  gf,
int  el,
int  v0,
int  v1,
int  v2 
)
protected

Definition at line 3906 of file mesh.cpp.

void mfem::Mesh::AddVertex ( const double *  x)

Definition at line 760 of file mesh.cpp.

void mfem::Mesh::AverageVertices ( int *  indexes,
int  n,
int  result 
)
protected

Averages the vertices with given indexes and save the result in vertices[result].

Definition at line 4988 of file mesh.cpp.

void mfem::Mesh::Bisection ( int  i,
const DSTable v_to_v,
int *  edge1,
int *  edge2,
int *  middle 
)
protected

Bisection. Element with index i is bisected.

Definition at line 5942 of file mesh.cpp.

void mfem::Mesh::Bisection ( int  i,
const DSTable v_to_v,
int *  middle 
)
protected

Bisection. Boundary element with index i is bisected.

Definition at line 6159 of file mesh.cpp.

void mfem::Mesh::BisectTetTrans ( DenseMatrix pointmat,
Tetrahedron tet,
int  child 
)
protected

Definition at line 6722 of file mesh.cpp.

void mfem::Mesh::BisectTriTrans ( DenseMatrix pointmat,
Triangle tri,
int  child 
)
protected

Definition at line 6698 of file mesh.cpp.

int * mfem::Mesh::CartesianPartitioning ( int  nxyz[])

Definition at line 4209 of file mesh.cpp.

void mfem::Mesh::CheckBdrElementOrientation ( bool  fix_it = true)

Check the orientation of the boundary elements.

Definition at line 3266 of file mesh.cpp.

void mfem::Mesh::CheckDisplacements ( const Vector displacements,
double &  tmax 
)

Definition at line 4819 of file mesh.cpp.

void mfem::Mesh::CheckElementOrientation ( bool  fix_it = true)

Check the orientation of the elements.

Definition at line 3098 of file mesh.cpp.

void mfem::Mesh::CheckPartitioning ( int *  partitioning)

Definition at line 4533 of file mesh.cpp.

void mfem::Mesh::DegreeElevate ( int  t)

Definition at line 2878 of file mesh.cpp.

void mfem::Mesh::DeleteCoarseTables ( )
protected

Delete the 'el_to_el', 'face_edge' and 'edge_vertex' tables. Usefull in refinement methods to destroy the coarse tables.

Definition at line 671 of file mesh.cpp.

void mfem::Mesh::DeleteTables ( )
protected

Definition at line 645 of file mesh.cpp.

int mfem::Mesh::Dimension ( ) const
inline

Definition at line 417 of file mesh.hpp.

void mfem::Mesh::DoNodeReorder ( DSTable old_v_to_v,
Table old_elem_vert 
)
protected

Definition at line 1046 of file mesh.cpp.

const Table & mfem::Mesh::ElementToEdgeTable ( ) const

Definition at line 3855 of file mesh.cpp.

const Table & mfem::Mesh::ElementToElementTable ( )

Definition at line 3783 of file mesh.cpp.

const Table & mfem::Mesh::ElementToFaceTable ( ) const

Definition at line 3848 of file mesh.cpp.

int mfem::Mesh::EulerNumber ( ) const
inline

Equals 1 + num_holes - num_loops.

Definition at line 411 of file mesh.hpp.

int mfem::Mesh::EulerNumber2D ( ) const
inline

Equals 1 - num_holes.

Definition at line 414 of file mesh.hpp.

bool mfem::Mesh::FaceIsInterior ( int  FaceNo) const
inline

Return true if the given face is interior.

Definition at line 570 of file mesh.hpp.

bool mfem::Mesh::FaceIsTrueInterior ( int  FaceNo) const
inlineprotected

For a serial Mesh, return true if the face is interior. For a parallel ParMesh return true if the face is interior or shared. In parallel, this method only works if the face neighbor data is exchanged.

Definition at line 259 of file mesh.hpp.

void mfem::Mesh::FinalizeHexMesh ( int  generate_edges = 0,
int  refine = 0,
bool  fix_orientation = true 
)

Definition at line 1292 of file mesh.cpp.

void mfem::Mesh::FinalizeQuadMesh ( int  generate_edges = 0,
int  refine = 0,
bool  fix_orientation = true 
)

Definition at line 916 of file mesh.cpp.

void mfem::Mesh::FinalizeTetMesh ( int  generate_edges = 0,
int  refine = 0,
bool  fix_orientation = true 
)

Definition at line 1254 of file mesh.cpp.

void mfem::Mesh::FinalizeTriMesh ( int  generate_edges = 0,
int  refine = 0,
bool  fix_orientation = true 
)

Definition at line 892 of file mesh.cpp.

void mfem::Mesh::FreeElement ( Element E)
protected

Definition at line 8206 of file mesh.cpp.

void mfem::Mesh::GeneralRefinement ( Array< Refinement > &  refinements,
int  nonconforming = -1,
int  nc_limit = 0 
)

Refine selected mesh elements. Refinement type can be specified for each element. The function can do conforming refinement of triangles and tetrahedra and non-conforming refinement (i.e., with hanging-nodes) of triangles, quadrilaterals and hexahedrons. If 'nonconforming' = -1, suitable refinement method is selected automatically (namely, conforming refinement for triangles). Use noncoforming = 0/1 to force the method. For nonconforming refinements, nc_limit optionally specifies the maximum level of hanging nodes (unlimited by default).

Definition at line 5894 of file mesh.cpp.

void mfem::Mesh::GeneralRefinement ( Array< int > &  el_to_refine,
int  nonconforming = -1,
int  nc_limit = 0 
)

Simplified version of GeneralRefinement taking a simple list of elements to refine, without refinement types.

Definition at line 5932 of file mesh.cpp.

void mfem::Mesh::GenerateBoundaryElements ( )

Definition at line 845 of file mesh.cpp.

void mfem::Mesh::GenerateFaces ( )
protected

Definition at line 3954 of file mesh.cpp.

int * mfem::Mesh::GeneratePartitioning ( int  nparts,
int  part_method = 1 
)

Definition at line 4252 of file mesh.cpp.

int mfem::Mesh::GetAttribute ( int  i) const
inline

Return the attribute of element i.

Definition at line 583 of file mesh.hpp.

int mfem::Mesh::GetBdrAttribute ( int  i) const
inline

Return the attribute of boundary element i.

Definition at line 586 of file mesh.hpp.

const Element* mfem::Mesh::GetBdrElement ( int  i) const
inline

Definition at line 428 of file mesh.hpp.

Element* mfem::Mesh::GetBdrElement ( int  i)
inline

Definition at line 430 of file mesh.hpp.

int mfem::Mesh::GetBdrElementBaseGeometry ( int  i) const
inline

Definition at line 439 of file mesh.hpp.

int mfem::Mesh::GetBdrElementEdgeIndex ( int  i) const

Return the edge index of boundary element i. (2D) return the face index of boundary element i. (3D)

Definition at line 3625 of file mesh.cpp.

void mfem::Mesh::GetBdrElementEdges ( int  i,
Array< int > &  edges,
Array< int > &  cor 
) const

Return the indices and the orientations of all edges of bdr element i.

Definition at line 3364 of file mesh.cpp.

void mfem::Mesh::GetBdrElementFace ( int  i,
int *  f,
int *  o 
) const

Return the index and the orientation of the face of bdr element i. (3D)

Definition at line 3560 of file mesh.cpp.

ElementTransformation * mfem::Mesh::GetBdrElementTransformation ( int  i)

Returns the transformation defining the i-th boundary element.

Definition at line 218 of file mesh.cpp.

void mfem::Mesh::GetBdrElementTransformation ( int  i,
IsoparametricTransformation ElTr 
)

Definition at line 224 of file mesh.cpp.

int mfem::Mesh::GetBdrElementType ( int  i) const

Returns the type of boundary element i.

Definition at line 3647 of file mesh.cpp.

void mfem::Mesh::GetBdrElementVertices ( int  i,
Array< int > &  dofs 
) const
inline

Returns the indices of the dofs of boundary element i.

Definition at line 447 of file mesh.hpp.

FaceElementTransformations * mfem::Mesh::GetBdrFaceTransformations ( int  BdrElemNo)

Definition at line 598 of file mesh.cpp.

void mfem::Mesh::GetBdrPointMatrix ( int  i,
DenseMatrix pointmat 
) const

Definition at line 3674 of file mesh.cpp.

int mfem::Mesh::GetBisectionHierarchy ( Element E)
protected

Definition at line 6502 of file mesh.cpp.

void mfem::Mesh::GetEdgeOrdering ( DSTable v_to_v,
Array< int > &  order 
)
protected

Definition at line 963 of file mesh.cpp.

void mfem::Mesh::GetEdgeTransformation ( int  i,
IsoparametricTransformation EdTr 
)

Returns the transformation defining the given edge element. The transformation is stored in a user-defined variable.

Definition at line 329 of file mesh.cpp.

ElementTransformation * mfem::Mesh::GetEdgeTransformation ( int  EdgeNo)

Returns the transformation defining the given face element.

Definition at line 367 of file mesh.cpp.

Table * mfem::Mesh::GetEdgeVertexTable ( ) const

Returns the edge-to-vertex Table (3D)

Definition at line 3451 of file mesh.cpp.

void mfem::Mesh::GetEdgeVertices ( int  i,
Array< int > &  vert 
) const

Returns the indices of the vertices of edge i.

Definition at line 3421 of file mesh.cpp.

const Element* mfem::Mesh::GetElement ( int  i) const
inline

Definition at line 424 of file mesh.hpp.

Element* mfem::Mesh::GetElement ( int  i)
inline

Definition at line 426 of file mesh.hpp.

void mfem::Mesh::GetElementArrayEdgeTable ( const Array< Element * > &  elem_array,
const DSTable v_to_v,
Table el_to_edge 
)
staticprotected

Definition at line 3701 of file mesh.cpp.

int mfem::Mesh::GetElementBaseGeometry ( int  i) const
inline

Definition at line 436 of file mesh.hpp.

void mfem::Mesh::GetElementColoring ( Array< int > &  colors,
int  el0 = 0 
)

Definition at line 7515 of file mesh.cpp.

void mfem::Mesh::GetElementEdges ( int  i,
Array< int > &  edges,
Array< int > &  cor 
) const

Return the indices and the orientations of all edges of element i.

Definition at line 3345 of file mesh.cpp.

void mfem::Mesh::GetElementFaces ( int  i,
Array< int > &  fcs,
Array< int > &  cor 
) const

Return the indices and the orientations of all faces of element i.

Definition at line 3534 of file mesh.cpp.

void mfem::Mesh::GetElementJacobian ( int  i,
DenseMatrix J 
)
protected

Compute the Jacobian of the transformation from the perfect reference element at the center of the element.

Definition at line 31 of file mesh.cpp.

double mfem::Mesh::GetElementSize ( int  i,
int  type = 0 
)

Get the size of the i-th element relative to the perfect reference element.

Definition at line 39 of file mesh.cpp.

double mfem::Mesh::GetElementSize ( int  i,
const Vector dir 
)

Definition at line 51 of file mesh.cpp.

int mfem::Mesh::GetElementToEdgeTable ( Table e_to_f,
Array< int > &  be_to_f 
)
protected

Return element to edge table and the indeces for the boundary edges. The entries in the table are ordered according to the order of the nodes in the elements. For example, if T is the element to edge table T(i, 0) gives the index of edge in element i that connects vertex 0 to vertex 1, etc. Returns the number of the edges.

Definition at line 3748 of file mesh.cpp.

STable3D * mfem::Mesh::GetElementToFaceTable ( int  ret_ftbl = 0)
protected

Definition at line 4053 of file mesh.cpp.

void mfem::Mesh::GetElementTransformation ( int  i,
IsoparametricTransformation ElTr 
)

Builds the transformation defining the i-th element in the user-defined variable.

Definition at line 157 of file mesh.cpp.

ElementTransformation * mfem::Mesh::GetElementTransformation ( int  i)

Returns the transformation defining the i-th element.

Definition at line 211 of file mesh.cpp.

void mfem::Mesh::GetElementTransformation ( int  i,
const Vector nodes,
IsoparametricTransformation ElTr 
)

Return the transformation defining the i-th element assuming the position of the vertices/nodes are given by 'nodes'.

Definition at line 181 of file mesh.cpp.

int mfem::Mesh::GetElementType ( int  i) const

Returns the type of element i.

Definition at line 3632 of file mesh.cpp.

void mfem::Mesh::GetElementVertices ( int  i,
Array< int > &  dofs 
) const
inline

Returns the indices of the dofs of element i.

Definition at line 443 of file mesh.hpp.

double mfem::Mesh::GetElementVolume ( int  i)

Definition at line 60 of file mesh.cpp.

const Element* mfem::Mesh::GetFace ( int  i) const
inline

Definition at line 432 of file mesh.hpp.

int mfem::Mesh::GetFaceBaseGeometry ( int  i) const

Definition at line 3589 of file mesh.cpp.

void mfem::Mesh::GetFaceEdges ( int  i,
Array< int > &  edges,
Array< int > &  o 
) const

Return the indices and the orientations of all edges of face i. Works for both 2D (face=edge) and 3D faces.

Definition at line 3393 of file mesh.cpp.

Table * mfem::Mesh::GetFaceEdgeTable ( ) const

Returns the face-to-edge Table (3D)

Definition at line 3429 of file mesh.cpp.

void mfem::Mesh::GetFaceElements ( int  Face,
int *  Elem1,
int *  Elem2 
)

Definition at line 615 of file mesh.cpp.

FaceElementTransformations * mfem::Mesh::GetFaceElementTransformations ( int  FaceNo,
int  mask = 31 
)

Returns (a poiter to a structure containing) the following data: 1) Elem1No - the index of the first element that contains this face this is the element that has the same outward unit normal vector as the face; 2) Elem2No - the index of the second element that contains this face this element has outward unit normal vector as the face multiplied with -1; 3) Elem1, Elem2 - pointers to the ElementTransformation's of the first and the second element respectively; 4) Face - pointer to the ElementTransformation of the face; 5) Loc1, Loc2 - IntegrationPointTransformation's mapping the face coordinate system to the element coordinate system (both in their reference elements). Used to transform IntegrationPoints from face to element. More formally, let: TL1, TL2 be the transformations represented by Loc1, Loc2, TE1, TE2 - the transformations represented by Eleme1, Elem2, TF - the transformation represented by Face, then TF(x) = TE1(TL1(x)) = TE2(TL2(x)) for all x in the reference face. 6) FaceGeom - the base geometry for the face. The mask specifies which fields in the structure to return: mask & 1 - Elem1, mask & 2 - Elem2, mask & 4 - Loc1, mask & 8 - Loc2, mask & 16 - Face.

Definition at line 499 of file mesh.cpp.

void mfem::Mesh::GetFaceInfos ( int  Face,
int *  Inf1,
int *  Inf2 
)

Definition at line 621 of file mesh.cpp.

STable3D * mfem::Mesh::GetFacesTable ( )
protected

Definition at line 4020 of file mesh.cpp.

Table * mfem::Mesh::GetFaceToElementTable ( ) const

Return the "face"-element Table. Here "face" refers to face (3D), edge (2D), or vertex (1D). The returned Table must be destroyed by the caller.

Definition at line 3506 of file mesh.cpp.

void mfem::Mesh::GetFaceTransformation ( int  i,
IsoparametricTransformation FTr 
)

Returns the transformation defining the given face element. The transformation is stored in a user-defined variable.

Definition at line 248 of file mesh.cpp.

ElementTransformation * mfem::Mesh::GetFaceTransformation ( int  FaceNo)

Returns the transformation defining the given face element.

Definition at line 323 of file mesh.cpp.

void mfem::Mesh::GetFaceVertices ( int  i,
Array< int > &  vert 
) const
inline

Returns the indices of the vertices of face i.

Definition at line 461 of file mesh.hpp.

int mfem::Mesh::GetFineElem ( int  i,
int  j 
)

For a given coarse element i and index j of one of its subelements return the global index of the fine element in the fine mesh.

Definition at line 6601 of file mesh.cpp.

int mfem::Mesh::GetFineElemPath ( int  i,
int  j 
)
protected

Definition at line 6769 of file mesh.cpp.

ElementTransformation * mfem::Mesh::GetFineElemTrans ( int  i,
int  j 
)

For a given coarse element i and index j of one of its subelements return the transformation that transforms the fine element into the coordinate system of the coarse element. Clear, isn't it? :-)

Definition at line 6802 of file mesh.cpp.

FaceElementTransformations* mfem::Mesh::GetInteriorFaceTransformations ( int  FaceNo)
inline

Definition at line 563 of file mesh.hpp.

double mfem::Mesh::GetLength ( int  i,
int  j 
) const
protected

Return the length of the segment from node i to node j.

Definition at line 3688 of file mesh.cpp.

void mfem::Mesh::GetLocalPtToSegTransformation ( IsoparametricTransformation Transf,
int  i 
)
protected

Used in GetFaceElementTransformations (...)

Definition at line 374 of file mesh.cpp.

void mfem::Mesh::GetLocalQuadToHexTransformation ( IsoparametricTransformation loc,
int  i 
)
protected

Used in GetFaceElementTransformations (...)

Definition at line 478 of file mesh.cpp.

void mfem::Mesh::GetLocalSegToQuadTransformation ( IsoparametricTransformation loc,
int  i 
)
protected

Definition at line 411 of file mesh.cpp.

void mfem::Mesh::GetLocalSegToTriTransformation ( IsoparametricTransformation loc,
int  i 
)
protected

Definition at line 388 of file mesh.cpp.

void mfem::Mesh::GetLocalTriToTetTransformation ( IsoparametricTransformation loc,
int  i 
)
protected

Used in GetFaceElementTransformations (...)

Definition at line 455 of file mesh.cpp.

int mfem::Mesh::GetNBE ( ) const
inline

Returns number of boundary elements.

Definition at line 399 of file mesh.hpp.

int mfem::Mesh::GetNE ( ) const
inline

Returns number of elements.

Definition at line 396 of file mesh.hpp.

int mfem::Mesh::GetNEdges ( ) const
inline

Return the number of edges.

Definition at line 402 of file mesh.hpp.

int mfem::Mesh::GetNFaces ( ) const
inline

Return the number of faces in a 3D mesh.

Definition at line 405 of file mesh.hpp.

const FiniteElementSpace * mfem::Mesh::GetNodalFESpace ( )

Return the FiniteElementSpace on which the current mesh nodes are defined or NULL if the mesh does not have nodes.

Definition at line 3078 of file mesh.cpp.

void mfem::Mesh::GetNode ( int  i,
double *  coord 
)

Definition at line 4909 of file mesh.cpp.

void mfem::Mesh::GetNodes ( Vector node_coord) const

Definition at line 4948 of file mesh.cpp.

GridFunction* mfem::Mesh::GetNodes ( )
inline

Return a pointer to the internal node GridFunction (may be NULL).

Definition at line 638 of file mesh.hpp.

void mfem::Mesh::GetNodes ( GridFunction nodes) const

Return the mesh nodes/vertices projected on the given GridFunction.

Definition at line 3054 of file mesh.cpp.

int mfem::Mesh::GetNumFaces ( ) const

Return the number of faces (3D), edges (2D) or vertices (1D).

Definition at line 3083 of file mesh.cpp.

int mfem::Mesh::GetNumFineElems ( int  i)

For a given coarse element i returns the number of subelements it is divided into.

Definition at line 6445 of file mesh.cpp.

int mfem::Mesh::GetNV ( ) const
inline

Returns number of vertices. Vertices are only at the corners of elements, where you would expect them in the lowest-order mesh.

Definition at line 393 of file mesh.hpp.

void mfem::Mesh::GetPointMatrix ( int  i,
DenseMatrix pointmat 
) const

Definition at line 3660 of file mesh.cpp.

int mfem::Mesh::GetQuadOrientation ( const int *  base,
const int *  test 
)
staticprotected

Returns the orientation of "test" relative to "base".

Definition at line 3230 of file mesh.cpp.

int mfem::Mesh::GetRefinementType ( int  i)

For a given coarse element i returns its refinement type.

Definition at line 6541 of file mesh.cpp.

int mfem::Mesh::GetState ( ) const
inline

Definition at line 694 of file mesh.hpp.

FiniteElement * mfem::Mesh::GetTransformationFEforElementType ( int  ElemType)
static

Definition at line 141 of file mesh.cpp.

int mfem::Mesh::GetTriOrientation ( const int *  base,
const int *  test 
)
staticprotected

Returns the orientation of "test" relative to "base".

Definition at line 3200 of file mesh.cpp.

const double* mfem::Mesh::GetVertex ( int  i) const
inline

Return pointer to vertex i's coordinates.

Definition at line 421 of file mesh.hpp.

double* mfem::Mesh::GetVertex ( int  i)
inline

Definition at line 422 of file mesh.hpp.

Table * mfem::Mesh::GetVertexToElementTable ( )

The returned Table must be destroyed by the caller.

Definition at line 3475 of file mesh.cpp.

void mfem::Mesh::GetVertexToVertexTable ( DSTable v_to_v) const
protected

Return vertex to vertex table. The connections stored in the table are from smaller to bigger vertex index, i.e. if i<j and (i, j) is in the table, then (j, i) is not stored.

Definition at line 3723 of file mesh.cpp.

void mfem::Mesh::GetVertices ( Vector vert_coord) const

Definition at line 4893 of file mesh.cpp.

void mfem::Mesh::GreenRefinement ( int  i,
const DSTable v_to_v,
int *  edge1,
int *  edge2,
int *  middle 
)
inlineprotected

Green refinement. Element with index i is refined. The default refinement for now is Bisection.

Definition at line 159 of file mesh.hpp.

void mfem::Mesh::HexUniformRefinement ( )
protectedvirtual

Refine hexahedral mesh.

Definition at line 5156 of file mesh.cpp.

void mfem::Mesh::Init ( )
protected

Definition at line 627 of file mesh.cpp.

void mfem::Mesh::InitMesh ( int  _Dim,
int  _spaceDim,
int  NVert,
int  NElem,
int  NBdrElem 
)
protected

Begin construction of a mesh.

Definition at line 742 of file mesh.cpp.

void mfem::Mesh::InitTables ( )
protected

Definition at line 639 of file mesh.cpp.

void mfem::Mesh::KnotInsert ( Array< KnotVector * > &  kv)

Definition at line 2853 of file mesh.cpp.

void mfem::Mesh::Load ( std::istream &  input,
int  generate_edges = 0,
int  refine = 1,
bool  fix_orientation = true 
)

Definition at line 1772 of file mesh.cpp.

void mfem::Mesh::LoadPatchTopo ( std::istream &  input,
Array< int > &  edge_to_knot 
)
protected

Read NURBS patch/macro-element mesh.

Definition at line 2951 of file mesh.cpp.

void mfem::Mesh::LocalRefinement ( const Array< int > &  marked_el,
int  type = 3 
)
protectedvirtual

This function is not public anymore. Use GeneralRefinement instead.

Reimplemented in mfem::ParMesh.

Definition at line 5371 of file mesh.cpp.

void mfem::Mesh::Make1D ( int  n,
double  sx = 1.0 
)
protected

Creates a 1D mesh for the interval [0,sx] divided into n equal intervals.

Definition at line 1621 of file mesh.cpp.

void mfem::Mesh::Make2D ( int  nx,
int  ny,
Element::Type  type,
int  generate_edges,
double  sx,
double  sy 
)
protected

Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny quadrilaterals if type = QUADRILATERAL or into 2*nx*ny triangles if type = TRIANGLE. If generate_edges = 0 (default) edges are not generated, if 1 edges are generated.

Definition at line 1473 of file mesh.cpp.

void mfem::Mesh::Make3D ( int  nx,
int  ny,
int  nz,
Element::Type  type,
int  generate_edges,
double  sx,
double  sy,
double  sz 
)
protected

Creates mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into nx*ny*nz hexahedrals if type=HEXAHEDRON or into 6*nx*ny*nz tetrahedrons if type=TETRAHEDRON. If generate_edges = 0 (default) edges are not generated, if 1 edges are generated.

Definition at line 1317 of file mesh.cpp.

void mfem::Mesh::MarkForRefinement ( )
protected

Definition at line 939 of file mesh.cpp.

void mfem::Mesh::MarkTetMeshForRefinement ( )
protected

Definition at line 987 of file mesh.cpp.

void mfem::Mesh::MarkTriMeshForRefinement ( )
protected

Definition at line 950 of file mesh.cpp.

int mfem::Mesh::MeshGenerator ( )
inline

Truegrid or NetGen?

Definition at line 389 of file mesh.hpp.

void mfem::Mesh::MesquiteSmooth ( const int  mesquite_option = 0)

Definition at line 1051 of file mesquite.cpp.

void mfem::Mesh::MoveNodes ( const Vector displacements)

Definition at line 4940 of file mesh.cpp.

void mfem::Mesh::MoveVertices ( const Vector displacements)

Definition at line 4886 of file mesh.cpp.

Element * mfem::Mesh::NewElement ( int  geom)

Definition at line 1670 of file mesh.cpp.

void mfem::Mesh::NewNodes ( GridFunction nodes,
bool  make_owner = false 
)

Replace the internal node GridFunction with the given GridFunction.

Definition at line 4964 of file mesh.cpp.

void mfem::Mesh::NonconformingRefinement ( const Array< Refinement > &  refinements,
int  nc_limit = 0 
)
protected

This function is not public anymore. Use GeneralRefinement instead.

Definition at line 5738 of file mesh.cpp.

void mfem::Mesh::NURBSUniformRefinement ( )
protectedvirtual

Refine NURBS mesh.

Definition at line 2868 of file mesh.cpp.

void mfem::Mesh::PrepareNodeReorder ( DSTable **  old_v_to_v,
Table **  old_elem_vert 
)
protected

Definition at line 1006 of file mesh.cpp.

void mfem::Mesh::Print ( std::ostream &  out = std::cout) const
virtual

Print the mesh to the given stream using the default MFEM mesh format.

Reimplemented in mfem::ParMesh.

Definition at line 7136 of file mesh.cpp.

void mfem::Mesh::PrintCharacteristics ( Vector Vh = NULL,
Vector Vk = NULL 
)

Definition at line 76 of file mesh.cpp.

void mfem::Mesh::PrintElement ( const Element el,
std::ostream &  out 
)
staticprotected

Definition at line 1727 of file mesh.cpp.

void mfem::Mesh::PrintElementsWithPartitioning ( int *  partitioning,
std::ostream &  out,
int  interior_faces = 0 
)

Definition at line 7678 of file mesh.cpp.

void mfem::Mesh::PrintElementWithoutAttr ( const Element el,
std::ostream &  out 
)
staticprotected

Definition at line 1705 of file mesh.cpp.

void mfem::Mesh::PrintSurfaces ( const Table Aface_face,
std::ostream &  out 
) const

Print set of disjoint surfaces:

If Aface_face(i,j) != 0, print face j as a boundary element with attribute i+1.

Definition at line 8000 of file mesh.cpp.

void mfem::Mesh::PrintTopo ( std::ostream &  out,
const Array< int > &  e_to_k 
) const
protected

Definition at line 7189 of file mesh.cpp.

void mfem::Mesh::PrintVTK ( std::ostream &  out)

Print the mesh in VTK format (linear and quadratic meshes only).

Definition at line 7225 of file mesh.cpp.

void mfem::Mesh::PrintVTK ( std::ostream &  out,
int  ref,
int  field_data = 0 
)

Print the mesh in VTK format. The parameter ref specifies an element subdivision number (useful for high order fields and curved meshes). If the optional field_data is set, we also add a FIELD section in the beginning of the file with additional dataset information.

Definition at line 7377 of file mesh.cpp.

void mfem::Mesh::PrintWithPartitioning ( int *  partitioning,
std::ostream &  out,
int  elem_attr = 0 
) const

Prints the mesh with bdr elements given by the boundary of the subdomains, so that the boundary of subdomain i has bdr attribute i+1.

Definition at line 7582 of file mesh.cpp.

void mfem::Mesh::PrintXG ( std::ostream &  out = std::cout) const
virtual

Print the mesh to the given stream using Netgen/Truegrid format.

Reimplemented in mfem::ParMesh.

Definition at line 6992 of file mesh.cpp.

void mfem::Mesh::QuadUniformRefinement ( )
protectedvirtual

Refine quadrilateral mesh.

Definition at line 5009 of file mesh.cpp.

Element * mfem::Mesh::ReadElement ( std::istream &  input)
protected

Definition at line 1715 of file mesh.cpp.

Element * mfem::Mesh::ReadElementWithoutAttr ( std::istream &  input)
protected

Definition at line 1690 of file mesh.cpp.

void mfem::Mesh::RedRefinement ( int  i,
const DSTable v_to_v,
int *  edge1,
int *  edge2,
int *  middle 
)
inlineprotected

Red refinement. Element with index i is refined. The default red refinement for now is Uniform.

Definition at line 153 of file mesh.hpp.

void mfem::Mesh::ReorientTetMesh ( )
virtual

This method modifies a tetrahedral mesh so that Nedelec spaces of order greater than 1 can be defined on the mesh. Specifically, we 1) rotate all tets in the mesh so that the vertices {v0, v1, v2, v3} satisfy: v0 < v1 < min(v2, v3). 2) rotate all boundary triangles so that the vertices {v0, v1, v2} satisfy: v0 < min(v1, v2).

Note: refinement does not work after a call to this method!

Reimplemented in mfem::ParMesh.

Definition at line 4118 of file mesh.cpp.

void mfem::Mesh::Rotate3 ( int &  a,
int &  b,
int &  c 
)
inlinestaticprotected

Definition at line 801 of file mesh.hpp.

void mfem::Mesh::ScaleElements ( double  sf)

Definition at line 8121 of file mesh.cpp.

void mfem::Mesh::ScaleSubdomains ( double  sf)

Definition at line 8059 of file mesh.cpp.

void mfem::Mesh::SetAttributes ( )

Definition at line 683 of file mesh.cpp.

void mfem::Mesh::SetMeshGen ( )
protected

Definition at line 1733 of file mesh.cpp.

void mfem::Mesh::SetNodalFESpace ( FiniteElementSpace nfes)

Replace the internal node GridFunction with a new GridFunction defined on the given FiniteElementSpace. The new node coordinates are projected (derived) from the current nodes/vertices.

Definition at line 3066 of file mesh.cpp.

void mfem::Mesh::SetNodalGridFunction ( GridFunction nodes,
bool  make_owner = false 
)

Replace the internal node GridFunction with the given GridFunction. The given GridFunction is updated with node coordinates projected (derived) from the current nodes/vertices.

Definition at line 3072 of file mesh.cpp.

void mfem::Mesh::SetNode ( int  i,
const double *  coord 
)

Definition at line 4924 of file mesh.cpp.

void mfem::Mesh::SetNodes ( const Vector node_coord)

Definition at line 4956 of file mesh.cpp.

void mfem::Mesh::SetState ( int  s)

Change the mesh state to NORMAL, TWO_LEVEL_COARSE, TWO_LEVEL_FINE.

Definition at line 6298 of file mesh.cpp.

void mfem::Mesh::SetVertices ( const Vector vert_coord)

Definition at line 4902 of file mesh.cpp.

void mfem::Mesh::ShiftL2R ( int &  a,
int &  b,
int &  c 
)
inlinestaticprotected

Definition at line 795 of file mesh.hpp.

int mfem::Mesh::SpaceDimension ( ) const
inline

Definition at line 418 of file mesh.hpp.

void mfem::Mesh::Swap ( Mesh other,
bool  non_geometry = false 
)
protected

Swaps internal data with another mesh. By default, non-geometry members like 'ncmesh' and 'NURBSExt' are only swapped when 'non_geometry' is set.

Definition at line 5832 of file mesh.cpp.

void mfem::Mesh::SwapNodes ( GridFunction *&  nodes,
int &  own_nodes_ 
)

Swap the internal node GridFunction pointer and onwership flag members with the given ones.

Definition at line 4978 of file mesh.cpp.

void mfem::Mesh::Transform ( void(*)(const Vector &, Vector &)  f)

Definition at line 8183 of file mesh.cpp.

void mfem::Mesh::UniformRefinement ( int  i,
const DSTable v_to_v,
int *  edge1,
int *  edge2,
int *  middle 
)
protected

Uniform Refinement. Element with index i is refined uniformly.

Definition at line 6225 of file mesh.cpp.

void mfem::Mesh::UniformRefinement ( )

Refine all mesh elements.

Definition at line 5874 of file mesh.cpp.

void mfem::Mesh::UpdateNodes ( )
protected

Update the nodes of a curved mesh after refinement.

Definition at line 5003 of file mesh.cpp.

void mfem::Mesh::UpdateNURBS ( )
protected

Definition at line 2896 of file mesh.cpp.

void mfem::Mesh::UseTwoLevelState ( int  use)
inline

Sets or clears the flag that indicates that mesh refinement methods should put the mesh in two-level state.

Definition at line 684 of file mesh.hpp.

Friends And Related Function Documentation

friend class NURBSExtension
friend

Definition at line 45 of file mesh.hpp.

friend class ParMesh
friend

Definition at line 43 of file mesh.hpp.

friend class Tetrahedron
friend

Definition at line 109 of file mesh.hpp.

Member Data Documentation

Array<int> mfem::Mesh::attributes

Definition at line 304 of file mesh.hpp.

Array<int> mfem::Mesh::bdr_attributes

Definition at line 305 of file mesh.hpp.

Array<int> mfem::Mesh::be_to_edge
protected

Definition at line 77 of file mesh.hpp.

Array<int> mfem::Mesh::be_to_face
protected

Definition at line 79 of file mesh.hpp.

Table* mfem::Mesh::bel_to_edge
protected

Definition at line 78 of file mesh.hpp.

MemAlloc<BisectedElement, 1024> mfem::Mesh::BEMemory
protected

Definition at line 112 of file mesh.hpp.

Array<Element *> mfem::Mesh::boundary
protected

Definition at line 68 of file mesh.hpp.

Table * mfem::Mesh::c_bel_to_edge
protected

Definition at line 83 of file mesh.hpp.

Table* mfem::Mesh::c_el_to_edge
protected

Definition at line 83 of file mesh.hpp.

Table* mfem::Mesh::c_el_to_face
protected

Definition at line 85 of file mesh.hpp.

int mfem::Mesh::c_NumOfBdrElements
protected

Definition at line 59 of file mesh.hpp.

int mfem::Mesh::c_NumOfEdges
protected

Definition at line 61 of file mesh.hpp.

int mfem::Mesh::c_NumOfElements
protected

Definition at line 59 of file mesh.hpp.

int mfem::Mesh::c_NumOfFaces
protected

Definition at line 61 of file mesh.hpp.

int mfem::Mesh::c_NumOfVertices
protected

Definition at line 59 of file mesh.hpp.

int mfem::Mesh::Dim
protected

Definition at line 48 of file mesh.hpp.

Table* mfem::Mesh::edge_vertex
mutableprotected

Definition at line 81 of file mesh.hpp.

IsoparametricTransformation mfem::Mesh::EdgeTransformation
protected

Definition at line 89 of file mesh.hpp.

Table* mfem::Mesh::el_to_edge
protected

Definition at line 74 of file mesh.hpp.

Table* mfem::Mesh::el_to_el
protected

Definition at line 76 of file mesh.hpp.

Table* mfem::Mesh::el_to_face
protected

Definition at line 75 of file mesh.hpp.

Array<Element *> mfem::Mesh::elements
protected

Definition at line 64 of file mesh.hpp.

Table * mfem::Mesh::f_bel_to_edge
protected

Definition at line 83 of file mesh.hpp.

Table * mfem::Mesh::f_el_to_edge
protected

Definition at line 83 of file mesh.hpp.

Table * mfem::Mesh::f_el_to_face
protected

Definition at line 85 of file mesh.hpp.

int mfem::Mesh::f_NumOfBdrElements
protected

Definition at line 60 of file mesh.hpp.

int mfem::Mesh::f_NumOfEdges
protected

Definition at line 62 of file mesh.hpp.

int mfem::Mesh::f_NumOfElements
protected

Definition at line 60 of file mesh.hpp.

int mfem::Mesh::f_NumOfFaces
protected

Definition at line 62 of file mesh.hpp.

int mfem::Mesh::f_NumOfVertices
protected

Definition at line 60 of file mesh.hpp.

Table* mfem::Mesh::face_edge
mutableprotected

Definition at line 80 of file mesh.hpp.

FaceElementTransformations mfem::Mesh::FaceElemTr
protected

Definition at line 90 of file mesh.hpp.

Array<Element *> mfem::Mesh::faces
protected

Definition at line 69 of file mesh.hpp.

Array<FaceInfo> mfem::Mesh::faces_info
protected

Definition at line 72 of file mesh.hpp.

IsoparametricTransformation mfem::Mesh::FaceTransformation
protected

Definition at line 89 of file mesh.hpp.

Array<int> mfem::Mesh::fc_be_to_edge
protected

Definition at line 84 of file mesh.hpp.

Array<FaceInfo> mfem::Mesh::fc_faces_info
protected

Definition at line 86 of file mesh.hpp.

const int mfem::Mesh::hex_faces
staticprotected
Initial value:
=
{{3, 2, 1, 0}, {0, 1, 5, 4},
{1, 2, 6, 5}, {2, 3, 7, 6},
{3, 0, 4, 7}, {4, 5, 6, 7}}

Definition at line 103 of file mesh.hpp.

int mfem::Mesh::meshgen
protected

Definition at line 57 of file mesh.hpp.

Mesh* mfem::Mesh::nc_coarse_level
protected

Definition at line 100 of file mesh.hpp.

NCMesh* mfem::Mesh::ncmesh

Definition at line 308 of file mesh.hpp.

GridFunction* mfem::Mesh::Nodes
protected

Definition at line 95 of file mesh.hpp.

int mfem::Mesh::NumOfBdrElements
protected

Definition at line 51 of file mesh.hpp.

int mfem::Mesh::NumOfEdges
protected

Definition at line 52 of file mesh.hpp.

int mfem::Mesh::NumOfElements
protected

Definition at line 51 of file mesh.hpp.

int mfem::Mesh::NumOfFaces
protected

Definition at line 52 of file mesh.hpp.

int mfem::Mesh::NumOfVertices
protected

Definition at line 51 of file mesh.hpp.

NURBSExtension* mfem::Mesh::NURBSext

Definition at line 307 of file mesh.hpp.

int mfem::Mesh::own_nodes
protected

Definition at line 96 of file mesh.hpp.

const int mfem::Mesh::quad_orientations
staticprotected
Initial value:
=
{{0, 1, 2, 3}, {0, 3, 2, 1},
{1, 2, 3, 0}, {1, 0, 3, 2},
{2, 3, 0, 1}, {2, 1, 0, 3},
{3, 0, 1, 2}, {3, 2, 1, 0}}

Definition at line 106 of file mesh.hpp.

int mfem::Mesh::spaceDim
protected

Definition at line 49 of file mesh.hpp.

int mfem::Mesh::State
protected

Definition at line 54 of file mesh.hpp.

const int mfem::Mesh::tet_faces
staticprotected
Initial value:
=
{{1, 2, 3}, {0, 3, 2},
{0, 1, 3}, {0, 2, 1}}

Definition at line 102 of file mesh.hpp.

MemAlloc<Tetrahedron, 1024> mfem::Mesh::TetMemory
protected

Definition at line 111 of file mesh.hpp.

IsoparametricTransformation mfem::Mesh::Transformation
protected

Definition at line 88 of file mesh.hpp.

IsoparametricTransformation mfem::Mesh::Transformation2
protected

Definition at line 88 of file mesh.hpp.

const int mfem::Mesh::tri_orientations
staticprotected
Initial value:
=
{{0, 1, 2}, {1, 0, 2},
{2, 0, 1}, {2, 1, 0},
{1, 2, 0}, {0, 2, 1}}

Definition at line 105 of file mesh.hpp.

Array<Vertex> mfem::Mesh::vertices
protected

Definition at line 67 of file mesh.hpp.

int mfem::Mesh::WantTwoLevelState
protected

Definition at line 54 of file mesh.hpp.


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