Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MeshElements/Triangle.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: Triangle.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 triangle object.
33 //
34 ////////////////////////////////////////////////////////////////////////////////
35 
37 #include <LocalRegions/TriExp.h>
39 
41 
42 using namespace std;
43 
44 namespace Nektar
45 {
46 namespace NekMeshUtils
47 {
48 
49 LibUtilities::ShapeType Triangle::m_type =
51  LibUtilities::eTriangle, Triangle::create, "Triangle");
52 
53 /**
54  * @brief Create a triangle element.
55  */
57  vector<NodeSharedPtr> pNodeList,
58  vector<int> pTagList)
59  : Element(pConf, GetNumNodes(pConf), pNodeList.size())
60 {
61  m_tag = "T";
62  m_dim = 2;
63  m_taglist = pTagList;
65  int n = m_conf.m_order - 1;
66 
67  // Create a map to relate edge nodes to a pair of vertices
68  // defining an edge. This is based on the ordering produced by
69  // gmsh.
70  map<pair<int, int>, int> edgeNodeMap;
71  map<pair<int, int>, int>::iterator it;
72  edgeNodeMap[pair<int, int>(1, 2)] = 4;
73  edgeNodeMap[pair<int, int>(2, 3)] = 4 + n;
74  edgeNodeMap[pair<int, int>(3, 1)] = 4 + 2 * n;
75 
76  // Add vertices. This logic will determine (in 2D) whether the
77  // element is clockwise (sum > 0) or counter-clockwise (sum < 0).
78  NekDouble sum = 0.0;
79  for (int i = 0; i < 3; ++i)
80  {
81  int o = (i + 1) % 3;
82  m_vertex.push_back(pNodeList[i]);
83  sum += (pNodeList[o]->m_x - pNodeList[i]->m_x) *
84  (pNodeList[o]->m_y + pNodeList[i]->m_y);
85  }
86 
87  // Create edges (with corresponding set of edge points)
88  for (it = edgeNodeMap.begin(); it != edgeNodeMap.end(); ++it)
89  {
90  vector<NodeSharedPtr> edgeNodes;
91  if (m_conf.m_order > 1)
92  {
93  for (int j = it->second; j < it->second + n; ++j)
94  {
95  edgeNodes.push_back(pNodeList[j - 1]);
96  }
97  }
98  m_edge.push_back(EdgeSharedPtr(new Edge(pNodeList[it->first.first - 1],
99  pNodeList[it->first.second - 1],
100  edgeNodes,
102  }
103 
104  if (pConf.m_reorient)
105  {
106  if (sum > 0.0)
107  {
108  reverse(m_edge.begin(), m_edge.end());
109  }
110  }
111 
112  if (m_conf.m_faceNodes)
113  {
114  m_volumeNodes.insert(m_volumeNodes.begin(),
115  pNodeList.begin() + 3 * m_conf.m_order,
116  pNodeList.end());
117  }
118 }
119 
121 {
125 
126  for (int i = 0; i < 3; ++i)
127  {
128  edges[i] = m_edge[i]->GetGeom(coordDim);
129  verts[i] = m_vertex[i]->GetGeom(coordDim);
130  }
131 
132  StdRegions::Orientation edgeorient[3] = {
133  SpatialDomains::SegGeom::GetEdgeOrientation(*edges[0], *edges[1]),
134  SpatialDomains::SegGeom::GetEdgeOrientation(*edges[1], *edges[2]),
135  SpatialDomains::SegGeom::GetEdgeOrientation(*edges[2], *edges[0])};
136 
138  m_id, verts, edges, edgeorient);
139 
140  return ret;
141 }
142 
144 {
145  int locVert = edgeId;
146  if (edge->m_n1 == m_vertex[locVert])
147  {
148  return StdRegions::eForwards;
149  }
150  else if (edge->m_n2 == m_vertex[locVert])
151  {
152  return StdRegions::eBackwards;
153  }
154  else
155  {
156  ASSERTL1(false, "Edge is not connected to this triangle.");
157  }
158 
160 }
161 
162 /**
163  * @brief Return the number of nodes defining a triangle.
164  */
166 {
167  int n = pConf.m_order;
168  if (!pConf.m_faceNodes)
169  return (n + 1) + 2 * (n - 1) + 1;
170  else
171  return (n + 1) * (n + 2) / 2;
172 }
173 
174 void Triangle::GetCurvedNodes(std::vector<NodeSharedPtr> &nodeList) const
175 {
176  int n = m_edge[0]->GetNodeCount();
177  nodeList.resize(n * (n + 1) / 2);
178 
179  // Populate nodelist
180  std::copy(m_vertex.begin(), m_vertex.end(), nodeList.begin());
181  for (int i = 0; i < 3; ++i)
182  {
183  std::copy(m_edge[i]->m_edgeNodes.begin(),
184  m_edge[i]->m_edgeNodes.end(),
185  nodeList.begin() + 3 + i * (n - 2));
186  if (m_edge[i]->m_n1 != m_vertex[i])
187  {
188  // If edge orientation is reversed relative to node ordering, we
189  // need to reverse order of nodes.
190  std::reverse(nodeList.begin() + 3 + i * (n - 2),
191  nodeList.begin() + 3 + (i + 1) * (n - 2));
192  }
193  }
194 
195  // Copy volume nodes.
196  std::copy(m_volumeNodes.begin(),
197  m_volumeNodes.end(),
198  nodeList.begin() + 3 * (n - 1));
199 }
200 
201 void Triangle::MakeOrder(int order,
204  int coordDim,
205  int &id,
206  bool justConfig)
207 
208 {
209  m_conf.m_order = order;
210  m_curveType = pType;
211  m_conf.m_volumeNodes = false;
212  m_volumeNodes.clear();
213 
214  // Triangles of order < 3 have no interior volume points.
215  if (order == 1 || order == 2)
216  {
217  m_conf.m_faceNodes = false;
218  return;
219  }
220 
221  m_conf.m_faceNodes = true;
222 
223  if (justConfig)
224  {
225  return;
226  }
227 
228  int nPoints = order + 1;
229  StdRegions::StdExpansionSharedPtr xmap = geom->GetXmap();
230 
231  Array<OneD, NekDouble> px, py;
232  LibUtilities::PointsKey pKey(nPoints, pType);
233  ASSERTL1(pKey.GetPointsDim() == 2, "Points distribution must be 2D");
234  LibUtilities::PointsManager()[pKey]->GetPoints(px, py);
235 
236  Array<OneD, Array<OneD, NekDouble> > phys(coordDim);
237 
238  for (int i = 0; i < coordDim; ++i)
239  {
240  phys[i] = Array<OneD, NekDouble>(xmap->GetTotPoints());
241  xmap->BwdTrans(geom->GetCoeffs(i), phys[i]);
242  }
243 
244  const int nTriPts = nPoints * (nPoints + 1) / 2;
245  const int nTriIntPts = (nPoints - 3) * (nPoints - 2) / 2;
246  m_volumeNodes.resize(nTriIntPts);
247 
248  for (int i = 3 + 3*(nPoints-2), cnt = 0; i < nTriPts; ++i, ++cnt)
249  {
251  xp[0] = px[i];
252  xp[1] = py[i];
253 
254  Array<OneD, NekDouble> x(3, 0.0);
255  for (int j = 0; j < coordDim; ++j)
256  {
257  x[j] = xmap->PhysEvaluate(xp, phys[j]);
258  }
259 
260  m_volumeNodes[cnt] = boost::shared_ptr<Node>(
261  new Node(id++, x[0], x[1], x[2]));
262  }
263 
264  m_conf.m_order = order;
265  m_conf.m_faceNodes = true;
266  m_conf.m_volumeNodes = false;
267 }
268 
270 {
271  Array<OneD, NekDouble> ret(3,0.0);
272 
273  ret[0] = (m_vertex[1]->m_y - m_vertex[0]->m_y) * (m_vertex[2]->m_z - m_vertex[0]->m_z) -
274  (m_vertex[1]->m_z - m_vertex[0]->m_z) * (m_vertex[2]->m_y - m_vertex[0]->m_y);
275  ret[1] = (m_vertex[1]->m_z - m_vertex[0]->m_z) * (m_vertex[2]->m_x - m_vertex[0]->m_x) -
276  (m_vertex[1]->m_x - m_vertex[0]->m_x) * (m_vertex[2]->m_z - m_vertex[0]->m_z);
277  ret[2] = (m_vertex[1]->m_x - m_vertex[0]->m_x) * (m_vertex[2]->m_y - m_vertex[0]->m_y) -
278  (m_vertex[1]->m_y - m_vertex[0]->m_y) * (m_vertex[2]->m_x - m_vertex[0]->m_x);
279 
280  NekDouble mt = ret[0] * ret[0] + ret[1] * ret[1] + ret[2] * ret[2];
281  mt = sqrt(mt);
282 
283  ret[0] /= mt;
284  ret[1] /= mt;
285  ret[2] /= mt;
286 
287  if(m_parentCAD)
288  {
289  //has cad so can orientate based on that
290  if(m_parentCAD->Orientation() == CADOrientation::eBackwards)
291  {
292  ret[0] *= -1.0;
293  ret[1] *= -1.0;
294  ret[2] *= -1.0;
295  }
296 
297  //by default normals point outwards so if we want inward for BLs
298  if(inward)
299  {
300  ret[0] *= -1.0;
301  ret[1] *= -1.0;
302  ret[2] *= -1.0;
303  }
304  }
305  return ret;
306 }
307 
308 }
309 }
virtual NEKMESHUTILS_EXPORT SpatialDomains::GeometrySharedPtr GetGeom(int coordDim)
Generate a Nektar++ geometry object for this element.
bool m_faceNodes
Denotes whether the element contains face nodes. For 2D elements, if this is true then the element co...
Definition: ElementConfig.h:80
CADObjectSharedPtr m_parentCAD
Definition: Element.h:372
Basic information about an element.
Definition: ElementConfig.h:50
virtual NEKMESHUTILS_EXPORT void MakeOrder(int order, SpatialDomains::GeometrySharedPtr geom, LibUtilities::PointsType pType, int coordDim, int &id, bool justConfig=false)
Insert interior (i.e. volume) points into this element to make the geometry an order order representa...
Represents an edge which joins two points.
Definition: Edge.h:58
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
STL namespace.
ElementFactory & GetElementFactory()
Definition: Element.cpp:47
ElmtConfig m_conf
Contains configuration of the element.
Definition: Element.h:380
static StdRegions::Orientation GetEdgeOrientation(const SegGeom &edge1, const SegGeom &edge2)
Get the orientation of edge1.
Definition: SegGeom.cpp:293
std::vector< int > m_taglist
List of integers specifying properties of the element.
Definition: Element.h:384
LibUtilities::PointsType m_edgeCurveType
Distribution of points in edges.
Definition: ElementConfig.h:92
unsigned int m_order
Order of the element.
Definition: ElementConfig.h:87
std::vector< NodeSharedPtr > m_vertex
List of element vertex nodes.
Definition: Element.h:386
unsigned int m_dim
Dimension of the element.
Definition: Element.h:378
bool m_volumeNodes
Denotes whether the element contains volume (i.e. interior) nodes. These are not supported by either ...
Definition: ElementConfig.h:85
boost::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:60
std::vector< EdgeSharedPtr > m_edge
List of element edges.
Definition: Element.h:388
unsigned int GetPointsDim() const
Definition: Points.h:149
PointsManagerT & PointsManager(void)
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:392
static NEKMESHUTILS_EXPORT unsigned int GetNumNodes(ElmtConfig pConf)
Return the number of nodes defining a triangle.
std::string m_tag
Tag character describing the element.
Definition: Element.h:382
boost::shared_ptr< Edge > EdgeSharedPtr
Shared pointer to an edge.
Definition: Edge.h:136
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
NEKMESHUTILS_EXPORT Array< OneD, NekDouble > Normal(bool inward=false)
returns the normal to the element
virtual NEKMESHUTILS_EXPORT StdRegions::Orientation GetEdgeOrient(int edgeId, EdgeSharedPtr edge)
Get the edge orientation of edge with respect to the local element, which lies at edge index edgeId...
unsigned int m_id
ID of the element.
Definition: Element.h:376
boost::shared_ptr< TriGeom > TriGeomSharedPtr
Definition: TriGeom.h:58
LibUtilities::PointsType m_curveType
Volume curve type.
Definition: Element.h:394
bool m_reorient
Denotes whether the element needs to be re-orientated for a spectral element framework.
Definition: ElementConfig.h:90
2D Evenly-spaced points on a Triangle
Definition: PointsType.h:72
boost::shared_ptr< StdExpansion > StdExpansionSharedPtr
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
boost::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:53
virtual NEKMESHUTILS_EXPORT void GetCurvedNodes(std::vector< NodeSharedPtr > &nodeList) const
get list of volume interior nodes
boost::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:60
Base class for element definitions.
Definition: Element.h:59
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215