MFEM v2.0
Classes | Public Types | Public Member Functions | Public Attributes | Protected Member Functions | Static Protected Member Functions | Protected Attributes | Friends
Mesh Class Reference

#include <mesh.hpp>

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

List of all members.

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)
void AddVertex (double *)
void AddElement (Element *elem)
void AddTri (int *vi, int attr=1)
void AddTriangle (int *vi, int attr=1)
void AddQuad (int *vi, int attr=1)
void AddTet (int *vi, int attr=1)
void AddHex (int *vi, int attr=1)
void AddBdrSegment (int *vi, int attr=1)
void AddBdrTriangle (int *vi, int attr=1)
void AddBdrQuad (int *vi, int attr=1)
void GenerateBoundaryElements ()
void FinalizeTriMesh (int generate_edges=0, int refine=0)
void FinalizeQuadMesh (int generate_edges=0, int refine=0)
void FinalizeTetMesh (int generate_edges=0, int refine=0)
void FinalizeHexMesh (int generate_edges=0, int refine=0)
void SetAttributes ()
 Mesh (int nx, int ny, Element::Type type, int generate_edges=0, double sx=1.0, double sy=1.0)
 Mesh (int n)
 Mesh (istream &input, int generate_edges=0, int refine=1)
 Mesh (Mesh *mesh_array[], int num_pieces)
 Create a disjoint mesh from the given mesh array.
void Load (istream &input, int generate_edges=0, int refine=1)
void SetNodalFESpace (FiniteElementSpace *nfes)
void SetNodalGridFunction (GridFunction *nodes)
const FiniteElementSpaceGetNodalFESpace ()
int MeshGenerator ()
 Truegrid or NetGen?
int GetNV () const
 Returns number of vertices.
int GetNE () const
 Returns number of elements.
int GetNBE () const
 Returns number of boundary elements.
int GetNEdges () const
 Return the number of edges.
int GetNFaces () const
int EulerNumber () const
 Equals 1 + num_holes - num_loops.
int EulerNumber2D () const
 Equals 1 - num_holes.
int Dimension () const
const double * GetVertex (int i) const
 Return pointer to vertex i's coordinates.
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.
void GetBdrElementVertices (int i, Array< int > &dofs) const
 Returns the indices of the dofs of boundary element i.
void GetElementEdges (int i, Array< int > &, Array< int > &) const
 Return the indices and the orientations of all edges of element i.
void GetBdrElementEdges (int i, Array< int > &, Array< int > &) const
 Return the indices and the orientations of all edges of bdr element i.
void GetFaceEdges (int i, Array< int > &, Array< int > &) const
 Return the indices and the orientations of all edges of face i.
void GetFaceVertices (int i, Array< int > &vert) const
 Returns the indices of the vertices of face i.
void GetEdgeVertices (int i, Array< int > &vert) const
 Returns the indices of the vertices of edge i.
TableGetFaceEdgeTable () const
 Returns the face-to-edge Table (3D)
TableGetEdgeVertexTable () const
 Returns the edge-to-vertex Table (3D)
void GetElementFaces (int i, Array< int > &, Array< int > &) const
 Return the indices and the orientations of all faces of element i.
void GetBdrElementFace (int i, int *, int *) const
 Return the index and the orientation of the face of bdr element i. (3D)
int GetBdrElementEdgeIndex (int i) const
int GetElementType (int i) const
 Returns the type of element i.
int GetBdrElementType (int i) const
 Returns the type of boundary element i.
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.
void GetElementTransformation (int i, const Vector &nodes, IsoparametricTransformation *ElTr)
ElementTransformationGetBdrElementTransformation (int i)
 Returns the transformation defining the i-th boundary element.
void GetFaceTransformation (int i, IsoparametricTransformation *FTr)
ElementTransformationGetFaceTransformation (int FaceNo)
 Returns the transformation defining the given face element.
FaceElementTransformationsGetFaceElementTransformations (int FaceNo, int mask=31)
FaceElementTransformationsGetInteriorFaceTransformations (int FaceNo)
FaceElementTransformationsGetBdrFaceTransformations (int BdrElemNo)
void GetFaceElements (int Face, int *Elem1, int *Elem2)
void GetFaceInfos (int Face, int *Inf1, int *Inf2)
void CheckElementOrientation ()
 Check the orientation of the elements.
void CheckBdrElementOrientation ()
 Check the orientation of the boundary elements.
int GetAttribute (int i) const
 Return the attribute of element i.
int GetBdrAttribute (int i) const
 Return the attribute of boundary element i.
const TableElementToElementTable ()
const TableElementToFaceTable () const
const TableElementToEdgeTable () const
TableGetVertexToElementTable ()
 The returned Table must be destroyed by the caller.
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 MoveNodes (const Vector &displacements)
void GetNodes (Vector &node_coord) const
void SetNodes (const Vector &node_coord)
GridFunctionGetNodes ()
 Return a pointer to the internal node grid function.
void NewNodes (GridFunction &nodes)
virtual void LocalRefinement (const Array< int > &marked_el, int type=3)
 Refine the marked elements.
void UniformRefinement ()
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.
int GetNumFineElems (int i)
int GetRefinementType (int i)
 For a given coarse element i returns its refinement type.
int GetFineElem (int i, int j)
ElementTransformationGetFineElemTrans (int i, int j)
virtual void PrintXG (ostream &out=cout) const
 Print the mesh to the given stream using Netgen/Truegrid format.
virtual void Print (ostream &out=cout) const
 Print the mesh to the given stream using the default MFEM mesh format.
void PrintVTK (ostream &out)
 Print the mesh in VTK format (linear and quadratic meshes only).
