Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Triangle.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: MeshElements.cpp
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 
37 #include <LocalRegions/TriExp.h>
39 
40 using namespace std;
41 
42 namespace Nektar
43 {
44 namespace NekMeshUtils
45 {
46 
47 LibUtilities::ShapeType Triangle::m_type =
49  LibUtilities::eTriangle, Triangle::create, "Triangle");
50 
51 /**
52  * @brief Create a triangle element.
53  */
55  vector<NodeSharedPtr> pNodeList,
56  vector<int> pTagList)
57  : Element(pConf, GetNumNodes(pConf), pNodeList.size())
58 {
59  m_tag = "T";
60  m_dim = 2;
61  m_taglist = pTagList;
63  int n = m_conf.m_order - 1;
64 
65  // Create a map to relate edge nodes to a pair of vertices
66  // defining an edge. This is based on the ordering produced by
67  // gmsh.
68  map<pair<int, int>, int> edgeNodeMap;
69  map<pair<int, int>, int>::iterator it;
70  edgeNodeMap[pair<int, int>(1, 2)] = 4;
71  edgeNodeMap[pair<int, int>(2, 3)] = 4 + n;
72  edgeNodeMap[pair<int, int>(3, 1)] = 4 + 2 * n;
73 
74  // Add vertices. This logic will determine (in 2D) whether the
75  // element is clockwise (sum > 0) or counter-clockwise (sum < 0).
76  NekDouble sum = 0.0;
77  for (int i = 0; i < 3; ++i)
78  {
79  int o = (i + 1) % 3;
80  m_vertex.push_back(pNodeList[i]);
81  sum += (pNodeList[o]->m_x - pNodeList[i]->m_x) *
82  (pNodeList[o]->m_y + pNodeList[i]->m_y);
83  }
84 
85  // Create edges (with corresponding set of edge points)
86  for (it = edgeNodeMap.begin(); it != edgeNodeMap.end(); ++it)
87  {
88  vector<NodeSharedPtr> edgeNodes;
89  if (m_conf.m_order > 1)
90  {
91  for (int j = it->second; j < it->second + n; ++j)
92  {
93  edgeNodes.push_back(pNodeList[j - 1]);
94  }
95  }
96  m_edge.push_back(EdgeSharedPtr(new Edge(pNodeList[it->first.first - 1],
97  pNodeList[it->first.second - 1],
98  edgeNodes,
100  }
101 
102  if (pConf.m_reorient)
103  {
104  if (sum > 0.0)
105  {
106  reverse(m_edge.begin(), m_edge.end());
107  }
108  }
109 
110  if (m_conf.m_faceNodes)
111  {
112  m_volumeNodes.insert(m_volumeNodes.begin(),
113  pNodeList.begin() + 3 * m_conf.m_order,
114  pNodeList.end());
115  }
116 }
117 
119 {
123 
124  for (int i = 0; i < 3; ++i)
125  {
126  edges[i] = m_edge[i]->GetGeom(coordDim);
127  verts[i] = m_vertex[i]->GetGeom(coordDim);
128  }
129 
130  StdRegions::Orientation edgeorient[3] = {
131  SpatialDomains::SegGeom::GetEdgeOrientation(*edges[0], *edges[1]),
132  SpatialDomains::SegGeom::GetEdgeOrientation(*edges[1], *edges[2]),
133  SpatialDomains::SegGeom::GetEdgeOrientation(*edges[2], *edges[0])};
134 
136  m_id, verts, edges, edgeorient);
137 
138  return ret;
139 }
140 
141 /**
142  * @brief Return the number of nodes defining a triangle.
143  */
145 {
146  int n = pConf.m_order;
147  if (!pConf.m_faceNodes)
148  return (n + 1) + 2 * (n - 1) + 1;
149  else
150  return (n + 1) * (n + 2) / 2;
151 }
152 
153 void Triangle::Complete(int order)
154 {
155  int i, j;
156 
157  // Create basis key for a nodal tetrahedron.
160  order + 1,
161  LibUtilities::PointsKey(order + 1,
165  order + 1,
166  LibUtilities::PointsKey(order + 1,
168 
169  // Create a standard nodal triangle in order to get the
170  // Vandermonde matrix to perform interpolation to nodal points.
174 
176  boost::dynamic_pointer_cast<SpatialDomains::TriGeom>(this->GetGeom(3));
177 
178  // Create basis key for a triangle.
181  order + 1,
182  LibUtilities::PointsKey(order + 1,
186  order + 1,
187  LibUtilities::PointsKey(order + 1,
189 
190  // Create a triangle.
193 
194  // Get coordinate array for tetrahedron.
195  int nqtot = tri->GetTotPoints();
196  Array<OneD, NekDouble> alloc(6 * nqtot);
197  Array<OneD, NekDouble> xi(alloc);
198  Array<OneD, NekDouble> yi(alloc + nqtot);
199  Array<OneD, NekDouble> zi(alloc + 2 * nqtot);
200  Array<OneD, NekDouble> xo(alloc + 3 * nqtot);
201  Array<OneD, NekDouble> yo(alloc + 4 * nqtot);
202  Array<OneD, NekDouble> zo(alloc + 5 * nqtot);
204 
205  tri->GetCoords(xi, yi, zi);
206 
207  for (i = 0; i < 3; ++i)
208  {
209  Array<OneD, NekDouble> coeffs(nodalTri->GetNcoeffs());
210  tri->FwdTrans(alloc + i * nqtot, coeffs);
211  // Apply Vandermonde matrix to project onto nodal space.
212  nodalTri->ModalToNodal(coeffs, tmp = alloc + (i + 3) * nqtot);
213  }
214 
215  // Now extract points from the co-ordinate arrays into the
216  // edge/face/volume nodes. First, extract edge-interior nodes.
217  for (i = 0; i < 3; ++i)
218  {
219  int pos = 3 + i * (order - 1);
220  m_edge[i]->m_edgeNodes.clear();
221  for (j = 0; j < order - 1; ++j)
222  {
223  m_edge[i]->m_edgeNodes.push_back(NodeSharedPtr(
224  new Node(0, xo[pos + j], yo[pos + j], zo[pos + j])));
225  }
226  }
227 
228  // Extract face-interior nodes.
229  int pos = 3 + 3 * (order - 1);
230  for (i = pos; i < (order + 1) * (order + 2) / 2; ++i)
231  {
232  m_volumeNodes.push_back(
233  NodeSharedPtr(new Node(0, xo[i], yo[i], zo[i])));
234  }
235 
236  m_conf.m_order = order;
237  m_conf.m_faceNodes = true;
238  m_conf.m_volumeNodes = true;
239 }
240 }
241 }
virtual NEKMESHUTILS_EXPORT SpatialDomains::GeometrySharedPtr GetGeom(int coordDim)
Generate a Nektar++ geometry object for this element.
Definition: Triangle.cpp:118
bool m_faceNodes
Denotes whether the element contains face nodes. For 2D elements, if this is true then the element co...
Definition: Element.h:89
Basic information about an element.
Definition: Element.h:58
Represents an edge which joins two points.
Definition: Edge.h:61
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
virtual NEKMESHUTILS_EXPORT void Complete(int order)
Complete this object.
Definition: Triangle.cpp:153
STL namespace.
ElementFactory & GetElementFactory()
Definition: Element.cpp:47
ElmtConfig m_conf
Contains configuration of the element.
Definition: Element.h:496
boost::shared_ptr< StdNodalTriExp > StdNodalTriExpSharedPtr
static StdRegions::Orientation GetEdgeOrientation(const SegGeom &edge1, const SegGeom &edge2)
Get the orientation of edge1.
Definition: SegGeom.cpp:293
Represents a point in the domain.
Definition: Node.h:60
std::vector< int > m_taglist
List of integers specifying properties of the element.
Definition: Element.h:500
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:57
LibUtilities::PointsType m_edgeCurveType
Distribution of points in edges.
Definition: Element.h:101
unsigned int m_order
Order of the element.
Definition: Element.h:96
Principle Orthogonal Functions .
Definition: BasisType.h:47
std::vector< NodeSharedPtr > m_vertex
List of element vertex nodes.
Definition: Element.h:502
unsigned int m_dim
Dimension of the element.
Definition: Element.h:494
bool m_volumeNodes
Denotes whether the element contains volume (i.e. interior) nodes. These are not supported by either ...
Definition: Element.h:94
boost::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:60
std::vector< EdgeSharedPtr > m_edge
List of element edges.
Definition: Element.h:504
boost::shared_ptr< Node > NodeSharedPtr
Definition: Node.h:50
Principle Orthogonal Functions .
Definition: BasisType.h:46
Defines a specification for a set of points.
Definition: Points.h:58
double NekDouble
std::vector< NodeSharedPtr > m_volumeNodes
List of element volume nodes.
Definition: Element.h:508
static NEKMESHUTILS_EXPORT unsigned int GetNumNodes(ElmtConfig pConf)
Return the number of nodes defining a triangle.
Definition: Triangle.cpp:144
std::string m_tag
Tag character describing the element.
Definition: Element.h:498
boost::shared_ptr< Edge > EdgeSharedPtr
Shared pointer to an edge.
Definition: Edge.h:196
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
unsigned int m_id
ID of the element.
Definition: Element.h:492
boost::shared_ptr< TriGeom > TriGeomSharedPtr
Definition: TriGeom.h:58
LibUtilities::PointsType m_curveType
Volume curve type.
Definition: Element.h:510
bool m_reorient
Denotes whether the element needs to be re-orientated for a spectral element framework.
Definition: Element.h:99
2D Evenly-spaced points on a Triangle
Definition: PointsType.h:70
boost::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:53
Describes the specification for a Basis.
Definition: Basis.h:50
boost::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:60
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:50
Base class for element definitions.
Definition: Element.h:115
boost::shared_ptr< TriExp > TriExpSharedPtr
Definition: TriExp.h:291
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215