Nektar++
SegGeom.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: SegGeom.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 <LibUtilities/Foundations/ManagerAccess.h> // for PointsManager, etc
42
44{
46{
48}
49
50SegGeom::SegGeom(int id, const int coordim, const PointGeomSharedPtr vertex[],
51 const CurveSharedPtr curve)
52 : Geometry1D(coordim)
53{
55 m_globalID = id;
57 m_curve = curve;
58
59 m_verts[0] = vertex[0];
60 m_verts[1] = vertex[1];
61}
62
64{
65 // From Geometry class
67
68 // info from EdgeComponent class
70 m_xmap = in.m_xmap;
71 SetUpCoeffs(m_xmap->GetNcoeffs());
72
73 // info from SegGeom class
75 m_verts[0] = in.m_verts[0];
76 m_verts[1] = in.m_verts[1];
77
78 m_state = in.m_state;
79}
80
82{
83 if (m_curve)
84 {
85 int npts = m_curve->m_points.size();
86 LibUtilities::PointsKey pkey(npts + 1,
90 }
91 else
92 {
97 }
98}
99
100/**
101 * \brief Generate a one dimensional space segment geometry where the vert[0]
102 * has the same x value and vert[1] is set to vert[0] plus the length of the
103 * original segment
104 **/
106{
108
109 // info about numbering
110 returnval->m_globalID = m_globalID;
111
112 // geometric information.
113 returnval->m_coordim = 1;
114 NekDouble x0 = (*m_verts[0])[0];
116 1, m_verts[0]->GetGlobalID(), x0, 0.0, 0.0);
117 vert0->SetGlobalID(vert0->GetGlobalID());
118 returnval->m_verts[0] = vert0;
119
120 // Get information to calculate length.
122 m_xmap->GetBase();
124 v.push_back(base[0]->GetPointsKey());
126
127 const Array<OneD, const NekDouble> jac = m_geomFactors->GetJac(v);
128
129 NekDouble len = 0.0;
130 if (jac.size() == 1)
131 {
132 len = jac[0] * 2.0;
133 }
134 else
135 {
136 Array<OneD, const NekDouble> w0 = base[0]->GetW();
137 len = 0.0;
138
139 for (int i = 0; i < jac.size(); ++i)
140 {
141 len += jac[i] * w0[i];
142 }
143 }
144 // Set up second vertex.
146 1, m_verts[1]->GetGlobalID(), x0 + len, 0.0, 0.0);
147 vert1->SetGlobalID(vert1->GetGlobalID());
148
149 returnval->m_verts[1] = vert1;
150
151 // at present just use previous m_xmap[0];
152 returnval->m_xmap = m_xmap;
153 returnval->SetUpCoeffs(m_xmap->GetNcoeffs());
154 returnval->m_state = eNotFilled;
155
156 return returnval;
157}
158
160{
161}
162
164{
166}
167
169 const Array<OneD, const NekDouble> &Lcoord)
170{
171 ASSERTL1(m_state == ePtsFilled, "Geometry is not in physical space");
172
173 Array<OneD, NekDouble> tmp(m_xmap->GetTotPoints());
174 m_xmap->BwdTrans(m_coeffs[i], tmp);
175
176 return m_xmap->PhysEvaluate(Lcoord, tmp);
177}
178
179/**
180 * @brief Get the orientation of @p edge1.
181 *
182 * If @p edge1 is connected to @p edge2 in the same direction as the points
183 * comprising @p edge1 then it is forward, otherwise it is backward.
184 *
185 * For example, assume @p edge1 is comprised of points 1 and 2, and @p edge2 is
186 * comprised of points 2 and 3, then @p edge1 is forward.
187 *
188 * If @p edge1 is comprised of points 2 and 1 and @p edge2 is comprised of
189 * points 3 and 2, then @p edge1 is backward.
190 *
191 * Since both edges are passed, it does not need any information from the
192 * EdgeComponent instance.
193 */
195 const SegGeom &edge2)
196{
198
199 if ((*edge1.GetVertex(0) == *edge2.GetVertex(0)) ||
200 (*edge1.GetVertex(0) == *edge2.GetVertex(1)))
201 {
202 // Backward direction. Vertex 0 is connected to edge 2.
203 returnval = StdRegions::eBackwards;
204 }
205 else if ((*edge1.GetVertex(1) != *edge2.GetVertex(0)) &&
206 (*edge1.GetVertex(1) != *edge2.GetVertex(1)))
207 {
208 // Not forward either, then we have a problem.
209 std::ostringstream errstrm;
210 errstrm << "Connected edges do not share a vertex. Edges ";
211 errstrm << edge1.GetGlobalID() << ", " << edge2.GetGlobalID();
212 ASSERTL0(false, errstrm.str());
213 }
214
215 return returnval;
216}
217
219{
220 if (!m_setupState)
221 {
223 }
224
226 {
229
230 if (m_xmap->GetBasisNumModes(0) != 2)
231 {
232 gType = eDeformed;
233 }
234
236 gType, m_coordim, m_xmap, m_coeffs);
238 }
239}
240
242{
243 if (m_state != ePtsFilled)
244 {
245 int i;
246
247 if (m_coordim > 0 && m_curve)
248 {
249 int npts = m_curve->m_points.size();
250 LibUtilities::PointsKey pkey(npts + 1,
252 Array<OneD, NekDouble> tmp(npts);
253
254 if (m_verts[0]->dist(*(m_curve->m_points[0])) >
256 {
257 std::string err =
258 "Vertex 0 is separated from first point by more than ";
259 std::stringstream strstrm;
260 strstrm << NekConstants::kVertexTheSameDouble << " in edge "
261 << m_globalID;
262 err += strstrm.str();
263 NEKERROR(ErrorUtil::ewarning, err.c_str());
264 }
265
266 if (m_verts[1]->dist(*(m_curve->m_points[npts - 1])) >
268 {
269 std::string err =
270 "Vertex 1 is separated from last point by more than ";
271 std::stringstream strstrm;
272 strstrm << NekConstants::kVertexTheSameDouble << " in edge "
273 << m_globalID;
274 err += strstrm.str();
275 NEKERROR(ErrorUtil::ewarning, err.c_str());
276 }
277
278 LibUtilities::PointsKey fkey(npts, m_curve->m_ptype);
280 LibUtilities::PointsManager()[fkey]->GetI(pkey);
281 NekVector<NekDouble> out(npts + 1);
282
283 for (int i = 0; i < m_coordim; ++i)
284 {
285 // Load up coordinate values into tmp
286 for (int j = 0; j < npts; ++j)
287 {
288 tmp[j] = (m_curve->m_points[j]->GetPtr())[i];
289 }
290
291 // Interpolate to GLL points
292 NekVector<NekDouble> in(npts, tmp, eWrapper);
293 out = (*I0) * in;
294
295 m_xmap->FwdTrans(out.GetPtr(), m_coeffs[i]);
296 }
297 }
298
299 for (i = 0; i < m_coordim; ++i)
300 {
301 m_coeffs[i][0] = (*m_verts[0])[i];
302 m_coeffs[i][1] = (*m_verts[1])[i];
303 }
304
306 }
307}
308
309void SegGeom::v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
310{
311 Geometry::v_Reset(curvedEdges, curvedFaces);
312 CurveMap::iterator it = curvedEdges.find(m_globalID);
313
314 if (it != curvedEdges.end())
315 {
316 m_curve = it->second;
317 }
318
319 SetUpXmap();
320 SetUpCoeffs(m_xmap->GetNcoeffs());
321}
322
324{
325 if (!m_setupState)
326 {
327 SetUpXmap();
328 SetUpCoeffs(m_xmap->GetNcoeffs());
329 m_setupState = true;
330 }
331}
332
334{
335 PointGeomSharedPtr returnval;
336
337 if (i >= 0 && i < kNverts)
338 {
339 returnval = m_verts[i];
340 }
341
342 return returnval;
343}
344
346{
347 return kNverts;
348}
349
352{
353 if (m_geomFactors->GetGtype() == eRegular)
354 {
355 xiOut = Array<OneD, NekDouble>(1, 0.0);
356
357 GetLocCoords(xs, xiOut);
358 ClampLocCoords(xiOut);
359
361 NekDouble tmp = 0;
362 for (int i = 0; i < m_coordim; ++i)
363 {
364 gloCoord[i] = GetCoord(i, xiOut);
365 tmp += (xs[i] - gloCoord[i]) * (xs[i] - gloCoord[i]);
366 }
367
368 return sqrt(tmp);
369 }
370 // If deformed edge then the inverse mapping is non-linear so need to
371 // numerically solve for the local coordinate
372 else if (m_geomFactors->GetGtype() == eDeformed)
373 {
374 Array<OneD, NekDouble> xi(1, 0.0);
375
376 // Armijo constants:
377 // https://en.wikipedia.org/wiki/Backtracking_line_search
378 const NekDouble c1 = 1e-4, c2 = 0.9;
379
380 int dim = GetCoordim();
381 int nq = m_xmap->GetTotPoints();
382
383 Array<OneD, Array<OneD, NekDouble>> x(dim), xder(dim), xder2(dim);
384 // Get x,y,z phys values from coefficients
385 for (int i = 0; i < dim; ++i)
386 {
387 x[i] = Array<OneD, NekDouble>(nq);
388 xder[i] = Array<OneD, NekDouble>(nq);
389 xder2[i] = Array<OneD, NekDouble>(nq);
390
391 m_xmap->BwdTrans(m_coeffs[i], x[i]);
392 }
393
394 NekDouble fx_prev = std::numeric_limits<NekDouble>::max();
395
396 // Minimisation loop (Quasi-newton method)
397 for (int i = 0; i < NekConstants::kNewtonIterations; ++i)
398 {
399 // Compute the objective function, f(x_k) and its derivatives
403 NekDouble fx = 0, fxp = 0, fxp2 = 0, xcDiff = 0;
404 for (int j = 0; j < dim; ++j)
405 {
406 xc[j] = m_xmap->PhysEvaluate(xi, x[j], xc_der[j], xc_der2[j]);
407
408 xcDiff = xc[j] - xs[j];
409 // Objective function is the distance to the search point
410 fx += xcDiff * xcDiff;
411 fxp += xc_der[j][0] * xcDiff;
412 fxp2 += xc_der2[j][0] * xcDiff + xc_der[j][0] * xc_der[j][0];
413 }
414
415 fxp *= 2;
416 fxp2 *= 2;
417
418 // Check for convergence
419 if (std::abs(fx - fx_prev) < 1e-12)
420 {
421 fx_prev = fx;
422 break;
423 }
424 else
425 {
426 fx_prev = fx;
427 }
428
429 NekDouble gamma = 1.0;
430 bool conv = false;
431
432 // Search direction: Newton's method
433 NekDouble pk = -fxp / fxp2;
434
435 // Perform backtracking line search
436 while (gamma > 1e-10)
437 {
438 Array<OneD, NekDouble> xi_pk(1);
439 xi_pk[0] = xi[0] + pk * gamma;
440
441 if (xi_pk[0] < -1.0 || xi_pk[0] > 1.0)
442 {
443 gamma /= 2.0;
444 continue;
445 }
446
447 Array<OneD, NekDouble> xc_pk(dim);
449 NekDouble fx_pk = 0, fxp_pk = 0, xc_pkDiff = 0;
450 for (int j = 0; j < dim; ++j)
451 {
452 xc_pk[j] = m_xmap->PhysEvaluate(xi_pk, x[j], xc_der_pk[j]);
453
454 xc_pkDiff = xc_pk[j] - xs[j];
455 fx_pk += xc_pkDiff * xc_pkDiff;
456 fxp_pk += xc_der_pk[j][0] * xc_pkDiff;
457 }
458
459 fxp_pk *= 2;
460
461 // Check Wolfe conditions using Armijo constants
462 // https://en.wikipedia.org/wiki/Wolfe_conditions
463 if ((fx_pk - (fx + c1 * gamma * pk * fxp)) <
464 std::numeric_limits<NekDouble>::epsilon() &&
465 (-pk * fxp_pk + c2 * pk * fxp) <
466 std::numeric_limits<NekDouble>::epsilon())
467 {
468 conv = true;
469 break;
470 }
471
472 gamma /= 2.0;
473 }
474
475 if (!conv)
476 {
477 break;
478 }
479
480 xi[0] += gamma * pk;
481 }
482
483 xiOut = xi;
484 return sqrt(fx_prev);
485 }
486 else
487 {
488 ASSERTL0(false, "Geometry type unknown")
489 }
490
491 return -1.0;
492}
493
494} // 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.
Array< OneD, DataType > & GetPtr()
Definition: NekVector.cpp:217
1D geometry information
Definition: Geometry1D.h:53
NekDouble GetCoord(const int i, const Array< OneD, const NekDouble > &Lcoord)
Given local collapsed coordinate Lcoord, return the value of physical coordinate in direction i.
Definition: Geometry.h:552
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
NekDouble GetLocCoords(const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
Determine the local collapsed coordinates that correspond to a given Cartesian coordinate for this ge...
Definition: Geometry.h:542
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
bool ClampLocCoords(Array< OneD, NekDouble > &locCoord, NekDouble tol=std::numeric_limits< NekDouble >::epsilon())
Clamp local coords to be within standard regions [-1, 1]^dim.
Definition: Geometry.cpp:550
int GetGlobalID(void) const
Get the ID of this object.
Definition: Geometry.h:324
int GetCoordim() const
Return the coordinate dimension of this object (i.e. the dimension of the space in which this object ...
Definition: Geometry.h:281
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
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:185
int m_coordim
Coordinate dimension of this geometry object.
Definition: Geometry.h:183
PointGeomSharedPtr v_GetVertex(const int i) const override
Definition: SegGeom.cpp:333
int v_GetNumVerts() const override
Get the number of vertices of this object.
Definition: SegGeom.cpp:345
virtual LibUtilities::ShapeType v_GetShapeType() const
Definition: SegGeom.cpp:163
void v_FillGeom() override
Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
Definition: SegGeom.cpp:241
SegGeomSharedPtr GenerateOneSpaceDimGeom(void)
Generate a one dimensional space segment geometry where the vert[0] has the same x value and vert[1] ...
Definition: SegGeom.cpp:105
CurveSharedPtr m_curve
Boolean indicating whether object owns the data.
Definition: SegGeom.h:93
void v_GenGeomFactors() override
Definition: SegGeom.cpp:218
SpatialDomains::PointGeomSharedPtr m_verts[kNverts]
Definition: SegGeom.h:76
static StdRegions::Orientation GetEdgeOrientation(const SegGeom &edge1, const SegGeom &edge2)
Get the orientation of edge1.
Definition: SegGeom.cpp:194
NekDouble v_FindDistance(const Array< OneD, const NekDouble > &xs, Array< OneD, NekDouble > &xi) override
Definition: SegGeom.cpp:350
static const int kNverts
Definition: SegGeom.h:73
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: SegGeom.cpp:168
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: SegGeom.cpp:309
PointsManagerT & PointsManager(void)
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:231
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48
static const NekDouble kVertexTheSameDouble
static const unsigned int kNewtonIterations
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
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:57
@ eNotFilled
Geometric information has not been generated.
@ ePtsFilled
Geometric information has been generated.
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:298
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294