void PrintVTK (ostream &out, int ref, int field_data=0)
void GetElementColoring (Array< int > &colors, int el0=0)
void PrintWithPartitioning (int *partitioning, ostream &out, int elem_attr=0) const
void PrintElementsWithPartitioning (int *partitioning, ostream &out, int interior_faces=0)
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)
virtual ~Mesh ()
 Destroys mesh.

Public Attributes

Array< int > attributes
Array< int > bdr_attributes
NURBSExtensionNURBSext

Protected Member Functions

ElementNewElement (int geom)
void Init ()
void InitTables ()
void DeleteTables ()
void DeleteCoarseTables ()
ElementReadElement (istream &)
double GetLength (int i, int j) const
 Return the length of the segment from node i to node j.
void GetElementJacobian (int i, DenseMatrix &J)
void MarkForRefinement ()
void MarkTriMeshForRefinement ()
void MarkTetMeshForRefinement ()
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.
virtual void QuadUniformRefinement ()
 Refine quadrilateral mesh.
virtual void HexUniformRefinement ()
 Refine hexahedral mesh.
virtual void NURBSUniformRefinement ()
 Refine NURBS mesh.
void LoadPatchTopo (istream &input, Array< int > &edge_to_knot)
 Read NURBS patch/macro-element mesh.
void UpdateNURBS ()
void PrintTopo (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)
FiniteElementGetTransformationFEforElementType (int)
void GetLocalSegToTriTransformation (IsoparametricTransformation &loc, int i)
 Used in GetFaceElementTransformations (...)
void GetLocalSegToQuadTransformation (IsoparametricTransformation &loc, int i)
void GetLocalTriToTetTransformation (IsoparametricTransformation &loc, int i)
 Used in GetFaceElementTransformations (...)
void GetLocalQuadToHexTransformation (IsoparametricTransformation &loc, int i)
 Used in GetFaceElementTransformations (...)
void GetVertexToVertexTable (DSTable &) const
int GetElementToEdgeTable (Table &, Array< int > &)
void AddSegmentFaceElement (int lf, int gf, int el, int v0, int v1)
 Used in GenerateFaces()
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)
void FreeElement (Element *E)
void GenerateFaces ()

Static Protected Member Functions

static void PrintElement (Element *, ostream &)
static int GetTriOrientation (const int *base, const int *test)
 Returns the orientation of "test" relative to "base".
static int GetQuadOrientation (const int *base, const int *test)
 Returns the orientation of "test" relative to "base".
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 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
FaceElementTransformations FaceElemTr
GridFunctionNodes
int own_nodes

Friends

class NURBSExtension

Detailed Description

Definition at line 25 of file mesh.hpp.


Member Enumeration Documentation

anonymous enum
Enumerator:
NORMAL 
TWO_LEVEL_COARSE 
TWO_LEVEL_FINE 

Definition at line 222 of file mesh.hpp.


Constructor & Destructor Documentation

Mesh::Mesh ( ) [inline]

Definition at line 229 of file mesh.hpp.

References Dim, Init(), InitTables(), and meshgen.

Mesh::Mesh ( int  _Dim,
int  NVert,
int  NElem,
int  NBdrElem = 0 
)
Mesh::Mesh ( int  nx,
int  ny,
Element::Type  type,
int  generate_edges = 0,
double  sx = 1.0,
double  sy = 1.0 
)

Creates mesh for the unit square, 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 864 of file mesh.cpp.

References Array< T >::Append(), attributes, bdr_attributes, be_to_edge, boundary, CheckBdrElementOrientation(), CheckElementOrientation(), Dim, el_to_edge, elements, GenerateFaces(), GetElementToEdgeTable(), Init(), InitTables(), MarkTriMeshForRefinement(), meshgen, NumOfBdrElements, NumOfEdges, NumOfElements, NumOfFaces, NumOfVertices, Element::QUADRILATERAL, Array< T >::SetSize(), Element::TRIANGLE, and vertices.

Mesh::Mesh ( int  n) [explicit]

Creates 1D mesh , divided into n equal intervals.

Definition at line 1012 of file mesh.cpp.

References boundary, Dim, elements, Init(), InitTables(), meshgen, NumOfBdrElements, NumOfEdges, NumOfElements, NumOfFaces, NumOfVertices, Array< T >::SetSize(), and vertices.

Mesh::Mesh ( istream &  input,
int  generate_edges = 0,
int  refine = 1 
)

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

Definition at line 1048 of file mesh.cpp.

References Init(), InitTables(), and Load().

Mesh::Mesh ( Mesh mesh_array[],
int  num_pieces 
)
Mesh::~Mesh ( ) [virtual]

Member Function Documentation

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

Definition at line 643 of file mesh.cpp.

References boundary, and NumOfBdrElements.

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

Definition at line 633 of file mesh.cpp.

References boundary, and NumOfBdrElements.

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

Definition at line 638 of file mesh.cpp.

References boundary, and NumOfBdrElements.

void Mesh::AddElement ( Element elem) [inline]

Definition at line 233 of file mesh.hpp.

References elements, and NumOfElements.

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

Definition at line 628 of file mesh.cpp.

References elements, and NumOfElements.

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

Definition at line 610 of file mesh.cpp.

References elements, and NumOfElements.

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

Definition at line 3129 of file mesh.cpp.

References faces, faces_info, GetQuadOrientation(), GetVertices(), and mfem_error().

Referenced by GenerateFaces().

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

Used in GenerateFaces()

Definition at line 3085 of file mesh.cpp.

References faces, faces_info, and mfem_error().

Referenced by GenerateFaces().

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

Definition at line 615 of file mesh.cpp.

References elements, NumOfElements, Element::SetAttribute(), and Tetrahedron::SetVertices().

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

