Nektar++
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.
304  NodeSharedPtr m_n1;
305  /// Second vertex node.
306  NodeSharedPtr m_n2;
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 == 1)
717  {
718  n += 2;
719  }
720  else if (m_dim == 2)
721  {
722  for (int i = 0; i < m_edge.size(); ++i)
723  {
724  n += m_edge[i]->GetNodeCount();
725  }
726  n -= m_vertex.size();
727  }
728  else
729  {
730  cerr << "Not supported." << endl;
731  exit(1);
732  }
733  return n;
734  }
735  /// Access the list of tags associated with this element.
736  std::vector<int> GetTagList() const {
737  return m_taglist;
738  }
739  /// Returns the number of vertices.
740  unsigned int GetVertexCount() const {
741  return m_vertex.size();
742  }
743  /// Returns the number of edges.
744  unsigned int GetEdgeCount() const {
745  return m_edge.size();
746  }
747  /// Returns the number of faces.
748  unsigned int GetFaceCount() const {
749  return m_face.size();
750  }
751  /// Change the ID of the element.
752  void SetId(unsigned int p) {
753  m_id = p;
754  }
755  /// Replace a vertex with another vertex object.
756  void SetVertex(unsigned int p, NodeSharedPtr pNew);
757  /// Replace an edge with another edge object.
758  void SetEdge(unsigned int p, EdgeSharedPtr pNew);
759  /// Replace a face with another face object.
760  void SetFace(unsigned int p, FaceSharedPtr pNew);
761  /// Set a correspondence between this element and an edge
762  /// (2D boundary element).
763  void SetEdgeLink(EdgeSharedPtr pLink) {
764  m_edgeLink = pLink;
765  }
766  /// Get correspondence between this element and an edge.
767  EdgeSharedPtr GetEdgeLink() {
768  return m_edgeLink;
769  }
770  /// Set a correspondence between this element and a face
771  /// (3D boundary element).
772  void SetFaceLink(FaceSharedPtr pLink) {
773  m_faceLink = pLink;
774  }
775  /// Get correspondence between this element and a face.
776  FaceSharedPtr GetFaceLink() {
777  return m_faceLink;
778  }
779  /// Set a correspondence between edge or face i and its
780  /// representative boundary element m->element[expDim-1][j].
781  void SetBoundaryLink(int i, int j) {
782  m_boundaryLinks[i] = j;
783  }
784  /// Get the location of the boundary face/edge i for this element.
785  int GetBoundaryLink(int i) {
786  std::map<int,int>::iterator it = m_boundaryLinks.find(i);
787  if (it == m_boundaryLinks.end())
788  {
789  return -1;
790  }
791  else
792  {
793  return it->second;
794  }
795  }
796  /// Set the list of tags associated with this element.
797  void SetTagList(const std::vector<int> &tags)
798  {
799  m_taglist = tags;
800  }
801  /// Generate a list of vertices (1D), edges (2D), or faces (3D).
802  virtual std::string GetXmlString() const
803  {
804  std::stringstream s;
805  switch (m_dim)
806  {
807  case 1:
808  for(int j=0; j< m_vertex.size(); ++j)
809  {
810  s << std::setw(5) << m_vertex[j]->m_id << " ";
811  }
812  break;
813  case 2:
814  for(int j=0; j< m_edge.size(); ++j)
815  {
816  s << std::setw(5) << m_edge[j]->m_id << " ";
817  }
818  break;
819  case 3:
820  for(int j=0; j< m_face.size(); ++j)
821  {
822  s << std::setw(5) << m_face[j]->m_id << " ";
823  }
824  break;
825  }
826  return s.str();
827  }
828 
829  /// Generates a string listing the coordinates of all nodes
830  /// associated with this element.
831  std::string GetXmlCurveString() const
832  {
833  // Temporary node list for reordering
834  std::vector<NodeSharedPtr> nodeList;
835 
836  // Node orderings are different for different elements.
837  // Triangle
838  if (m_vertex.size() == 2)
839  {
840  nodeList.push_back(m_vertex[0]);
841  for (int i = 0; i < m_volumeNodes.size(); ++i)
842  {
843  nodeList.push_back(m_volumeNodes[i]);
844  }
845  nodeList.push_back(m_vertex[1]);
846  }
847  else if (m_vertex.size() == 3)
848  {
849  int n = m_edge[0]->GetNodeCount();
850  nodeList.resize(n*(n+1)/2);
851 
852  // Populate nodelist
853  std::copy(m_vertex.begin(), m_vertex.end(), nodeList.begin());
854  for (int i = 0; i < 3; ++i)
855  {
856  std::copy(m_edge[i]->m_edgeNodes.begin(),
857  m_edge[i]->m_edgeNodes.end(),
858  nodeList.begin() + 3 + i*(n-2));
859  if (m_edge[i]->m_n1 != m_vertex[i])
860  {
861  // If edge orientation is reversed relative to node
862  // ordering, we need to reverse order of nodes.
863  std::reverse(nodeList.begin() + 3 + i*(n-2),
864  nodeList.begin() + 3 + (i+1)*(n-2));
865  }
866  }
867 
868  // Copy volume nodes.
869  std::copy(m_volumeNodes.begin(), m_volumeNodes.end(),
870  nodeList.begin() + 3*(n-1));
871  }
872  // Quad
873  else if (m_dim == 2 && m_vertex.size() == 4)
874  {
875  int n = m_edge[0]->GetNodeCount();
876  nodeList.resize(n*n);
877 
878  // Write vertices
879  nodeList[0] = m_vertex[0];
880  nodeList[n-1] = m_vertex[1];
881  nodeList[n*n-1] = m_vertex[2];
882  nodeList[n*(n-1)] = m_vertex[3];
883 
884  // Write edge-interior
885  int skips[4][2] = {{0,1}, {n-1,n}, {n*n-1,-1}, {n*(n-1),-n}};
886  for (int i = 0; i < 4; ++i)
887  {
888  bool reverseEdge = m_edge[i]->m_n1 == m_vertex[i];
889 
890  if (!reverseEdge)
891  {
892  for (int j = 1; j < n-1; ++j)
893  {
894  nodeList[skips[i][0] + j*skips[i][1]] =
895  m_edge[i]->m_edgeNodes[n-2-j];
896  }
897  }
898  else
899  {
900  for (int j = 1; j < n-1; ++j)
901  {
902  nodeList[skips[i][0] + j*skips[i][1]] =
903  m_edge[i]->m_edgeNodes[j-1];
904  }
905  }
906  }
907 
908  // Write interior
909  for (int i = 1; i < n-1; ++i)
910  {
911  for (int j = 1; j < n-1; ++j)
912  {
913  nodeList[i*n+j] = m_volumeNodes[(i-1)*(n-2)+(j-1)];
914  }
915  }
916  }
917  else
918  {
919  cerr << "GetXmlCurveString for a " << m_vertex.size()
920  << "-vertex element is not yet implemented." << endl;
921  }
922 
923  // Finally generate the XML string corresponding to our new
924  // node reordering.
925  std::stringstream s;
926  std::string str;
927  for (int k = 0; k < nodeList.size(); ++k)
928  {
929  s << std::scientific << std::setprecision(8) << " "
930  << nodeList[k]->m_x << " " << nodeList[k]->m_y
931  << " " << nodeList[k]->m_z << " ";
932 
933  }
934  return s.str();
935  }
936 
937  /// Generate a Nektar++ geometry object for this element.
939  {
940  ASSERTL0(false, "This function should be implemented at a shape level.");
941  return boost::shared_ptr<SpatialDomains::Geometry>();
942  }
943  int GetMaxOrder();
944  /// Complete this object.
945  virtual void Complete(int order)
946  {
947  ASSERTL0(false, "This function should be implemented at a shape level.");
948  }
949  void Print()
950  {
951  int i, j;
952  for (i = 0; i < m_vertex.size(); ++i)
953  {
954  cout << m_vertex[i]->m_x << " " << m_vertex[i]->m_y << " " << m_vertex[i]->m_z << endl;
955  }
956  for (i = 0; i < m_edge.size(); ++i)
957  {
958  for (j = 0; j < m_edge[i]->m_edgeNodes.size(); ++j)
959  {
960  NodeSharedPtr n = m_edge[i]->m_edgeNodes[j];
961  cout << n->m_x << " " << n->m_y << " " << n->m_z << endl;
962  }
963  }
964  for (i = 0; i < m_face.size(); ++i)
965  {
966  for (j = 0; j < m_face[i]->m_faceNodes.size(); ++j)
967  {
968  NodeSharedPtr n = m_face[i]->m_faceNodes[j];
969  cout << n->m_x << " " << n->m_y << " " << n->m_z << endl;
970  }
971  }
972  }
973 
974  protected:
975  /// ID of the element.
976  unsigned int m_id;
977  /// Dimension of the element.
978  unsigned int m_dim;
979  /// Contains configuration of the element.
981  /// Tag character describing the element.
982  std::string m_tag;
983  /// List of integers specifying properties of the element.
984  std::vector<int> m_taglist;
985  /// List of element vertex nodes.
986  std::vector<NodeSharedPtr> m_vertex;
987  /// List of element edges.
988  std::vector<EdgeSharedPtr> m_edge;
989  /// List of element faces.
990  std::vector<FaceSharedPtr> m_face;
991  /// List of element volume nodes.
992  std::vector<NodeSharedPtr> m_volumeNodes;
993  /// Volume curve type
995  /// Pointer to the corresponding edge if element is a 2D boundary.
996  EdgeSharedPtr m_edgeLink;
997  /// Pointer to the corresponding face if element is a 3D boundary.
998  FaceSharedPtr m_faceLink;
999  /// Array mapping faces/edges to the location of the appropriate
1000  /// boundary elements in m->element.
1001  std::map<int,int> m_boundaryLinks;
1002  /// Nektar++ geometry object for this element.
1004  };
1005  /// Container for elements; key is expansion dimension, value is
1006  /// vector of elements of that dimension.
1007  typedef std::map<unsigned int, std::vector<ElementSharedPtr> > ElementMap;
1008  /// Element factory definition.
1010  ElmtConfig, std::vector<NodeSharedPtr>, std::vector<int> > ElementFactory;
1011  ElementFactory& GetElementFactory();
1012 
1013  /// Define element ordering based on ID.
1015  {
1016  typedef boost::shared_ptr<Element> pT;
1017  bool operator()(const pT a, const pT b) const
1018  {
1019  // check for 0
1020  if (a.get() == 0)
1021  {
1022  // if b is also 0, then they are equal, hence a is not
1023  // less than b
1024  return b.get() != 0;
1025  }
1026  else if (b.get() == 0)
1027  {
1028  return false;
1029  }
1030  else
1031  {
1032  return a->GetId() < b->GetId();
1033  }
1034  }
1035  };
1036 
1037 
1038  /**
1039  * @brief A composite is a collection of elements.
1040  *
1041  * All elements should be of the same type, i.e. have the same tag.
1042  */
1043  class Composite {
1044  public:
1045  Composite() : m_reorder(true) {}
1046 
1047  /// Generate the list of IDs of elements within this composite.
1048  std::string GetXmlString(bool doSort=true);
1049 
1050  /// ID of composite.
1051  unsigned int m_id;
1052  /// Element type tag.
1053  std::string m_tag;
1054  /// boundary label
1055  std::string m_label;
1056  /// Determines whether items can be reordered.
1058  /// List of elements in this composite.
1059  std::vector<ElementSharedPtr> m_items;
1060  };
1061  /// Shared pointer to a composite.
1062  typedef boost::shared_ptr<Composite> CompositeSharedPtr;
1063  /// Container of composites; key is the composite id, value is the
1064  /// composite.
1065  typedef std::map<unsigned int, CompositeSharedPtr> CompositeMap;
1066 
1067  /**
1068  * Enumeration of condition types (Dirichlet, Neumann, etc).
1069  */
1071  {
1078  };
1079 
1080  /**
1081  * @brief Defines a boundary condition.
1082  *
1083  * A boundary condition is defined by its type (e.g. Dirichlet), the
1084  * field it applies to, the value imposed on this field and the
1085  * composite which the boundary condition is applied to.
1086  */
1087  struct Condition
1088  {
1089  Condition() : type(), field(), value(), m_composite() {}
1090  std::vector<ConditionType> type;
1091  std::vector<std::string> field;
1092  std::vector<std::string> value;
1093  std::vector<int> m_composite;
1094  };
1095 
1096  typedef boost::shared_ptr<Condition> ConditionSharedPtr;
1097  typedef std::map<int,ConditionSharedPtr> ConditionMap;
1098 
1099  bool operator==(ConditionSharedPtr const &c1, ConditionSharedPtr const &c2);
1100 
1101  class Mesh
1102  {
1103  public:
1104  Mesh() : m_verbose(false) {}
1105 
1106  /// Verbose flag
1108  /// Dimension of the expansion.
1109  unsigned int m_expDim;
1110  /// Dimension of the space in which the mesh is defined.
1111  unsigned int m_spaceDim;
1112  /// List of mesh nodes.
1113  std::vector<NodeSharedPtr> m_node;
1114  /// Set of element vertices.
1115  NodeSet m_vertexSet;
1116  /// Set of element edges.
1117  EdgeSet m_edgeSet;
1118  /// Set of element faces.
1119  FaceSet m_faceSet;
1120  /// Map for elements.
1121  ElementMap m_element;
1122  /// Map for composites.
1123  CompositeMap m_composite;
1124  /// Boundary conditions maps tag to condition.
1125  ConditionMap m_condition;
1126  /// List of fields names.
1127  std::vector<std::string> m_fields;
1128  /// Map of vertex normals.
1129  boost::unordered_map<int, Node> m_vertexNormals;
1130  /// Set of all pairs of element ID and edge/face number on which to
1131  /// apply spherigon surface smoothing.
1132  set<pair<int,int> > m_spherigonSurfs;
1133  /// List of face labels for composite annotation
1134  map<int,string> m_faceLabels;
1135 
1136  /// Returns the total number of elements in the mesh with
1137  /// dimension expDim.
1138  unsigned int GetNumElements();
1139  /// Returns the total number of elements in the mesh with
1140  /// dimension < expDim.
1141  unsigned int GetNumBndryElements();
1142  /// Returns the total number of entities in the mesh.
1143  unsigned int GetNumEntities();
1144 
1145  };
1146  /// Shared pointer to a mesh.
1147  typedef boost::shared_ptr<Mesh> MeshSharedPtr;
1148 
1149 
1150  /**
1151  * @brief A 0-dimensional vertex.
1152  */
1153  class Point : public Element {
1154  public:
1155  /// Creates an instance of this class
1156  static ElementSharedPtr create(
1157  ElmtConfig pConf,
1158  std::vector<NodeSharedPtr> pNodeList,
1159  std::vector<int> pTagList)
1160  {
1161  return boost::shared_ptr<Element>(
1162  new Point(pConf, pNodeList, pTagList));
1163  }
1164  /// Element type
1165  static LibUtilities::ShapeType m_type;
1166 
1167  Point(ElmtConfig pConf,
1168  std::vector<NodeSharedPtr> pNodeList,
1169  std::vector<int> pTagList);
1170  Point(const Point& pSrc);
1171  virtual ~Point() {}
1172 
1173  static unsigned int GetNumNodes(ElmtConfig pConf);
1174  };
1175 
1176 
1177  /**
1178  * @brief A 1-dimensional line between two vertex nodes.
1179  */
1180  class Line : public Element {
1181  public:
1182  /// Creates an instance of this class
1183  static ElementSharedPtr create(
1184  ElmtConfig pConf,
1185  std::vector<NodeSharedPtr> pNodeList,
1186  std::vector<int> pTagList)
1187  {
1188  return boost::shared_ptr<Element>(
1189  new Line(pConf, pNodeList, pTagList));
1190  }
1191  /// Element type
1192  static LibUtilities::ShapeType m_type;
1193 
1194  Line(ElmtConfig pConf,
1195  std::vector<NodeSharedPtr> pNodeList,
1196  std::vector<int> pTagList);
1197  Line(const Point& pSrc);
1198  virtual ~Line() {}
1199 
1200  virtual SpatialDomains::GeometrySharedPtr GetGeom(int coordDim);
1201 
1202  static unsigned int GetNumNodes(ElmtConfig pConf);
1203  };
1204 
1205  /**
1206  * @brief A lightweight struct for dealing with high-order triangle
1207  * alignment.
1208  *
1209  * The logic underlying these routines is taken from the original Nektar
1210  * code.
1211  */
1212  template<typename T>
1213  struct HOTriangle
1214  {
1215  HOTriangle(vector<int> pVertId,
1216  vector<T> pSurfVerts) :
1217  vertId(pVertId), surfVerts(pSurfVerts) {}
1218  HOTriangle(vector<int> pVertId) : vertId(pVertId) {}
1219 
1220  /// The triangle vertex IDs
1221  vector<int> vertId;
1222 
1223  /// The triangle surface vertices -- templated so that this can
1224  /// either be nodes or IDs.
1225  vector<T> surfVerts;
1226 
1227  /**
1228  * @brief Rotates the triangle of data points inside #surfVerts
1229  * counter-clockwise nrot times.
1230  *
1231  * @param nrot Number of times to rotate triangle.
1232  */
1233  void Rotate(int nrot)
1234  {
1235  int n, i, j, cnt;
1236  int np = ((int)sqrt(8.0*surfVerts.size()+1.0)-1)/2;
1237  vector<T> tmp(np*np);
1238 
1239  for (n = 0; n < nrot; ++n)
1240  {
1241  for (cnt = i = 0; i < np; ++i)
1242  {
1243  for (j = 0; j < np-i; ++j, cnt++)
1244  {
1245  tmp[i*np+j] = surfVerts[cnt];
1246  }
1247  }
1248  for (cnt = i = 0; i < np; ++i)
1249  {
1250  for (j = 0; j < np-i; ++j,cnt++)
1251  {
1252  surfVerts[cnt] = tmp[(np-1-i-j)*np+i];
1253  }
1254  }
1255  }
1256  }
1257 
1258  /**
1259  * @brief Reflect data points inside #surfVerts.
1260  *
1261  * This applies a mapping essentially doing the following
1262  * reordering:
1263  *
1264  * 9 9
1265  * 7 8 -> 8 7
1266  * 4 5 6 6 5 4
1267  * 0 1 2 3 3 2 1 0
1268  */
1269  void Reflect()
1270  {
1271  int i, j, cnt;
1272  int np = ((int)sqrt(8.0*surfVerts.size()+1.0)-1)/2;
1273  vector<T> tmp(np*np);
1274 
1275  for (cnt = i = 0; i < np; ++i)
1276  {
1277  for (j = 0; j < np-i; ++j,cnt++)
1278  {
1279  tmp[i*np+np-i-1-j] = surfVerts[cnt];
1280  }
1281  }
1282 
1283  for (cnt = i = 0; i < np; ++i)
1284  {
1285  for(j = 0; j < np-i; ++j,cnt++)
1286  {
1287  surfVerts[cnt] = tmp[i*np+j];
1288  }
1289  }
1290  }
1291 
1292  /**
1293  * @brief Align this surface to a given vertex ID.
1294  */
1295  void Align(vector<int> vertId)
1296  {
1297  if (vertId[0] == this->vertId[0])
1298  {
1299  if (vertId[1] == this->vertId[1] ||
1300  vertId[1] == this->vertId[2])
1301  {
1302  if (vertId[1] == this->vertId[2])
1303  {
1304  Rotate(1);
1305  Reflect();
1306  }
1307  }
1308  }
1309  else if (vertId[0] == this->vertId[1])
1310  {
1311  if (vertId[1] == this->vertId[0] ||
1312  vertId[1] == this->vertId[2])
1313  {
1314  if (vertId[1] == this->vertId[0])
1315  {
1316  Reflect();
1317  }
1318  else
1319  {
1320  Rotate(2);
1321  }
1322  }
1323  }
1324  else if (vertId[0] == this->vertId[2])
1325  {
1326  if (vertId[1] == this->vertId[0] ||
1327  vertId[1] == this->vertId[1])
1328  {
1329  if (vertId[1] == this->vertId[1])
1330  {
1331  Rotate(2);
1332  Reflect();
1333  }
1334  else
1335  {
1336  Rotate(1);
1337  }
1338  }
1339  }
1340  }
1341  };
1342 
1343  /**
1344  * @brief A 2-dimensional three-sided element.
1345  */
1346  class Triangle : public Element {
1347  public:
1348  /// Creates an instance of this class
1349  static ElementSharedPtr create(
1350  ElmtConfig pConf,
1351  std::vector<NodeSharedPtr> pNodeList,
1352  std::vector<int> pTagList)
1353  {
1354  ElementSharedPtr e = boost::shared_ptr<Element>(
1355  new Triangle(pConf, pNodeList, pTagList));
1356  vector<EdgeSharedPtr> m_edges = e->GetEdgeList();
1357  for (int i = 0; i < m_edges.size(); ++i)
1358  {
1359  m_edges[i]->m_elLink.push_back(pair<ElementSharedPtr, int>(e,i));
1360  }
1361  return e;
1362  }
1363  /// Element type
1364  static LibUtilities::ShapeType m_type;
1365 
1366  Triangle(ElmtConfig pConf,
1367  std::vector<NodeSharedPtr> pNodeList,
1368  std::vector<int> pTagList);
1369  Triangle(const Triangle& pSrc);
1370  virtual ~Triangle() {}
1371 
1372  virtual SpatialDomains::GeometrySharedPtr GetGeom(int coordDim);
1373  virtual void Complete(int order);
1374 
1375  static unsigned int GetNumNodes(ElmtConfig pConf);
1376  };
1377 
1378 
1379  /**
1380  * @brief A 2-dimensional four-sided element.
1381  */
1382  class Quadrilateral : public Element {
1383  public:
1384  /// Creates an instance of this class
1385  static ElementSharedPtr create(
1386  ElmtConfig pConf,
1387  std::vector<NodeSharedPtr> pNodeList,
1388  std::vector<int> pTagList)
1389  {
1390  ElementSharedPtr e = boost::shared_ptr<Element>(
1391  new Quadrilateral(pConf, pNodeList, pTagList));
1392  vector<EdgeSharedPtr> m_edges = e->GetEdgeList();
1393  for (int i = 0; i < m_edges.size(); ++i)
1394  {
1395  m_edges[i]->m_elLink.push_back(pair<ElementSharedPtr, int>(e,i));
1396  }
1397  return e;
1398  }
1399  /// Element type
1400  static LibUtilities::ShapeType m_type;
1401 
1402  Quadrilateral(ElmtConfig pConf,
1403  std::vector<NodeSharedPtr> pNodeList,
1404  std::vector<int> pTagList);
1405  Quadrilateral(const Quadrilateral& pSrc);
1406  virtual ~Quadrilateral() {}
1407 
1408  virtual SpatialDomains::GeometrySharedPtr GetGeom(int coordDim);
1409  virtual void Complete(int order);
1410 
1411  static unsigned int GetNumNodes(ElmtConfig pConf);
1412  };
1413 
1414 
1415  /**
1416  * @brief A 3-dimensional four-faced element.
1417  */
1418  class Tetrahedron : public Element {
1419  public:
1420  /// Creates an instance of this class
1421  static ElementSharedPtr create(
1422  ElmtConfig pConf,
1423  std::vector<NodeSharedPtr> pNodeList,
1424  std::vector<int> pTagList)
1425  {
1426  ElementSharedPtr e = boost::shared_ptr<Element>(
1427  new Tetrahedron(pConf, pNodeList, pTagList));
1428  vector<FaceSharedPtr> faces = e->GetFaceList();
1429  for (int i = 0; i < faces.size(); ++i)
1430  {
1431  faces[i]->m_elLink.push_back(pair<ElementSharedPtr, int>(e,i));
1432  }
1433  return e;
1434  }
1435  /// Element type
1436  static LibUtilities::ShapeType m_type;
1437 
1438  Tetrahedron(ElmtConfig pConf,
1439  std::vector<NodeSharedPtr> pNodeList,
1440  std::vector<int> pTagList);
1441  Tetrahedron(const Tetrahedron& pSrc);
1442  virtual ~Tetrahedron() {}
1443 
1444  virtual SpatialDomains::GeometrySharedPtr GetGeom(int coordDim);
1445  virtual void Complete(int order);
1446 
1447  static unsigned int GetNumNodes(ElmtConfig pConf);
1448 
1449  int m_orientationMap[4];
1450  int m_origVertMap[4];
1451 
1452  protected:
1453  void OrientTet();
1454  };
1455 
1456 
1457  /**
1458  * @brief A 3-dimensional square-based pyramidic element
1459  */
1460  class Pyramid : public Element {
1461  public:
1462  /// Creates an instance of this class
1463  static ElementSharedPtr create(
1464  ElmtConfig pConf,
1465  std::vector<NodeSharedPtr> pNodeList,
1466  std::vector<int> pTagList)
1467  {
1468  ElementSharedPtr e = boost::shared_ptr<Element>(
1469  new Pyramid(pConf, pNodeList, pTagList));
1470  vector<FaceSharedPtr> faces = e->GetFaceList();
1471  for (int i = 0; i < faces.size(); ++i)
1472  {
1473  faces[i]->m_elLink.push_back(pair<ElementSharedPtr, int>(e,i));
1474  }
1475  return e;
1476  }
1477  /// Element type
1478  static LibUtilities::ShapeType type;
1479 
1480  Pyramid(ElmtConfig pConf,
1481  std::vector<NodeSharedPtr> pNodeList,
1482  std::vector<int> pTagList);
1483  Pyramid(const Pyramid& pSrc);
1484  virtual ~Pyramid() {}
1485 
1486  virtual SpatialDomains::GeometrySharedPtr GetGeom(int coordDim);
1487  static unsigned int GetNumNodes(ElmtConfig pConf);
1488 
1489  /**
1490  * Orientation of pyramid.
1491  */
1492  int orientationMap[5];
1493  };
1494 
1495  /**
1496  * @brief A 3-dimensional five-faced element (2 triangles, 3
1497  * quadrilaterals).
1498  */
1499  class Prism : public Element {
1500  public:
1501  /// Creates an instance of this class
1502  static ElementSharedPtr create(
1503  ElmtConfig pConf,
1504  std::vector<NodeSharedPtr> pNodeList,
1505  std::vector<int> pTagList)
1506  {
1507  ElementSharedPtr e = boost::shared_ptr<Element>(
1508  new Prism(pConf, pNodeList, pTagList));
1509  vector<FaceSharedPtr> faces = e->GetFaceList();
1510  for (int i = 0; i < faces.size(); ++i)
1511  {
1512  faces[i]->m_elLink.push_back(pair<ElementSharedPtr, int>(e,i));
1513  }
1514  return e;
1515  }
1516  /// Element type
1517  static LibUtilities::ShapeType m_type;
1518 
1519  Prism(ElmtConfig pConf,
1520  std::vector<NodeSharedPtr> pNodeList,
1521  std::vector<int> pTagList);
1522  Prism(const Prism& pSrc);
1523  virtual ~Prism() {}
1524 
1525  virtual SpatialDomains::GeometrySharedPtr GetGeom(int coordDim);
1526  virtual void Complete(int order);
1527 
1528  static unsigned int GetNumNodes(ElmtConfig pConf);
1529 
1530  /**
1531  * Orientation of prism; unchanged = 0; clockwise = 1;
1532  * counter-clockwise = 2. This is set by OrientPrism.
1533  */
1534  unsigned int m_orientation;
1535 
1536  protected:
1537  void OrientPrism();
1538  };
1539 
1540 
1541  /**
1542  * @brief A 3-dimensional six-faced element.
1543  */
1544  class Hexahedron : public Element {
1545  public:
1546  /// Creates an instance of this class
1547  static ElementSharedPtr create(
1548  ElmtConfig pConf,
1549  std::vector<NodeSharedPtr> pNodeList,
1550  std::vector<int> pTagList)
1551  {
1552  ElementSharedPtr e = boost::shared_ptr<Element>(
1553  new Hexahedron(pConf, pNodeList, pTagList));
1554  vector<FaceSharedPtr> faces = e->GetFaceList();
1555  for (int i = 0; i < faces.size(); ++i)
1556  {
1557  faces[i]->m_elLink.push_back(pair<ElementSharedPtr, int>(e,i));
1558  }
1559  return e;
1560  }
1561  /// Element type
1562  static LibUtilities::ShapeType m_type;
1563 
1564  Hexahedron(ElmtConfig pConf,
1565  std::vector<NodeSharedPtr> pNodeList,
1566  std::vector<int> pTagList);
1567  Hexahedron(const Hexahedron& pSrc);
1568  virtual ~Hexahedron() {}
1569 
1570  virtual SpatialDomains::GeometrySharedPtr GetGeom(int coordDim);
1571 
1572  static unsigned int GetNumNodes(ElmtConfig pConf);
1573  };
1574  }
1575 }
1576 
1577 #endif
static ElementSharedPtr create(ElmtConfig pConf, std::vector< NodeSharedPtr > pNodeList, std::vector< int > pTagList)
Creates an instance of this class.
vector< pair< ElementSharedPtr, int > > m_elLink
Element(s) which are linked to this face.
Definition: MeshElements.h:545
vector< int > vertId
The triangle vertex IDs.
Define element ordering based on ID.
Edge(NodeSharedPtr pVertex1, NodeSharedPtr pVertex2, std::vector< NodeSharedPtr > pEdgeNodes, LibUtilities::PointsType pCurveType)
Creates a new edge.
Definition: MeshElements.h:230
std::vector< NodeSharedPtr > m_node
List of mesh nodes.
static LibUtilities::ShapeType m_type
Element type.
EdgeSet m_edgeSet
Set of element edges.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
bool m_volumeNodes
Denotes whether the element contains volume (i.e. interior) nodes. These are not supported by either ...
Definition: MeshElements.h:627
std::ostream & operator<<(std::ostream &os, const ModuleKey &rhs)
unsigned int GetVertexCount() const
Returns the number of vertices.
Definition: MeshElements.h:740
SpatialDomains::PointGeomSharedPtr m_geom
Definition: MeshElements.h:192
LibUtilities::PointsType m_curveType
Volume curve type.
Definition: MeshElements.h:994
map< int, string > m_faceLabels
List of face labels for composite annotation.
void SetVolumeNodes(std::vector< NodeSharedPtr > &nodes)
Definition: MeshElements.h:702
NodeSharedPtr m_n2
Second vertex node.
Definition: MeshElements.h:306
A 3-dimensional square-based pyramidic element.
Node operator*(const Node &pSrc) const
Definition: MeshElements.h:119
std::vector< NodeSharedPtr > GetVolumeNodes() const
Access the list of volume nodes.
Definition: MeshElements.h:699
void Rotate(int nrot)
Rotates the triangle of data points inside surfVerts counter-clockwise nrot times.
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
Represents an edge which joins two points.
Definition: MeshElements.h:227
Node operator*(const NekDouble &alpha) const
Definition: MeshElements.h:124
static LibUtilities::ShapeType type
Element type.
A composite is a collection of elements.
std::vector< std::string > field
ConditionMap m_condition
Boundary conditions maps tag to condition.
std::vector< NodeSharedPtr > m_volumeNodes
List of element volume nodes.
Definition: MeshElements.h:992
Face(const Face &pSrc)
Copy an existing face.
Definition: MeshElements.h:366
std::vector< int > m_taglist
List of integers specifying properties of the element.
Definition: MeshElements.h:984
std::map< unsigned int, CompositeSharedPtr > CompositeMap
Container of composites; key is the composite id, value is the composite.
EdgeSharedPtr m_edgeLink
Pointer to the corresponding edge if element is a 2D boundary.
Definition: MeshElements.h:996
SpatialDomains::Geometry2DSharedPtr GetGeom(int coordDim)
Generate either SpatialDomains::TriGeom or SpatialDomains::QuadGeom for this element.
Definition: MeshElements.h:501
Node curl(const Node &pSrc) const
Definition: MeshElements.h:166
SpatialDomains::SegGeomSharedPtr m_geom
Definition: MeshElements.h:315
boost::shared_ptr< Edge > EdgeSharedPtr
Shared pointer to an edge.
Definition: MeshElements.h:318
unsigned int m_id
ID of edge.
Definition: MeshElements.h:302
std::vector< FaceSharedPtr > m_face
List of element faces.
Definition: MeshElements.h:990
ElmtConfig GetConf() const
Returns the configuration of the element.
Definition: MeshElements.h:665
EdgeSharedPtr GetEdge(unsigned int i) const
Access an edge.
Definition: MeshElements.h:679
boost::unordered_set< FaceSharedPtr, FaceHash > FaceSet
Definition: MeshElements.h:574
std::vector< ElementSharedPtr > m_items
List of elements in this composite.
Nektar::LibUtilities::NekFactory< LibUtilities::ShapeType, Element, ElmtConfig, std::vector< NodeSharedPtr >, std::vector< int > > ElementFactory
Element factory definition.
LibUtilities::PointsType m_curveType
Distribution of points in this face.
Definition: MeshElements.h:543
boost::shared_ptr< Face > FaceSharedPtr
Shared pointer to a face.
Definition: MeshElements.h:550
NodeSharedPtr GetVertex(unsigned int i) const
Access a vertex node.
Definition: MeshElements.h:675
static ElementSharedPtr create(ElmtConfig pConf, std::vector< NodeSharedPtr > pNodeList, std::vector< int > pTagList)
Creates an instance of this class.
EdgeSharedPtr GetEdgeLink()
Get correspondence between this element and an edge.
Definition: MeshElements.h:767
vector< T > surfVerts
The triangle surface vertices – templated so that this can either be nodes or IDs.
unsigned int GetNodeCount() const
Returns the total number of nodes (vertices, edge nodes and face nodes and volume nodes)...
Definition: MeshElements.h:713
NekDouble abs2() const
Definition: MeshElements.h:155
void SetTagList(const std::vector< int > &tags)
Set the list of tags associated with this element.
Definition: MeshElements.h:797
FaceSet m_faceSet
Set of element faces.
Node operator+(const Node &pSrc) const
Definition: MeshElements.h:109
boost::shared_ptr< Node > NodeSharedPtr
Shared pointer to a Node.
Definition: MeshElements.h:195
unsigned int GetEdgeCount() const
Returns the number of edges.
Definition: MeshElements.h:744
std::vector< NodeSharedPtr > m_edgeNodes
List of control nodes between the first and second vertices.
Definition: MeshElements.h:308
void SetEdgeLink(EdgeSharedPtr pLink)
Set a correspondence between this element and an edge (2D boundary element).
Definition: MeshElements.h:763
static ElementSharedPtr create(ElmtConfig pConf, std::vector< NodeSharedPtr > pNodeList, std::vector< int > pTagList)
Creates an instance of this class.
unsigned int m_spaceDim
Dimension of the space in which the mesh is defined.
A 2-dimensional three-sided element.
static StdRegions::Orientation GetEdgeOrientation(const SegGeom &edge1, const SegGeom &edge2)
Get the orientation of edge1.
Definition: SegGeom.cpp:293
std::vector< NodeSharedPtr > GetVertexList() const
Access the list of vertex nodes.
Definition: MeshElements.h:687
boost::shared_ptr< Curve > CurveSharedPtr
Definition: Curve.hpp:62
vector< pair< ElementSharedPtr, int > > m_elLink
Element(s) which are linked to this edge.
Definition: MeshElements.h:312
static LibUtilities::ShapeType m_type
Element type.
boost::shared_ptr< Composite > CompositeSharedPtr
Shared pointer to a composite.
SpatialDomains::GeometrySharedPtr m_geom
Nektar++ geometry object for this element.
NekDouble m_x
X-coordinate.
Definition: MeshElements.h:185
boost::shared_ptr< Condition > ConditionSharedPtr
boost::shared_ptr< Element > ElementSharedPtr
Shared pointer to an element.
Definition: MeshElements.h:63
bool operator==(Face &pSrc)
Equality is defined by matching all vertices.
Definition: MeshElements.h:373
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:63
std::vector< NodeSharedPtr > m_vertexList
List of vertex nodes.
Definition: MeshElements.h:537
static ElementSharedPtr create(ElmtConfig pConf, std::vector< NodeSharedPtr > pNodeList, std::vector< int > pTagList)
Creates an instance of this class.
NodeSet m_vertexSet
Set of element vertices.
bool m_reorient
Denotes whether the element needs to be re-orientated for a spectral element framework.
Definition: MeshElements.h:632
void operator+=(const Node &pSrc)
Definition: MeshElements.h:134
A 3-dimensional five-faced element (2 triangles, 3 quadrilaterals).
boost::unordered_set< NodeSharedPtr, NodeHash > NodeSet
Definition: MeshElements.h:218
HOTriangle(vector< int > pVertId, vector< T > pSurfVerts)
bool m_verbose
Verbose flag.
void operator/=(const NekDouble &alpha)
Definition: MeshElements.h:148
A 3-dimensional six-faced element.
LibUtilities::ShapeType m_e
Element type (e.g. triangle, quad, etc).
Definition: MeshElements.h:618
unsigned int m_order
Order of the element.
Definition: MeshElements.h:629
std::map< unsigned int, std::vector< ElementSharedPtr > > ElementMap
Container for elements; key is expansion dimension, value is vector of elements of that dimension...
string GetXmlString(char tag, vector< unsigned int > &ids)
std::map< int, int > m_boundaryLinks
Array mapping faces/edges to the location of the appropriate boundary elements in m->element...
Node operator/(const NekDouble &alpha) const
Definition: MeshElements.h:129
HOTriangle(vector< int > pVertId)
boost::shared_ptr< Element > pT
virtual SpatialDomains::GeometrySharedPtr GetGeom(int coordDim)
Generate a Nektar++ geometry object for this element.
Definition: MeshElements.h:938
CompositeMap m_composite
Map for composites.
std::vector< ConditionType > type
void SetId(unsigned int p)
Change the ID of the element.
Definition: MeshElements.h:752
boost::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:60
ElmtConfig(ElmtConfig const &p)
Definition: MeshElements.h:604
Face(std::vector< NodeSharedPtr > pVertexList, std::vector< NodeSharedPtr > pFaceNodes, std::vector< EdgeSharedPtr > pEdgeList, LibUtilities::PointsType pCurveType)
Create a new face.
Definition: MeshElements.h:355
unsigned int m_id
ID of the element.
Definition: MeshElements.h:976
std::size_t operator()(NodeSharedPtr const &p) const
Definition: MeshElements.h:209
std::string m_tag
Element type tag.
static LibUtilities::ShapeType m_type
Element type.
LibUtilities::PointsType m_faceCurveType
Distribution of points in faces.
Definition: MeshElements.h:636
A 1-dimensional line between two vertex nodes.
bool m_faceNodes
Denotes whether the element contains face nodes. For 2D elements, if this is true then the element co...
Definition: MeshElements.h:622
unsigned int GetFaceCount() const
Returns the number of faces.
Definition: MeshElements.h:748
boost::shared_ptr< Mesh > MeshSharedPtr
Shared pointer to a mesh.
SpatialDomains::PointGeomSharedPtr GetGeom(int coordDim)
Generate a SpatialDomains::PointGeom for this node.
Definition: MeshElements.h:173
static LibUtilities::ShapeType m_type
Element type.
Defines a hash function for nodes.
Definition: MeshElements.h:207
int m_id
ID of node.
Definition: MeshElements.h:183
NekDouble m_y
Y-coordinate.
Definition: MeshElements.h:187
ElmtConfig(LibUtilities::ShapeType pE, unsigned int pOrder, bool pFn, bool pVn, bool pReorient=true, LibUtilities::PointsType pECt=LibUtilities::ePolyEvenlySpaced, LibUtilities::PointsType pFCt=LibUtilities::ePolyEvenlySpaced)
Definition: MeshElements.h:585
std::vector< std::string > m_fields
List of fields names.
double NekDouble
bool operator<(NodeSharedPtr const &p1, NodeSharedPtr const &p2)
Defines ordering between two NodeSharedPtr objects.
static LibUtilities::ShapeType m_type
Element type.
Basic information about an element.
Definition: MeshElements.h:583
unsigned int GetId() const
Returns the ID of the element (or associated edge or face for boundary elements). ...
Definition: MeshElements.h:655
static ElementSharedPtr create(ElmtConfig pConf, std::vector< NodeSharedPtr > pNodeList, std::vector< int > pTagList)
Creates an instance of this class.
unsigned int GetNodeCount() const
Returns the total number of nodes defining the edge.
Definition: MeshElements.h:242
ElementMap m_element
Map for elements.
bool operator()(const pT a, const pT b) const
SpatialDomains::Geometry2DSharedPtr m_geom
Definition: MeshElements.h:547
Node operator-(const Node &pSrc) const
Definition: MeshElements.h:114
std::string m_tag
Tag character describing the element.
Definition: MeshElements.h:982
set< pair< int, int > > m_spherigonSurfs
Set of all pairs of element ID and edge/face number on which to apply spherigon surface smoothing...
static LibUtilities::ShapeType m_type
Element type.
FaceSharedPtr GetFace(unsigned int i) const
Access a face.
Definition: MeshElements.h:683
unsigned int m_dim
Dimension of the element.
Definition: MeshElements.h:978
Represents a point in the domain.
Definition: MeshElements.h:74
unsigned int m_expDim
Dimension of the expansion.
unsigned int m_id
ID of the face.
Definition: MeshElements.h:535
Node(const Node &pSrc)
Copy an existing node.
Definition: MeshElements.h:80
std::vector< NodeSharedPtr > m_faceNodes
List of face-interior nodes defining the shape of the face.
Definition: MeshElements.h:541
A 0-dimensional vertex.
static ElementSharedPtr create(ElmtConfig pConf, std::vector< NodeSharedPtr > pNodeList, std::vector< int > pTagList)
Creates an instance of this class.
void SetCurveType(LibUtilities::PointsType cT)
Definition: MeshElements.h:708
std::vector< EdgeSharedPtr > GetEdgeList() const
Access the list of edges.
Definition: MeshElements.h:691
void SetFaceLink(FaceSharedPtr pLink)
Set a correspondence between this element and a face (3D boundary element).
Definition: MeshElements.h:772
2D Nodal Fekete Points on a Triangle
Definition: PointsType.h:69
LibUtilities::PointsType m_curveType
Distributions of points along edge.
Definition: MeshElements.h:310
LibUtilities::PointsType m_edgeCurveType
Distribution of points in edges.
Definition: MeshElements.h:634
boost::shared_ptr< Geometry2D > Geometry2DSharedPtr
Definition: Geometry2D.h:59
void operator*=(const NekDouble &alpha)
Definition: MeshElements.h:141
virtual void Complete(int order)
Complete this object.
Definition: MeshElements.h:945
boost::unordered_set< EdgeSharedPtr, EdgeHash > EdgeSet
Definition: MeshElements.h:342
static LibUtilities::ShapeType m_type
Element type.
std::vector< int > m_composite
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
bool operator==(TriFaceIDs const &p1, TriFaceIDs const &p2)
LibUtilities::PointsType GetCurveType() const
Definition: MeshElements.h:705
std::vector< FaceSharedPtr > GetFaceList() const
Access the list of faces.
Definition: MeshElements.h:695
std::map< int, ConditionSharedPtr > ConditionMap
FaceSharedPtr m_faceLink
Pointer to the corresponding face if element is a 3D boundary.
Definition: MeshElements.h:998
NekDouble dot(const Node &pSrc) const
Definition: MeshElements.h:160
std::size_t operator()(EdgeSharedPtr const &p) const
Definition: MeshElements.h:332
std::vector< NodeSharedPtr > m_vertex
List of element vertex nodes.
Definition: MeshElements.h:986
A 2-dimensional four-sided element.
NekDouble m_z
Z-coordinate.
Definition: MeshElements.h:189
std::string m_label
boundary label
boost::unordered_map< int, Node > m_vertexNormals
Map of vertex normals.
Node(int pId, NekDouble pX, NekDouble pY, NekDouble pZ)
Create a new node at a specified coordinate.
Definition: MeshElements.h:77
virtual std::string GetXmlString() const
Generate a list of vertices (1D), edges (2D), or faces (3D).
Definition: MeshElements.h:802
std::vector< EdgeSharedPtr > m_edgeList
List of corresponding edges.
Definition: MeshElements.h:539
void SetBoundaryLink(int i, int j)
Set a correspondence between edge or face i and its representative boundary element m->element[expDim...
Definition: MeshElements.h:781
void Reflect()
Reflect data points inside surfVerts.
Defines a boundary condition.
FaceSharedPtr GetFaceLink()
Get correspondence between this element and a face.
Definition: MeshElements.h:776
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:312
std::string GetXmlCurveString() const
Generates a string listing the coordinates of all nodes associated with this face.
Definition: MeshElements.h:402
SpatialDomains::SegGeomSharedPtr GetGeom(int coordDim)
Generate a SpatialDomains::SegGeom object for this edge.
Definition: MeshElements.h:266
std::vector< std::string > value
static ElementSharedPtr create(ElmtConfig pConf, std::vector< NodeSharedPtr > pNodeList, std::vector< int > pTagList)
Creates an instance of this class.
unsigned int GetDim() const
Returns the expansion dimension of the element.
Definition: MeshElements.h:661
bool m_reorder
Determines whether items can be reordered.
2D Evenly-spaced points on a Triangle
Definition: PointsType.h:70
std::vector< int > GetTagList() const
Access the list of tags associated with this element.
Definition: MeshElements.h:736
Edge(const Edge &pSrc)
Copies an existing edge.
Definition: MeshElements.h:236
Represents a face comprised of three or more edges.
Definition: MeshElements.h:352
void Align(vector< int > vertId)
Align this surface to a given vertex ID.
std::string GetTag() const
Returns the tag which defines the element shape.
Definition: MeshElements.h:669
int GetBoundaryLink(int i)
Get the location of the boundary face/edge i for this element.
Definition: MeshElements.h:785
std::string GetXmlCurveString() const
Creates a Nektar++ string listing the coordinates of all the nodes.
Definition: MeshElements.h:249
Defines a hash function for edges.
Definition: MeshElements.h:330
unsigned int GetNodeCount() const
Returns the total number of nodes (vertices, edge nodes and face nodes).
Definition: MeshElements.h:389
NodeSharedPtr m_n1
First vertex node.
Definition: MeshElements.h:304
std::size_t operator()(FaceSharedPtr const &p) const
Definition: MeshElements.h:557
static ElementSharedPtr create(ElmtConfig pConf, std::vector< NodeSharedPtr > pNodeList, std::vector< int > pTagList)
Creates an instance of this class.
A lightweight struct for dealing with high-order triangle alignment.
A 3-dimensional four-faced element.
std::vector< EdgeSharedPtr > m_edge
List of element edges.
Definition: MeshElements.h:988
boost::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:53
std::string GetXmlCurveString() const
Generates a string listing the coordinates of all nodes associated with this element.
Definition: MeshElements.h:831
bool operator==(const Node &pSrc)
Define node equality based on coordinate.
Definition: MeshElements.h:104
ElementFactory & GetElementFactory()
boost::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:60
2D Nodal Electrostatic Points on a Triangle
Definition: PointsType.h:68
unsigned int m_id
ID of composite.
Base class for element definitions.
Definition: MeshElements.h:647
Provides a generic Factory class.
Definition: NekFactory.hpp:116
void SetID(int pId)
Reset the local id;.
Definition: MeshElements.h:87
bool operator<(const Node &pSrc)
Define node ordering based on ID.
Definition: MeshElements.h:99
int GetID(void)
Get the local id;.
Definition: MeshElements.h:93
ElmtConfig m_conf
Contains configuration of the element.
Definition: MeshElements.h:980