Nektar++
Loading...
Searching...
No Matches
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
44
46{
47
53
58
59QuadGeom::QuadGeom(const int id, std::array<SegGeom *, kNedges> edges,
60 const CurveSharedPtr curve)
61 : Geometry2D(edges[0]->GetVertex(0)->GetCoordim(), curve)
62{
63 int j;
64
66 m_globalID = id;
67
68 /// Copy the edge pointers
69 m_edges = edges;
70
71 for (j = 0; j < kNedges; ++j)
72 {
73 m_eorient[j] =
74 SegGeom::GetEdgeOrientation(*edges[j], *edges[(j + 1) % kNedges]);
75 m_verts[j] =
76 edges[j]->GetVertex(m_eorient[j] == StdRegions::eForwards ? 0 : 1);
77 }
78
79 for (j = 2; j < kNedges; ++j)
80 {
84 }
85
86 m_coordim = edges[0]->GetVertex(0)->GetCoordim();
87 ASSERTL0(m_coordim > 1, "Cannot call function with dim == 1");
88}
89
91{
92 // From Geometry
95
96 // From QuadGeom
97 m_verts = in.m_verts;
98 m_edges = in.m_edges;
99 for (int i = 0; i < kNedges; i++)
100 {
101 m_eorient[i] = in.m_eorient[i];
102 }
103}
104
106{
107 int nc = 1, d0 = m_manifold[0], d1 = m_manifold[1];
108 if (0 == m_edgeNormal.size())
109 {
112 x[0] = Array<OneD, NekDouble>(3);
113 x[1] = Array<OneD, NekDouble>(3);
114 m_verts[0]->GetCoords(x[0]);
115 int i0 = 1, i1 = 0, direction = 1;
116 for (size_t i = 0; i < m_verts.size(); ++i)
117 {
118 i0 ^= 1;
119 i1 ^= 1;
120 m_verts[(i + 1) % m_verts.size()]->GetCoords(x[i1]);
121 if (m_edges[i]->GetXmap()->GetBasis(0)->GetNumModes() > 2)
122 {
123 continue;
124 }
126 m_edgeNormal[i][0] = x[i0][d1] - x[i1][d1];
127 m_edgeNormal[i][1] = x[i1][d0] - x[i0][d0];
128 }
129 if (m_coordim == 3)
130 {
131 for (size_t i = 0; i < m_verts.size(); ++i)
132 {
133 if (m_edgeNormal[i].size() == 2)
134 {
135 m_verts[i]->GetCoords(x[0]);
136 m_verts[(i + 2) % m_verts.size()]->GetCoords(x[1]);
137 if (m_edgeNormal[i][0] * (x[1][d0] - x[0][d0]) <
138 m_edgeNormal[i][1] * (x[0][d1] - x[1][d1]))
139 {
140 direction = -1;
141 }
142 break;
143 }
144 }
145 }
146 if (direction == -1)
147 {
148 for (size_t i = 0; i < m_verts.size(); ++i)
149 {
150 if (m_edgeNormal[i].size() == 2)
151 {
152 m_edgeNormal[i][0] = -m_edgeNormal[i][0];
153 m_edgeNormal[i][1] = -m_edgeNormal[i][1];
154 }
155 }
156 }
157 }
158
159 Array<OneD, NekDouble> vertex(3);
160 for (size_t i = 0; i < m_verts.size(); ++i)
161 {
162 int i1 = (i + 1) % m_verts.size();
163 if (m_verts[i]->GetGlobalID() < m_verts[i1]->GetGlobalID())
164 {
165 m_verts[i]->GetCoords(vertex);
166 }
167 else
168 {
169 m_verts[i1]->GetCoords(vertex);
170 }
171 if (m_edgeNormal[i].size() == 0)
172 {
173 nc = 0; // not sure
174 continue;
175 }
176 if (m_edgeNormal[i][0] * (gloCoord[d0] - vertex[d0]) <
177 m_edgeNormal[i][1] * (vertex[d1] - gloCoord[d1]))
178 {
179 return -1; // outside
180 }
181 }
182 // 3D manifold needs to check the distance
183 if (m_coordim == 3)
184 {
185 nc = 0;
186 }
187 // nc: 1 (side element), 0 (maybe inside), -1 (outside)
188 return nc;
189}
190
192{
193 int order0 = std::max(m_edges[0]->GetXmap()->GetBasis(0)->GetNumModes(),
194 m_edges[2]->GetXmap()->GetBasis(0)->GetNumModes());
195 int order1 = std::max(m_edges[1]->GetXmap()->GetBasis(0)->GetNumModes(),
196 m_edges[3]->GetXmap()->GetBasis(0)->GetNumModes());
197
198 std::array<LibUtilities::BasisKey, 2> basis = {
201 LibUtilities::PointsKey(order0 + 1,
205 LibUtilities::PointsKey(order1 + 1,
207
208 m_xmap = GetStdQuadFactory().CreateInstance(basis);
209}
210
212 const Array<OneD, const NekDouble> &Lcoord)
213{
214 ASSERTL1(m_state == ePtsFilled, "Geometry is not in physical space");
215
216 Array<OneD, NekDouble> tmp(m_xmap->GetTotPoints());
217 m_xmap->BwdTrans(m_coeffs[i], tmp);
218
219 return m_xmap->PhysEvaluate(Lcoord, tmp);
220}
221
223 const QuadGeom &face2,
224 bool doRot, int dir,
225 NekDouble angle,
226 NekDouble tol)
227{
228 return GetFaceOrientation(face1.m_verts, face2.m_verts, doRot, dir, angle,
229 tol);
230}
231
232/**
233 * Calculate the orientation of face2 to face1 (note this is
234 * not face1 to face2!).
235 */
237 std::array<PointGeom *, 4> face1, std::array<PointGeom *, 4> face2,
238 bool doRot, int dir, NekDouble angle, NekDouble tol)
239{
240 int i, j, vmap[4] = {-1, -1, -1, -1};
241
242 if (doRot)
243 {
244 PointGeom rotPt;
245
246 for (i = 0; i < 4; ++i)
247 {
248 rotPt.Rotate((*face1[i]), dir, angle);
249 for (j = 0; j < 4; ++j)
250 {
251 if (rotPt.dist(*face2[j]) < tol)
252 {
253 vmap[j] = i;
254 break;
255 }
256 }
257 }
258 }
259 else
260 {
261
262 NekDouble x, y, z, x1, y1, z1, cx = 0.0, cy = 0.0, cz = 0.0;
263
264 // For periodic faces, we calculate the vector between the centre
265 // points of the two faces. (For connected faces this will be
266 // zero). We can then use this to determine alignment later in the
267 // algorithm.
268 for (i = 0; i < 4; ++i)
269 {
270 cx += (*face2[i])(0) - (*face1[i])(0);
271 cy += (*face2[i])(1) - (*face1[i])(1);
272 cz += (*face2[i])(2) - (*face1[i])(2);
273 }
274 cx /= 4;
275 cy /= 4;
276 cz /= 4;
277
278 // Now construct a mapping which takes us from the vertices of one
279 // face to the other. That is, vertex j of face2 corresponds to
280 // vertex vmap[j] of face1.
281 for (i = 0; i < 4; ++i)
282 {
283 x = (*face1[i])(0);
284 y = (*face1[i])(1);
285 z = (*face1[i])(2);
286 for (j = 0; j < 4; ++j)
287 {
288 x1 = (*face2[j])(0) - cx;
289 y1 = (*face2[j])(1) - cy;
290 z1 = (*face2[j])(2) - cz;
291 if (sqrt((x1 - x) * (x1 - x) + (y1 - y) * (y1 - y) +
292 (z1 - z) * (z1 - z)) < 1e-8)
293 {
294 vmap[j] = i;
295 break;
296 }
297 }
298 }
299 }
300
301 // Use the mapping to determine the eight alignment options between
302 // faces.
303 if (vmap[1] == (vmap[0] + 1) % 4)
304 {
305 switch (vmap[0])
306 {
307 case 0:
309 break;
310 case 1:
312 break;
313 case 2:
315 break;
316 case 3:
318 break;
319 }
320 }
321 else
322 {
323 switch (vmap[0])
324 {
325 case 0:
327 break;
328 case 1:
330 break;
331 case 2:
333 break;
334 case 3:
336 break;
337 }
338 }
339
340 ASSERTL0(false, "unable to determine face orientation");
342}
343
344/**
345 * Set up GeoFac for this geometry using Coord quadrature distribution
346 */
348{
349 if (!m_setupState)
350 {
352 }
353
355 {
356 GeomType Gtype = eRegular;
357
359
360 // We will first check whether we have a regular or deformed
361 // geometry. We will define regular as those cases where the
362 // Jacobian and the metric terms of the derivative are constants
363 // (i.e. not coordinate dependent)
364
365 // Check to see if expansions are linear
366 // If not linear => deformed geometry
367 m_straightEdge = 1;
368 if ((m_xmap->GetBasisNumModes(0) != 2) ||
369 (m_xmap->GetBasisNumModes(1) != 2))
370 {
371 Gtype = eDeformed;
372 m_straightEdge = 0;
373 }
374
375 // For linear expansions, the mapping from standard to local
376 // element is given by the relation:
377 // x_i = 0.25 * [ ( x_i^A + x_i^B + x_i^C + x_i^D) +
378 // (-x_i^A + x_i^B + x_i^C - x_i^D)*xi_1 +
379 // (-x_i^A - x_i^B + x_i^C + x_i^D)*xi_2 +
380 // ( x_i^A - x_i^B + x_i^C - x_i^D)*xi_1*xi_2 ]
381 //
382 // The jacobian of the transformation and the metric terms
383 // dxi_i/dx_j, involve only terms of the form dx_i/dxi_j (both
384 // for coordim == 2 or 3). Inspecting the formula above, it can
385 // be appreciated that the derivatives dx_i/dxi_j will be
386 // constant, if the coefficient of the non-linear term is zero.
387 //
388 // That is why for regular geometry, we require
389 //
390 // x_i^A - x_i^B + x_i^C - x_i^D = 0
391 //
392 // or equivalently
393 //
394 // x_i^A - x_i^B = x_i^D - x_i^C
395 //
396 // This corresponds to quadrilaterals which are paralellograms.
398 m_manifold[0] = 0;
399 m_manifold[1] = 1;
400 if (m_coordim == 3)
401 {
402 PointGeom e01, e21, norm;
403 e01.Sub(*m_verts[0], *m_verts[1]);
404 e21.Sub(*m_verts[3], *m_verts[1]);
405 norm.Mult(e01, e21);
406 int tmpi = 0;
407 double tmp = std::fabs(norm[0]);
408 if (tmp < fabs(norm[1]))
409 {
410 tmp = fabs(norm[1]);
411 tmpi = 1;
412 }
413 if (tmp < fabs(norm[2]))
414 {
415 tmpi = 2;
416 }
417 m_manifold[0] = (tmpi + 1) % 3;
418 m_manifold[1] = (tmpi + 2) % 3;
419 m_manifold[2] = (tmpi + 3) % 3;
420 }
421
422 if (Gtype == eRegular)
423 {
425 for (int i = 0; i < m_verts.size(); ++i)
426 {
427 verts[i] = Array<OneD, NekDouble>(3);
428 m_verts[i]->GetCoords(verts[i]);
429 }
430 // a00 + a01 xi1 + a02 xi2 + a03 xi1 xi2
431 // a10 + a11 xi1 + a12 xi2 + a03 xi1 xi2
433 for (int i = 0; i < 2; i++)
434 {
435 unsigned int d = m_manifold[i];
437 // Karniadakis, Sherwin 2005, Appendix D
438 NekDouble A = verts[0][d];
439 NekDouble B = verts[1][d];
440 NekDouble D = verts[2][d];
441 NekDouble C = verts[3][d];
442 m_isoParameter[i][0] = 0.25 * (A + B + C + D); // 1
443 m_isoParameter[i][1] = 0.25 * (-A + B - C + D); // xi1
444 m_isoParameter[i][2] = 0.25 * (-A - B + C + D); // xi2
445 m_isoParameter[i][3] = 0.25 * (A - B - C + D); // xi1*xi2
446 NekDouble tmp =
447 fabs(m_isoParameter[i][1]) + fabs(m_isoParameter[i][2]);
448 if (fabs(m_isoParameter[i][3]) >
450 {
451 Gtype = eDeformed;
452 }
453 }
454 }
455
456 if (Gtype == eRegular)
457 {
459 }
460 else if (m_straightEdge)
461 {
463 }
464
466 Gtype, m_coordim, m_xmap, m_coeffs);
468 }
469}
470
471/**
472 * Note verts and edges are listed according to anticlockwise
473 * convention but points in _coeffs have to be in array format from
474 * left to right.
475 */
477{
478 // check to see if geometry structure is already filled
479 if (m_state != ePtsFilled)
480 {
481 int i, j, k;
482 int nEdgeCoeffs;
483
484 if (m_curve)
485 {
486 int npts = m_curve->m_points.size();
487 int nEdgePts = (int)sqrt(static_cast<NekDouble>(npts));
488 Array<OneD, NekDouble> tmp(npts);
489 Array<OneD, NekDouble> tmp2(m_xmap->GetTotPoints());
490 LibUtilities::PointsKey curveKey(nEdgePts, m_curve->m_ptype);
491
492 // Sanity checks:
493 // - Curved faces should have square number of points;
494 // - Each edge should have sqrt(npts) points.
495 ASSERTL0(nEdgePts * nEdgePts == npts,
496 "NUMPOINTS should be a square number in"
497 " quadrilteral " +
498 std::to_string(m_globalID));
499
500 for (i = 0; i < kNedges; ++i)
501 {
502 ASSERTL0(m_edges[i]->GetXmap()->GetNcoeffs() == nEdgePts,
503 "Number of edge points does not correspond to "
504 "number of face points in quadrilateral " +
505 std::to_string(m_globalID));
506 }
507
508 for (i = 0; i < m_coordim; ++i)
509 {
510 for (j = 0; j < npts; ++j)
511 {
512 tmp[j] = (m_curve->m_points[j]->GetPtr())[i];
513 }
514
515 // Interpolate m_curve points to GLL points
516 LibUtilities::Interp2D(curveKey, curveKey, tmp,
517 m_xmap->GetBasis(0)->GetPointsKey(),
518 m_xmap->GetBasis(1)->GetPointsKey(),
519 tmp2);
520
521 // Forwards transform to get coefficient space.
522 m_xmap->FwdTrans(tmp2, m_coeffs[i]);
523 }
524 }
525
526 // Now fill in edges.
528 Array<OneD, int> signArray;
529
530 for (i = 0; i < kNedges; i++)
531 {
532 m_edges[i]->FillGeom();
533 m_xmap->GetTraceToElementMap(i, mapArray, signArray, m_eorient[i]);
534
535 nEdgeCoeffs = m_edges[i]->GetXmap()->GetNcoeffs();
536
537 for (j = 0; j < m_coordim; j++)
538 {
539 for (k = 0; k < nEdgeCoeffs; k++)
540 {
541 m_coeffs[j][mapArray[k]] =
542 signArray[k] * (m_edges[i]->GetCoeffs(j))[k];
543 }
544 }
545 }
546
548 }
549}
550
552{
553 int i0, i1, j1, j2;
554 if (fabs(m_isoParameter[0][3]) >= fabs(m_isoParameter[1][3]))
555 {
556 i0 = 0;
557 i1 = 1;
558 }
559 else
560 {
561 i1 = 0;
562 i0 = 1;
563 m_straightEdge |= 2;
564 }
565 NekDouble gamma = m_isoParameter[i1][3] / m_isoParameter[i0][3];
566 std::vector<NekDouble> c(3);
567 for (int i = 0; i < 3; ++i)
568 {
569 c[i] = m_isoParameter[i1][i] - gamma * m_isoParameter[i0][i];
570 }
571 if (fabs(c[1]) >= fabs(c[2]))
572 {
573 j1 = 1;
574 j2 = 2;
575 }
576 else
577 {
578 j1 = 2;
579 j2 = 1;
580 m_straightEdge |= 4;
581 }
582 NekDouble beta = c[j2] / c[j1];
583 if (i0 == 1)
584 {
586 }
587 if (j1 == 2)
588 {
589 NekDouble temp = m_isoParameter[0][j1];
590 m_isoParameter[0][j1] = m_isoParameter[0][j2];
591 m_isoParameter[0][j2] = temp;
592 }
593 m_isoParameter[0][2] -= m_isoParameter[0][1] * beta;
594 m_isoParameter[1][0] = c[0];
595 m_isoParameter[1][1] = 1. / c[j1];
596 m_isoParameter[1][2] = beta;
597 m_isoParameter[1][3] = gamma;
598}
599
600int QuadGeom::v_GetDir(const int i, [[maybe_unused]] const int j) const
601{
602 return i % 2;
603}
604
605void QuadGeom::v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
606{
607 Geometry::v_Reset(curvedEdges, curvedFaces);
608 CurveMap::iterator it = curvedFaces.find(m_globalID);
609
610 if (it != curvedFaces.end())
611 {
612 m_curve = it->second;
613 }
614
615 for (int i = 0; i < 4; ++i)
616 {
617 m_edges[i]->Reset(curvedEdges, curvedFaces);
618 }
619
620 SetUpXmap();
621 SetUpCoeffs(m_xmap->GetNcoeffs());
622}
623
625{
626 if (!m_setupState)
627 {
628 for (int i = 0; i < 4; ++i)
629 {
630 m_edges[i]->Setup();
631 }
632 SetUpXmap();
633 SetUpCoeffs(m_xmap->GetNcoeffs());
634 m_setupState = true;
635 }
636}
637
638} // namespace Nektar::SpatialDomains
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Describes the specification for a Basis.
Definition Basis.h:45
Defines a specification for a set of points.
Definition Points.h:50
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
2D geometry information
Definition Geometry2D.h:50
Array< OneD, Array< OneD, NekDouble > > m_edgeNormal
Definition Geometry2D.h:66
bool m_setupState
Wether or not the setup routines have been run.
Definition Geometry.h:193
GeomState m_state
Enumeration to dictate whether coefficients are filled.
Definition Geometry.h:191
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition Geometry.h:704
Array< OneD, Array< OneD, NekDouble > > m_isoParameter
Definition Geometry.h:204
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:372
int GetGlobalID(void) const
Get the ID of this object.
Definition Geometry.h:322
LibUtilities::ShapeType m_shapeType
Type of shape.
Definition Geometry.h:197
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition Geometry.h:201
GeomState m_geomFactorsState
State of the geometric factors.
Definition Geometry.h:187
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition Geometry.h:189
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition Geometry.h:439
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition Geometry.h:185
int m_coordim
Coordinate dimension of this geometry object.
Definition Geometry.h:183
void Sub(PointGeom &a, PointGeom &b)
void Mult(PointGeom &a, PointGeom &b)
_this = a x b
void Rotate(PointGeom &a, int dir, NekDouble angle)
_this = rotation of a by angle 'angle' around axis dir
NekDouble dist(PointGeom &a)
return distance between this and input a
std::array< PointGeom *, kNverts > m_verts
Definition QuadGeom.h:117
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:605
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:222
std::array< StdRegions::Orientation, kNedges > m_eorient
Definition QuadGeom.h:119
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:211
std::array< SegGeom *, kNedges > m_edges
Definition QuadGeom.h:118
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:600
int v_AllLeftCheck(const Array< OneD, const NekDouble > &gloCoord) override
Definition QuadGeom.cpp:105
static StdRegions::Orientation GetEdgeOrientation(const SegGeom &edge1, const SegGeom &edge2)
Get the orientation of edge1.
Definition SegGeom.cpp:209
A simple factory for Xmap objects that is based on the element type, the basis and quadrature selecti...
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:101
@ beta
Gauss Radau pinned at x=-1,.
Definition PointsType.h:59
@ 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::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.
XmapFactory< StdRegions::StdQuadExp, 2 > & GetStdQuadFactory()
Definition QuadGeom.cpp:48
@ ePtsFilled
Geometric information has been generated.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition Vmath.hpp:825
scalarT< T > sqrt(scalarT< T > in)
Definition scalar.hpp:290