Definition at line 600 of file mesh.cpp.

References elements, and NumOfElements.

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

Definition at line 605 of file mesh.cpp.

References elements, and NumOfElements.

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

Definition at line 3106 of file mesh.cpp.

References faces, faces_info, GetTriOrientation(), GetVertices(), and mfem_error().

Referenced by GenerateFaces().

void Mesh::AddVertex ( double *  x)

Definition at line 591 of file mesh.cpp.

References Dim, NumOfVertices, and vertices.

void 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 4100 of file mesh.cpp.

References Dim, and vertices.

Referenced by HexUniformRefinement(), and QuadUniformRefinement().

void Mesh::Bisection ( int  i,
const DSTable v_to_v,
int *  edge1,
int *  edge2,
int *  middle 
) [protected]
void Mesh::Bisection ( int  i,
const DSTable v_to_v,
int *  middle 
) [protected]
void Mesh::BisectTetTrans ( DenseMatrix pointmat,
Tetrahedron tet,
int  child 
) [protected]

Definition at line 5611 of file mesh.cpp.

References Tetrahedron::ParseRefinementFlag().

Referenced by GetFineElemTrans().

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

Definition at line 5587 of file mesh.cpp.

Referenced by GetFineElemTrans().

int * Mesh::CartesianPartitioning ( int  nxyz[])
void Mesh::CheckBdrElementOrientation ( )
void Mesh::CheckDisplacements ( const Vector displacements,
double &  tmax 
)
void Mesh::CheckElementOrientation ( )
void Mesh::CheckPartitioning ( int *  partitioning)
void Mesh::DegreeElevate ( int  t)
void 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 506 of file mesh.cpp.

References edge_vertex, el_to_el, and face_edge.

Referenced by HexUniformRefinement(), LocalRefinement(), QuadUniformRefinement(), and ReorientTetMesh().

void Mesh::DeleteTables ( ) [protected]

Definition at line 483 of file mesh.cpp.

References bel_to_edge, Dim, edge_vertex, el_to_edge, el_to_el, el_to_face, face_edge, and InitTables().

Referenced by Load(), and ~Mesh().

int Mesh::Dimension ( ) const [inline]
const Table & Mesh::ElementToEdgeTable ( ) const

Definition at line 3078 of file mesh.cpp.

References el_to_edge, and mfem_error().

Referenced by ElementToElementTable(), and PrintElementsWithPartitioning().

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

Definition at line 3071 of file mesh.cpp.

References el_to_face, and mfem_error().

int Mesh::EulerNumber ( ) const [inline]

Equals 1 + num_holes - num_loops.

Definition at line 296 of file mesh.hpp.

References NumOfEdges, NumOfElements, NumOfFaces, and NumOfVertices.

Referenced by PrintCharacteristics().

int Mesh::EulerNumber2D ( ) const [inline]

Equals 1 - num_holes.

Definition at line 299 of file mesh.hpp.

References NumOfEdges, NumOfElements, and NumOfVertices.

Referenced by PrintCharacteristics().

void Mesh::FinalizeHexMesh ( int  generate_edges = 0,
int  refine = 0 
)
void Mesh::FinalizeQuadMesh ( int  generate_edges = 0,
int  refine = 0 
)
void Mesh::FinalizeTetMesh ( int  generate_edges = 0,
int  refine = 0 
)
void Mesh::FinalizeTriMesh ( int  generate_edges = 0,
int  refine = 0 
)
void Mesh::FreeElement ( Element E) [protected]
void Mesh::GenerateBoundaryElements ( )
void Mesh::GenerateFaces ( ) [protected]
int * Mesh::GeneratePartitioning ( int  nparts,
int  part_method = 1 
)
int Mesh::GetAttribute ( int  i) const [inline]
int Mesh::GetBdrAttribute ( int  i) const [inline]

Return the attribute of boundary element i.

Definition at line 443 of file mesh.hpp.

References boundary.

Referenced by NURBSExtension::Get2DBdrElementTopo(), NURBSExtension::Get3DBdrElementTopo(), GetBdrElementTransformation(), and SetAttributes().

const Element* Mesh::GetBdrElement ( int  i) const [inline]

Definition at line 312 of file mesh.hpp.

References boundary.

Referenced by Mesh().

Element* Mesh::GetBdrElement ( int  i) [inline]

Definition at line 314 of file mesh.hpp.

References boundary.

int Mesh::GetBdrElementBaseGeometry ( int  i) const [inline]

Definition at line 323 of file mesh.hpp.

References boundary.

Referenced by FiniteElementSpace::GetBdrElementDofs(), and FiniteElementSpace::GetBE().

int 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 2867 of file mesh.cpp.

References be_to_edge, be_to_face, and Dim.

void Mesh::GetBdrElementEdges ( int  i,
Array< int > &  edges,
Array< int > &  cor 
) const
void 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 2805 of file mesh.cpp.

References be_to_face, boundary, faces, GetBdrElementType(), GetQuadOrientation(), GetTriOrientation(), mfem_error(), Element::QUADRILATERAL, State, Element::TRIANGLE, and TWO_LEVEL_COARSE.

Referenced by FiniteElementSpace::GetBdrElementDofs(), and NURBSPatchMap::GetBdrPatchKnotVectors().

ElementTransformation * Mesh::GetBdrElementTransformation ( int  i)
int Mesh::GetBdrElementType ( int  i) const

Returns the type of boundary element i.

Definition at line 2889 of file mesh.cpp.

References Element::BISECTED, boundary, Element::GetType(), and Element::QUADRISECTED.

Referenced by GetBdrElementFace(), GetBdrElementTransformation(), GetElementToFaceTable(), and ReorientTetMesh().

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

