Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Quadrilateral.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 
36 #include <LocalRegions/QuadExp.h>
38 
39 using namespace std;
40 
41 namespace Nektar
42 {
43 namespace NekMeshUtils
44 {
45 
46 LibUtilities::ShapeType Quadrilateral::m_type =
48  LibUtilities::eQuadrilateral, Quadrilateral::create, "Quadrilateral");
49 
50 /**
51  * @brief Create a quadrilateral element.
52  */
53 Quadrilateral::Quadrilateral(ElmtConfig pConf,
54  vector<NodeSharedPtr> pNodeList,
55  vector<int> pTagList)
56  : Element(pConf, GetNumNodes(pConf), pNodeList.size())
57 {
58  m_tag = "Q";
59  m_dim = 2;
60  m_taglist = pTagList;
61  int n = m_conf.m_order - 1;
62 
63  // Create a map to relate edge nodes to a pair of vertices
64  // defining an edge. This is based on the ordering produced by
65  // gmsh.
66  map<pair<int, int>, int> edgeNodeMap;
67  map<pair<int, int>, int>::iterator it;
68  edgeNodeMap[pair<int, int>(1, 2)] = 5;
69  edgeNodeMap[pair<int, int>(2, 3)] = 5 + n;
70  edgeNodeMap[pair<int, int>(3, 4)] = 5 + 2 * n;
71  edgeNodeMap[pair<int, int>(4, 1)] = 5 + 3 * n;
72 
73  // Add vertices. This logic will determine (in 2D) whether the
74  // element is clockwise (sum > 0) or counter-clockwise (sum < 0).
75  NekDouble sum = 0.0;
76  for (int i = 0; i < 4; ++i)
77  {
78  int o = (i + 1) % 4;
79  m_vertex.push_back(pNodeList[i]);
80  sum += (pNodeList[o]->m_x - pNodeList[i]->m_x) *
81  (pNodeList[o]->m_y + pNodeList[i]->m_y);
82  }
83 
84  // Create edges (with corresponding set of edge points)
85  for (it = edgeNodeMap.begin(); it != edgeNodeMap.end(); ++it)
86  {
87  vector<NodeSharedPtr> edgeNodes;
88  if (m_conf.m_order > 1)
89  {
90  for (int j = it->second; j < it->second + n; ++j)
91  {
92  edgeNodes.push_back(pNodeList[j - 1]);
93  }
94  }
95  m_edge.push_back(EdgeSharedPtr(new Edge(pNodeList[it->first.first - 1],
96  pNodeList[it->first.second - 1],
97  edgeNodes,
99  }
100 
101  if (pConf.m_reorient)
102  {
103  if (sum > 0.0)
104  {
105  reverse(m_edge.begin(), m_edge.end());
106  }
107  }
108 
109  if (m_conf.m_faceNodes)
110  {
111  m_volumeNodes.insert(m_volumeNodes.begin(),
112  pNodeList.begin() + 4 * m_conf.m_order,
113  pNodeList.end());
114  }
115 }
116 
117 void Quadrilateral::Complete(int order)
118 {
121  order + 1,
122  LibUtilities::PointsKey(order + 1,
124 
126  boost::dynamic_pointer_cast<SpatialDomains::QuadGeom>(this->GetGeom(3));
127 
128  // Create a quad.
131 
132  // Get coordinate array for quadrilateral.
133  int nqtot = quad->GetTotPoints();
134  Array<OneD, NekDouble> alloc(3 * nqtot);
135  Array<OneD, NekDouble> x(alloc);
136  Array<OneD, NekDouble> y(alloc + 1 * nqtot);
137  Array<OneD, NekDouble> z(alloc + 2 * nqtot);
138 
139  quad->GetCoords(x, y, z);
140 
141  // Now extract points from the co-ordinate arrays into the edge
142  // and face nodes. First, extract edge-interior nodes.
143  int edgeMap[4][2] = {{0, 1},
144  {order, order + 1},
145  {nqtot - 1, -1},
146  {order * (order + 1), -order - 1}};
147 
148  for (int i = 0; i < 4; ++i)
149  {
150  int pos = edgeMap[i][0] + edgeMap[i][1];
151  m_edge[i]->m_edgeNodes.clear();
152  for (int j = 1; j < order; ++j, pos += edgeMap[i][1])
153  {
154  m_edge[i]->m_edgeNodes.push_back(
155  NodeSharedPtr(new Node(0, x[pos], y[pos], z[pos])));
156  }
157  }
158 
159  // Extract face-interior nodes.
160  m_volumeNodes.clear();
161  for (int i = 1; i < order; ++i)
162  {
163  int pos = i * (order + 1);
164  for (int j = 1; j < order; ++j)
165  {
166  m_volumeNodes.push_back(
167  NodeSharedPtr(new Node(0, x[pos + j], y[pos + j], z[pos + j])));
168  }
169  }
170 
171  m_conf.m_order = order;
172  m_conf.m_faceNodes = true;
173  m_conf.m_volumeNodes = true;
174 }
175 
177 {
181 
182  for (int i = 0; i < 4; ++i)
183  {
184  edges[i] = m_edge[i]->GetGeom(coordDim);
185  verts[i] = m_vertex[i]->GetGeom(coordDim);
186  }
187 
188  StdRegions::Orientation edgeorient[4] = {
189  SpatialDomains::SegGeom::GetEdgeOrientation(*edges[0], *edges[1]),
190  SpatialDomains::SegGeom::GetEdgeOrientation(*edges[1], *edges[2]),
191  SpatialDomains::SegGeom::GetEdgeOrientation(*edges[2], *edges[3]),
192  SpatialDomains::SegGeom::GetEdgeOrientation(*edges[3], *edges[0])};
193 
195  m_id, verts, edges, edgeorient);
196 
197  return ret;
198 }
199 
200 /**
201  * @brief Return the number of nodes defining a quadrilateral.
202  */
204 {
205  int n = pConf.m_order;
206  if (!pConf.m_faceNodes)
207  return 4 * n;
208  else
209  return (n + 1) * (n + 1);
210 }
211 }
212 }
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.
static NEKMESHUTILS_EXPORT unsigned int GetNumNodes(ElmtConfig pConf)
Return the number of nodes defining a quadrilateral.
virtual NEKMESHUTILS_EXPORT SpatialDomains::GeometrySharedPtr GetGeom(int coordDim)
Generate a Nektar++ geometry object for this element.
boost::shared_ptr< QuadGeom > QuadGeomSharedPtr
Definition: QuadGeom.h:54
STL namespace.
ElementFactory & GetElementFactory()
Definition: Element.cpp:47
ElmtConfig m_conf
Contains configuration of the element.
Definition: Element.h:493
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:497
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
std::vector< NodeSharedPtr > m_vertex
List of element vertex nodes.
Definition: Element.h:499
unsigned int m_dim
Dimension of the element.
Definition: Element.h:491
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:501
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
boost::shared_ptr< QuadExp > QuadExpSharedPtr
Definition: QuadExp.h:293
double NekDouble
std::vector< NodeSharedPtr > m_volumeNodes
List of element volume nodes.
Definition: Element.h:505
std::string m_tag
Tag character describing the element.
Definition: Element.h:495
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
virtual NEKMESHUTILS_EXPORT void Complete(int order)
Complete this object.
unsigned int m_id
ID of the element.
Definition: Element.h:489
bool m_reorient
Denotes whether the element needs to be re-orientated for a spectral element framework.
Definition: Element.h:99
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
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215