Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
MeshElements.h
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: MeshElements.h
4 //
5 // For more information, please see: http://www.nektar.info/
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Mesh manipulation objects.
33 //
34 ////////////////////////////////////////////////////////////////////////////////
35 
36 #ifndef UTILITY_PREPROCESSING_MESHCONVERT_MESHELEMENTS
37 #define UTILITY_PREPROCESSING_MESHCONVERT_MESHELEMENTS
38 
39 #include <vector>
40 #include <string>
41 #include <iostream>
42 #include <iomanip>
43 
44 #include <boost/shared_ptr.hpp>
45 #include <boost/unordered_set.hpp>
46 #include <boost/unordered_map.hpp>
47 
51 
52 #include <SpatialDomains/SegGeom.h>
53 #include <SpatialDomains/TriGeom.h>
55 #include <SpatialDomains/Curve.hpp>
57 
58 namespace Nektar
59 {
60  namespace Utilities
61  {
62  // Forwards declaration for Element class.
63  class Element;
64  /// Shared pointer to an element.
65  typedef boost::shared_ptr<Element> ElementSharedPtr;
66 
67  /**
68  * @brief Represents a point in the domain.
69  *
70  * Such points may either be element vertices, or simply control
71  * points on high-order edges/faces, although this information is not
72  * contained within this class.
73  */
74  class Node {
75  public:
76  /// Create a new node at a specified coordinate.
77  Node(int pId, NekDouble pX, NekDouble pY, NekDouble pZ)
78  : m_id(pId), m_x(pX), m_y(pY), m_z(pZ), m_geom() {}
79  /// Copy an existing node.
80  Node(const Node& pSrc)
81  : m_id(pSrc.m_id), m_x(pSrc.m_x), m_y(pSrc.m_y),
82  m_z(pSrc.m_z), m_geom() {}
83  Node() : m_id(0), m_x(0.0), m_y(0.0), m_z(0.0), m_geom() {}
84  ~Node() {}
85 
86  /// Reset the local id;
87  void SetID(int pId)
88  {
89  m_id = pId;
90  }
91 
92  /// Get the local id;
93  int GetID(void)
94  {
95  return m_id;
96  }
97 
98  /// Define node ordering based on ID.
99  bool operator<(const Node& pSrc)
100  {
101  return (m_id < pSrc.m_id);
102  }
103  /// Define node equality based on coordinate.
104  bool operator==(const Node& pSrc)
105  {
106  return m_x == pSrc.m_x && m_y == pSrc.m_y && m_z == pSrc.m_z;
107  }
108 
109  Node operator+(const Node &pSrc) const
110  {
111  return Node(m_id, m_x+pSrc.m_x, m_y+pSrc.m_y, m_z+pSrc.m_z);
112  }
113 
114  Node operator-(const Node &pSrc) const
115  {
116  return Node(m_id, m_x-pSrc.m_x, m_y-pSrc.m_y, m_z-pSrc.m_z);
117  }
118 
119  Node operator*(const Node &pSrc) const
120  {
121  return Node(m_id, m_x*pSrc.m_x, m_y*pSrc.m_y, m_z*pSrc.m_z);
122  }
123 
124  Node operator*(const NekDouble &alpha) const
125  {
126  return Node(m_id, alpha*m_x, alpha*m_y, alpha*m_z);
127  }
128 
129  Node operator/(const NekDouble &alpha) const
130  {
131  return Node(m_id, m_x/alpha, m_y/alpha, m_z/alpha);
132  }
133 
134  void operator+=(const Node &pSrc)
135  {
136  m_x += pSrc.m_x;
137  m_y += pSrc.m_y;
138  m_z += pSrc.m_z;
139  }
140 
141  void operator*=(const NekDouble &alpha)
142  {
143  m_x *= alpha;
144  m_y *= alpha;
145  m_z *= alpha;
146  }
147 
148  void operator/=(const NekDouble &alpha)
149  {
150  m_x /= alpha;
151  m_y /= alpha;
152  m_z /= alpha;
153  }
154 
155  NekDouble abs2() const
156  {
157  return m_x*m_x+m_y*m_y+m_z*m_z;
158  }
159 
160  NekDouble dot(const Node &pSrc) const
161  {
162  return m_x*pSrc.m_x + m_y*pSrc.m_y + m_z*pSrc.m_z;
163  }
164 
165 
166  Node curl(const Node &pSrc) const
167  {
168  return Node(m_id, m_y*pSrc.m_z - m_z*pSrc.m_y,
169  m_z*pSrc.m_x-m_x*pSrc.m_z, m_x*pSrc.m_y-m_y*pSrc.m_x);
170  }
171 
172  /// Generate a %SpatialDomains::PointGeom for this node.
174  {
178 
179  return ret;
180  }
181 
182  /// ID of node.
183  int m_id;
184  /// X-coordinate.
186  /// Y-coordinate.
188  /// Z-coordinate.
190 
191  private:
193  };
194  /// Shared pointer to a Node.
195  typedef boost::shared_ptr<Node> NodeSharedPtr;
196 
197  bool operator==(NodeSharedPtr const &p1, NodeSharedPtr const &p2);
198  bool operator< (NodeSharedPtr const &p1, NodeSharedPtr const &p2);
199  std::ostream &operator<<(std::ostream &os, const NodeSharedPtr &n);
200 
201  /**
202  * @brief Defines a hash function for nodes.
203  *
204  * The hash of a node is straight-forward; a combination of the x, y,
205  * and z co-ordinates in this order.
206  */
207  struct NodeHash : std::unary_function<NodeSharedPtr, std::size_t>
208  {
209  std::size_t operator()(NodeSharedPtr const& p) const
210  {
211  std::size_t seed = 0;
212  boost::hash_combine(seed, p->m_x);
213  boost::hash_combine(seed, p->m_y);
214  boost::hash_combine(seed, p->m_z);
215  return seed;
216  }
217  };
218  typedef boost::unordered_set<NodeSharedPtr, NodeHash> NodeSet;
219 
220 
221  /**
222  * @brief Represents an edge which joins two points.
223  *
224  * An edge is defined by two nodes (vertices) and, for high-order edges,
225  * a set of control nodes defining the shape of the edge.
226  */
227  class Edge {
228  public:
229  /// Creates a new edge.
230  Edge(NodeSharedPtr pVertex1, NodeSharedPtr pVertex2,
231  std::vector<NodeSharedPtr> pEdgeNodes,
232  LibUtilities::PointsType pCurveType)
233  : m_n1(pVertex1), m_n2(pVertex2), m_edgeNodes(pEdgeNodes),
234  m_curveType(pCurveType), m_geom() {}
235  /// Copies an existing edge.
236  Edge(const Edge& pSrc)
237  : m_n1(pSrc.m_n1), m_n2(pSrc.m_n2), m_edgeNodes(pSrc.m_edgeNodes),
238  m_curveType(pSrc.m_curveType), m_geom(pSrc.m_geom) {}
239  ~Edge() {}
240 
241  /// Returns the total number of nodes defining the edge.
242  unsigned int GetNodeCount() const
243  {
244  return m_edgeNodes.size() + 2;
245  }
246 
247  /// Creates a Nektar++ string listing the coordinates of all the
248  /// nodes.
249  std::string GetXmlCurveString() const
250  {
251  std::stringstream s;
252  std::string str;
253  s << std::scientific << std::setprecision(8) << " "
254  << m_n1->m_x << " " << m_n1->m_y << " " << m_n1->m_z << " ";
255  for (int k = 0; k < m_edgeNodes.size(); ++k) {
256  s << std::scientific << std::setprecision(8) << " "
257  << m_edgeNodes[k]->m_x << " " << m_edgeNodes[k]->m_y
258  << " " << m_edgeNodes[k]->m_z << " ";
259  }
260  s << std::scientific << std::setprecision(8) << " "
261  << m_n2->m_x << " " << m_n2->m_y << " " << m_n2->m_z;
262  return s.str();
263  }
264 
265  /// Generate a SpatialDomains::SegGeom object for this edge.
267  {
268  // Create edge vertices.
271 
272  p[0] = m_n1->GetGeom(coordDim);
273  p[1] = m_n2->GetGeom(coordDim);
274 
275  // Create a curve if high-order information exists.
276  if (m_edgeNodes.size() > 0)
277  {
281 
282  c->m_points.push_back(p[0]);
283  for (int i = 0; i < m_edgeNodes.size(); ++i)
284  {
285  c->m_points.push_back(m_edgeNodes[i]->GetGeom(coordDim));
286  }
287  c->m_points.push_back(p[1]);
288 
290  AllocateSharedPtr(m_id, coordDim, p, c);
291  }
292  else
293  {
295  AllocateSharedPtr(m_id, coordDim, p);
296  }
297 
298  return ret;
299  }
300 
301  /// ID of edge.
302  unsigned int m_id;
303  /// First vertex node.
305  /// Second vertex node.
307  /// List of control nodes between the first and second vertices.
308  std::vector<NodeSharedPtr> m_edgeNodes;
309  /// Distributions of points along edge.
311  /// Element(s) which are linked to this edge.
312  vector<pair<ElementSharedPtr, int> > m_elLink;
313 
314  private:
316  };
317  /// Shared pointer to an edge.
318  typedef boost::shared_ptr<Edge> EdgeSharedPtr;
319 
320  bool operator==(EdgeSharedPtr const &p1, EdgeSharedPtr const &p2);
321  bool operator< (EdgeSharedPtr const &p1, EdgeSharedPtr const &p2);
322 
323  /**
324  * @brief Defines a hash function for edges.
325  *
326  * The hash of an edge is defined using the IDs of the two nodes which
327  * define it. First the minimum ID is hashed, then the maximum
328  * ID, which takes the two possible orientations into account.
329  */
330  struct EdgeHash : std::unary_function<EdgeSharedPtr, std::size_t>
331  {
332  std::size_t operator()(EdgeSharedPtr const& p) const
333  {
334  std::size_t seed = 0;
335  unsigned int id1 = p->m_n1->m_id;
336  unsigned int id2 = p->m_n2->m_id;
337  boost::hash_combine(seed, id1 < id2 ? id1 : id2);
338  boost::hash_combine(seed, id2 < id1 ? id1 : id2);
339  return seed;
340  }
341  };
342  typedef boost::unordered_set<EdgeSharedPtr, EdgeHash> EdgeSet;
343 
344 
345  /**
346  * @brief Represents a face comprised of three or more edges.
347  *
348  * A face is defined by a list of vertices, a list of edges joining
349  * these vertices, and a list of control nodes within the interior of
350  * the face, defining the shape of the face.
351  */
352  class Face {
353  public:
354  /// Create a new face.
355  Face(std::vector<NodeSharedPtr> pVertexList,
356  std::vector<NodeSharedPtr> pFaceNodes,
357  std::vector<EdgeSharedPtr> pEdgeList,
358  LibUtilities::PointsType pCurveType)
359  : m_vertexList(pVertexList),
360  m_edgeList (pEdgeList),
361  m_faceNodes (pFaceNodes),
362  m_curveType (pCurveType),
363  m_geom () {}
364 
365  /// Copy an existing face.
366  Face(const Face& pSrc)
369  m_geom (pSrc.m_geom) {}
370  ~Face() {}
371 
372  /// Equality is defined by matching all vertices.
373  bool operator==(Face& pSrc)
374  {
376  for (it1 = m_vertexList.begin(); it1 != m_vertexList.end(); ++it1)
377  {
378  if (find(pSrc.m_vertexList.begin(), pSrc.m_vertexList.end(), *it1)
379  == pSrc.m_vertexList.end())
380  {
381  return false;
382  }
383  }
384  return true;
385  }
386 
387  /// Returns the total number of nodes (vertices, edge nodes and
388  /// face nodes).
389  unsigned int GetNodeCount() const
390  {
391  unsigned int n = m_faceNodes.size();
392  for (int i = 0; i < m_edgeList.size(); ++i)
393  {
394  n += m_edgeList[i]->GetNodeCount();
395  }
396  n -= m_vertexList.size();
397  return n;
398  }
399 
400  /// Generates a string listing the coordinates of all nodes
401  /// associated with this face.
402  std::string GetXmlCurveString() const
403  {
404  std::stringstream s;
405  std::string str;
406 
407  // Treat 2D point distributions differently to 3D.
411  {
412  vector<NodeSharedPtr> tmp;
413  int n = m_edgeList[0]->GetNodeCount();
414 
415  tmp.insert(tmp.end(), m_vertexList.begin(), m_vertexList.end());
416  for (int k = 0; k < m_edgeList.size(); ++k)
417  {
418  tmp.insert(tmp.end(), m_edgeList[k]->m_edgeNodes.begin(),
419  m_edgeList[k]->m_edgeNodes.end());
420  if (m_edgeList[k]->m_n1 != m_vertexList[k])
421  {
422  // If edge orientation is reversed relative to node
423  // ordering, we need to reverse order of nodes.
424  std::reverse(tmp.begin() + 3 + k*(n-2),
425  tmp.begin() + 3 + (k+1)*(n-2));
426  }
427  }
428  tmp.insert(tmp.end(), m_faceNodes.begin(), m_faceNodes.end());
429 
430  for (int k = 0; k < tmp.size(); ++k) {
431  s << std::scientific << std::setprecision(8) << " "
432  << tmp[k]->m_x << " " << tmp[k]->m_y
433  << " " << tmp[k]->m_z << " ";
434  }
435 
436  return s.str();
437  }
438  else
439  {
440  // Write out in 2D tensor product order.
441  ASSERTL0(m_vertexList.size() == 4,
442  "Face nodes of tensor product only supported "
443  "for quadrilaterals.");
444 
445  int n = (int)sqrt((NekDouble)GetNodeCount());
446  vector<NodeSharedPtr> tmp(n*n);
447 
448  ASSERTL0(n*n == GetNodeCount(), "Wrong number of modes?");
449 
450  // Write vertices
451  tmp[0] = m_vertexList[0];
452  tmp[n-1] = m_vertexList[1];
453  tmp[n*n-1] = m_vertexList[2];
454  tmp[n*(n-1)] = m_vertexList[3];
455 
456  // Write edge-interior
457  int skips[4][2] = {{0,1}, {n-1,n}, {n*n-1,-1}, {n*(n-1),-n}};
458  for (int i = 0; i < 4; ++i)
459  {
460  bool reverseEdge = m_edgeList[i]->m_n1 == m_vertexList[i];
461 
462  if (!reverseEdge)
463  {
464  for (int j = 1; j < n-1; ++j)
465  {
466  tmp[skips[i][0] + j*skips[i][1]] =
467  m_edgeList[i]->m_edgeNodes[n-2-j];
468  }
469  }
470  else
471  {
472  for (int j = 1; j < n-1; ++j)
473  {
474  tmp[skips[i][0] + j*skips[i][1]] =
475  m_edgeList[i]->m_edgeNodes[j-1];
476  }
477  }
478  }
479 
480  // Write interior
481  for (int i = 1; i < n-1; ++i)
482  {
483  for (int j = 1; j < n-1; ++j)
484  {
485  tmp[i*n+j] = m_faceNodes[(i-1)*(n-2)+(j-1)];
486  }
487  }
488 
489  for (int k = 0; k < tmp.size(); ++k) {
490  s << std::scientific << std::setprecision(8) << " "
491  << tmp[k]->m_x << " " << tmp[k]->m_y
492  << " " << tmp[k]->m_z << " ";
493  }
494 
495  return s.str();
496  }
497  }
498 
499  /// Generate either SpatialDomains::TriGeom or
500  /// SpatialDomains::QuadGeom for this element.
502  {
503  int nEdge = m_edgeList.size();
504 
507  StdRegions::Orientation edgeo[4];
508 
509  for (int i = 0; i < nEdge; ++i)
510  {
511  edges[i] = m_edgeList[i]->GetGeom(coordDim);
512  }
513 
514  for (int i = 0; i < nEdge; ++i)
515  {
517  *edges[i], *edges[(i+1) % nEdge]);
518  }
519 
520  if (nEdge == 3)
521  {
523  AllocateSharedPtr(m_id, edges, edgeo);
524  }
525  else
526  {
528  AllocateSharedPtr(m_id, edges, edgeo);
529  }
530 
531  return ret;
532  }
533 
534  /// ID of the face.
535  unsigned int m_id;
536  /// List of vertex nodes.
537  std::vector<NodeSharedPtr> m_vertexList;
538  /// List of corresponding edges.
539  std::vector<EdgeSharedPtr> m_edgeList;
540  /// List of face-interior nodes defining the shape of the face.
541  std::vector<NodeSharedPtr> m_faceNodes;
542  /// Distribution of points in this face.
544  /// Element(s) which are linked to this face.
545  vector<pair<ElementSharedPtr, int> > m_elLink;
546 
548  };
549  /// Shared pointer to a face.
550  typedef boost::shared_ptr<Face> FaceSharedPtr;
551 
552  bool operator==(FaceSharedPtr const &p1, FaceSharedPtr const &p2);
553  bool operator< (FaceSharedPtr const &p1, FaceSharedPtr const &p2);
554 
555  struct FaceHash : std::unary_function<FaceSharedPtr, std::size_t>
556  {
557  std::size_t operator()(FaceSharedPtr const& p) const
558  {
559  unsigned int nVert = p->m_vertexList.size();
560  std::size_t seed = 0;
561  std::vector<unsigned int> ids(nVert);
562 
563  for (int i = 0; i < nVert; ++i)
564  {
565  ids[i] = p->m_vertexList[i]->m_id;
566  }
567 
568  std::sort(ids.begin(), ids.end());
569  boost::hash_range(seed, ids.begin(), ids.end());
570 
571  return seed;
572  }
573  };
574  typedef boost::unordered_set<FaceSharedPtr, FaceHash> FaceSet;
575 
576 
577  /**
578  * @brief Basic information about an element.
579  *
580  * ElmtConfig contains four member variables which denote the
581  * properties of an element when it is created.
582  */
583  struct ElmtConfig
584  {
586  unsigned int pOrder,
587  bool pFn,
588  bool pVn,
589  bool pReorient = true,
594  : m_e (pE),
595  m_faceNodes (pFn),
596  m_volumeNodes (pVn),
597  m_order (pOrder),
598  m_reorient (pReorient),
599  m_edgeCurveType(pECt),
600  m_faceCurveType(pFCt)
601  {
602  }
603 
604  ElmtConfig(ElmtConfig const &p) :
605  m_e (p.m_e),
607  m_volumeNodes (p.m_volumeNodes),
608  m_order (p.m_order),
609  m_reorient (p.m_reorient),
610  m_edgeCurveType(p.m_edgeCurveType),
611  m_faceCurveType(p.m_faceCurveType)
612  {
613  }
614 
616 
617  /// Element type (e.g. triangle, quad, etc).
619  /// Denotes whether the element contains face nodes. For 2D
620  /// elements, if this is true then the element contains interior
621  /// nodes.
623  /// Denotes whether the element contains volume (i.e. interior)
624  /// nodes. These are not supported by either the mesh converter or
625  /// Nektar++ but are included for completeness and are required
626  /// for some output modules (e.g. Gmsh).
628  /// Order of the element.
629  unsigned int m_order;
630  /// Denotes whether the element needs to be re-orientated for a
631  /// spectral element framework.
633  /// Distribution of points in edges.
635  /// Distribution of points in faces.
637  };
638 
639 
640  /**
641  * @brief Base class for element definitions.
642  *
643  * An element is defined by a list of vertices, edges and faces
644  * (depending on the dimension of the problem). This base class
645  * provides the underlying structure.
646  */
647  class Element {
648  public:
649  Element(ElmtConfig pConf,
650  unsigned int pNumNodes,
651  unsigned int pGotNodes);
652 
653  /// Returns the ID of the element (or associated edge or face for
654  /// boundary elements).
655  unsigned int GetId() const {
656  if (m_faceLink.get() != 0) return m_faceLink->m_id;
657  if (m_edgeLink.get() != 0) return m_edgeLink->m_id;
658  return m_id;
659  }
660  /// Returns the expansion dimension of the element.
661  unsigned int GetDim() const {
662  return m_dim;
663  }
664  /// Returns the configuration of the element.
665  ElmtConfig GetConf() const {
666  return m_conf;
667  }
668  /// Returns the tag which defines the element shape.
669  std::string GetTag() const {
670  if (m_faceLink.get() != 0) return "F";
671  if (m_edgeLink.get() != 0) return "E";
672  return m_tag;
673  }
674  /// Access a vertex node.
675  NodeSharedPtr GetVertex(unsigned int i) const {
676  return m_vertex[i];
677  }
678  /// Access an edge.
679  EdgeSharedPtr GetEdge(unsigned int i) const {
680  return m_edge[i];
681  }
682  /// Access a face.
683  FaceSharedPtr GetFace(unsigned int i) const {
684  return m_face[i];
685  }
686  /// Access the list of vertex nodes.
687  std::vector<NodeSharedPtr> GetVertexList() const {
688  return m_vertex;
689  }
690  /// Access the list of edges.
691  std::vector<EdgeSharedPtr> GetEdgeList() const {
692  return m_edge;
693  }
694  /// Access the list of faces.
695  std::vector<FaceSharedPtr> GetFaceList() const {
696  return m_face;
697  }
698  /// Access the list of volume nodes.
699  std::vector<NodeSharedPtr> GetVolumeNodes() const {
700  return m_volumeNodes;
701  }
702  void SetVolumeNodes(std::vector<NodeSharedPtr> &nodes) {
703  m_volumeNodes = nodes;
704  }
706  return m_curveType;
707  }
709  m_curveType = cT;
710  }
711  /// Returns the total number of nodes (vertices, edge nodes and
712  /// face nodes and volume nodes).
713  unsigned int GetNodeCount() const
714  {
715  unsigned int n = m_volumeNodes.size();
716  if (m_dim == 2)
717  {
718  for (int i = 0; i < m_edge.size(); ++i)
719  {
720  n += m_edge[i]->GetNodeCount();
721  }
722  n -= m_vertex.size();
723  }
724  else
725  {
726  cerr << "Not supported." << endl;
727  exit(1);
728  }
729  return n;
730  }
731  /// Access the list of tags associated with this element.
732  std::vector<int> GetTagList() const {
733  return m_taglist;
734  }
735  /// Returns the number of vertices.
736  unsigned int GetVertexCount() const {
737  return m_vertex.size();
738  }
739  /// Returns the number of edges.
740  unsigned int GetEdgeCount() const {
741  return m_edge.size();
742  }
743  /// Returns the number of faces.
744  unsigned int GetFaceCount() const {
745  return m_face.size();
746  }
747  /// Change the ID of the element.
748  void SetId(unsigned int p) {
749  m_id = p;
750  }
751  /// Replace a vertex with another vertex object.
752  void SetVertex(unsigned int p, NodeSharedPtr pNew);
753  /// Replace an edge with another edge object.
754  void SetEdge(unsigned int p, EdgeSharedPtr pNew);
755  /// Replace a face with another face object.
756  void SetFace(unsigned int p, FaceSharedPtr pNew);
757  /// Set a correspondence between this element and an edge
758  /// (2D boundary element).
760  m_edgeLink = pLink;
761  }
762  /// Get correspondence between this element and an edge.
764  return m_edgeLink;
765  }
766  /// Set a correspondence between this element and a face
767  /// (3D boundary element).
768  void SetFaceLink(FaceSharedPtr pLink) {
769  m_faceLink = pLink;
770  }
771  /// Get correspondence between this element and a face.
772  FaceSharedPtr GetFaceLink() {
773  return m_faceLink;
774  }
775  /// Set a correspondence between edge or face i and its
776  /// representative boundary element m->element[expDim-1][j].
777  void SetBoundaryLink(int i, int j) {
778  m_boundaryLinks[i] = j;
779  }
780  /// Get the location of the boundary face/edge i for this element.
781  int GetBoundaryLink(int i) {
782  std::map<int,int>::iterator it = m_boundaryLinks.find(i);
783  if (it == m_boundaryLinks.end())
784  {
785  return -1;
786  }
787  else
788  {
789  return it->second;
790  }
791  }
792  /// Set the list of tags associated with this element.
793  void SetTagList(const std::vector<int> &tags)
794  {
795  m_taglist = tags;
796  }
797  /// Generate a list of vertices (1D), edges (2D), or faces (3D).
798  virtual std::string GetXmlString() const
799  {
800  std::stringstream s;
801  switch (m_dim)
802  {
803  case 1:
804  for(int j=0; j< m_vertex.size(); ++j)
805  {
806  s << std::setw(5) << m_vertex[j]->m_id << " ";
807  }
808  break;
809  case 2:
810  for(int j=0; j< m_edge.size(); ++j)
811  {
812  s << std::setw(5) << m_edge[j]->m_id << " ";
813  }
814  break;
815  case 3:
816  for(int j=0; j< m_face.size(); ++j)
817  {
818  s << std::setw(5) << m_face[j]->m_id << " ";
819  }
820  break;
821  }
822  return s.str();
823  }
824 
825  /// Generates a string listing the coordinates of all nodes
826  /// associated with this element.
827  std::string GetXmlCurveString() const
828  {
829  // Temporary node list for reordering
830  std::vector<NodeSharedPtr> nodeList;
831 
832  // Node orderings are different for different elements.
833  // Triangle
834  if (m_vertex.size() == 3)
835  {
836  int n = m_edge[0]->GetNodeCount();
837  nodeList.resize(n*(n+1)/2);
838 
839  // Populate nodelist
840  std::copy(m_vertex.begin(), m_vertex.end(), nodeList.begin());
841  for (int i = 0; i < 3; ++i)
842  {
843  std::copy(m_edge[i]->m_edgeNodes.begin(),
844  m_edge[i]->m_edgeNodes.end(),
845  nodeList.begin() + 3 + i*(n-2));
846  if (m_edge[i]->m_n1 != m_vertex[i])
847  {
848  // If edge orientation is reversed relative to node
849  // ordering, we need to reverse order of nodes.
850  std::reverse(nodeList.begin() + 3 + i*(n-2),
851  nodeList.begin() + 3 + (i+1)*(n-2));
852  }
853  }
854 
855  // Copy volume nodes.
856  std::copy(m_volumeNodes.begin(), m_volumeNodes.end(),
857  nodeList.begin() + 3*(n-1));
858  }
859  // Quad
860  else if (m_dim == 2 && m_vertex.size() == 4)
861  {
862  int n = m_edge[0]->GetNodeCount();
863  nodeList.resize(n*n);
864 
865  // Write vertices
866  nodeList[0] = m_vertex[0];
867  nodeList[n-1] = m_vertex[1];
868  nodeList[n*n-1] = m_vertex[2];
869  nodeList[n*(n-1)] = m_vertex[3];
870 
871  // Write edge-interior
872  int skips[4][2] = {{0,1}, {n-1,n}, {n*n-1,-1}, {n*(n-1),-n}};
873  for (int i = 0; i < 4; ++i)
874  {
875  bool reverseEdge = m_edge[i]->m_n1 == m_vertex[i];
876 
877  if (!reverseEdge)
878  {
879  for (int j = 1; j < n-1; ++j)
880  {
881  nodeList[skips[i][0] + j*skips[i][1]] =
882  m_edge[i]->m_edgeNodes[n-2-j];
883  }
884  }
885  else
886  {
887  for (int j = 1; j < n-1; ++j)
888  {
889  nodeList[skips[i][0] + j*skips[i][1]] =
890  m_edge[i]->m_edgeNodes[j-1];
891  }
892  }
893  }
894 
895  // Write interior
896  for (int i = 1; i < n-1; ++i)
897  {
898  for (int j = 1; j < n-1; ++j)
899  {
900  nodeList[i*n+j] = m_volumeNodes[(i-1)*(n-2)+(j-1)];
901  }
902  }
903  }
904  else
905  {
906  cerr << "GetXmlCurveString for a " << m_vertex.size()
907  << "-vertex element is not yet implemented." << endl;
908  }
909 
910  // Finally generate the XML string corresponding to our new
911  // node reordering.
912  std::stringstream s;
913  std::string str;
914  for (int k = 0; k < nodeList.size(); ++k)
915  {
916  s << std::scientific << std::setprecision(8) << " "
917  << nodeList[k]->m_x << " " << nodeList[k]->m_y
918  << " " << nodeList[k]->m_z << " ";
919 
920  }
921  return s.str();
922  }
923 
924  /// Generate a Nektar++ geometry object for this element.
926  {
927  ASSERTL0(false, "This function should be implemented at a shape level.");
928  return boost::shared_ptr<SpatialDomains::Geometry>();
929  }
930  int GetMaxOrder();
931  /// Complete this object.
932  virtual void Complete(int order)
933  {
934  ASSERTL0(false, "This function should be implemented at a shape level.");
935  }
936  void Print()
937  {
938  int i, j;
939  for (i = 0; i < m_vertex.size(); ++i)
940  {
941  cout << m_vertex[i]->m_x << " " << m_vertex[i]->m_y << " " << m_vertex[i]->m_z << endl;
942  }
943  for (i = 0; i < m_edge.size(); ++i)
944  {
945  for (j = 0; j < m_edge[i]->m_edgeNodes.size(); ++j)
946  {
947  NodeSharedPtr n = m_edge[i]->m_edgeNodes[j];
948  cout << n->m_x << " " << n->m_y << " " << n->m_z << endl;
949  }
950  }
951  for (i = 0; i < m_face.size(); ++i)
952  {
953  for (j = 0; j < m_face[i]->m_faceNodes.size(); ++j)
954  {
955  NodeSharedPtr n = m_face[i]->m_faceNodes[j];
956  cout << n->m_x << " " << n->m_y << " " << n->m_z << endl;
957  }
958  }
959  }
960 
961  protected:
962  /// ID of the element.
963  unsigned int m_id;
964  /// Dimension of the element.
965  unsigned int m_dim;
966  /// Contains configuration of the element.
968  /// Tag character describing the element.
969  std::string m_tag;
970  /// List of integers specifying properties of the element.
971  std::vector<int> m_taglist;
972  /// List of element vertex nodes.
973  std::vector<NodeSharedPtr> m_vertex;
974  /// List of element edges.
975  std::vector<EdgeSharedPtr> m_edge;
976  /// List of element faces.
977  std::vector<FaceSharedPtr> m_face;
978  /// List of element volume nodes.
979  std::vector<NodeSharedPtr> m_volumeNodes;
980  /// Volume curve type
982  /// Pointer to the corresponding edge if element is a 2D boundary.
984  /// Pointer to the corresponding face if element is a 3D boundary.
985  FaceSharedPtr m_faceLink;
986  /// Array mapping faces/edges to the location of the appropriate
987  /// boundary elements in m->element.
988  std::map<int,int> m_boundaryLinks;
989  /// Nektar++ geometry object for this element.
991  };
992  /// Container for elements; key is expansion dimension, value is
993  /// vector of elements of that dimension.
994  typedef std::map<unsigned int, std::vector<ElementSharedPtr> > ElementMap;
995  /// Element factory definition.
997  ElmtConfig, std::vector<NodeSharedPtr>, std::vector<int> > ElementFactory;
998  ElementFactory& GetElementFactory();
999 
1000  /// Define element ordering based on ID.
1002  {
1003  typedef boost::shared_ptr<Element> pT;
1004  bool operator()(const pT a, const pT b) const
1005  {
1006  // check for 0
1007  if (a.get() == 0)
1008  {
1009  // if b is also 0, then they are equal, hence a is not
1010  // less than b
1011  return b.get() != 0;
1012  }
1013  else if (b.get() == 0)
1014  {
1015  return false;
1016  }
1017  else
1018  {
1019  return a->GetId() < b->GetId();
1020  }
1021  }
1022  };
1023 
1024 
1025  /**
1026  * @brief A composite is a collection of elements.
1027  *
1028  * All elements should be of the same type, i.e. have the same tag.
1029  */
1030  class Composite {
1031  public:
1032  Composite() : m_reorder(true) {}
1033 
1034  /// Generate the list of IDs of elements within this composite.
1035  std::string GetXmlString(bool doSort=true);
1036 
1037  /// ID of composite.
1038  unsigned int m_id;
1039  /// Element type tag.
1040  std::string m_tag;
1041  /// Determines whether items can be reordered.
1043  /// List of elements in this composite.
1044  std::vector<ElementSharedPtr> m_items;
1045  };
1046  /// Shared pointer to a composite.
1047  typedef boost::shared_ptr<Composite> CompositeSharedPtr;
1048  /// Container of composites; key is the composite id, value is the
1049  /// composite.
1050  typedef std::map<unsigned int, CompositeSharedPtr> CompositeMap;
1051 
1052  /**
1053  * Enumeration of condition types (Dirichlet, Neumann, etc).
1054  */
1056  {
1063  };
1064 
1065  /**
1066  * @brief Defines a boundary condition.
1067  *
1068  * A boundary condition is defined by its type (e.g. Dirichlet), the
1069  * field it applies to, the value imposed on this field and the
1070  * composite which the boundary condition is applied to.
1071  */
1072  struct Condition
1073  {
1074  Condition() : type(), field(), value(), m_composite() {}
1075  std::vector<ConditionType> type;
1076  std::vector<std::string> field;
1077  std::vector<std::string> value;
1078  std::vector<int> m_composite;
1079  };
1080 
1081  typedef boost::shared_ptr<Condition> ConditionSharedPtr;
1082  typedef std::map<int,ConditionSharedPtr> ConditionMap;
1083 
1084  bool operator==(ConditionSharedPtr const &c1, ConditionSharedPtr const &c2);
1085 
1086  class Mesh
1087  {
1088  public:
1089  Mesh() : m_verbose(false) {}
1090 
1091  /// Verbose flag
1093  /// Dimension of the expansion.
1094  unsigned int m_expDim;
1095  /// Dimension of the space in which the mesh is defined.
1096  unsigned int m_spaceDim;
1097  /// List of mesh nodes.
1098  std::vector<NodeSharedPtr> m_node;
1099  /// Set of element vertices.
1101  /// Set of element edges.
1103  /// Set of element faces.
1104  FaceSet m_faceSet;
1105  /// Map for elements.
1106  ElementMap m_element;
1107  /// Map for composites.
1108  CompositeMap m_composite;
1109  /// Boundary conditions maps tag to condition.
1110  ConditionMap m_condition;
1111  /// List of fields names.
1112  std::vector<std::string> m_fields;
1113  /// Map of vertex normals.
1114  boost::unordered_map<int, Node> m_vertexNormals;
1115  /// Set of all pairs of element ID and edge/face number on which to
1116  /// apply spherigon surface smoothing.
1117  set<pair<int,int> > m_spherigonSurfs;
1118  /// Returns the total number of elements in the mesh with
1119  /// dimension expDim.
1120  unsigned int GetNumElements();
1121  /// Returns the total number of elements in the mesh with
1122  /// dimension < expDim.
1123  unsigned int GetNumBndryElements();
1124  /// Returns the total number of entities in the mesh.
1125  unsigned int GetNumEntities();
1126  };
1127  /// Shared pointer to a mesh.
1128  typedef boost::shared_ptr<Mesh> MeshSharedPtr;
1129 
1130 
1131  /**
1132  * @brief A 0-dimensional vertex.
1133  */
1134  class Point : public Element {
1135  public:
1136  /// Creates an instance of this class
1137  static ElementSharedPtr create(
1138  ElmtConfig pConf,
1139  std::vector<NodeSharedPtr> pNodeList,
1140  std::vector<int> pTagList)
1141  {
1142  return boost::shared_ptr<Element>(
1143  new Point(pConf, pNodeList, pTagList));
1144  }
1145  /// Element type
1146  static LibUtilities::ShapeType m_type;
1147 
1148  Point(ElmtConfig pConf,
1149  std::vector<NodeSharedPtr> pNodeList,
1150  std::vector<int> pTagList);
1151  Point(const Point& pSrc);
1152  virtual ~Point() {}
1153 
1154  static unsigned int GetNumNodes(ElmtConfig pConf);
1155  };
1156 
1157 
1158  /**
1159  * @brief A 1-dimensional line between two vertex nodes.
1160  */
1161  class Line : public Element {
1162  public:
1163  /// Creates an instance of this class
1164  static ElementSharedPtr create(
1165  ElmtConfig pConf,
1166  std::vector<NodeSharedPtr> pNodeList,
1167  std::vector<int> pTagList)
1168  {
1169  return boost::shared_ptr<Element>(
1170  new Line(pConf, pNodeList, pTagList));
1171  }
1172  /// Element type
1173  static LibUtilities::ShapeType m_type;
1174 
1175  Line(ElmtConfig pConf,
1176  std::vector<NodeSharedPtr> pNodeList,
1177  std::vector<int> pTagList);
1178  Line(const Point& pSrc);
1179  virtual ~Line() {}
1180 
1181  virtual SpatialDomains::GeometrySharedPtr GetGeom(int coordDim);
1182 
1183  static unsigned int GetNumNodes(ElmtConfig pConf);
1184  };
1185 
1186  /**
1187  * @brief A lightweight struct for dealing with high-order triangle
1188  * alignment.
1189  *
1190  * The logic underlying these routines is taken from the original Nektar
1191  * code.
1192  */
1193  template<typename T>
1194  struct HOTriangle
1195  {
1196  HOTriangle(vector<int> pVertId,
1197  vector<T> pSurfVerts) :
1198  vertId(pVertId), surfVerts(pSurfVerts) {}
1199  HOTriangle(vector<int> pVertId) : vertId(pVertId) {}
1200 
1201  /// The triangle vertex IDs
1202  vector<int> vertId;
1203 
1204  /// The triangle surface vertices -- templated so that this can
1205  /// either be nodes or IDs.
1206  vector<T> surfVerts;
1207 
1208  /**
1209  * @brief Rotates the triangle of data points inside #surfVerts
1210  * counter-clockwise nrot times.
1211  *
1212  * @param nrot Number of times to rotate triangle.
1213  */
1214  void Rotate(int nrot)
1215  {
1216  int n, i, j, cnt;
1217  int np = ((int)sqrt(8.0*surfVerts.size()+1.0)-1)/2;
1218  vector<T> tmp(np*np);
1219 
1220  for (n = 0; n < nrot; ++n)
1221  {
1222  for (cnt = i = 0; i < np; ++i)
1223  {
1224  for (j = 0; j < np-i; ++j, cnt++)
1225  {
1226  tmp[i*np+j] = surfVerts[cnt];
1227  }
1228  }
1229  for (cnt = i = 0; i < np; ++i)
1230  {
1231  for (j = 0; j < np-i; ++j,cnt++)
1232  {
1233  surfVerts[cnt] = tmp[(np-1-i-j)*np+i];
1234  }
1235  }
1236  }
1237  }
1238 
1239  /**
1240  * @brief Reflect data points inside #surfVerts.
1241  *
1242  * This applies a mapping essentially doing the following
1243  * reordering:
1244  *
1245  * 9 9
1246  * 7 8 -> 8 7
1247  * 4 5 6 6 5 4
1248  * 0 1 2 3 3 2 1 0
1249  */
1250  void Reflect()
1251  {
1252  int i, j, cnt;
1253  int np = ((int)sqrt(8.0*surfVerts.size()+1.0)-1)/2;
1254  vector<T> tmp(np*np);
1255 
1256  for (cnt = i = 0; i < np; ++i)
1257  {
1258  for (j = 0; j < np-i; ++j,cnt++)
1259  {
1260  tmp[i*np+np-i-1-j] = surfVerts[cnt];
1261  }
1262  }
1263 
1264  for (cnt = i = 0; i < np; ++i)
1265  {
1266  for(j = 0; j < np-i; ++j,cnt++)
1267  {
1268  surfVerts[cnt] = tmp[i*np+j];
1269  }
1270  }
1271  }
1272 
1273  /**
1274  * @brief Align this surface to a given vertex ID.
1275  */
1276  void Align(vector<int> vertId)
1277  {
1278  if (vertId[0] == this->vertId[0])
1279  {
1280  if (vertId[1] == this->vertId[1] ||
1281  vertId[1] == this->vertId[2])
1282  {
1283  if (vertId[1] == this->vertId[2])
1284  {
1285  Rotate(1);
1286  Reflect();
1287  }
1288  }
1289  }
1290  else if (vertId[0] == this->vertId[1])
1291  {
1292  if (vertId[1] == this->vertId[0] ||
1293  vertId[1] == this->vertId[2])
1294  {
1295  if (vertId[1] == this->vertId[0])
1296  {
1297  Reflect();
1298  }
1299  else
1300  {
1301  Rotate(2);
1302  }
1303  }
1304  }
1305  else if (vertId[0] == this->vertId[2])
1306  {
1307  if (vertId[1] == this->vertId[0] ||
1308  vertId[1] == this->vertId[1])
1309  {
1310  if (vertId[1] == this->vertId[1])
1311  {
1312  Rotate(2);
1313  Reflect();
1314  }
1315  else
1316  {
1317  Rotate(1);
1318  }
1319  }
1320  }
1321  }
1322  };
1323 
1324  /**
1325  * @brief A 2-dimensional three-sided element.
1326  */
1327  class Triangle : public Element {
1328  public:
1329  /// Creates an instance of this class
1330  static ElementSharedPtr create(
1331  ElmtConfig pConf,
1332  std::vector<NodeSharedPtr> pNodeList,
1333  std::vector<int> pTagList)
1334  {
1335  ElementSharedPtr e = boost::shared_ptr<Element>(
1336  new Triangle(pConf, pNodeList, pTagList));
1337  vector<EdgeSharedPtr> m_edges = e->GetEdgeList();
1338  for (int i = 0; i < m_edges.size(); ++i)
1339  {
1340  m_edges[i]->m_elLink.push_back(pair<ElementSharedPtr, int>(e,i));
1341  }
1342  return e;
1343  }
1344  /// Element type
1345  static LibUtilities::ShapeType m_type;
1346 
1347  Triangle(ElmtConfig pConf,
1348  std::vector<NodeSharedPtr> pNodeList,
1349  std::vector<int> pTagList);
1350  Triangle(const Triangle& pSrc);
1351  virtual ~Triangle() {}
1352 
1353  virtual SpatialDomains::GeometrySharedPtr GetGeom(int coordDim);
1354  virtual void Complete(int order);
1355 
1356  static unsigned int GetNumNodes(ElmtConfig pConf);
1357  };
1358 
1359 
1360  /**
1361  * @brief A 2-dimensional four-sided element.
1362  */
1363  class Quadrilateral : public Element {
1364  public:
1365  /// Creates an instance of this class
1366  static ElementSharedPtr create(
1367  ElmtConfig pConf,
1368  std::vector<NodeSharedPtr> pNodeList,
1369  std::vector<int> pTagList)
1370  {
1371  ElementSharedPtr e = boost::shared_ptr<Element>(
1372  new Quadrilateral(pConf, pNodeList, pTagList));
1373  vector<EdgeSharedPtr> m_edges = e->GetEdgeList();
1374  for (int i = 0; i < m_edges.size(); ++i)
1375  {
1376  m_edges[i]->m_elLink.push_back(pair<ElementSharedPtr, int>(e,i));
1377  }
1378  return e;
1379  }
1380  /// Element type
1381  static LibUtilities::ShapeType m_type;
1382 
1383  Quadrilateral(ElmtConfig pConf,
1384  std::vector<NodeSharedPtr> pNodeList,
1385  std::vector<int> pTagList);
1386  Quadrilateral(const Quadrilateral& pSrc);
1387  virtual ~Quadrilateral() {}
1388 
1389  virtual SpatialDomains::GeometrySharedPtr GetGeom(int coordDim);
1390  virtual void Complete(int order);
1391 
1392  static unsigned int GetNumNodes(ElmtConfig pConf);
1393  };
1394 
1395 
1396  /**
1397  * @brief A 3-dimensional four-faced element.
1398  */
1399  class Tetrahedron : public Element {
1400  public:
1401  /// Creates an instance of this class
1402  static ElementSharedPtr create(
1403  ElmtConfig pConf,
1404  std::vector<NodeSharedPtr> pNodeList,
1405  std::vector<int> pTagList)
1406  {
1407  ElementSharedPtr e = boost::shared_ptr<Element>(
1408  new Tetrahedron(pConf, pNodeList, pTagList));
1409  vector<FaceSharedPtr> faces = e->GetFaceList();
1410  for (int i = 0; i < faces.size(); ++i)
1411  {
1412  faces[i]->m_elLink.push_back(pair<ElementSharedPtr, int>(e,i));
1413  }
1414  return e;
1415  }
1416  /// Element type
1417  static LibUtilities::ShapeType m_type;
1418 
1419  Tetrahedron(ElmtConfig pConf,
1420  std::vector<NodeSharedPtr> pNodeList,
1421  std::vector<int> pTagList);
1422  Tetrahedron(const Tetrahedron& pSrc);
1423  virtual ~Tetrahedron() {}
1424 
1425  virtual SpatialDomains::GeometrySharedPtr GetGeom(int coordDim);
1426  virtual void Complete(int order);
1427 
1428  static unsigned int GetNumNodes(ElmtConfig pConf);
1429 
1430  int m_orientationMap[4];
1431  int m_origVertMap[4];
1432 
1433  protected:
1434  void OrientTet();
1435  };
1436 
1437 
1438  /**
1439  * @brief A 3-dimensional square-based pyramidic element
1440  */
1441  class Pyramid : public Element {
1442  public:
1443  /// Creates an instance of this class
1444  static ElementSharedPtr create(
1445  ElmtConfig pConf,
1446  std::vector<NodeSharedPtr> pNodeList,
1447  std::vector<int> pTagList)
1448  {
1449  ElementSharedPtr e = boost::shared_ptr<Element>(
1450  new Pyramid(pConf, pNodeList, pTagList));
1451  vector<FaceSharedPtr> faces = e->GetFaceList();
1452  for (int i = 0; i < faces.size(); ++i)
1453  {
1454  faces[i]->m_elLink.push_back(pair<ElementSharedPtr, int>(e,i));
1455  }
1456  return e;
1457  }
1458  /// Element type
1459  static LibUtilities::ShapeType type;
1460 
1461  Pyramid(ElmtConfig pConf,
1462  std::vector<NodeSharedPtr> pNodeList,
1463  std::vector<int> pTagList);
1464  Pyramid(const Pyramid& pSrc);
1465  virtual ~Pyramid() {}
1466 
1467  virtual SpatialDomains::GeometrySharedPtr GetGeom(int coordDim);
1468  static unsigned int GetNumNodes(ElmtConfig pConf);
1469 
1470  /**
1471  * Orientation of pyramid.
1472  */
1473  int orientationMap[5];
1474  };
1475 
1476  /**
1477  * @brief A 3-dimensional five-faced element (2 triangles, 3
1478  * quadrilaterals).
1479  */
1480  class Prism : public Element {
1481  public:
1482  /// Creates an instance of this class
1483  static ElementSharedPtr create(
1484  ElmtConfig pConf,
1485  std::vector<NodeSharedPtr> pNodeList,
1486  std::vector<int> pTagList)
1487  {
1488  ElementSharedPtr e = boost::shared_ptr<Element>(
1489  new Prism(pConf, pNodeList, pTagList));
1490  vector<FaceSharedPtr> faces = e->GetFaceList();
1491  for (int i = 0; i < faces.size(); ++i)
1492  {
1493  faces[i]->m_elLink.push_back(pair<ElementSharedPtr, int>(e,i));
1494  }
1495  return e;
1496  }
1497  /// Element type
1498  static LibUtilities::ShapeType m_type;
1499 
1500  Prism(ElmtConfig pConf,
1501  std::vector<NodeSharedPtr> pNodeList,
1502  std::vector<int> pTagList);
1503  Prism(const Prism& pSrc);
1504  virtual ~Prism() {}
1505 
1506  virtual SpatialDomains::GeometrySharedPtr GetGeom(int coordDim);
1507  virtual void Complete(int order);
1508 
1509  static unsigned int GetNumNodes(ElmtConfig pConf);
1510 
1511  /**
1512  * Orientation of prism; unchanged = 0; clockwise = 1;
1513  * counter-clockwise = 2. This is set by OrientPrism.
1514  */
1515  unsigned int m_orientation;
1516 
1517  protected:
1518  void OrientPrism();
1519  };
1520 
1521 
1522  /**
1523  * @brief A 3-dimensional six-faced element.
1524  */
1525  class Hexahedron : public Element {
1526  public:
1527  /// Creates an instance of this class
1528  static ElementSharedPtr create(
1529  ElmtConfig pConf,
1530  std::vector<NodeSharedPtr> pNodeList,
1531  std::vector<int> pTagList)
1532  {
1533  ElementSharedPtr e = boost::shared_ptr<Element>(
1534  new Hexahedron(pConf, pNodeList, pTagList));
1535  vector<FaceSharedPtr> faces = e->GetFaceList();
1536  for (int i = 0; i < faces.size(); ++i)
1537  {
1538  faces[i]->m_elLink.push_back(pair<ElementSharedPtr, int>(e,i));
1539  }
1540  return e;
1541  }
1542  /// Element type
1543  static LibUtilities::ShapeType m_type;
1544 
1545  Hexahedron(ElmtConfig pConf,
1546  std::vector<NodeSharedPtr> pNodeList,
1547  std::vector<int> pTagList);
1548  Hexahedron(const Hexahedron& pSrc);
1549  virtual ~Hexahedron() {}
1550 
1551  virtual SpatialDomains::GeometrySharedPtr GetGeom(int coordDim);
1552 
1553  static unsigned int GetNumNodes(ElmtConfig pConf);
1554  };
1555  }
1556 }
1557 
1558 #endif