Returns the indices of the dofs of boundary element i.

Definition at line 331 of file mesh.hpp.

References boundary.

Referenced by FiniteElementSpace::GetBdrElementDofs(), and NURBSPatchMap::GetBdrPatchKnotVectors().

FaceElementTransformations * Mesh::GetBdrFaceTransformations ( int  BdrElemNo)
void Mesh::GetBdrPointMatrix ( int  i,
DenseMatrix pointmat 
) const

Definition at line 2916 of file mesh.cpp.

References boundary, Dim, DenseMatrix::SetSize(), and vertices.

Referenced by GetBdrElementTransformation().

int Mesh::GetBisectionHierarchy ( Element E) [protected]

Definition at line 5405 of file mesh.cpp.

References Element::BISECTED, elements, and Element::GetType().

Referenced by GetRefinementType().

Table * Mesh::GetEdgeVertexTable ( ) const

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

Definition at line 2724 of file mesh.cpp.

References edge_vertex, Table::Finalize(), GetVertexToVertexTable(), DSTable::NumberOfEntries(), NumOfVertices, and Table::Push().

Referenced by GetEdgeVertices().

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

Returns the indices of the vertices of edge i.

Definition at line 2688 of file mesh.cpp.

References Dim, edge_vertex, faces, GetEdgeVertexTable(), Table::GetRow(), NumOfEdges, and Array< T >::Size().

Referenced by FiniteElementSpace::GetEdgeDofs(), and PrintElementsWithPartitioning().

Element* Mesh::GetElement ( int  i) [inline]

Definition at line 310 of file mesh.hpp.

References elements.

const Element* Mesh::GetElement ( int  i) const [inline]

Definition at line 308 of file mesh.hpp.

References elements.

Referenced by Mesh().

void Mesh::GetElementArrayEdgeTable ( const Array< Element * > &  elem_array,
const DSTable v_to_v,
Table el_to_edge 
) [static, protected]
int Mesh::GetElementBaseGeometry ( int  i) const [inline]
void Mesh::GetElementColoring ( Array< int > &  colors,
int  el0 = 0 
)

Definition at line 6369 of file mesh.cpp.

References el_to_el, ElementToElementTable(), Table::GetI(), Table::GetJ(), GetNE(), and Array< T >::SetSize().

Referenced by PrintVTK().

void Mesh::GetElementEdges ( int  i,
Array< int > &  edges,
Array< int > &  cor 
) const
void 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 2779 of file mesh.cpp.

References el_to_face, faces_info, Table::GetRow(), mfem_error(), Array< T >::SetSize(), and Array< T >::Size().

Referenced by FiniteElementSpace::GetElementDofs(), and NURBSPatchMap::GetPatchKnotVectors().

void 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 25 of file mesh.cpp.

References Geometries, Geometry::GetCenter(), GetElementBaseGeometry(), GetElementTransformation(), ElementTransformation::Jacobian(), Geometry::JacToPerfJac(), and ElementTransformation::SetIntPoint().

Referenced by GetElementSize(), and PrintCharacteristics().

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

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

Definition at line 33 of file mesh.cpp.

References DenseMatrix::CalcSingularvalue(), DenseMatrix::Det(), Dim, and GetElementJacobian().

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

Definition at line 45 of file mesh.cpp.

References Dim, GetElementJacobian(), and DenseMatrix::MultTranspose().

int 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 2990 of file mesh.cpp.

References bel_to_edge, boundary, Dim, elements, GetElementArrayEdgeTable(), GetVertexToVertexTable(), mfem_error(), DSTable::NumberOfEntries(), NumOfBdrElements, NumOfVertices, and Array< T >::SetSize().

Referenced by FinalizeHexMesh(), FinalizeQuadMesh(), FinalizeTetMesh(), FinalizeTriMesh(), HexUniformRefinement(), Load(), LoadPatchTopo(), LocalRefinement(), Mesh(), QuadUniformRefinement(), ReorientTetMesh(), and UpdateNURBS().

STable3D * Mesh::GetElementToFaceTable ( int  ret_ftbl = 0) [protected]
void Mesh::GetElementTransformation ( int  i,
IsoparametricTransformation ElTr 
)
ElementTransformation * Mesh::GetElementTransformation ( int  i)

Returns the transformation defining the i-th element.

Definition at line 201 of file mesh.cpp.

References GetElementTransformation(), and Transformation.

void Mesh::GetElementTransformation ( int  i,
const Vector nodes,
IsoparametricTransformation ElTr 
)
int Mesh::GetElementType ( int  i) const
void Mesh::GetElementVertices ( int  i,
Array< int > &  dofs 
) const [inline]

Returns the indices of the dofs of element i.

Definition at line 327 of file mesh.hpp.

References elements.

Referenced by FiniteElementSpace::GetElementDofs(), NURBSPatchMap::GetPatchKnotVectors(), GridFunction::ProjectCoefficient(), ScaleElements(), and ScaleSubdomains().

double Mesh::GetElementVolume ( int  i)
const Element* Mesh::GetFace ( int  i) const [inline]

Definition at line 316 of file mesh.hpp.

References faces.

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

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

Definition at line 2669 of file mesh.cpp.

References Dim, face_edge, faces, GetFaceEdgeTable(), Table::GetRow(), and Array< T >::SetSize().

Referenced by NURBSExtension::GenerateOffsets(), and FiniteElementSpace::GetFaceDofs().

Table * Mesh::GetFaceEdgeTable ( ) const

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

Definition at line 2702 of file mesh.cpp.

References Dim, face_edge, faces, GetElementArrayEdgeTable(), GetVertexToVertexTable(), mfem_error(), NumOfFaces, NumOfVertices, and Array< T >::Size().

