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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description:
33 //
34 //
35 ////////////////////////////////////////////////////////////////////////////////
36 
37 #include <SpatialDomains/TriGeom.h>
41 
42 #include <SpatialDomains/SegGeom.h>
44 
45 
46 namespace Nektar
47 {
48  namespace SpatialDomains
49  {
50  /**
51  *
52  */
54  {
56  }
57 
58  /**
59  *
60  */
61  TriGeom::TriGeom(int id, const int coordim):
62  Geometry2D(coordim), m_fid(id)
63  {
64  m_globalID = m_fid;
65  SetUpXmap();
66  SetUpCoeffs(m_xmap->GetNcoeffs());
67  }
68 
69  /**
70  *
71  */
72  TriGeom::TriGeom(const int id,
73  const PointGeomSharedPtr verts[],
74  const SegGeomSharedPtr edges[],
75  const StdRegions::Orientation eorient[]):
76  Geometry2D(verts[0]->GetCoordim()),
77  m_fid(id)
78  {
79  m_globalID = m_fid;
81 
82  /// Copy the vert shared pointers.
83  m_verts.insert(m_verts.begin(), verts, verts+TriGeom::kNverts);
84 
85  /// Copy the edge shared pointers.
86  m_edges.insert(m_edges.begin(), edges, edges+TriGeom::kNedges);
87 
88  for (int j=0; j<kNedges; ++j)
89  {
90  m_eorient[j] = eorient[j];
91  }
92 
93  m_coordim = verts[0]->GetCoordim();
94  ASSERTL0(m_coordim > 1,
95  "Cannot call function with dim == 1");
96 
97  SetUpXmap();
98  SetUpCoeffs(m_xmap->GetNcoeffs());
99  }
100 
101 
102  /**
103  *
104  */
105  TriGeom::TriGeom(const int id, const SegGeomSharedPtr edges[],
106  const StdRegions::Orientation eorient[]):
107  Geometry2D(edges[0]->GetVertex(0)->GetCoordim()),
108  m_fid(id)
109  {
110  m_globalID = m_fid;
112 
113  /// Copy the edge shared pointers.
114  m_edges.insert(m_edges.begin(), edges, edges+TriGeom::kNedges);
115 
116  for(int j=0; j <kNedges; ++j)
117  {
118  if(eorient[j] == StdRegions::eForwards)
119  {
120  m_verts.push_back(edges[j]->GetVertex(0));
121  }
122  else
123  {
124  m_verts.push_back(edges[j]->GetVertex(1));
125  }
126  }
127 
128  for (int j=0; j<kNedges; ++j)
129  {
130  m_eorient[j] = eorient[j];
131  }
132 
133  m_coordim = edges[0]->GetVertex(0)->GetCoordim();
134  ASSERTL0(m_coordim > 1,"Cannot call function with dim == 1");
135 
136  SetUpXmap();
137  SetUpCoeffs(m_xmap->GetNcoeffs());
138  }
139 
140 
141  /**
142  *
143  */
145  const int id,
146  const SegGeomSharedPtr edges[],
147  const StdRegions::Orientation eorient[],
148  const CurveSharedPtr &curve)
149  : Geometry2D(edges[0]->GetVertex(0)->GetCoordim()),
150  m_fid(id),
151  m_curve(curve)
152  {
153  m_globalID = m_fid;
155 
156  /// Copy the edge shared pointers.
157  m_edges.insert(m_edges.begin(), edges, edges+TriGeom::kNedges);
158 
159  for(int j=0; j <kNedges; ++j)
160  {
161  if(eorient[j] == StdRegions::eForwards)
162  {
163  m_verts.push_back(edges[j]->GetVertex(0));
164  }
165  else
166  {
167  m_verts.push_back(edges[j]->GetVertex(1));
168  }
169  }
170 
171  for (int j=0; j<kNedges; ++j)
172  {
173  m_eorient[j] = eorient[j];
174  }
175 
176  m_coordim = edges[0]->GetVertex(0)->GetCoordim();
177  ASSERTL0(m_coordim > 1,"Cannot call function with dim == 1");
178 
179  SetUpXmap();
180  SetUpCoeffs(m_xmap->GetNcoeffs());
181  }
182 
183  /**
184  *
185  */
187  {
188  // From Geomtry
190 
191  // From TriFaceComponent
192  m_fid = in.m_fid;
193  m_globalID = m_fid;
194  m_ownVerts = in.m_ownVerts;
195  std::list<CompToElmt>::const_iterator def;
196  for(def = in.m_elmtMap.begin(); def != in.m_elmtMap.end(); def++)
197  {
198  m_elmtMap.push_back(*def);
199  }
200 
201  // From TriGeom
202  m_verts = in.m_verts;
203  m_edges = in.m_edges;
204  for (int i = 0; i < kNedges; i++)
205  {
206  m_eorient[i] = in.m_eorient[i];
207  }
208  m_ownData = in.m_ownData;
209  }
210 
211 
212  /**
213  *
214  */
216  {
217  }
218 
219 
220  /**
221  * Given local collapsed coordinate Lcoord return the value of physical
222  * coordinate in direction i
223  */
225  const Array<OneD, const NekDouble> &Lcoord)
226  {
228  "Geometry is not in physical space");
229 
230  Array<OneD, NekDouble> tmp(m_xmap->GetTotPoints());
231  m_xmap->BwdTrans(m_coeffs[i], tmp);
232 
233  return m_xmap->PhysEvaluate(Lcoord, tmp);
234  }
235 
237  const TriGeom &face1,
238  const TriGeom &face2)
239  {
240  return GetFaceOrientation(face1.m_verts, face2.m_verts);
241  }
242 
244  const PointGeomVector &face1,
245  const PointGeomVector &face2)
246  {
247  int i, j, vmap[3] = {-1,-1,-1};
248  NekDouble x, y, z, x1, y1, z1, cx = 0.0, cy = 0.0, cz = 0.0;
249 
250  // For periodic faces, we calculate the vector between the centre
251  // points of the two faces. (For connected faces this will be
252  // zero). We can then use this to determine alignment later in the
253  // algorithm.
254  for (i = 0; i < 3; ++i)
255  {
256  cx += (*face2[i])(0) - (*face1[i])(0);
257  cy += (*face2[i])(1) - (*face1[i])(1);
258  cz += (*face2[i])(2) - (*face1[i])(2);
259  }
260  cx /= 3;
261  cy /= 3;
262  cz /= 3;
263 
264  // Now construct a mapping which takes us from the vertices of one
265  // face to the other. That is, vertex j of face2 corresponds to
266  // vertex vmap[j] of face1.
267  for (i = 0; i < 3; ++i)
268  {
269  x = (*face1[i])(0);
270  y = (*face1[i])(1);
271  z = (*face1[i])(2);
272  for (j = 0; j < 3; ++j)
273  {
274  x1 = (*face2[j])(0)-cx;
275  y1 = (*face2[j])(1)-cy;
276  z1 = (*face2[j])(2)-cz;
277  if (sqrt((x1-x)*(x1-x)+(y1-y)*(y1-y)+(z1-z)*(z1-z)) < 1e-8)
278  {
279  vmap[j] = i;
280  break;
281  }
282  }
283  }
284 
285  if (vmap[1] == (vmap[0]+1) % 3)
286  {
287  switch (vmap[0])
288  {
289  case 0:
291  break;
292  case 1:
294  break;
295  case 2:
297  break;
298  }
299  }
300  else
301  {
302  switch (vmap[0])
303  {
304  case 0:
306  break;
307  case 1:
309  break;
310  case 2:
312  break;
313  }
314  }
315 
316  ASSERTL0(false, "Unable to determine triangle orientation");
318  }
319 
320 
321  /**
322  *
323  */
324  void TriGeom::v_AddElmtConnected(int gvo_id, int locid)
325  {
326  CompToElmt ee(gvo_id,locid);
327  m_elmtMap.push_back(ee);
328  }
329 
330 
331  /**
332  *
333  */
335  {
336  return int(m_elmtMap.size());
337  }
338 
339 
340  /**
341  *
342  */
343  bool TriGeom::v_IsElmtConnected(int gvo_id, int locid) const
344  {
345  std::list<CompToElmt>::const_iterator def;
346  CompToElmt ee(gvo_id,locid);
347 
348  def = find(m_elmtMap.begin(),m_elmtMap.end(),ee);
349 
350  // Found the element connectivity object in the list
351  return (def != m_elmtMap.end());
352  }
353 
354 
355  /**
356  *
357  */
358  int TriGeom::v_GetFid() const
359  {
360  return m_fid;
361  }
362 
363 
364  /**
365  *
366  */
368  {
369  return m_coordim;
370  }
371 
372 
373  /**
374  *
375  */
377  {
378  return m_xmap->GetBasis(i);
379  }
380 
381 
382  /**
383  *
384  */
386  {
387  ASSERTL1(i <= 2, "edge is out of range");
388 
389  if(i == 0)
390  {
391  return m_xmap->GetBasis(0);
392  }
393  else
394  {
395  return m_xmap->GetBasis(1);
396  }
397  }
398 
399 
400  /**
401  *
402  */
404  {
405  return GetCoord(i,Lcoord);
406  }
407 
408 
409  /**
410  * Set up GeoFac for this geometry using Coord quadrature distribution
411  */
413  {
415  {
416  GeomType Gtype = eRegular;
417 
419 
420  // check to see if expansions are linear
421  for(int i = 0; i < m_coordim; ++i)
422  {
423  if(m_xmap->GetBasisNumModes(0) != 2 ||
424  m_xmap->GetBasisNumModes(1) != 2)
425  {
426  Gtype = eDeformed;
427  }
428  }
429 
431  Gtype, m_coordim, m_xmap, m_coeffs);
432 
434  }
435  }
436 
437  /**
438  *
439  */
441  {
442  m_ownData = true;
443  }
444 
445 
446  /**
447  * Note verts and edges are listed according to anticlockwise
448  * convention but points in _coeffs have to be in array format from
449  * left to right.
450  */
452  {
453  // check to see if geometry structure is already filled
454  if(m_state != ePtsFilled)
455  {
456  int i,j,k;
457  int nEdgeCoeffs = m_xmap->GetEdgeNcoeffs(0);
458 
459  if (m_curve)
460  {
461  int pdim = LibUtilities::PointsManager()[
462  LibUtilities::PointsKey(2, m_curve->m_ptype)]
463  ->GetPointsDim();
464 
465  // Deal with 2D points type separately
466  // (e.g. electrostatic or Fekete points) to 1D tensor
467  // product.
468  if (pdim == 2)
469  {
470  int N = m_curve->m_points.size();
471  int nEdgePts = (
472  -1+(int)sqrt(static_cast<NekDouble>(8*N+1)))/2;
473 
474  ASSERTL0(nEdgePts*(nEdgePts+1)/2 == N,
475  "NUMPOINTS should be a triangle number for"
476  " triangle "
477  + boost::lexical_cast<string>(m_globalID));
478 
479  for (i = 0; i < kNedges; ++i)
480  {
481  ASSERTL0(
482  m_edges[i]->GetXmap()->GetNcoeffs() == nEdgePts,
483  "Number of edge points does not correspond "
484  "to number of face points in triangle "
485  + boost::lexical_cast<string>(m_globalID));
486  }
487 
488  const LibUtilities::PointsKey P0(
490  const LibUtilities::PointsKey P1(
492  const LibUtilities::BasisKey T0(
493  LibUtilities::eOrtho_A, nEdgePts, P0);
494  const LibUtilities::BasisKey T1(
495  LibUtilities::eOrtho_B, nEdgePts, P1);
497  max(nEdgePts*nEdgePts, m_xmap->GetTotPoints()));
498  Array<OneD, NekDouble> tmp(nEdgePts*nEdgePts);
499 
500  for(i = 0; i < m_coordim; ++i)
501  {
502  // Create a StdNodalTriExp.
505  ::AllocateSharedPtr(T0, T1, m_curve->m_ptype);
506 
507  for (j = 0; j < N; ++j)
508  {
509  phys[j] = (m_curve->m_points[j]->GetPtr())[i];
510  }
511 
512  t->BwdTrans(phys, tmp);
513 
514  // Interpolate points to standard region.
516  P0, P1, tmp,
517  m_xmap->GetBasis(0)->GetPointsKey(),
518  m_xmap->GetBasis(1)->GetPointsKey(),
519  phys);
520 
521  // Forwards transform to get coefficient space.
522  m_xmap->FwdTrans(phys, m_coeffs[i]);
523  }
524  }
525  else if (pdim == 1)
526  {
527  int npts = m_curve->m_points.size();
528  int nEdgePts = (int)sqrt(static_cast<NekDouble>(npts));
529  Array<OneD, NekDouble> tmp (npts);
530  Array<OneD, NekDouble> phys(m_xmap->GetTotPoints());
531  LibUtilities::PointsKey curveKey(
532  nEdgePts, m_curve->m_ptype);
533 
534  // Sanity checks:
535  // - Curved faces should have square number of points;
536  // - Each edge should have sqrt(npts) points.
537  ASSERTL0(nEdgePts * nEdgePts == npts,
538  "NUMPOINTS should be a square number for"
539  " triangle "
540  + boost::lexical_cast<string>(m_globalID));
541 
542  for (i = 0; i < kNedges; ++i)
543  {
544  ASSERTL0(
545  m_edges[i]->GetXmap()->GetNcoeffs() == nEdgePts,
546  "Number of edge points does not correspond to "
547  "number of face points in triangle "
548  + boost::lexical_cast<string>(m_globalID));
549  }
550 
551  for (i = 0; i < m_coordim; ++i)
552  {
553  for (j = 0; j < npts; ++j)
554  {
555  tmp[j] = (m_curve->m_points[j]->GetPtr())[i];
556  }
557 
558  // Interpolate curve points to standard triangle
559  // points.
561  curveKey, curveKey, tmp,
562  m_xmap->GetBasis(0)->GetPointsKey(),
563  m_xmap->GetBasis(1)->GetPointsKey(),
564  phys);
565 
566  // Forwards transform to get coefficient space.
567  m_xmap->FwdTrans(phys, m_coeffs[i]);
568  }
569  }
570  else
571  {
572  ASSERTL0(false, "Only 1D/2D points distributions "
573  "supported.");
574  }
575  }
576 
577  Array<OneD, unsigned int> mapArray (nEdgeCoeffs);
578  Array<OneD, int> signArray(nEdgeCoeffs);
579 
580  for(i = 0; i < kNedges; i++)
581  {
582  m_edges[i]->FillGeom();
583  m_xmap->GetEdgeToElementMap(i,m_eorient[i],mapArray,signArray);
584 
585  nEdgeCoeffs = m_edges[i]->GetXmap()->GetNcoeffs();
586 
587  for(j = 0 ; j < m_coordim; j++)
588  {
589  for(k = 0; k < nEdgeCoeffs; k++)
590  {
591  m_coeffs[j][mapArray[k]] =
592  signArray[k] * m_edges[i]->GetCoeffs(j)[k];
593  }
594  }
595  }
596 
598  }
599  }
600 
601 
602  /**
603  *
604  */
606  Array<OneD,NekDouble> &Lcoords)
607  {
608  NekDouble resid = 0.0;
610 
611  // calculate local coordinate for coord
612  if(GetMetricInfo()->GetGtype() == eRegular)
613  {
614  NekDouble coords2 = (m_coordim == 3)? coords[2]: 0.0;
615  PointGeom dv1, dv2, norm, orth1, orth2;
616  PointGeom xin(m_coordim,0,coords[0],coords[1],coords2);
617 
618  // Calculate edge vectors from 0-1 and 0-2 edges.
619  dv1.Sub(*m_verts[1],*m_verts[0]);
620  dv2.Sub(*m_verts[2],*m_verts[0]);
621 
622  // Obtain normal to plane in which dv1 and dv2 lie
623  norm.Mult(dv1,dv2);
624 
625  // Obtain vector which are proportional to normal of dv1 and dv2.
626  orth1.Mult(norm,dv1);
627  orth2.Mult(norm,dv2);
628 
629  // Start with vector of desired points minus vertex_0
630  xin -= *m_verts[0];
631 
632  // Calculate length using L/|dv1| = (x-v0).n1/(dv1.n1) for coordiante 1
633  // Then rescale to [-1,1].
634  Lcoords[0] = xin.dot(orth2)/dv1.dot(orth2);
635  Lcoords[0] = 2*Lcoords[0]-1;
636  Lcoords[1] = xin.dot(orth1)/dv2.dot(orth1);
637  Lcoords[1] = 2*Lcoords[1]-1;
638  }
639  else
640  {
641  // Determine nearest point of coords to values in m_xmap
642  int npts = m_xmap->GetTotPoints();
643  Array<OneD, NekDouble> ptsx(npts), ptsy(npts);
644  Array<OneD, NekDouble> tmpx(npts), tmpy(npts);
645 
646  m_xmap->BwdTrans(m_coeffs[0], ptsx);
647  m_xmap->BwdTrans(m_coeffs[1], ptsy);
648 
649  const Array<OneD, const NekDouble> za = m_xmap->GetPoints(0);
650  const Array<OneD, const NekDouble> zb = m_xmap->GetPoints(1);
651 
652  //guess the first local coords based on nearest point
653  Vmath::Sadd(npts, -coords[0], ptsx,1,tmpx,1);
654  Vmath::Sadd(npts, -coords[1], ptsy,1,tmpy,1);
655  Vmath::Vmul (npts, tmpx,1,tmpx,1,tmpx,1);
656  Vmath::Vvtvp(npts, tmpy,1,tmpy,1,tmpx,1,tmpx,1);
657 
658  int min_i = Vmath::Imin(npts,tmpx,1);
659 
660  Lcoords[0] = za[min_i%za.num_elements()];
661  Lcoords[1] = zb[min_i/za.num_elements()];
662 
663  // recover cartesian coordinate from collapsed coordinate.
664  Lcoords[0] = (1.0+Lcoords[0])*(1.0-Lcoords[1])/2 -1.0;
665 
666  // Perform newton iteration to find local coordinates
667  NewtonIterationForLocCoord(coords, ptsx, ptsy, Lcoords,resid);
668  }
669  return resid;
670  }
671 
672 
673  /**
674  *
675  */
676  int TriGeom::v_GetEid(int i) const
677  {
678  ASSERTL2((i >=0) && (i <= 2),"Edge id must be between 0 and 2");
679  return m_edges[i]->GetEid();
680  }
681 
682 
683  /**
684  *
685  */
686  int TriGeom::v_GetVid(int i) const
687  {
688  ASSERTL2((i >=0) && (i <= 2),"Vertex id must be between 0 and 2");
689  return m_verts[i]->GetVid();
690  }
691 
692 
693  /**
694  *
695  */
697  {
698  ASSERTL2((i >=0) && (i <= 2),"Vertex id must be between 0 and 2");
699  return m_verts[i];
700  }
701 
702 
703  /**
704  *
705  */
707  {
708  ASSERTL2((i >=0) && (i <= 2),"Edge id must be between 0 and 3");
709  return m_edges[i];
710  }
711 
712 
713  /**
714  *
715  */
717  {
718  ASSERTL2((i >=0) && (i <= 2),"Edge id must be between 0 and 2");
719  return m_eorient[i];
720  }
721 
722 
723  /**
724  *
725  */
727  {
728  ASSERTL2((i >=0) && (i <= 3),"Edge id must be between 0 and 3");
729  if(i < 2)
730  {
731  return m_eorient[i];
732  }
733  else
734  {
736  {
737  return StdRegions::eBackwards;
738  }
739  else
740  {
741  return StdRegions::eForwards;
742  }
743  }
744  }
745 
746 
747  /**
748  *
749  */
751  {
752  int returnval = -1;
753 
754  SegGeomVector::iterator edgeIter;
755  int i;
756 
757  for (i=0,edgeIter = m_edges.begin(); edgeIter != m_edges.end(); ++edgeIter,++i)
758  {
759  if (*edgeIter == edge)
760  {
761  returnval = i;
762  break;
763  }
764  }
765 
766  return returnval;
767  }
768 
769 
770  /**
771  *
772  */
774  {
775  return kNverts;
776  }
777 
778 
779  /**
780  *
781  */
783  {
784  return kNedges;
785  }
786 
787 
788  /**
789  * @brief Determines if a point specified in global coordinates is
790  * located within this tetrahedral geometry.
791  */
793  const Array<OneD, const NekDouble> &gloCoord, NekDouble tol)
794  {
795  Array<OneD,NekDouble> locCoord(GetCoordim(),0.0);
796  return v_ContainsPoint(gloCoord,locCoord,tol);
797 
798  }
799 
800  /**
801  *
802  */
804  Array<OneD, NekDouble> &stdCoord,
805  NekDouble tol)
806  {
807  NekDouble resid;
808  return v_ContainsPoint(gloCoord,stdCoord,tol,resid);
809  }
810 
812  Array<OneD, NekDouble> &stdCoord,
813  NekDouble tol,
814  NekDouble &resid)
815  {
816  ASSERTL1(gloCoord.num_elements() >= 2,
817  "Two dimensional geometry expects at least two coordinates.");
818 
819  resid = GetLocCoords(gloCoord, stdCoord);
820  if (stdCoord[0] >= -(1+tol) && stdCoord[1] >= -(1+tol)
821  && stdCoord[0] + stdCoord[1] <= tol)
822  {
823  return true;
824  }
825  return false;
826  }
827 
829  CurveMap &curvedEdges,
830  CurveMap &curvedFaces)
831  {
832  Geometry::v_Reset(curvedEdges, curvedFaces);
833  CurveMap::iterator it = curvedFaces.find(m_globalID);
834 
835  if (it != curvedFaces.end())
836  {
837  m_curve = it->second;
838  }
839 
840  for (int i = 0; i < 3; ++i)
841  {
842  m_edges[i]->Reset(curvedEdges, curvedFaces);
843  }
844 
845  SetUpXmap();
846  SetUpCoeffs(m_xmap->GetNcoeffs());
847  }
848 
850  {
851  int order0 = m_edges[0]->GetBasis(0)->GetNumModes();
852  int points0 = m_edges[0]->GetBasis(0)->GetNumPoints();
853  int order1 = max(order0,
854  max(m_edges[1]->GetBasis(0)->GetNumModes(),
855  m_edges[2]->GetBasis(0)->GetNumModes()));
856  int points1 = max(points0,
857  max(m_edges[1]->GetBasis(0)->GetNumPoints(),
858  m_edges[2]->GetBasis(0)->GetNumPoints()));
859 
860 
861  const LibUtilities::BasisKey B0(
864  const LibUtilities::BasisKey B1(
867 
869  ::AllocateSharedPtr(B0,B1);
870  }
871  }; //end of namespace
872 }; //end of namespace
StdRegions::StdExpansionSharedPtr m_xmap
Definition: Geometry.h:172
virtual int v_NumElmtConnected() const
Definition: TriGeom.cpp:334
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
virtual int v_GetFid() const
Definition: TriGeom.cpp:358
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
std::vector< PointGeomSharedPtr > PointGeomVector
Definition: Geometry3D.h:61
GeomFactorsSharedPtr m_geomFactors
Definition: Geometry.h:170
static const int kNverts
Definition: TriGeom.h:99
Structure holding graphvertexobject id and local element facet id.
NekDouble dot(PointGeom &a)
Definition: PointGeom.cpp:192
int Imin(int n, const T *x, const int incx)
Return the index of the minimum element in x.
Definition: Vmath.cpp:824
virtual void v_GenGeomFactors()
Definition: TriGeom.cpp:412
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:428
Principle Modified Functions .
Definition: BasisType.h:49
void Sub(PointGeom &a, PointGeom &b)
Definition: PointGeom.cpp:168
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state and remove allocated GeomFactors.
Definition: Geometry.cpp:307
StdRegions::StdExpansionSharedPtr GetXmap() const
Definition: Geometry.h:383
boost::shared_ptr< StdNodalTriExp > StdNodalTriExpSharedPtr
virtual void v_FillGeom()
Put all quadrature information into edge structure.
Definition: TriGeom.cpp:451
boost::shared_ptr< Curve > CurveSharedPtr
Definition: Curve.hpp:62
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:57
StdRegions::Orientation m_eorient[kNedges]
Definition: TriGeom.h:113
const LibUtilities::BasisSharedPtr GetBasis(const int i)
Return the j-th basis of the i-th co-ordinate dimension.
Definition: Geometry.h:475
Principle Orthogonal Functions .
Definition: BasisType.h:47
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:116
virtual PointGeomSharedPtr v_GetVertex(int i) const
Definition: TriGeom.cpp:696
virtual StdRegions::Orientation v_GetEorient(const int i) const
Definition: TriGeom.cpp:716
boost::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:60
Principle Modified Functions .
Definition: BasisType.h:50
std::list< CompToElmt > m_elmtMap
Definition: TriGeom.h:116
virtual const LibUtilities::BasisSharedPtr v_GetEdgeBasis(const int i)
Definition: TriGeom.cpp:385
void NewtonIterationForLocCoord(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &ptsx, const Array< OneD, const NekDouble > &ptsy, Array< OneD, NekDouble > &Lcoords, NekDouble &resid)
Definition: Geometry2D.cpp:100
virtual NekDouble v_GetLocCoords(const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
Definition: TriGeom.cpp:605
static std::string npts
Definition: InputFld.cpp:43
virtual void v_AddElmtConnected(int gvo_id, int locid)
Definition: TriGeom.cpp:324
PointsManagerT & PointsManager(void)
Principle Orthogonal Functions .
Definition: BasisType.h:46
Defines a specification for a set of points.
Definition: Points.h:58
virtual int v_GetVid(int i) const
Definition: TriGeom.cpp:686
double NekDouble
virtual int v_GetCoordim() const
Definition: TriGeom.cpp:367
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.cpp:301
virtual NekDouble v_GetCoord(const int i, const Array< OneD, const NekDouble > &Lcoord)
Definition: TriGeom.cpp:403
PointGeomSharedPtr GetVertex(int i) const
Definition: Geometry.h:348
static const int kNedges
Get the orientation of face1.
Definition: TriGeom.h:98
virtual const LibUtilities::BasisSharedPtr v_GetBasis(const int i)
Definition: TriGeom.cpp:376
void Mult(PointGeom &a, PointGeom &b)
Definition: PointGeom.cpp:177
NekDouble GetCoord(const int i, const Array< OneD, const NekDouble > &Lcoord)
Definition: TriGeom.cpp:224
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
virtual StdRegions::Orientation v_GetCartesianEorient(const int i) const
Definition: TriGeom.cpp:726
virtual int v_WhichEdge(SegGeomSharedPtr edge)
Return the edge number of the given edge.
Definition: TriGeom.cpp:750
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Definition: Geometry.h:180
2D geometry information
Definition: Geometry2D.h:65
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:187
Geometry is straight-sided with constant geometric factors.
virtual int v_GetEid(int i) const
Definition: TriGeom.cpp:676
boost::shared_ptr< Geometry1D > Geometry1DSharedPtr
Definition: Geometry1D.h:48
NekDouble GetLocCoords(const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
Definition: Geometry.h:450
LibUtilities::ShapeType m_shapeType
Definition: Geometry.h:177
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:312
static StdRegions::Orientation GetFaceOrientation(const TriGeom &face1, const TriGeom &face2)
Definition: TriGeom.cpp:236
virtual bool v_ContainsPoint(const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0)
Determines if a point specified in global coordinates is located within this tetrahedral geometry...
Definition: TriGeom.cpp:792
boost::unordered_map< int, CurveSharedPtr > CurveMap
Definition: Curve.hpp:63
virtual int v_GetNumEdges() const
Definition: TriGeom.cpp:782
Geometric information has been generated.
void SetUpCoeffs(const int nCoeffs)
Initialise the m_coeffs array.
Definition: Geometry.h:484
GeomFactorsSharedPtr GetMetricInfo()
Definition: Geometry.h:299
virtual int v_GetNumVerts() const
Definition: TriGeom.cpp:773
GeomType
Indicates the type of element geometry.
boost::shared_ptr< Basis > BasisSharedPtr
virtual bool v_IsElmtConnected(int gvo_id, int locid) const
Definition: TriGeom.cpp:343
GeomState m_state
enum identifier to determine if quad points are filled
Definition: Geometry.h:175
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:165
Geometry is curved or has non-constant factors.
int m_coordim
coordinate dimension
Definition: Geometry.h:169
Describes the specification for a Basis.
Definition: Basis.h:50
virtual const Geometry1DSharedPtr v_GetEdge(int i) const
Definition: TriGeom.cpp:706
boost::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:60
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:50
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state and remove allocated GeomFactors.
Definition: TriGeom.cpp:828