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
46namespace Nektar
47{
48namespace SpatialDomains
49{
50
52{
54}
55
56QuadGeom::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 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 m_straightEdge = true;
283 if ((m_xmap->GetBasisNumModes(0) != 2) ||
284 (m_xmap->GetBasisNumModes(1) != 2))
285 {
286 Gtype = eDeformed;
287 m_straightEdge = false;
288 }
289
290 // For linear expansions, the mapping from standard to local
291 // element is given by the relation:
292 // x_i = 0.25 * [ ( x_i^A + x_i^B + x_i^C + x_i^D) +
293 // (-x_i^A + x_i^B + x_i^C - x_i^D)*xi_1 +
294 // (-x_i^A - x_i^B + x_i^C + x_i^D)*xi_2 +
295 // ( x_i^A - x_i^B + x_i^C - x_i^D)*xi_1*xi_2 ]
296 //
297 // The jacobian of the transformation and the metric terms
298 // dxi_i/dx_j, involve only terms of the form dx_i/dxi_j (both
299 // for coordim == 2 or 3). Inspecting the formula above, it can
300 // be appreciated that the derivatives dx_i/dxi_j will be
301 // constant, if the coefficient of the non-linear term is zero.
302 //
303 // That is why for regular geometry, we require
304 //
305 // x_i^A - x_i^B + x_i^C - x_i^D = 0
306 //
307 // or equivalently
308 //
309 // x_i^A - x_i^B = x_i^D - x_i^C
310 //
311 // This corresponds to quadrilaterals which are paralellograms.
313 m_manifold[0] = 0;
314 m_manifold[1] = 1;
315 if (m_coordim == 3)
316 {
317 PointGeom e01, e21, norm;
318 e01.Sub(*m_verts[0], *m_verts[1]);
319 e21.Sub(*m_verts[3], *m_verts[1]);
320 norm.Mult(e01, e21);
321 int tmpi = 0;
322 double tmp = std::fabs(norm[0]);
323 if (tmp < fabs(norm[1]))
324 {
325 tmp = fabs(norm[1]);
326 tmpi = 1;
327 }
328 if (tmp < fabs(norm[2]))
329 {
330 tmpi = 2;
331 }
332 m_manifold[0] = (tmpi + 1) % 3;
333 m_manifold[1] = (tmpi + 2) % 3;
334 }
335
336 if (Gtype == eRegular)
337 {
339 for (int i = 0; i < m_verts.size(); ++i)
340 {
341 verts[i] = Array<OneD, NekDouble>(3);
342 m_verts[i]->GetCoords(verts[i]);
343 }
344 // a00 + a01 xi1 + a02 xi2 + a03 xi1 xi2
345 // a10 + a11 xi1 + a12 xi2 + a03 xi1 xi2
347 for (int i = 0; i < 2; i++)
348 {
349 unsigned int d = m_manifold[i];
351 // Karniadakis, Sherwin 2005, Appendix D
352 NekDouble A = verts[0][d];
353 NekDouble B = verts[1][d];
354 NekDouble D = verts[2][d];
355 NekDouble C = verts[3][d];
356 m_isoParameter[i][0] = 0.25 * (A + B + C + D); // 1
357 m_isoParameter[i][1] = 0.25 * (-A + B - C + D); // xi1
358 m_isoParameter[i][2] = 0.25 * (-A - B + C + D); // xi2
359 m_isoParameter[i][3] = 0.25 * (A - B - C + D); // xi1*xi2
360 NekDouble tmp =
361 fabs(m_isoParameter[i][1]) + fabs(m_isoParameter[i][2]);
362 if (fabs(m_isoParameter[i][3]) >
364 {
365 Gtype = eDeformed;
366 }
367 }
368 }
369
370 if (Gtype == eRegular)
371 {
373 }
374
376 Gtype, m_coordim, m_xmap, m_coeffs);
378 }
379}
380
381/**
382 * Note verts and edges are listed according to anticlockwise
383 * convention but points in _coeffs have to be in array format from
384 * left to right.
385 */
387{
388 // check to see if geometry structure is already filled
389 if (m_state != ePtsFilled)
390 {
391 int i, j, k;
392 int nEdgeCoeffs;
393
394 if (m_curve)
395 {
396 int npts = m_curve->m_points.size();
397 int nEdgePts = (int)sqrt(static_cast<NekDouble>(npts));
398 Array<OneD, NekDouble> tmp(npts);
399 Array<OneD, NekDouble> tmp2(m_xmap->GetTotPoints());
400 LibUtilities::PointsKey curveKey(nEdgePts, m_curve->m_ptype);
401
402 // Sanity checks:
403 // - Curved faces should have square number of points;
404 // - Each edge should have sqrt(npts) points.
405 ASSERTL0(nEdgePts * nEdgePts == npts,
406 "NUMPOINTS should be a square number in"
407 " quadrilteral " +
408 boost::lexical_cast<string>(m_globalID));
409
410 for (i = 0; i < kNedges; ++i)
411 {
412 ASSERTL0(m_edges[i]->GetXmap()->GetNcoeffs() == nEdgePts,
413 "Number of edge points does not correspond to "
414 "number of face points in quadrilateral " +
415 boost::lexical_cast<string>(m_globalID));
416 }
417
418 for (i = 0; i < m_coordim; ++i)
419 {
420 for (j = 0; j < npts; ++j)
421 {
422 tmp[j] = (m_curve->m_points[j]->GetPtr())[i];
423 }
424
425 // Interpolate m_curve points to GLL points
426 LibUtilities::Interp2D(curveKey, curveKey, tmp,
427 m_xmap->GetBasis(0)->GetPointsKey(),
428 m_xmap->GetBasis(1)->GetPointsKey(),
429 tmp2);
430
431 // Forwards transform to get coefficient space.
432 m_xmap->FwdTrans(tmp2, m_coeffs[i]);
433 }
434 }
435
436 // Now fill in edges.
438 Array<OneD, int> signArray;
439
440 for (i = 0; i < kNedges; i++)
441 {
442 m_edges[i]->FillGeom();
443 m_xmap->GetTraceToElementMap(i, mapArray, signArray, m_eorient[i]);
444
445 nEdgeCoeffs = m_edges[i]->GetXmap()->GetNcoeffs();
446
447 for (j = 0; j < m_coordim; j++)
448 {
449 for (k = 0; k < nEdgeCoeffs; k++)
450 {
451 m_coeffs[j][mapArray[k]] =
452 signArray[k] * (m_edges[i]->GetCoeffs(j))[k];
453 }
454 }
455 }
456
458 }
459}
460
461int QuadGeom::v_GetDir(const int i, const int j) const
462{
463 boost::ignore_unused(j); // required in 3D shapes
464
465 return i % 2;
466}
467
468void QuadGeom::v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
469{
470 Geometry::v_Reset(curvedEdges, curvedFaces);
471 CurveMap::iterator it = curvedFaces.find(m_globalID);
472
473 if (it != curvedFaces.end())
474 {
475 m_curve = it->second;
476 }
477
478 for (int i = 0; i < 4; ++i)
479 {
480 m_edges[i]->Reset(curvedEdges, curvedFaces);
481 }
482
483 SetUpXmap();
484 SetUpCoeffs(m_xmap->GetNcoeffs());
485}
486
488{
489 if (!m_setupState)
490 {
491 for (int i = 0; i < 4; ++i)
492 {
493 m_edges[i]->Setup();
494 }
495 SetUpXmap();
496 SetUpCoeffs(m_xmap->GetNcoeffs());
497 m_setupState = true;
498 }
499}
500
501} // namespace SpatialDomains
502} // 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:47
Defines a specification for a set of points.
Definition: Points.h:55
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
2D geometry information
Definition: Geometry2D.h:69
Array< OneD, int > m_manifold
Definition: Geometry2D.h:86
std::vector< StdRegions::Orientation > m_eorient
Definition: Geometry2D.h:84
virtual void v_CalculateInverseIsoParam() override
Definition: Geometry2D.cpp:625
bool m_setupState
Wether or not the setup routines have been run.
Definition: Geometry.h:196
PointGeomSharedPtr GetVertex(int i) const
Returns vertex i of this object.
Definition: Geometry.h:358
GeomState m_state
Enumeration to dictate whether coefficients are filled.
Definition: Geometry.h:194
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition: Geometry.h:690
Array< OneD, Array< OneD, NekDouble > > m_isoParameter
Definition: Geometry.h:207
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:390
LibUtilities::ShapeType m_shapeType
Type of shape.
Definition: Geometry.h:200
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:204
GeomState m_geomFactorsState
State of the geometric factors.
Definition: Geometry.h:190
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:192
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.h:436
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:188
int m_coordim
Coordinate dimension of this geometry object.
Definition: Geometry.h:186
void Sub(PointGeom &a, PointGeom &b)
Definition: PointGeom.cpp:124
void Mult(PointGeom &a, PointGeom &b)
_this = a x b
Definition: PointGeom.cpp:133
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:468
virtual void v_Setup() override
Definition: QuadGeom.cpp:487
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:386
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:461
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:103
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:53
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50
static const NekDouble kNekZeroTol
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.
std::vector< double > z(NPUPPER)
std::vector< double > d(NPUPPER *NPUPPER)
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