Referenced by GetFaceEdges().

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

Definition at line 455 of file mesh.cpp.

References faces_info.

FaceElementTransformations * 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. 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 355 of file mesh.cpp.

References FaceElementTransformations::Elem1, FaceElementTransformations::Elem1No, FaceElementTransformations::Elem2, FaceElementTransformations::Elem2No, FaceElementTransformations::Face, FaceElemTr, FaceElementTransformations::FaceGeom, faces, faces_info, GetElementTransformation(), GetElementType(), GetFaceTransformation(), GetLocalQuadToHexTransformation(), GetLocalSegToQuadTransformation(), GetLocalSegToTriTransformation(), GetLocalTriToTetTransformation(), FaceElementTransformations::Loc1, FaceElementTransformations::Loc2, mfem_error(), NURBSext, Element::QUADRILATERAL, Element::SEGMENT, IntegrationPointTransformation::Transf, Transformation, Transformation2, and Element::TRIANGLE.

Referenced by GridFunction::ComputeH1Error(), GetBdrFaceTransformations(), GridFunction::GetFaceValues(), GridFunction::GetFaceVectorValues(), and GetInteriorFaceTransformations().

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

Definition at line 461 of file mesh.cpp.

References faces_info.

STable3D * Mesh::GetFacesTable ( ) [protected]
ElementTransformation * Mesh::GetFaceTransformation ( int  FaceNo)

Returns the transformation defining the given face element.

Definition at line 248 of file mesh.cpp.

References FaceTransformation, and GetFaceTransformation().

void 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 234 of file mesh.cpp.

References ElementTransformation::Attribute, Dim, ElementTransformation::ElementNo, faces, IsoparametricTransformation::GetPointMat(), GetTransformationFEforElementType(), IsoparametricTransformation::SetFE(), DenseMatrix::SetSize(), and vertices.

Referenced by GetFaceElementTransformations(), and GetFaceTransformation().

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

Returns the indices of the vertices of face i.

Definition at line 344 of file mesh.hpp.

References faces, and GetVertices().

Referenced by FiniteElementSpace::GetFaceDofs().

int 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 5499 of file mesh.cpp.

References Element::BISECTED, QuadrisectedElement::Child2, QuadrisectedElement::Child3, QuadrisectedElement::Child4, RefinedElement::CoarseElem, Dim, elements, RefinedElement::FirstChild, Element::GetType(), Element::OCTASECTED, Element::QUADRILATERAL, Element::QUADRISECTED, and BisectedElement::SecondChild.

int Mesh::GetFineElemPath ( int  i,
int  j 
) [protected]
ElementTransformation * Mesh::GetFineElemTrans ( int  i,
int  j 
)
FaceElementTransformations* Mesh::GetInteriorFaceTransformations ( int  FaceNo) [inline]

Definition at line 425 of file mesh.hpp.

References faces_info, and GetFaceElementTransformations().

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

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

Definition at line 2930 of file mesh.cpp.

References Dim, and vertices.

Referenced by MarkTetMeshForRefinement().

void Mesh::GetLocalQuadToHexTransformation ( IsoparametricTransformation loc,
int  i 
) [protected]
void Mesh::GetLocalSegToQuadTransformation ( IsoparametricTransformation loc,
int  i 
) [protected]
void Mesh::GetLocalSegToTriTransformation ( IsoparametricTransformation loc,
int  i 
) [protected]
void Mesh::GetLocalTriToTetTransformation ( IsoparametricTransformation loc,
int  i 
) [protected]
int Mesh::GetNBE ( ) const [inline]

Returns number of boundary elements.

Definition at line 288 of file mesh.hpp.

References NumOfBdrElements.

Referenced by NURBSExtension::GetNBP(), Mesh(), PrintCharacteristics(), and SetAttributes().

int Mesh::GetNE ( ) const [inline]
int Mesh::GetNEdges ( ) const [inline]
int Mesh::GetNFaces ( ) const [inline]
const FiniteElementSpace * Mesh::GetNodalFESpace ( )

Definition at line 2391 of file mesh.cpp.

References GridFunction::FESpace(), and Nodes.

void Mesh::GetNodes ( Vector node_coord) const

Definition at line 4077 of file mesh.cpp.

References GetVertices(), and Nodes.

Referenced by main(), and Mesh().

GridFunction* Mesh::GetNodes ( ) [inline]

Return a pointer to the internal node grid function.

Definition at line 481 of file mesh.hpp.

References Nodes.

Referenced by Mesh().

int Mesh::GetNumFineElems ( int  i)

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

Definition at line 5351 of file mesh.cpp.

References Element::BISECTED, Dim, elements, RefinedElement::FirstChild, Element::GetType(), Element::OCTASECTED, Element::QUADRISECTED, and BisectedElement::SecondChild.

int Mesh::GetNV ( ) const [inline]
void Mesh::GetPointMatrix ( int  i,
DenseMatrix pointmat 
) const
int Mesh::GetQuadOrientation ( const int *  base,
const int *  test 
) [static, protected]

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

Definition at line 2490 of file mesh.cpp.

References mfem_error().

Referenced by AddQuadFaceElement(), CheckBdrElementOrientation(), GetBdrElementFace(), and Load().

int Mesh::GetRefinementType ( int  i)
FiniteElement * Mesh::GetTransformationFEforElementType ( int  ElemType) [protected]
int Mesh::GetTriOrientation ( const int *  base,
const int *  test 
) [static, protected]

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

Definition at line 2457 of file mesh.cpp.

References mfem_error().

Referenced by AddTriangleFaceElement(), CheckBdrElementOrientation(), GetBdrElementFace(), and Load().

const double* Mesh::GetVertex ( int  i) const [inline]

