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 <SpatialDomains/Curve.hpp>
41 #include <SpatialDomains/SegGeom.h>
42 #include <StdRegions/StdQuadExp.h>
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, const SegGeomSharedPtr edges[],
57  const CurveSharedPtr curve)
58  : Geometry2D(edges[0]->GetVertex(0)->GetCoordim(), curve)
59 {
60  int j;
61 
63  m_globalID = id;
64 
65  /// Copy the edge shared pointers.
66  m_edges.insert(m_edges.begin(), edges, edges + QuadGeom::kNedges);
67  m_eorient.resize(kNedges);
68 
69  for (j = 0; j < kNedges; ++j)
70  {
71  m_eorient[j] =
72  SegGeom::GetEdgeOrientation(*edges[j], *edges[(j + 1) % kNedges]);
73  m_verts.push_back(
74  edges[j]->GetVertex(m_eorient[j] == StdRegions::eForwards ? 0 : 1));
75  }
76 
77  for (j = 2; j < kNedges; ++j)
78  {
82  }
83 
84  m_coordim = edges[0]->GetVertex(0)->GetCoordim();
85  ASSERTL0(m_coordim > 1, "Cannot call function with dim == 1");
86 }
87 
89 {
90  // From Geometry
93 
94  // From QuadGeom
95  m_verts = in.m_verts;
96  m_edges = in.m_edges;
97  for (int i = 0; i < kNedges; i++)
98  {
99  m_eorient[i] = in.m_eorient[i];
100  }
101 }
102 
104 {
105 }
106 
108 {
109  int order0 = max(m_edges[0]->GetXmap()->GetBasis(0)->GetNumModes(),
110  m_edges[2]->GetXmap()->GetBasis(0)->GetNumModes());
111  int order1 = max(m_edges[1]->GetXmap()->GetBasis(0)->GetNumModes(),
112  m_edges[3]->GetXmap()->GetBasis(0)->GetNumModes());
113 
114  const LibUtilities::BasisKey B0(
116  LibUtilities::PointsKey(order0 + 1,
118  const LibUtilities::BasisKey B1(
120  LibUtilities::PointsKey(order1 + 1,
122 
124 }
125 
127  const Array<OneD, const NekDouble> &Lcoord)
128 {
129  ASSERTL1(m_state == ePtsFilled, "Geometry is not in physical space");
130 
131  Array<OneD, NekDouble> tmp(m_xmap->GetTotPoints());
132  m_xmap->BwdTrans(m_coeffs[i], tmp);
133 
134  return m_xmap->PhysEvaluate(Lcoord, tmp);
135 }
136 
138  const QuadGeom &face2,
139  bool doRot, int dir,
140  NekDouble angle,
141  NekDouble tol)
142 {
143  return GetFaceOrientation(face1.m_verts, face2.m_verts, doRot, dir, angle,
144  tol);
145 }
146 
147 /**
148  * Calculate the orientation of face2 to face1 (note this is
149  * not face1 to face2!).
150  */
152  const PointGeomVector &face1, const PointGeomVector &face2, bool doRot,
153  int dir, NekDouble angle, NekDouble tol)
154 {
155  int i, j, vmap[4] = {-1, -1, -1, -1};
156 
157  if (doRot)
158  {
159  PointGeom rotPt;
160 
161  for (i = 0; i < 4; ++i)
162  {
163  rotPt.Rotate((*face1[i]), dir, angle);
164  for (j = 0; j < 4; ++j)
165  {
166  if (rotPt.dist(*face2[j]) < tol)
167  {
168  vmap[j] = i;
169  break;
170  }
171  }
172  }
173  }
174  else
175  {
176 
177  NekDouble x, y, z, x1, y1, z1, cx = 0.0, cy = 0.0, cz = 0.0;
178 
179  // For periodic faces, we calculate the vector between the centre
180  // points of the two faces. (For connected faces this will be
181  // zero). We can then use this to determine alignment later in the
182  // algorithm.
183  for (i = 0; i < 4; ++i)
184  {
185  cx += (*face2[i])(0) - (*face1[i])(0);
186  cy += (*face2[i])(1) - (*face1[i])(1);
187  cz += (*face2[i])(2) - (*face1[i])(2);
188  }
189  cx /= 4;
190  cy /= 4;
191  cz /= 4;
192 
193  // Now construct a mapping which takes us from the vertices of one
194  // face to the other. That is, vertex j of face2 corresponds to
195  // vertex vmap[j] of face1.
196  for (i = 0; i < 4; ++i)
197  {
198  x = (*face1[i])(0);
199  y = (*face1[i])(1);
200  z = (*face1[i])(2);
201  for (j = 0; j < 4; ++j)
202  {
203  x1 = (*face2[j])(0) - cx;
204  y1 = (*face2[j])(1) - cy;
205  z1 = (*face2[j])(2) - cz;
206  if (sqrt((x1 - x) * (x1 - x) + (y1 - y) * (y1 - y) +
207  (z1 - z) * (z1 - z)) < 1e-8)
208  {
209  vmap[j] = i;
210  break;
211  }
212  }
213  }
214  }
215 
216  // Use the mapping to determine the eight alignment options between
217  // faces.
218  if (vmap[1] == (vmap[0] + 1) % 4)
219  {
220  switch (vmap[0])
221  {
222  case 0:
224  break;
225  case 1:
227  break;
228  case 2:
230  break;
231  case 3:
233  break;
234  }
235  }
236  else
237  {
238  switch (vmap[0])
239  {
240  case 0:
242  break;
243  case 1:
245  break;
246  case 2:
248  break;
249  case 3:
251  break;
252  }
253  }
254 
255  ASSERTL0(false, "unable to determine face orientation");
257 }
258 
259 /**
260  * Set up GeoFac for this geometry using Coord quadrature distribution
261  */
263 {
264  if (!m_setupState)
265  {
267  }
268 
270  {
271  int i;
272  GeomType Gtype = eRegular;
273 
275 
276  // We will first check whether we have a regular or deformed
277  // geometry. We will define regular as those cases where the
278  // Jacobian and the metric terms of the derivative are constants
279  // (i.e. not coordinate dependent)
280 
281  // Check to see if expansions are linear
282  // If not linear => deformed geometry
283  for (i = 0; i < m_coordim; ++i)
284  {
285  if ((m_xmap->GetBasisNumModes(0) != 2) ||
286  (m_xmap->GetBasisNumModes(1) != 2))
287  {
288  Gtype = eDeformed;
289  }
290  }
291 
292  // For linear expansions, the mapping from standard to local
293  // element is given by the relation:
294  // x_i = 0.25 * [ ( x_i^A + x_i^B + x_i^C + x_i^D) +
295  // (-x_i^A + x_i^B + x_i^C - x_i^D)*xi_1 +
296  // (-x_i^A - x_i^B + x_i^C + x_i^D)*xi_2 +
297  // ( x_i^A - x_i^B + x_i^C - x_i^D)*xi_1*xi_2 ]
298  //
299  // The jacobian of the transformation and the metric terms
300  // dxi_i/dx_j, involve only terms of the form dx_i/dxi_j (both
301  // for coordim == 2 or 3). Inspecting the formula above, it can
302  // be appreciated that the derivatives dx_i/dxi_j will be
303  // constant, if the coefficient of the non-linear term is zero.
304  //
305  // That is why for regular geometry, we require
306  //
307  // x_i^A - x_i^B + x_i^C - x_i^D = 0
308  //
309  // or equivalently
310  //
311  // x_i^A - x_i^B = x_i^D - x_i^C
312  //
313  // This corresponds to quadrilaterals which are paralellograms.
314  if (Gtype == eRegular)
315  {
316  for (i = 0; i < m_coordim; i++)
317  {
318  if (fabs((*m_verts[0])(i) - (*m_verts[1])(i) +
319  (*m_verts[2])(i) - (*m_verts[3])(i)) >
321  {
322  Gtype = eDeformed;
323  break;
324  }
325  }
326  }
327 
329  Gtype, m_coordim, m_xmap, m_coeffs);
331  }
332 }
333 
334 /**
335  * Note verts and edges are listed according to anticlockwise
336  * convention but points in _coeffs have to be in array format from
337  * left to right.
338  */
340 {
341  // check to see if geometry structure is already filled
342  if (m_state != ePtsFilled)
343  {
344  int i, j, k;
345  int nEdgeCoeffs;
346 
347  if (m_curve)
348  {
349  int npts = m_curve->m_points.size();
350  int nEdgePts = (int)sqrt(static_cast<NekDouble>(npts));
351  Array<OneD, NekDouble> tmp(npts);
352  Array<OneD, NekDouble> tmp2(m_xmap->GetTotPoints());
353  LibUtilities::PointsKey curveKey(nEdgePts, m_curve->m_ptype);
354 
355  // Sanity checks:
356  // - Curved faces should have square number of points;
357  // - Each edge should have sqrt(npts) points.
358  ASSERTL0(nEdgePts * nEdgePts == npts,
359  "NUMPOINTS should be a square number in"
360  " quadrilteral " +
361  boost::lexical_cast<string>(m_globalID));
362 
363  for (i = 0; i < kNedges; ++i)
364  {
365  ASSERTL0(m_edges[i]->GetXmap()->GetNcoeffs() == nEdgePts,
366  "Number of edge points does not correspond to "
367  "number of face points in quadrilateral " +
368  boost::lexical_cast<string>(m_globalID));
369  }
370 
371  for (i = 0; i < m_coordim; ++i)
372  {
373  for (j = 0; j < npts; ++j)
374  {
375  tmp[j] = (m_curve->m_points[j]->GetPtr())[i];
376  }
377 
378  // Interpolate m_curve points to GLL points
379  LibUtilities::Interp2D(curveKey, curveKey, tmp,
380  m_xmap->GetBasis(0)->GetPointsKey(),
381  m_xmap->GetBasis(1)->GetPointsKey(),
382  tmp2);
383 
384  // Forwards transform to get coefficient space.
385  m_xmap->FwdTrans(tmp2, m_coeffs[i]);
386  }
387  }
388 
389  // Now fill in edges.
390  Array<OneD, unsigned int> mapArray;
391  Array<OneD, int> signArray;
392 
393  for (i = 0; i < kNedges; i++)
394  {
395  m_edges[i]->FillGeom();
396  m_xmap->GetTraceToElementMap(i, mapArray, signArray, m_eorient[i]);
397 
398  nEdgeCoeffs = m_edges[i]->GetXmap()->GetNcoeffs();
399 
400  for (j = 0; j < m_coordim; j++)
401  {
402  for (k = 0; k < nEdgeCoeffs; k++)
403  {
404  m_coeffs[j][mapArray[k]] =
405  signArray[k] * (m_edges[i]->GetCoeffs(j))[k];
406  }
407  }
408  }
409 
411  }
412 }
413 
414 int QuadGeom::v_GetDir(const int i, const int j) const
415 {
416  boost::ignore_unused(j); // required in 3D shapes
417 
418  return i % 2;
419 }
420 
421 void QuadGeom::v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
422 {
423  Geometry::v_Reset(curvedEdges, curvedFaces);
424  CurveMap::iterator it = curvedFaces.find(m_globalID);
425 
426  if (it != curvedFaces.end())
427  {
428  m_curve = it->second;
429  }
430 
431  for (int i = 0; i < 4; ++i)
432  {
433  m_edges[i]->Reset(curvedEdges, curvedFaces);
434  }
435 
436  SetUpXmap();
437  SetUpCoeffs(m_xmap->GetNcoeffs());
438 }
439 
441 {
442  if (!m_setupState)
443  {
444  for (int i = 0; i < 4; ++i)
445  {
446  m_edges[i]->Setup();
447  }
448  SetUpXmap();
449  SetUpCoeffs(m_xmap->GetNcoeffs());
450  m_setupState = true;
451  }
452 }
453 
454 } // namespace SpatialDomains
455 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
Describes the specification for a Basis.
Definition: Basis.h:50
Defines a specification for a set of points.
Definition: Points.h:59
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:84
bool m_setupState
Wether or not the setup routines have been run.
Definition: Geometry.h:194
PointGeomSharedPtr GetVertex(int i) const
Returns vertex i of this object.
Definition: Geometry.h:351
GeomState m_state
Enumeration to dictate whether coefficients are filled.
Definition: Geometry.h:192
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition: Geometry.h:683
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:376
LibUtilities::ShapeType m_shapeType
Type of shape.
Definition: Geometry.h:198
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:202
GeomState m_geomFactorsState
State of the geometric factors.
Definition: Geometry.h:188
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:190
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.h:429
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:186
int m_coordim
Coordinate dimension of this geometry object.
Definition: Geometry.h:184
void Rotate(PointGeom &a, int dir, NekDouble angle)
_this = rotation of a by angle 'angle' around axis dir
Definition: PointGeom.cpp:142
NekDouble dist(PointGeom &a)
return distance between this and input a
Definition: PointGeom.cpp:180
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces) override
Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated Geo...
Definition: QuadGeom.cpp:421
virtual void v_Setup() override
Definition: QuadGeom.cpp:440
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:137
virtual NekDouble v_GetCoord(const int i, const Array< OneD, const NekDouble > &Lcoord) override
Given local collapsed coordinate Lcoord, return the value of physical coordinate in direction i.
Definition: QuadGeom.cpp:126
virtual void v_GenGeomFactors() override
Definition: QuadGeom.cpp:262
static const int kNedges
Definition: QuadGeom.h:76
virtual void v_FillGeom() override
Definition: QuadGeom.cpp:339
virtual int v_GetDir(const int faceidx, const int facedir) const override
Returns the element coordinate direction corresponding to a given face coordinate direction.
Definition: QuadGeom.cpp:414
static StdRegions::Orientation GetEdgeOrientation(const SegGeom &edge1, const SegGeom &edge2)
Get the orientation of edge1.
Definition: SegGeom.cpp:196
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:106
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:53
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50
static const NekDouble kGeomFactorsTol
std::vector< PointGeomSharedPtr > PointGeomVector
Definition: Geometry2D.h:64
std::shared_ptr< Curve > CurveSharedPtr
Definition: Curve.hpp:60
std::unordered_map< int, CurveSharedPtr > CurveMap
Definition: Curve.hpp:61
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:2
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294