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
43
44using namespace std;
45
47{
48
50{
52}
53
54QuadGeom::QuadGeom(const int id, const SegGeomSharedPtr edges[],
55 const CurveSharedPtr curve)
56 : Geometry2D(edges[0]->GetVertex(0)->GetCoordim(), curve)
57{
58 int j;
59
61 m_globalID = id;
62
63 /// Copy the edge shared pointers.
64 m_edges.insert(m_edges.begin(), edges, edges + QuadGeom::kNedges);
65 m_eorient.resize(kNedges);
66
67 for (j = 0; j < kNedges; ++j)
68 {
69 m_eorient[j] =
70 SegGeom::GetEdgeOrientation(*edges[j], *edges[(j + 1) % kNedges]);
71 m_verts.push_back(
72 edges[j]->GetVertex(m_eorient[j] == StdRegions::eForwards ? 0 : 1));
73 }
74
75 for (j = 2; j < kNedges; ++j)
76 {
80 }
81
82 m_coordim = edges[0]->GetVertex(0)->GetCoordim();
83 ASSERTL0(m_coordim > 1, "Cannot call function with dim == 1");
84}
85
87{
88 // From Geometry
91
92 // From QuadGeom
93 m_verts = in.m_verts;
94 m_edges = in.m_edges;
95 for (int i = 0; i < kNedges; i++)
96 {
97 m_eorient[i] = in.m_eorient[i];
98 }
99}
100
102{
103}
104
106{
107 int order0 = max(m_edges[0]->GetXmap()->GetBasis(0)->GetNumModes(),
108 m_edges[2]->GetXmap()->GetBasis(0)->GetNumModes());
109 int order1 = max(m_edges[1]->GetXmap()->GetBasis(0)->GetNumModes(),
110 m_edges[3]->GetXmap()->GetBasis(0)->GetNumModes());
111
112 const LibUtilities::BasisKey B0(
114 LibUtilities::PointsKey(order0 + 1,
116 const LibUtilities::BasisKey B1(
118 LibUtilities::PointsKey(order1 + 1,
120
122}
123
125 const Array<OneD, const NekDouble> &Lcoord)
126{
127 ASSERTL1(m_state == ePtsFilled, "Geometry is not in physical space");
128
129 Array<OneD, NekDouble> tmp(m_xmap->GetTotPoints());
130 m_xmap->BwdTrans(m_coeffs[i], tmp);
131
132 return m_xmap->PhysEvaluate(Lcoord, tmp);
133}
134
136 const QuadGeom &face2,
137 bool doRot, int dir,
138 NekDouble angle,
139 NekDouble tol)
140{
141 return GetFaceOrientation(face1.m_verts, face2.m_verts, doRot, dir, angle,
142 tol);
143}
144
145/**
146 * Calculate the orientation of face2 to face1 (note this is
147 * not face1 to face2!).
148 */
150 const PointGeomVector &face1, const PointGeomVector &face2, bool doRot,
151 int dir, NekDouble angle, NekDouble tol)
152{
153 int i, j, vmap[4] = {-1, -1, -1, -1};
154
155 if (doRot)
156 {
157 PointGeom rotPt;
158
159 for (i = 0; i < 4; ++i)
160 {
161 rotPt.Rotate((*face1[i]), dir, angle);
162 for (j = 0; j < 4; ++j)
163 {
164 if (rotPt.dist(*face2[j]) < tol)
165 {
166 vmap[j] = i;
167 break;
168 }
169 }
170 }
171 }
172 else
173 {
174
175 NekDouble x, y, z, x1, y1, z1, cx = 0.0, cy = 0.0, cz = 0.0;
176
177 // For periodic faces, we calculate the vector between the centre
178 // points of the two faces. (For connected faces this will be
179 // zero). We can then use this to determine alignment later in the
180 // algorithm.
181 for (i = 0; i < 4; ++i)
182 {
183 cx += (*face2[i])(0) - (*face1[i])(0);
184 cy += (*face2[i])(1) - (*face1[i])(1);
185 cz += (*face2[i])(2) - (*face1[i])(2);
186 }
187 cx /= 4;
188 cy /= 4;
189 cz /= 4;
190
191 // Now construct a mapping which takes us from the vertices of one
192 // face to the other. That is, vertex j of face2 corresponds to
193 // vertex vmap[j] of face1.
194 for (i = 0; i < 4; ++i)
195 {
196 x = (*face1[i])(0);
197 y = (*face1[i])(1);
198 z = (*face1[i])(2);
199 for (j = 0; j < 4; ++j)
200 {
201 x1 = (*face2[j])(0) - cx;
202 y1 = (*face2[j])(1) - cy;
203 z1 = (*face2[j])(2) - cz;
204 if (sqrt((x1 - x) * (x1 - x) + (y1 - y) * (y1 - y) +
205 (z1 - z) * (z1 - z)) < 1e-8)
206 {
207 vmap[j] = i;
208 break;
209 }
210 }
211 }
212 }
213
214 // Use the mapping to determine the eight alignment options between
215 // faces.
216 if (vmap[1] == (vmap[0] + 1) % 4)
217 {
218 switch (vmap[0])
219 {
220 case 0:
222 break;
223 case 1:
225 break;
226 case 2:
228 break;
229 case 3:
231 break;
232 }
233 }
234 else
235 {
236 switch (vmap[0])
237 {
238 case 0:
240 break;
241 case 1:
243 break;
244 case 2:
246 break;
247 case 3:
249 break;
250 }
251 }
252
253 ASSERTL0(false, "unable to determine face orientation");
255}
256
257/**
258 * Set up GeoFac for this geometry using Coord quadrature distribution
259 */
261{
262 if (!m_setupState)
263 {
265 }
266
268 {
269 GeomType Gtype = eRegular;
270
272
273 // We will first check whether we have a regular or deformed
274 // geometry. We will define regular as those cases where the
275 // Jacobian and the metric terms of the derivative are constants
276 // (i.e. not coordinate dependent)
277
278 // Check to see if expansions are linear
279 // If not linear => deformed geometry
280 m_straightEdge = 1;
281 if ((m_xmap->GetBasisNumModes(0) != 2) ||
282 (m_xmap->GetBasisNumModes(1) != 2))
283 {
284 Gtype = eDeformed;
285 m_straightEdge = 0;
286 }
287
288 // For linear expansions, the mapping from standard to local
289 // element is given by the relation:
290 // x_i = 0.25 * [ ( x_i^A + x_i^B + x_i^C + x_i^D) +
291 // (-x_i^A + x_i^B + x_i^C - x_i^D)*xi_1 +
292 // (-x_i^A - x_i^B + x_i^C + x_i^D)*xi_2 +
293 // ( x_i^A - x_i^B + x_i^C - x_i^D)*xi_1*xi_2 ]
294 //
295 // The jacobian of the transformation and the metric terms
296 // dxi_i/dx_j, involve only terms of the form dx_i/dxi_j (both
297 // for coordim == 2 or 3). Inspecting the formula above, it can
298 // be appreciated that the derivatives dx_i/dxi_j will be
299 // constant, if the coefficient of the non-linear term is zero.
300 //
301 // That is why for regular geometry, we require
302 //
303 // x_i^A - x_i^B + x_i^C - x_i^D = 0
304 //
305 // or equivalently
306 //
307 // x_i^A - x_i^B = x_i^D - x_i^C
308 //
309 // This corresponds to quadrilaterals which are paralellograms.
311 m_manifold[0] = 0;
312 m_manifold[1] = 1;
313 if (m_coordim == 3)
314 {
315 PointGeom e01, e21, norm;
316 e01.Sub(*m_verts[0], *m_verts[1]);
317 e21.Sub(*m_verts[3], *m_verts[1]);
318 norm.Mult(e01, e21);
319 int tmpi = 0;
320 double tmp = std::fabs(norm[0]);
321 if (tmp < fabs(norm[1]))
322 {
323 tmp = fabs(norm[1]);
324 tmpi = 1;
325 }
326 if (tmp < fabs(norm[2]))
327 {
328 tmpi = 2;
329 }
330 m_manifold[0] = (tmpi + 1) % 3;
331 m_manifold[1] = (tmpi + 2) % 3;
332 }
333
334 if (Gtype == eRegular)
335 {
337 for (int i = 0; i < m_verts.size(); ++i)
338 {
339 verts[i] = Array<OneD, NekDouble>(3);
340 m_verts[i]->GetCoords(verts[i]);
341 }
342 // a00 + a01 xi1 + a02 xi2 + a03 xi1 xi2
343 // a10 + a11 xi1 + a12 xi2 + a03 xi1 xi2
345 for (int i = 0; i < 2; i++)
346 {
347 unsigned int d = m_manifold[i];
349 // Karniadakis, Sherwin 2005, Appendix D
350 NekDouble A = verts[0][d];
351 NekDouble B = verts[1][d];
352 NekDouble D = verts[2][d];
353 NekDouble C = verts[3][d];
354 m_isoParameter[i][0] = 0.25 * (A + B + C + D); // 1
355 m_isoParameter[i][1] = 0.25 * (-A + B - C + D); // xi1
356 m_isoParameter[i][2] = 0.25 * (-A - B + C + D); // xi2
357 m_isoParameter[i][3] = 0.25 * (A - B - C + D); // xi1*xi2
358 NekDouble tmp =
359 fabs(m_isoParameter[i][1]) + fabs(m_isoParameter[i][2]);
360 if (fabs(m_isoParameter[i][3]) >
362 {
363 Gtype = eDeformed;
364 }
365 }
366 }
367
368 if (Gtype == eRegular)
369 {
371 }
372 else if (m_straightEdge)
373 {
375 }
376
378 Gtype, m_coordim, m_xmap, m_coeffs);
380 }
381}
382
383/**
384 * Note verts and edges are listed according to anticlockwise
385 * convention but points in _coeffs have to be in array format from
386 * left to right.
387 */
389{
390 // check to see if geometry structure is already filled
391 if (m_state != ePtsFilled)
392 {
393 int i, j, k;
394 int nEdgeCoeffs;
395
396 if (m_curve)
397 {
398 int npts = m_curve->m_points.size();
399 int nEdgePts = (int)sqrt(static_cast<NekDouble>(npts));
400 Array<OneD, NekDouble> tmp(npts);
401 Array<OneD, NekDouble> tmp2(m_xmap->GetTotPoints());
402 LibUtilities::PointsKey curveKey(nEdgePts, m_curve->m_ptype);
403
404 // Sanity checks:
405 // - Curved faces should have square number of points;
406 // - Each edge should have sqrt(npts) points.
407 ASSERTL0(nEdgePts * nEdgePts == npts,
408 "NUMPOINTS should be a square number in"
409 " quadrilteral " +
410 boost::lexical_cast<string>(m_globalID));
411
412 for (i = 0; i < kNedges; ++i)
413 {
414 ASSERTL0(m_edges[i]->GetXmap()->GetNcoeffs() == nEdgePts,
415 "Number of edge points does not correspond to "
416 "number of face points in quadrilateral " +
417 boost::lexical_cast<string>(m_globalID));
418 }
419
420 for (i = 0; i < m_coordim; ++i)
421 {
422 for (j = 0; j < npts; ++j)
423 {
424 tmp[j] = (m_curve->m_points[j]->GetPtr())[i];
425 }
426
427 // Interpolate m_curve points to GLL points
428 LibUtilities::Interp2D(curveKey, curveKey, tmp,
429 m_xmap->GetBasis(0)->GetPointsKey(),
430 m_xmap->GetBasis(1)->GetPointsKey(),
431 tmp2);
432
433 // Forwards transform to get coefficient space.
434 m_xmap->FwdTrans(tmp2, m_coeffs[i]);
435 }
436 }
437
438 // Now fill in edges.
440 Array<OneD, int> signArray;
441
442 for (i = 0; i < kNedges; i++)
443 {
444 m_edges[i]->FillGeom();
445 m_xmap->GetTraceToElementMap(i, mapArray, signArray, m_eorient[i]);
446
447 nEdgeCoeffs = m_edges[i]->GetXmap()->GetNcoeffs();
448
449 for (j = 0; j < m_coordim; j++)
450 {
451 for (k = 0; k < nEdgeCoeffs; k++)
452 {
453 m_coeffs[j][mapArray[k]] =
454 signArray[k] * (m_edges[i]->GetCoeffs(j))[k];
455 }
456 }
457 }
458
460 }
461}
462
464{
465 int i0, i1, j1, j2;
466 if (fabs(m_isoParameter[0][3]) >= fabs(m_isoParameter[1][3]))
467 {
468 i0 = 0;
469 i1 = 1;
470 }
471 else
472 {
473 i1 = 0;
474 i0 = 1;
475 m_straightEdge |= 2;
476 }
477 NekDouble gamma = m_isoParameter[i1][3] / m_isoParameter[i0][3];
478 std::vector<NekDouble> c(3);
479 for (int i = 0; i < 3; ++i)
480 {
481 c[i] = m_isoParameter[i1][i] - gamma * m_isoParameter[i0][i];
482 }
483 if (fabs(c[1]) >= fabs(c[2]))
484 {
485 j1 = 1;
486 j2 = 2;
487 }
488 else
489 {
490 j1 = 2;
491 j2 = 1;
492 m_straightEdge |= 4;
493 }
494 NekDouble beta = c[j2] / c[j1];
495 if (i0 == 1)
496 {
498 }
499 if (j1 == 2)
500 {
501 NekDouble temp = m_isoParameter[0][j1];
502 m_isoParameter[0][j1] = m_isoParameter[0][j2];
503 m_isoParameter[0][j2] = temp;
504 }
505 m_isoParameter[0][2] -= m_isoParameter[0][1] * beta;
506 m_isoParameter[1][0] = c[0];
507 m_isoParameter[1][1] = 1. / c[j1];
508 m_isoParameter[1][2] = beta;
509 m_isoParameter[1][3] = gamma;
510}
511
512int QuadGeom::v_GetDir(const int i, [[maybe_unused]] const int j) const
513{
514 return i % 2;
515}
516
517void QuadGeom::v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
518{
519 Geometry::v_Reset(curvedEdges, curvedFaces);
520 CurveMap::iterator it = curvedFaces.find(m_globalID);
521
522 if (it != curvedFaces.end())
523 {
524 m_curve = it->second;
525 }
526
527 for (int i = 0; i < 4; ++i)
528 {
529 m_edges[i]->Reset(curvedEdges, curvedFaces);
530 }
531
532 SetUpXmap();
533 SetUpCoeffs(m_xmap->GetNcoeffs());
534}
535
537{
538 if (!m_setupState)
539 {
540 for (int i = 0; i < 4; ++i)
541 {
542 m_edges[i]->Setup();
543 }
544 SetUpXmap();
545 SetUpCoeffs(m_xmap->GetNcoeffs());
546 m_setupState = true;
547 }
548}
549
550} // namespace Nektar::SpatialDomains
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
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:66
Array< OneD, int > m_manifold
Definition: Geometry2D.h:83
std::vector< StdRegions::Orientation > m_eorient
Definition: Geometry2D.h:81
void v_CalculateInverseIsoParam() override
Definition: Geometry2D.cpp:575
bool m_setupState
Wether or not the setup routines have been run.
Definition: Geometry.h:193
PointGeomSharedPtr GetVertex(int i) const
Returns vertex i of this object.
Definition: Geometry.h:355
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:687
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:384
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:433
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)
Definition: PointGeom.cpp:122
void Mult(PointGeom &a, PointGeom &b)
_this = a x b
Definition: PointGeom.cpp:131
void Rotate(PointGeom &a, int dir, NekDouble angle)
_this = rotation of a by angle 'angle' around axis dir
Definition: PointGeom.cpp:140
NekDouble dist(PointGeom &a)
return distance between this and input a
Definition: PointGeom.cpp:178
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:517
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:135
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:124
void v_GenGeomFactors() override
Definition: QuadGeom.cpp:260
static const int kNedges
Definition: QuadGeom.h:74
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:512
static StdRegions::Orientation GetEdgeOrientation(const SegGeom &edge1, const SegGeom &edge2)
Get the orientation of edge1.
Definition: SegGeom.cpp:194
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::vector< PointGeomSharedPtr > PointGeomVector
Definition: Geometry2D.h:61
std::shared_ptr< Curve > CurveSharedPtr
Definition: Curve.hpp:58
std::unordered_map< int, CurveSharedPtr > CurveMap
Definition: Curve.hpp:59
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:59
@ ePtsFilled
Geometric information has been generated.
std::vector< double > z(NPUPPER)
std::vector< double > d(NPUPPER *NPUPPER)
double NekDouble
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:294