Return pointer to vertex i's coordinates.

Definition at line 305 of file mesh.hpp.

References vertices.

Referenced by Mesh(), and GridFunction::ProjectCoefficient().

double* Mesh::GetVertex ( int  i) [inline]

Definition at line 306 of file mesh.hpp.

References vertices.

Table * Mesh::GetVertexToElementTable ( )

The returned Table must be destroyed by the caller.

Definition at line 2748 of file mesh.cpp.

References Table::AddAColumnInRow(), Table::AddConnection(), elements, Table::MakeI(), Table::MakeJ(), NumOfElements, NumOfVertices, and Table::ShiftUpI().

void 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 2965 of file mesh.cpp.

References edge_vertex, elements, Table::GetRow(), NumOfElements, DSTable::Push(), and Table::Size().

Referenced by GetEdgeVertexTable(), GetElementToEdgeTable(), GetFaceEdgeTable(), Load(), LocalRefinement(), and MarkTetMeshForRefinement().

void Mesh::GetVertices ( Vector vert_coord) const
void Mesh::GreenRefinement ( int  i,
const DSTable v_to_v,
int *  edge1,
int *  edge2,
int *  middle 
) [inline, protected]

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

Definition at line 121 of file mesh.hpp.

References Bisection().

Referenced by LocalRefinement().

void Mesh::HexUniformRefinement ( ) [protected, virtual]
void Mesh::Init ( ) [protected]
void Mesh::InitTables ( ) [protected]

Definition at line 477 of file mesh.cpp.

References bel_to_edge, edge_vertex, el_to_edge, el_to_el, el_to_face, and face_edge.

Referenced by DeleteTables(), Load(), LoadPatchTopo(), and Mesh().

void Mesh::KnotInsert ( Array< KnotVector * > &  kv)
void Mesh::Load ( istream &  input,
int  generate_edges = 0,
int  refine = 1 
)

Definition at line 1120 of file mesh.cpp.

References Table::AddColumnsInRow(), Table::AddConnections(), be_to_edge, be_to_face, boundary, Ordering::byVDIM, c_el_to_edge, CheckBdrElementOrientation(), CheckElementOrientation(), Geometry::CUBE, Array< T >::DeleteAll(), DeleteTables(), Vector::Destroy(), Dim, NURBSExtension::Dimension(), FiniteElementCollection::DofForGeometry(), FiniteElementCollection::DofOrderForOrientation(), FiniteElementSpace::DofsToVDofs(), el_to_edge, elements, faces, faces_info, FiniteElementSpace::FEColl(), GridFunction::FESpace(), FreeElement(), GenerateBoundaryElements(), GenerateFaces(), NURBSExtension::GetBdrElementTopo(), GetElementBaseGeometry(), FiniteElementSpace::GetElementDofs(), GetElementToEdgeTable(), GetElementToFaceTable(), NURBSExtension::GetElementTopo(), GetFaceBaseGeometry(), NURBSExtension::GetNBE(), NURBSExtension::GetNE(), GridFunction::GetNodalValues(), NURBSExtension::GetNV(), NURBSExtension::GetOrder(), GetQuadOrientation(), Table::GetRow(), GetTriOrientation(), GetVertexToVertexTable(), NURBSExtension::HavePatches(), Element::HEXAHEDRON, InitTables(), Table::MakeI(), Table::MakeJ(), GridFunction::MakeOwner(), MarkForRefinement(), meshgen, mfem_error(), Nodes, NumOfBdrElements, NumOfEdges, NumOfElements, NumOfFaces, NumOfVertices, NURBSext, NURBSExtension, own_nodes, Geometry::POINT, Element::QUADRILATERAL, ReadElement(), Table::RowSize(), Geometry::SEGMENT, Element::SEGMENT, Element::SetAttribute(), SetAttributes(), NURBSExtension::SetCoordsFromPatches(), Vector::SetSize(), Array< T >::SetSize(), Tetrahedron::SetVertices(), Table::ShiftUpI(), Vector::Size(), Array< T >::Size(), skip_comment_lines(), Geometry::SQUARE, Element::TETRAHEDRON, Geometry::TETRAHEDRON, Element::TRIANGLE, Geometry::TRIANGLE, GridFunction::VectorDim(), and vertices.

Referenced by Mesh().

void Mesh::LoadPatchTopo ( istream &  input,
Array< int > &  edge_to_knot 
) [protected]
void Mesh::LocalRefinement ( const Array< int > &  marked_el,
int  type = 3 
) [virtual]
void Mesh::MarkForRefinement ( ) [protected]

Definition at line 740 of file mesh.cpp.

References Dim, MarkTetMeshForRefinement(), MarkTriMeshForRefinement(), and meshgen.

Referenced by Load().

void Mesh::MarkTetMeshForRefinement ( ) [protected]
void Mesh::MarkTriMeshForRefinement ( ) [protected]

Definition at line 751 of file mesh.cpp.

References elements, GetPointMatrix(), NumOfElements, and Element::TRIANGLE.

Referenced by FinalizeTriMesh(), MarkForRefinement(), and Mesh().

int Mesh::MeshGenerator ( ) [inline]

Truegrid or NetGen?

Definition at line 279 of file mesh.hpp.

References meshgen.

Referenced by Mesh().

void Mesh::MoveNodes ( const Vector displacements)

Definition at line 4069 of file mesh.cpp.

References MoveVertices(), and Nodes.

void Mesh::MoveVertices ( const Vector displacements)

Definition at line 4046 of file mesh.cpp.

References Dim, Array< T >::Size(), and vertices.

Referenced by MoveNodes().

Element * Mesh::NewElement ( int  geom) [protected]
void Mesh::NewNodes ( GridFunction nodes)

