Nektar++
TriGeom.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: TriGeom.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
40
43
44using namespace std;
45
47{
48
50{
52}
53
54TriGeom::TriGeom(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 + TriGeom::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
78
79 m_coordim = edges[0]->GetVertex(0)->GetCoordim();
80 ASSERTL0(m_coordim > 1, "Cannot call function with dim == 1");
81}
82
84{
85 // From Geometry
87
88 // From TriFaceComponent
90
91 // From TriGeom
92 m_verts = in.m_verts;
93 m_edges = in.m_edges;
94 for (int i = 0; i < kNedges; i++)
95 {
96 m_eorient[i] = in.m_eorient[i];
97 }
98}
99
101{
102}
103
105 const Array<OneD, const NekDouble> &Lcoord)
106{
107 ASSERTL1(m_state == ePtsFilled, "Geometry is not in physical space");
108
109 Array<OneD, NekDouble> tmp(m_xmap->GetTotPoints());
110 m_xmap->BwdTrans(m_coeffs[i], tmp);
111
112 return m_xmap->PhysEvaluate(Lcoord, tmp);
113}
114
116 const TriGeom &face2,
117 bool doRot, int dir,
118 NekDouble angle,
119 NekDouble tol)
120{
121 return GetFaceOrientation(face1.m_verts, face2.m_verts, doRot, dir, angle,
122 tol);
123}
124
126 const PointGeomVector &face1, const PointGeomVector &face2, bool doRot,
127 int dir, NekDouble angle, NekDouble tol)
128{
129 int i, j, vmap[3] = {-1, -1, -1};
130
131 if (doRot)
132 {
133 PointGeom rotPt;
134
135 for (i = 0; i < 3; ++i)
136 {
137 rotPt.Rotate((*face1[i]), dir, angle);
138 for (j = 0; j < 3; ++j)
139 {
140 if (rotPt.dist(*face2[j]) < tol)
141 {
142 vmap[j] = i;
143 break;
144 }
145 }
146 }
147 }
148 else
149 {
150
151 NekDouble x, y, z, x1, y1, z1, cx = 0.0, cy = 0.0, cz = 0.0;
152
153 // For periodic faces, we calculate the vector between the centre
154 // points of the two faces. (For connected faces this will be
155 // zero). We can then use this to determine alignment later in the
156 // algorithm.
157 for (i = 0; i < 3; ++i)
158 {
159 cx += (*face2[i])(0) - (*face1[i])(0);
160 cy += (*face2[i])(1) - (*face1[i])(1);
161 cz += (*face2[i])(2) - (*face1[i])(2);
162 }
163 cx /= 3;
164 cy /= 3;
165 cz /= 3;
166
167 // Now construct a mapping which takes us from the vertices of one
168 // face to the other. That is, vertex j of face2 corresponds to
169 // vertex vmap[j] of face1.
170 for (i = 0; i < 3; ++i)
171 {
172 x = (*face1[i])(0);
173 y = (*face1[i])(1);
174 z = (*face1[i])(2);
175 for (j = 0; j < 3; ++j)
176 {
177 x1 = (*face2[j])(0) - cx;
178 y1 = (*face2[j])(1) - cy;
179 z1 = (*face2[j])(2) - cz;
180 if (sqrt((x1 - x) * (x1 - x) + (y1 - y) * (y1 - y) +
181 (z1 - z) * (z1 - z)) < 1e-8)
182 {
183 vmap[j] = i;
184 break;
185 }
186 }
187 }
188 }
189
190 if (vmap[1] == (vmap[0] + 1) % 3)
191 {
192 switch (vmap[0])
193 {
194 case 0:
196 break;
197 case 1:
199 break;
200 case 2:
202 break;
203 }
204 }
205 else
206 {
207 switch (vmap[0])
208 {
209 case 0:
211 break;
212 case 1:
214 break;
215 case 2:
217 break;
218 }
219 }
220
221 ASSERTL0(false, "Unable to determine triangle orientation");
223}
224
226{
227 if (!m_setupState)
228 {
230 }
231
233 {
234 GeomType Gtype = eRegular;
235
237
238 // check to see if expansions are linear
239 m_straightEdge = 1;
240 if (m_xmap->GetBasisNumModes(0) != 2 ||
241 m_xmap->GetBasisNumModes(1) != 2)
242 {
243 Gtype = eDeformed;
244 m_straightEdge = 0;
245 }
246
248 m_manifold[0] = 0;
249 m_manifold[1] = 1;
250 if (m_coordim == 3)
251 {
252 PointGeom e01, e21, norm;
253 e01.Sub(*m_verts[0], *m_verts[1]);
254 e21.Sub(*m_verts[2], *m_verts[1]);
255 norm.Mult(e01, e21);
256 int tmpi = 0;
257 double tmp = std::fabs(norm[0]);
258 if (tmp < fabs(norm[1]))
259 {
260 tmp = fabs(norm[1]);
261 tmpi = 1;
262 }
263 if (tmp < fabs(norm[2]))
264 {
265 tmpi = 2;
266 }
267 m_manifold[0] = (tmpi + 1) % 3;
268 m_manifold[1] = (tmpi + 2) % 3;
269 }
270 if (Gtype == eRegular)
271 {
273 for (int i = 0; i < m_verts.size(); ++i)
274 {
275 verts[i] = Array<OneD, NekDouble>(3);
276 m_verts[i]->GetCoords(verts[i]);
277 }
278 // a00 + a01 xi1 + a02 xi2
279 // a10 + a11 xi1 + a12 xi2
281 for (int i = 0; i < 2; ++i)
282 {
283 unsigned int d = m_manifold[i];
285 NekDouble A = verts[0][d];
286 NekDouble B = verts[1][d];
287 NekDouble C = verts[2][d];
288 m_isoParameter[i][0] = 0.5 * (B + C); // 1
289 m_isoParameter[i][1] = 0.5 * (-A + B); // xi1
290 m_isoParameter[i][2] = 0.5 * (-A + C); // xi2
291 }
293 }
294
296 Gtype, m_coordim, m_xmap, m_coeffs);
297
299 }
300}
301
302/**
303 * Note verts and edges are listed according to anticlockwise
304 * convention but points in _coeffs have to be in array format from
305 * left to right.
306 */
308{
309 // check to see if geometry structure is already filled
310 if (m_state == ePtsFilled)
311 {
312 return;
313 }
314
315 int i, j, k;
316 int nEdgeCoeffs = m_xmap->GetTraceNcoeffs(0);
317
318 if (m_curve)
319 {
321 2, m_curve->m_ptype)]
322 ->GetPointsDim();
323
324 // Deal with 2D points type separately
325 // (e.g. electrostatic or Fekete points) to 1D tensor
326 // product.
327 if (pdim == 2)
328 {
329 int N = m_curve->m_points.size();
330 int nEdgePts =
331 (-1 + (int)sqrt(static_cast<NekDouble>(8 * N + 1))) / 2;
332
333 ASSERTL0(nEdgePts * (nEdgePts + 1) / 2 == N,
334 "NUMPOINTS should be a triangle number for"
335 " triangle curved face " +
336 boost::lexical_cast<string>(m_globalID));
337
338 // Sanity check 1: are curved vertices consistent with
339 // triangle vertices?
340 for (i = 0; i < 3; ++i)
341 {
342 NekDouble dist = m_verts[i]->dist(*(m_curve->m_points[i]));
344 {
345 std::stringstream ss;
346 ss << "Curved vertex " << i << " of triangle " << m_globalID
347 << " is separated from expansion vertex by"
348 << " more than " << NekConstants::kVertexTheSameDouble
349 << " (dist = " << dist << ")";
350 NEKERROR(ErrorUtil::ewarning, ss.str().c_str());
351 }
352 }
353
354 // Sanity check 2: are curved edges from the face curvature
355 // consistent with curved edges?
356 for (i = 0; i < kNedges; ++i)
357 {
358 CurveSharedPtr edgeCurve = m_edges[i]->GetCurve();
359
360 ASSERTL0(edgeCurve->m_points.size() == nEdgePts,
361 "Number of edge points does not correspond "
362 "to number of face points in triangle " +
363 boost::lexical_cast<string>(m_globalID));
364
365 const int offset = 3 + i * (nEdgePts - 2);
366 NekDouble maxDist = 0.0;
367
368 // Account for different ordering of nodal coordinates
369 // vs. Cartesian ordering of element.
371
372 if (i == 2)
373 {
374 orient = orient == StdRegions::eForwards
377 }
378
379 if (orient == StdRegions::eForwards)
380 {
381 for (j = 0; j < nEdgePts - 2; ++j)
382 {
383 NekDouble dist = m_curve->m_points[offset + j]->dist(
384 *(edgeCurve->m_points[j + 1]));
385 maxDist = dist > maxDist ? dist : maxDist;
386 }
387 }
388 else
389 {
390 for (j = 0; j < nEdgePts - 2; ++j)
391 {
392 NekDouble dist = m_curve->m_points[offset + j]->dist(
393 *(edgeCurve->m_points[nEdgePts - 2 - j]));
394 maxDist = dist > maxDist ? dist : maxDist;
395 }
396 }
397
399 {
400 std::stringstream ss;
401 ss << "Curved edge " << i << " of triangle " << m_globalID
402 << " has a point separated from edge interior"
403 << " points by more than "
405 << " (maxdist = " << maxDist << ")";
406 NEKERROR(ErrorUtil::ewarning, ss.str().c_str());
407 }
408 }
409
413 nEdgePts, LibUtilities::eGaussRadauMAlpha1Beta0);
415 P0);
417 P1);
419 max(nEdgePts * nEdgePts, m_xmap->GetTotPoints()));
420 Array<OneD, NekDouble> tmp(nEdgePts * nEdgePts);
421
422 for (i = 0; i < m_coordim; ++i)
423 {
424 // Create a StdNodalTriExp.
427 AllocateSharedPtr(T0, T1, m_curve->m_ptype);
428
429 for (j = 0; j < N; ++j)
430 {
431 phys[j] = (m_curve->m_points[j]->GetPtr())[i];
432 }
433
434 t->BwdTrans(phys, tmp);
435
436 // Interpolate points to standard region.
438 P0, P1, tmp, m_xmap->GetBasis(0)->GetPointsKey(),
439 m_xmap->GetBasis(1)->GetPointsKey(), phys);
440
441 // Forwards transform to get coefficient space.
442 m_xmap->FwdTrans(phys, m_coeffs[i]);
443 }
444 }
445 else if (pdim == 1)
446 {
447 int npts = m_curve->m_points.size();
448 int nEdgePts = (int)sqrt(static_cast<NekDouble>(npts));
449 Array<OneD, NekDouble> tmp(npts);
450 Array<OneD, NekDouble> phys(m_xmap->GetTotPoints());
451 LibUtilities::PointsKey curveKey(nEdgePts, m_curve->m_ptype);
452
453 // Sanity checks:
454 // - Curved faces should have square number of points;
455 // - Each edge should have sqrt(npts) points.
456 ASSERTL0(nEdgePts * nEdgePts == npts,
457 "NUMPOINTS should be a square number for"
458 " triangle " +
459 boost::lexical_cast<string>(m_globalID));
460
461 for (i = 0; i < kNedges; ++i)
462 {
463 ASSERTL0(m_edges[i]->GetXmap()->GetNcoeffs() == nEdgePts,
464 "Number of edge points does not correspond to "
465 "number of face points in triangle " +
466 boost::lexical_cast<string>(m_globalID));
467 }
468
469 for (i = 0; i < m_coordim; ++i)
470 {
471 for (j = 0; j < npts; ++j)
472 {
473 tmp[j] = (m_curve->m_points[j]->GetPtr())[i];
474 }
475
476 // Interpolate curve points to standard triangle
477 // points.
478 LibUtilities::Interp2D(curveKey, curveKey, tmp,
479 m_xmap->GetBasis(0)->GetPointsKey(),
480 m_xmap->GetBasis(1)->GetPointsKey(),
481 phys);
482
483 // Forwards transform to get coefficient space.
484 m_xmap->FwdTrans(phys, m_coeffs[i]);
485 }
486 }
487 else
488 {
489 ASSERTL0(false, "Only 1D/2D points distributions "
490 "supported.");
491 }
492 }
493
494 Array<OneD, unsigned int> mapArray(nEdgeCoeffs);
495 Array<OneD, int> signArray(nEdgeCoeffs);
496
497 for (i = 0; i < kNedges; i++)
498 {
499 m_edges[i]->FillGeom();
500 m_xmap->GetTraceToElementMap(i, mapArray, signArray, m_eorient[i]);
501
502 nEdgeCoeffs = m_edges[i]->GetXmap()->GetNcoeffs();
503
504 for (j = 0; j < m_coordim; j++)
505 {
506 for (k = 0; k < nEdgeCoeffs; k++)
507 {
508 m_coeffs[j][mapArray[k]] =
509 signArray[k] * m_edges[i]->GetCoeffs(j)[k];
510 }
511 }
512 }
513
515}
516
517int TriGeom::v_GetDir(const int i, [[maybe_unused]] const int j) const
518{
519 return i == 0 ? 0 : 1;
520}
521
522void TriGeom::v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
523{
524 Geometry::v_Reset(curvedEdges, curvedFaces);
525 CurveMap::iterator it = curvedFaces.find(m_globalID);
526
527 if (it != curvedFaces.end())
528 {
529 m_curve = it->second;
530 }
531
532 for (int i = 0; i < 3; ++i)
533 {
534 m_edges[i]->Reset(curvedEdges, curvedFaces);
535 }
536
537 SetUpXmap();
538 SetUpCoeffs(m_xmap->GetNcoeffs());
539}
540
542{
543 if (!m_setupState)
544 {
545 for (int i = 0; i < 3; ++i)
546 {
547 m_edges[i]->Setup();
548 }
549 SetUpXmap();
550 SetUpCoeffs(m_xmap->GetNcoeffs());
551 m_setupState = true;
552 }
553}
554
556{
557 int order0 = m_edges[0]->GetXmap()->GetBasis(0)->GetNumModes();
558 int order1 =
559 max(order0, max(m_edges[1]->GetXmap()->GetBasis(0)->GetNumModes(),
560 m_edges[2]->GetXmap()->GetBasis(0)->GetNumModes()));
561
562 const LibUtilities::BasisKey B0(
564 LibUtilities::PointsKey(order0 + 1,
566 const LibUtilities::BasisKey B1(
568 LibUtilities::PointsKey(order1, LibUtilities::eGaussRadauMAlpha1Beta0));
569
571}
572
573} // namespace Nektar::SpatialDomains
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
#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
static StdRegions::Orientation GetEdgeOrientation(const SegGeom &edge1, const SegGeom &edge2)
Get the orientation of edge1.
Definition: SegGeom.cpp:194
static const int kNedges
Get the orientation of face1.
Definition: TriGeom.h:70
void v_GenGeomFactors() override
Definition: TriGeom.cpp:225
int v_GetDir(const int faceidx, const int facedir) const override
Returns the element coordinate direction corresponding to a given face coordinate direction.
Definition: TriGeom.cpp:517
static StdRegions::Orientation GetFaceOrientation(const TriGeom &face1, const TriGeom &face2, bool doRot, int dir, NekDouble angle, NekDouble tol)
Definition: TriGeom.cpp:115
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: TriGeom.cpp:104
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: TriGeom.cpp:522
PointsManagerT & PointsManager(void)
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
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:49
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:42
@ eOrtho_B
Principle Orthogonal Functions .
Definition: BasisType.h:44
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48
static const NekDouble kVertexTheSameDouble
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::shared_ptr< StdNodalTriExp > StdNodalTriExpSharedPtr
std::vector< double > z(NPUPPER)
std::vector< double > d(NPUPPER *NPUPPER)
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294