Nektar++
QuadGeom.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: QuadGeom.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 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description:
32 //
33 //
34 ////////////////////////////////////////////////////////////////////////////////
35 
38 
39 #include <StdRegions/StdQuadExp.h>
40 #include <SpatialDomains/SegGeom.h>
41 #include <SpatialDomains/Curve.hpp>
43 
44 using namespace std;
45 
46 namespace Nektar
47 {
48 namespace SpatialDomains
49 {
50 
51 QuadGeom::QuadGeom()
52 {
53  m_shapeType = LibUtilities::eQuadrilateral;
54 }
55 
56 QuadGeom::QuadGeom(const int id,
57  const SegGeomSharedPtr edges[],
58  const CurveSharedPtr curve)
59  : Geometry2D(edges[0]->GetVertex(0)->GetCoordim(), curve)
60 {
61  int j;
62 
64  m_globalID = id;
65 
66  /// Copy the edge shared pointers.
67  m_edges.insert(m_edges.begin(), edges, edges + QuadGeom::kNedges);
68  m_eorient.resize(kNedges);
69 
70  for (j = 0; j < kNedges; ++j)
71  {
72  m_eorient[j] =
73  SegGeom::GetEdgeOrientation(*edges[j], *edges[(j + 1) % kNedges]);
74  m_verts.push_back(
75  edges[j]->GetVertex(m_eorient[j] == StdRegions::eForwards ? 0 : 1));
76  }
77 
78  for (j = 2; j < kNedges; ++j)
79  {
82  }
83 
84  m_coordim = edges[0]->GetVertex(0)->GetCoordim();
85  ASSERTL0(m_coordim > 1, "Cannot call function with dim == 1");
86 }
87 
89  : Geometry2D(in)
90 {
91  // From Geometry
94 
95  // From QuadGeom
96  m_verts = in.m_verts;
97  m_edges = in.m_edges;
98  for (int i = 0; i < kNedges; i++)
99  {
100  m_eorient[i] = in.m_eorient[i];
101  }
102 }
103 
105 {
106 }
107 
109 {
110  int order0 = max(m_edges[0]->GetXmap()->GetBasis(0)->GetNumModes(),
111  m_edges[2]->GetXmap()->GetBasis(0)->GetNumModes());
112  int order1 = max(m_edges[1]->GetXmap()->GetBasis(0)->GetNumModes(),
113  m_edges[3]->GetXmap()->GetBasis(0)->GetNumModes());
114 
115  const LibUtilities::BasisKey B0(
117  order0,
119  const LibUtilities::BasisKey B1(
121  order1,
123 
125 }
126 
128  const Array<OneD, const NekDouble> &Lcoord)
129 {
130  ASSERTL1(m_state == ePtsFilled, "Geometry is not in physical space");
131 
132  Array<OneD, NekDouble> tmp(m_xmap->GetTotPoints());
133  m_xmap->BwdTrans(m_coeffs[i], tmp);
134 
135  return m_xmap->PhysEvaluate(Lcoord, tmp);
136 }
137 
139  const QuadGeom &face2,
140  bool doRot, int dir, NekDouble angle, NekDouble tol)
141 {
142  return GetFaceOrientation(face1.m_verts, face2.m_verts,
143  doRot, dir, angle, tol);
144 }
145 
146 /**
147  * Calculate the orientation of face2 to face1 (note this is
148  * not face1 to face2!).
149  */
151  const PointGeomVector &face1, const PointGeomVector &face2,
152  bool doRot, int dir, NekDouble angle, NekDouble tol)
153 {
154  int i, j, vmap[4] = {-1, -1, -1, -1};
155 
156  if(doRot)
157  {
158  PointGeom rotPt;
159 
160  for (i = 0; i < 4; ++i)
161  {
162  rotPt.Rotate((*face1[i]), dir, angle);
163  for (j = 0; j < 4; ++j)
164  {
165  if (rotPt.dist(*face2[j]) < tol)
166  {
167  vmap[j] = i;
168  break;
169  }
170  }
171  }
172  }
173  else
174  {
175 
176  NekDouble x, y, z, x1, y1, z1, cx = 0.0, cy = 0.0, cz = 0.0;
177 
178  // For periodic faces, we calculate the vector between the centre
179  // points of the two faces. (For connected faces this will be
180  // zero). We can then use this to determine alignment later in the
181  // algorithm.
182  for (i = 0; i < 4; ++i)
183  {
184  cx += (*face2[i])(0) - (*face1[i])(0);
185  cy += (*face2[i])(1) - (*face1[i])(1);
186  cz += (*face2[i])(2) - (*face1[i])(2);
187  }
188  cx /= 4;
189  cy /= 4;
190  cz /= 4;
191 
192  // Now construct a mapping which takes us from the vertices of one
193  // face to the other. That is, vertex j of face2 corresponds to
194  // vertex vmap[j] of face1.
195  for (i = 0; i < 4; ++i)
196  {
197  x = (*face1[i])(0);
198  y = (*face1[i])(1);
199  z = (*face1[i])(2);
200  for (j = 0; j < 4; ++j)
201  {
202  x1 = (*face2[j])(0) - cx;
203  y1 = (*face2[j])(1) - cy;
204  z1 = (*face2[j])(2) - cz;
205  if (sqrt((x1 - x) * (x1 - x) + (y1 - y) * (y1 - y) +
206  (z1 - z) * (z1 - z)) < 1e-8)
207  {
208  vmap[j] = i;
209  break;
210  }
211  }
212  }
213  }
214 
215  // Use the mapping to determine the eight alignment options between
216  // faces.
217  if (vmap[1] == (vmap[0] + 1) % 4)
218  {
219  switch (vmap[0])
220  {
221  case 0:
223  break;
224  case 1:
226  break;
227  case 2:
229  break;
230  case 3:
232  break;
233  }
234  }
235  else
236  {
237  switch (vmap[0])
238  {
239  case 0:
241  break;
242  case 1:
244  break;
245  case 2:
247  break;
248  case 3:
250  break;
251  }
252  }
253 
254  ASSERTL0(false, "unable to determine face orientation");
256 }
257 
258 /**
259  * Set up GeoFac for this geometry using Coord quadrature distribution
260  */
262 {
263  if(!m_setupState)
264  {
266  }
267 
269  {
270  int i;
271  GeomType Gtype = eRegular;
272 
274 
275  // We will first check whether we have a regular or deformed
276  // geometry. We will define regular as those cases where the
277  // Jacobian and the metric terms of the derivative are constants
278  // (i.e. not coordinate dependent)
279 
280  // Check to see if expansions are linear
281  // If not linear => deformed geometry
282  for (i = 0; i < m_coordim; ++i)
283  {
284  if ((m_xmap->GetBasisNumModes(0) != 2) ||
285  (m_xmap->GetBasisNumModes(1) != 2))
286  {
287  Gtype = eDeformed;
288  }
289  }
290 
291  // For linear expansions, the mapping from standard to local
292  // element is given by the relation:
293  // x_i = 0.25 * [ ( x_i^A + x_i^B + x_i^C + x_i^D) +
294  // (-x_i^A + x_i^B + x_i^C - x_i^D)*xi_1 +
295  // (-x_i^A - x_i^B + x_i^C + x_i^D)*xi_2 +
296  // ( x_i^A - x_i^B + x_i^C - x_i^D)*xi_1*xi_2 ]
297  //
298  // The jacobian of the transformation and the metric terms
299  // dxi_i/dx_j, involve only terms of the form dx_i/dxi_j (both
300  // for coordim == 2 or 3). Inspecting the formula above, it can
301  // be appreciated that the derivatives dx_i/dxi_j will be
302  // constant, if the coefficient of the non-linear term is zero.
303  //
304  // That is why for regular geometry, we require
305  //
306  // x_i^A - x_i^B + x_i^C - x_i^D = 0
307  //
308  // or equivalently
309  //
310  // x_i^A - x_i^B = x_i^D - x_i^C
311  //
312  // This corresponds to quadrilaterals which are paralellograms.
313  if (Gtype == eRegular)
314  {
315  for (i = 0; i < m_coordim; i++)
316  {
317  if (fabs((*m_verts[0])(i) - (*m_verts[1])(i) +
318  (*m_verts[2])(i) - (*m_verts[3])(i)) >
320  {
321  Gtype = eDeformed;
322  break;
323  }
324  }
325  }
326 
328  Gtype, m_coordim, m_xmap, m_coeffs);
330  }
331 }
332 
333 /**
334  * Note verts and edges are listed according to anticlockwise
335  * convention but points in _coeffs have to be in array format from
336  * left to right.
337  */
339 {
340  // check to see if geometry structure is already filled
341  if (m_state != ePtsFilled)
342  {
343  int i, j, k;
344  int nEdgeCoeffs;
345 
346  if (m_curve)
347  {
348  int npts = m_curve->m_points.size();
349  int nEdgePts = (int)sqrt(static_cast<NekDouble>(npts));
350  Array<OneD, NekDouble> tmp(npts);
351  Array<OneD, NekDouble> tmp2(m_xmap->GetTotPoints());
352  LibUtilities::PointsKey curveKey(nEdgePts, m_curve->m_ptype);
353 
354  // Sanity checks:
355  // - Curved faces should have square number of points;
356  // - Each edge should have sqrt(npts) points.
357  ASSERTL0(nEdgePts * nEdgePts == npts,
358  "NUMPOINTS should be a square number in"
359  " quadrilteral " +
360  boost::lexical_cast<string>(m_globalID));
361 
362  for (i = 0; i < kNedges; ++i)
363  {
364  ASSERTL0(m_edges[i]->GetXmap()->GetNcoeffs() == nEdgePts,
365  "Number of edge points does not correspond to "
366  "number of face points in quadrilateral " +
367  boost::lexical_cast<string>(m_globalID));
368  }
369 
370  for (i = 0; i < m_coordim; ++i)
371  {
372  for (j = 0; j < npts; ++j)
373  {
374  tmp[j] = (m_curve->m_points[j]->GetPtr())[i];
375  }
376 
377  // Interpolate m_curve points to GLL points
378  LibUtilities::Interp2D(curveKey,
379  curveKey,
380  tmp,
381  m_xmap->GetBasis(0)->GetPointsKey(),
382  m_xmap->GetBasis(1)->GetPointsKey(),
383  tmp2);
384 
385  // Forwards transform to get coefficient space.
386  m_xmap->FwdTrans(tmp2, m_coeffs[i]);
387  }
388  }
389 
390  // Now fill in edges.
391  Array<OneD, unsigned int> mapArray;
392  Array<OneD, int> signArray;
393 
394  for (i = 0; i < kNedges; i++)
395  {
396  m_edges[i]->FillGeom();
397  m_xmap->GetTraceToElementMap(i, mapArray, signArray, m_eorient[i]);
398 
399  nEdgeCoeffs = m_edges[i]->GetXmap()->GetNcoeffs();
400 
401  for (j = 0; j < m_coordim; j++)
402  {
403  for (k = 0; k < nEdgeCoeffs; k++)
404  {
405  m_coeffs[j][mapArray[k]] =
406  signArray[k] * (m_edges[i]->GetCoeffs(j))[k];
407  }
408  }
409  }
410 
412  }
413 }
414 
415 int QuadGeom::v_GetDir(const int i, const int j) const
416 {
417  boost::ignore_unused(j); // required in 3D shapes
418 
419  return i%2;
420 }
421 
422 void QuadGeom::v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
423 {
424  Geometry::v_Reset(curvedEdges, curvedFaces);
425  CurveMap::iterator it = curvedFaces.find(m_globalID);
426 
427  if (it != curvedFaces.end())
428  {
429  m_curve = it->second;
430  }
431 
432  for (int i = 0; i < 4; ++i)
433  {
434  m_edges[i]->Reset(curvedEdges, curvedFaces);
435  }
436 
437  SetUpXmap();
438  SetUpCoeffs(m_xmap->GetNcoeffs());
439 }
440 
442 {
443  if(!m_setupState)
444  {
445  for (int i = 0; i < 4; ++i)
446  {
447  m_edges[i]->Setup();
448  }
449  SetUpXmap();
450  SetUpCoeffs(m_xmap->GetNcoeffs());
451  m_setupState = true;
452  }
453 }
454 
455 } // end of namespace
456 } // end of namespace
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
Describes the specification for a Basis.
Definition: Basis.h:50
Defines a specification for a set of points.
Definition: Points.h:60
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
2D geometry information
Definition: Geometry2D.h:69
std::vector< StdRegions::Orientation > m_eorient
Definition: Geometry2D.h:80
bool m_setupState
Wether or not the setup routines have been run.
Definition: Geometry.h:195
PointGeomSharedPtr GetVertex(int i) const
Returns vertex i of this object.
Definition: Geometry.h:348
GeomState m_state
Enumeration to dictate whether coefficients are filled.
Definition: Geometry.h:193
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition: Geometry.h:658
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated Geo...
Definition: Geometry.cpp:355
LibUtilities::ShapeType m_shapeType
Type of shape.
Definition: Geometry.h:199
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:203
GeomState m_geomFactorsState
State of the geometric factors.
Definition: Geometry.h:189
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:191
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.h:426
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:187
int m_coordim
Coordinate dimension of this geometry object.
Definition: Geometry.h:185
void Rotate(PointGeom &a, int dir, NekDouble angle)
_this = rotation of a by angle 'angle' around axis dir
Definition: PointGeom.cpp:143
NekDouble dist(PointGeom &a)
return distance between this and input a
Definition: PointGeom.cpp:181
virtual NekDouble v_GetCoord(const int i, const Array< OneD, const NekDouble > &Lcoord)
Given local collapsed coordinate Lcoord, return the value of physical coordinate in direction i.
Definition: QuadGeom.cpp:127
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated Geo...
Definition: QuadGeom.cpp:422
static StdRegions::Orientation GetFaceOrientation(const QuadGeom &face1, const QuadGeom &face2, bool doRot=false, int dir=0, NekDouble angle=0.0, NekDouble tol=1e-8)
Get the orientation of face1.
Definition: QuadGeom.cpp:138
static const int kNedges
Definition: QuadGeom.h:78
virtual int v_GetDir(const int faceidx, const int facedir) const
Returns the element coordinate direction corresponding to a given face coordinate direction.
Definition: QuadGeom.cpp:415
static StdRegions::Orientation GetEdgeOrientation(const SegGeom &edge1, const SegGeom &edge2)
Get the orientation of edge1.
Definition: SegGeom.cpp:200
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis,...
Definition: Interp.cpp:115
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48
static const NekDouble kNekZeroTol
std::vector< PointGeomSharedPtr > PointGeomVector
Definition: Geometry2D.h:64
std::shared_ptr< Curve > CurveSharedPtr
Definition: Curve.hpp:61
std::unordered_map< int, CurveSharedPtr > CurveMap
Definition: Curve.hpp:62
GeomType
Indicates the type of element geometry.
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eDeformed
Geometry is curved or has non-constant factors.
std::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:62
@ ePtsFilled
Geometric information has been generated.
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267