Definition at line 4093 of file mesh.cpp.

References Nodes, and own_nodes.

void Mesh::NURBSUniformRefinement ( ) [protected, virtual]

Refine NURBS mesh.

Reimplemented in ParMesh.

Definition at line 2193 of file mesh.cpp.

References NURBSExtension::ConvertToPatches(), Nodes, NURBSext, NURBSExtension::UniformRefinement(), and UpdateNURBS().

Referenced by UniformRefinement().

void Mesh::Print ( ostream &  out = cout) const [virtual]

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

Reimplemented in ParMesh.

Definition at line 6002 of file mesh.cpp.

References boundary, Dim, elements, Nodes, NumOfBdrElements, NumOfElements, NumOfVertices, NURBSext, NURBSExtension::Print(), PrintElement(), GridFunction::Save(), and vertices.

Referenced by main().

void Mesh::PrintCharacteristics ( Vector Vh = NULL,
Vector Vk = NULL 
)
void Mesh::PrintElement ( Element el,
ostream &  out 
) [static, protected]
void Mesh::PrintElementsWithPartitioning ( int *  partitioning,
ostream &  out,
int  interior_faces = 0 
)
void Mesh::PrintTopo ( ostream &  out,
const Array< int > &  e_to_k 
) const [protected]
void Mesh::PrintVTK ( 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 6238 of file mesh.cpp.

References attributes, Geometry::CUBE, Geometries, GetAttribute(), GetElementBaseGeometry(), GetElementColoring(), GetElementTransformation(), GetNE(), IntegrationRule::GetNPoints(), Geometry::GetVertices(), GlobGeometryRefiner, DenseMatrix::Height(), RefinedGeometry::RefGeoms, GeometryRefiner::Refine(), RefinedGeometry::RefPts, Array< T >::Size(), Geometry::SQUARE, Geometry::TETRAHEDRON, Geometry::TRIANGLE, and DenseMatrix::Width().

void Mesh::PrintVTK ( ostream &  out)
void Mesh::PrintWithPartitioning ( int *  partitioning,
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.

Definition at line 6436 of file mesh.cpp.

References Dim, elements, faces, faces_info, GetAttribute(), Nodes, NumOfElements, NumOfFaces, NumOfVertices, GridFunction::Save(), and vertices.

void Mesh::PrintXG ( ostream &  out = cout) const [virtual]

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

Reimplemented in ParMesh.

Definition at line 5859 of file mesh.cpp.

References boundary, Dim, elements, meshgen, mfem_error(), Nodes, NumOfBdrElements, NumOfElements, NumOfVertices, GridFunction::Save(), Array< T >::Size(), and vertices.

void Mesh::QuadUniformRefinement ( ) [protected, virtual]
Element * Mesh::ReadElement ( istream &  input) [protected]

Definition at line 1074 of file mesh.cpp.

References Element::GetNVertices(), Element::GetVertices(), NewElement(), and Element::SetAttribute().

Referenced by Load(), and LoadPatchTopo().

void Mesh::RedRefinement ( int  i,
const DSTable v_to_v,
int *  edge1,
int *  edge2,
int *  middle 
) [inline, protected]

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

Definition at line 115 of file mesh.hpp.

References UniformRefinement().

Referenced by LocalRefinement().

void 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 ParMesh.

Definition at line 3291 of file mesh.cpp.

References be_to_edge, boundary, DeleteCoarseTables(), Dim, el_to_edge, elements, GenerateFaces(), GetBdrElementType(), GetElementToEdgeTable(), GetElementToFaceTable(), GetElementType(), meshgen, NumOfBdrElements, NumOfEdges, NumOfElements, Rotate3(), ShiftL2R(), Element::TETRAHEDRON, and Element::TRIANGLE.

void Mesh::Rotate3 ( int &  a,
int &  b,
int &  c 
) [inline, static, protected]

Definition at line 576 of file mesh.hpp.

References ShiftL2R().

Referenced by ReorientTetMesh().

void Mesh::ScaleElements ( double  sf)
void Mesh::ScaleSubdomains ( double  sf)
void Mesh::SetAttributes ( )
void Mesh::SetNodalFESpace ( FiniteElementSpace nfes)

Definition at line 2365 of file mesh.cpp.

References Dim, Nodes, own_nodes, GridFunction::ProjectCoefficient(), and XYZ_VectorFunction().

Referenced by main().

void Mesh::SetNodalGridFunction ( GridFunction nodes)
void Mesh::SetNodes ( const Vector node_coord)

Definition at line 4085 of file mesh.cpp.

References Nodes, and SetVertices().

void Mesh::SetState ( int  s)
void Mesh::SetVertices ( const Vector vert_coord)

Definition at line 4062 of file mesh.cpp.

References Dim, Array< T >::Size(), and vertices.

Referenced by SetNodes().

void Mesh::ShiftL2R ( int &  a,
int &  b,
int &  c 
) [inline, static, protected]

Definition at line 570 of file mesh.hpp.

Referenced by ReorientTetMesh(), and Rotate3().

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

Definition at line 6977 of file mesh.cpp.

References Dim, GridFunction::FESpace(), Nodes, Vector::SetData(), Array< T >::Size(), and vertices.

void Mesh::UniformRefinement ( int  i,
const DSTable v_to_v,
int *  edge1,
int *  edge2,
int *  middle 
) [protected]
void Mesh::UniformRefinement ( )
void Mesh::UpdateNodes ( ) [protected]
void Mesh::UpdateNURBS ( ) [protected]
void 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 496 of file mesh.hpp.

References NORMAL, SetState(), State, and WantTwoLevelState.

Referenced by HexUniformRefinement(), LocalRefinement(), and QuadUniformRefinement().


Friends And Related Function Documentation

friend class NURBSExtension [friend]

Definition at line 30 of file mesh.hpp.

Referenced by Load(), and Mesh().


Member Data Documentation

Definition at line 224 of file mesh.hpp.

Referenced by main(), Mesh(), PrintVTK(), ScaleSubdomains(), and SetAttributes().

Definition at line 225 of file mesh.hpp.

Referenced by main(), Mesh(), and SetAttributes().

Array<int> Mesh::be_to_edge [protected]
Array<int> Mesh::be_to_face [protected]
Table* Mesh::bel_to_edge [protected]
Array<Element *> Mesh::boundary [protected]
Table * Mesh::c_bel_to_edge [protected]

Definition at line 65 of file mesh.hpp.

Referenced by HexUniformRefinement(), LocalRefinement(), and SetState().

Table* Mesh::c_el_to_edge [protected]

Definition at line 65 of file mesh.hpp.

Referenced by HexUniformRefinement(), Load(), LocalRefinement(), QuadUniformRefinement(), and SetState().

Table* Mesh::c_el_to_face [protected]

Definition at line 67 of file mesh.hpp.

Referenced by HexUniformRefinement(), LocalRefinement(), and SetState().

int Mesh::c_NumOfBdrElements [protected]

Definition at line 43 of file mesh.hpp.

Referenced by HexUniformRefinement(), LocalRefinement(), QuadUniformRefinement(), and SetState().

int Mesh::c_NumOfEdges [protected]

Definition at line 45 of file mesh.hpp.

Referenced by HexUniformRefinement(), LocalRefinement(), QuadUniformRefinement(), and SetState().

int Mesh::c_NumOfElements [protected]

Definition at line 43 of file mesh.hpp.

Referenced by HexUniformRefinement(), LocalRefinement(), QuadUniformRefinement(), and SetState().

int Mesh::c_NumOfFaces [protected]

Definition at line 45 of file mesh.hpp.

Referenced by HexUniformRefinement(), LocalRefinement(), and SetState().

int Mesh::c_NumOfVertices [protected]

Definition at line 43 of file mesh.hpp.

Referenced by HexUniformRefinement(), LocalRefinement(), QuadUniformRefinement(), and SetState().

int Mesh::Dim [protected]
Table* Mesh::edge_vertex [mutable, protected]
Table* Mesh::el_to_edge [protected]
Table* Mesh::el_to_el [protected]
Table* Mesh::el_to_face [protected]
Array<Element *> Mesh::elements [protected]
Table * Mesh::f_bel_to_edge [protected]

Definition at line 65 of file mesh.hpp.

Referenced by HexUniformRefinement(), LocalRefinement(), and SetState().

Table * Mesh::f_el_to_edge [protected]

Definition at line 65 of file mesh.hpp.

Referenced by HexUniformRefinement(), LocalRefinement(), QuadUniformRefinement(), and SetState().

Table * Mesh::f_el_to_face [protected]

Definition at line 67 of file mesh.hpp.

Referenced by HexUniformRefinement(), LocalRefinement(), and SetState().

int Mesh::f_NumOfBdrElements [protected]

Definition at line 44 of file mesh.hpp.

Referenced by HexUniformRefinement(), LocalRefinement(), QuadUniformRefinement(), and SetState().

int Mesh::f_NumOfEdges [protected]

Definition at line 46 of file mesh.hpp.

Referenced by HexUniformRefinement(), LocalRefinement(), QuadUniformRefinement(), and SetState().

int Mesh::f_NumOfElements [protected]

Definition at line 44 of file mesh.hpp.

Referenced by HexUniformRefinement(), LocalRefinement(), QuadUniformRefinement(), and SetState().

int Mesh::f_NumOfFaces [protected]

Definition at line 46 of file mesh.hpp.

Referenced by HexUniformRefinement(), LocalRefinement(), and SetState().

int Mesh::f_NumOfVertices [protected]

Definition at line 44 of file mesh.hpp.

Referenced by HexUniformRefinement(), LocalRefinement(), QuadUniformRefinement(), and SetState().

Table* Mesh::face_edge [mutable, protected]

Definition at line 62 of file mesh.hpp.

Referenced by DeleteCoarseTables(), DeleteTables(), GetFaceEdges(), GetFaceEdgeTable(), and InitTables().

Definition at line 72 of file mesh.hpp.

Referenced by GetFaceElementTransformations().

Array<Element *> Mesh::faces [protected]

Definition at line 71 of file mesh.hpp.

Referenced by GetBdrElementTransformation(), and GetFaceTransformation().

Array<int> Mesh::fc_be_to_edge [protected]

Definition at line 66 of file mesh.hpp.

Referenced by LocalRefinement(), QuadUniformRefinement(), and SetState().

Definition at line 68 of file mesh.hpp.

Referenced by HexUniformRefinement(), LocalRefinement(), and SetState().

int Mesh::meshgen [protected]
GridFunction* Mesh::Nodes [protected]
int Mesh::NumOfBdrElements [protected]
int Mesh::NumOfEdges [protected]
int Mesh::NumOfElements [protected]
int Mesh::NumOfFaces [protected]
int Mesh::NumOfVertices [protected]
int Mesh::own_nodes [protected]

Definition at line 75 of file mesh.hpp.

Referenced by Init(), Load(), Mesh(), NewNodes(), SetNodalFESpace(), SetNodalGridFunction(), and ~Mesh().

int Mesh::State [protected]

Definition at line 70 of file mesh.hpp.

Referenced by GetFaceElementTransformations().

Array<Vertex> Mesh::vertices [protected]
int Mesh::WantTwoLevelState [protected]

The documentation for this class was generated from the following files:
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Friends Defines