Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 using namespace std;
46 
47 namespace Nektar
48 {
49  namespace SpatialDomains
50  {
51  /**
52  *
53  */
54  TriGeom::TriGeom()
55  {
56  m_shapeType = LibUtilities::eTriangle;
57  }
58 
59  /**
60  *
61  */
62  TriGeom::TriGeom(int id, const int coordim):
63  Geometry2D(coordim), m_fid(id)
64  {
65  m_globalID = m_fid;
66  SetUpXmap();
67  SetUpCoeffs(m_xmap->GetNcoeffs());
68  }
69 
70  /**
71  *
72  */
73  TriGeom::TriGeom(const int id,
74  const PointGeomSharedPtr verts[],
75  const SegGeomSharedPtr edges[],
76  const StdRegions::Orientation eorient[]):
77  Geometry2D(verts[0]->GetCoordim()),
78  m_fid(id)
79  {
80  m_globalID = m_fid;
82 
83  /// Copy the vert shared pointers.
84  m_verts.insert(m_verts.begin(), verts, verts+TriGeom::kNverts);
85 
86  /// Copy the edge shared pointers.
87  m_edges.insert(m_edges.begin(), edges, edges+TriGeom::kNedges);
88 
89  for (int j=0; j<kNedges; ++j)
90  {
91  m_eorient[j] = eorient[j];
92  }
93 
94  m_coordim = verts[0]->GetCoordim();
95  ASSERTL0(m_coordim > 1,
96  "Cannot call function with dim == 1");
97 
98  SetUpXmap();
99  SetUpCoeffs(m_xmap->GetNcoeffs());
100  }
101 
102 
103  /**
104  *
105  */
106  TriGeom::TriGeom(const int id, const SegGeomSharedPtr edges[],
107  const StdRegions::Orientation eorient[]):
108  Geometry2D(edges[0]->GetVertex(0)->GetCoordim()),
109  m_fid(id)
110  {
111  m_globalID = m_fid;
113 
114  /// Copy the edge shared pointers.
115  m_edges.insert(m_edges.begin(), edges, edges+TriGeom::kNedges);
116 
117  for(int j=0; j <kNedges; ++j)
118  {
119  if(eorient[j] == StdRegions::eForwards)
120  {
121  m_verts.push_back(edges[j]->GetVertex(0));
122  }
123  else
124  {
125  m_verts.push_back(edges[j]->GetVertex(1));
126  }
127  }
128 
129  for (int j=0; j<kNedges; ++j)
130  {
131  m_eorient[j] = eorient[j];
132  }
133 
134  m_coordim = edges[0]->GetVertex(0)->GetCoordim();
135  ASSERTL0(m_coordim > 1,"Cannot call function with dim == 1");
136 
137  SetUpXmap();
138  SetUpCoeffs(m_xmap->GetNcoeffs());
139  }
140 
141 
142  /**
143  *
144  */
146  const int id,
147  const SegGeomSharedPtr edges[],
148  const StdRegions::Orientation eorient[],
149  const CurveSharedPtr &curve)
150  : Geometry2D(edges[0]->GetVertex(0)->GetCoordim()),
151  m_fid(id),
152  m_curve(curve)
153  {
154  m_globalID = m_fid;
156 
157  /// Copy the edge shared pointers.
158  m_edges.insert(m_edges.begin(), edges, edges+TriGeom::kNedges);
159 
160  for(int j=0; j <kNedges; ++j)
161  {
162  if(eorient[j] == StdRegions::eForwards)
163  {
164  m_verts.push_back(edges[j]->GetVertex(0));
165  }
166  else
167  {
168  m_verts.push_back(edges[j]->GetVertex(1));
169  }
170  }
171 
172  for (int j=0; j<kNedges; ++j)
173  {
174  m_eorient[j] = eorient[j];
175  }
176 
177  m_coordim = edges[0]->GetVertex(0)->GetCoordim();
178  ASSERTL0(m_coordim > 1,"Cannot call function with dim == 1");
179 
180  SetUpXmap();
181  SetUpCoeffs(m_xmap->GetNcoeffs());
182  }
183 
184  /**
185  *
186  */
188  {
189  // From Geomtry
191 
192  // From TriFaceComponent
193  m_fid = in.m_fid;
194  m_globalID = m_fid;
195  m_ownVerts = in.m_ownVerts;
196  std::list<CompToElmt>::const_iterator def;
197  for(def = in.m_elmtMap.begin(); def != in.m_elmtMap.end(); def++)
198  {
199  m_elmtMap.push_back(*def);
200  }
201 
202  // From TriGeom
203  m_verts = in.m_verts;
204  m_edges = in.m_edges;
205  for (int i = 0; i < kNedges; i++)
206  {
207  m_eorient[i] = in.m_eorient[i];
208  }
209  m_ownData = in.m_ownData;
210  }
211 
212 
213  /**
214  *
215  */
217  {
218  }
219 
220 
221  /**
222  * Given local collapsed coordinate Lcoord return the value of physical
223  * coordinate in direction i
224  */
226  const Array<OneD, const NekDouble> &Lcoord)
227  {
229  "Geometry is not in physical space");
230 
231  Array<OneD, NekDouble> tmp(m_xmap->GetTotPoints());
232  m_xmap->BwdTrans(m_coeffs[i], tmp);
233 
234  return m_xmap->PhysEvaluate(Lcoord, tmp);
235  }
236 
238  const TriGeom &face1,
239  const TriGeom &face2)
240  {
241  return GetFaceOrientation(face1.m_verts, face2.m_verts);
242  }
243 
245  const PointGeomVector &face1,
246  const PointGeomVector &face2)
247  {
248  int i, j, vmap[3] = {-1,-1,-1};
249  NekDouble x, y, z, x1, y1, z1, cx = 0.0, cy = 0.0, cz = 0.0;
250 
251  // For periodic faces, we calculate the vector between the centre
252  // points of the two faces. (For connected faces this will be
253  // zero). We can then use this to determine alignment later in the
254  // algorithm.
255  for (i = 0; i < 3; ++i)
256  {
257  cx += (*face2[i])(0) - (*face1[i])(0);
258  cy += (*face2[i])(1) - (*face1[i])(1);
259  cz += (*face2[i])(2) - (*face1[i])(2);
260  }
261  cx /= 3;
262  cy /= 3;
263  cz /= 3;
264 
265  // Now construct a mapping which takes us from the vertices of one
266  // face to the other. That is, vertex j of face2 corresponds to
267  // vertex vmap[j] of face1.
268  for (i = 0; i < 3; ++i)
269  {
270  x = (*face1[i])(0);
271  y = (*face1[i])(1);
272  z = (*face1[i])(2);
273  for (j = 0; j < 3; ++j)
274  {
275  x1 = (*face2[j])(0)-cx;
276  y1 = (*face2[j])(1)-cy;
277  z1 = (*face2[j])(2)-cz;
278  if (sqrt((x1-x)*(x1-x)+(y1-y)*(y1-y)+(z1-z)*(z1-z)) < 1e-8)
279  {
280  vmap[j] = i;
281  break;
282  }
283  }
284  }
285 
286  if (vmap[1] == (vmap[0]+1) % 3)
287  {
288  switch (vmap[0])
289  {
290  case 0:
292  break;
293  case 1:
295  break;
296  case 2:
298  break;
299  }
300  }
301  else
302  {
303  switch (vmap[0])
304  {
305  case 0:
307  break;
308  case 1:
310  break;
311  case 2:
313  break;
314  }
315  }
316 
317  ASSERTL0(false, "Unable to determine triangle orientation");
319  }
320 
321 
322  /**
323  *
324  */
325  void TriGeom::v_AddElmtConnected(int gvo_id, int locid)
326  {
327  CompToElmt ee(gvo_id,locid);
328  m_elmtMap.push_back(ee);
329  }
330 
331 
332  /**
333  *
334  */
336  {
337  return int(m_elmtMap.size());
338  }
339 
340 
341  /**
342  *
343  */
344  bool TriGeom::v_IsElmtConnected(int gvo_id, int locid) const
345  {
346  std::list<CompToElmt>::const_iterator def;
347  CompToElmt ee(gvo_id,locid);
348 
349  def = find(m_elmtMap.begin(),m_elmtMap.end(),ee);
350 
351  // Found the element connectivity object in the list
352  return (def != m_elmtMap.end());
353  }
354 
355 
356  /**
357  *
358  */
359  int TriGeom::v_GetFid() const
360  {
361  return m_fid;
362  }
363 
364 
365  /**
366  *
367  */
369  {
370  return m_coordim;
371  }
372 
373 
374  /**
375  *
376  */
378  {
379  return m_xmap->GetBasis(i);
380  }
381 
382 
383  /**
384  *
385  */
387  {
388  ASSERTL1(i <= 2, "edge is out of range");
389 
390  if(i == 0)
391  {
392  return m_xmap->GetBasis(0);
393  }
394  else
395  {
396  return m_xmap->GetBasis(1);
397  }
398  }
399 
400 
401  /**
402  *
403  */
405  {
406  return GetCoord(i,Lcoord);
407  }
408 
409 
410  /**
411  * Set up GeoFac for this geometry using Coord quadrature distribution
412  */
414  {
416  {
417  GeomType Gtype = eRegular;
418 
420 
421  // check to see if expansions are linear
422  for(int i = 0; i < m_coordim; ++i)
423  {
424  if(m_xmap->GetBasisNumModes(0) != 2 ||
425  m_xmap->GetBasisNumModes(1) != 2)
426  {
427  Gtype = eDeformed;
428  }
429  }
430 
432  Gtype, m_coordim, m_xmap, m_coeffs);
433 
435  }
436  }
437 
438  /**
439  *
440  */
442  {
443  m_ownData = true;
444  }
445 
446 
447  /**
448  * Note verts and edges are listed according to anticlockwise
449  * convention but points in _coeffs have to be in array format from
450  * left to right.
451  */
453  {
454  // check to see if geometry structure is already filled
455  if(m_state != ePtsFilled)
456  {
457  int i,j,k;
458  int nEdgeCoeffs = m_xmap->GetEdgeNcoeffs(0);
459 
460  if (m_curve)
461  {
462  int pdim = LibUtilities::PointsManager()[
463  LibUtilities::PointsKey(2, m_curve->m_ptype)]
464  ->GetPointsDim();
465 
466  // Deal with 2D points type separately
467  // (e.g. electrostatic or Fekete points) to 1D tensor
468  // product.
469  if (pdim == 2)
470  {
471  int N = m_curve->m_points.size();
472  int nEdgePts = (
473  -1+(int)sqrt(static_cast<NekDouble>(8*N+1)))/2;
474 
475  ASSERTL0(nEdgePts*(nEdgePts+1)/2 == N,
476  "NUMPOINTS should be a triangle number for"
477  " triangle "
478  + boost::lexical_cast<string>(m_globalID));
479 
480  for (i = 0; i < kNedges; ++i)
481  {
482  ASSERTL0(
483  m_edges[i]->GetXmap()->GetNcoeffs() == nEdgePts,
484  "Number of edge points does not correspond "
485  "to number of face points in triangle "
486  + boost::lexical_cast<string>(m_globalID));
487  }
488 
489  const LibUtilities::PointsKey P0(
491  const LibUtilities::PointsKey P1(
493  const LibUtilities::BasisKey T0(
494  LibUtilities::eOrtho_A, nEdgePts, P0);
495  const LibUtilities::BasisKey T1(
496  LibUtilities::eOrtho_B, nEdgePts, P1);
498  max(nEdgePts*nEdgePts, m_xmap->GetTotPoints()));
499  Array<OneD, NekDouble> tmp(nEdgePts*nEdgePts);
500 
501  for(i = 0; i < m_coordim; ++i)
502  {
503  // Create a StdNodalTriExp.
506  ::AllocateSharedPtr(T0, T1, m_curve->m_ptype);
507 
508  for (j = 0; j < N; ++j)
509  {
510  phys[j] = (m_curve->m_points[j]->GetPtr())[i];
511  }
512 
513  t->BwdTrans(phys, tmp);
514 
515  // Interpolate points to standard region.
517  P0, P1, tmp,
518  m_xmap->GetBasis(0)->GetPointsKey(),
519  m_xmap->GetBasis(1)->GetPointsKey(),
520  phys);
521 
522  // Forwards transform to get coefficient space.
523  m_xmap->FwdTrans(phys, m_coeffs[i]);
524  }
525  }
526  else if (pdim == 1)
527  {
528  int npts = m_curve->m_points.size();
529  int nEdgePts = (int)sqrt(static_cast<NekDouble>(npts));
530  Array<OneD, NekDouble> tmp (npts);
531  Array<OneD, NekDouble> phys(m_xmap->GetTotPoints());
532  LibUtilities::PointsKey curveKey(
533  nEdgePts, m_curve->m_ptype);
534 
535  // Sanity checks:
536  // - Curved faces should have square number of points;
537  // - Each edge should have sqrt(npts) points.
538  ASSERTL0(nEdgePts * nEdgePts == npts,
539  "NUMPOINTS should be a square number for"
540  " triangle "
541  + boost::lexical_cast<string>(m_globalID));
542 
543  for (i = 0; i < kNedges; ++i)
544  {
545  ASSERTL0(
546  m_edges[i]->GetXmap()->GetNcoeffs() == nEdgePts,
547  "Number of edge points does not correspond to "
548  "number of face points in triangle "
549  + boost::lexical_cast<string>(m_globalID));
550  }
551 
552  for (i = 0; i < m_coordim; ++i)
553  {
554  for (j = 0; j < npts; ++j)
555  {
556  tmp[j] = (m_curve->m_points[j]->GetPtr())[i];
557  }
558 
559  // Interpolate curve points to standard triangle
560  // points.
562  curveKey, curveKey, tmp,
563  m_xmap->GetBasis(0)->GetPointsKey(),
564  m_xmap->GetBasis(1)->GetPointsKey(),
565  phys);
566 
567  // Forwards transform to get coefficient space.
568  m_xmap->FwdTrans(phys, m_coeffs[i]);
569  }
570  }
571  else
572  {
573  ASSERTL0(false, "Only 1D/2D points distributions "
574  "supported.");
575  }
576  }
577 
578  Array<OneD, unsigned int> mapArray (nEdgeCoeffs);
579  Array<OneD, int> signArray(nEdgeCoeffs);
580 
581  for(i = 0; i < kNedges; i++)
582  {
583  m_edges[i]->FillGeom();
584  m_xmap->GetEdgeToElementMap(i,m_eorient[i],mapArray,signArray);
585 
586  nEdgeCoeffs = m_edges[i]->GetXmap()->GetNcoeffs();
587 
588  for(j = 0 ; j < m_coordim; j++)
589  {
590  for(k = 0; k < nEdgeCoeffs; k++)
591  {
592  m_coeffs[j][mapArray[k]] =
593  signArray[k] * m_edges[i]->GetCoeffs(j)[k];
594  }
595  }
596  }
597 
599  }
600  }
601 
602 
603  /**
604  *
605  */
607  Array<OneD,NekDouble> &Lcoords)
608  {
609  NekDouble resid = 0.0;
611 
612  // calculate local coordinate for coord
613  if(GetMetricInfo()->GetGtype() == eRegular)
614  {
615  NekDouble coords2 = (m_coordim == 3)? coords[2]: 0.0;
616  PointGeom dv1, dv2, norm, orth1, orth2;
617  PointGeom xin(m_coordim,0,coords[0],coords[1],coords2);
618 
619  // Calculate edge vectors from 0-1 and 0-2 edges.
620  dv1.Sub(*m_verts[1],*m_verts[0]);
621  dv2.Sub(*m_verts[2],*m_verts[0]);
622 
623  // Obtain normal to plane in which dv1 and dv2 lie
624  norm.Mult(dv1,dv2);
625 
626  // Obtain vector which are proportional to normal of dv1 and dv2.
627  orth1.Mult(norm,dv1);
628  orth2.Mult(norm,dv2);
629 
630  // Start with vector of desired points minus vertex_0
631  xin -= *m_verts[0];
632 
633  // Calculate length using L/|dv1| = (x-v0).n1/(dv1.n1) for coordiante 1
634  // Then rescale to [-1,1].
635  Lcoords[0] = xin.dot(orth2)/dv1.dot(orth2);
636  Lcoords[0] = 2*Lcoords[0]-1;
637  Lcoords[1] = xin.dot(orth1)/dv2.dot(orth1);
638  Lcoords[1] = 2*Lcoords[1]-1;
639  }
640  else
641  {
642  // Determine nearest point of coords to values in m_xmap
643  int npts = m_xmap->GetTotPoints();
644  Array<OneD, NekDouble> ptsx(npts), ptsy(npts);
645  Array<OneD, NekDouble> tmpx(npts), tmpy(npts);
646 
647  m_xmap->BwdTrans(m_coeffs[0], ptsx);
648  m_xmap->BwdTrans(m_coeffs[1], ptsy);
649 
650  const Array<OneD, const NekDouble> za = m_xmap->GetPoints(0);
651  const Array<OneD, const NekDouble> zb = m_xmap->GetPoints(1);
652 
653  //guess the first local coords based on nearest point
654  Vmath::Sadd(npts, -coords[0], ptsx,1,tmpx,1);
655  Vmath::Sadd(npts, -coords[1], ptsy,1,tmpy,1);
656  Vmath::Vmul (npts, tmpx,1,tmpx,1,tmpx,1);
657  Vmath::Vvtvp(npts, tmpy,1,tmpy,1,tmpx,1,tmpx,1);
658 
659  int min_i = Vmath::Imin(npts,tmpx,1);
660 
661  Lcoords[0] = za[min_i%za.num_elements()];
662  Lcoords[1] = zb[min_i/za.num_elements()];
663 
664  // recover cartesian coordinate from collapsed coordinate.
665  Lcoords[0] = (1.0+Lcoords[0])*(1.0-Lcoords[1])/2 -1.0;
666 
667  // Perform newton iteration to find local coordinates
668  NewtonIterationForLocCoord(coords, ptsx, ptsy, Lcoords,resid);
669  }
670  return resid;
671  }
672 
673 
674  /**
675  *
676  */
677  int TriGeom::v_GetEid(int i) const
678  {
679  ASSERTL2((i >=0) && (i <= 2),"Edge id must be between 0 and 2");
680  return m_edges[i]->GetEid();
681  }
682 
683 
684  /**
685  *
686  */
687  int TriGeom::v_GetVid(int i) const
688  {
689  ASSERTL2((i >=0) && (i <= 2),"Vertex id must be between 0 and 2");
690  return m_verts[i]->GetVid();
691  }
692 
693 
694  /**
695  *
696  */
698  {
699  ASSERTL2((i >=0) && (i <= 2),"Vertex id must be between 0 and 2");
700  return m_verts[i];
701  }
702 
703 
704  /**
705  *
706  */
708  {
709  ASSERTL2((i >=0) && (i <= 2),"Edge id must be between 0 and 3");
710  return m_edges[i];
711  }
712 
713 
714  /**
715  *
716  */
718  {
719  ASSERTL2((i >=0) && (i <= 2),"Edge id must be between 0 and 2");
720  return m_eorient[i];
721  }
722 
723 
724  /**
725  *
726  */
728  {
729  ASSERTL2((i >=0) && (i <= 3),"Edge id must be between 0 and 3");
730  if(i < 2)
731  {
732  return m_eorient[i];
733  }
734  else
735  {
737  {
738  return StdRegions::eBackwards;
739  }
740  else
741  {
742  return StdRegions::eForwards;
743  }
744  }
745  }
746 
747 
748  /**
749  *
750  */
752  {
753  int returnval = -1;
754 
755  SegGeomVector::iterator edgeIter;
756  int i;
757 
758  for (i=0,edgeIter = m_edges.begin(); edgeIter != m_edges.end(); ++edgeIter,++i)
759  {
760  if (*edgeIter == edge)
761  {
762  returnval = i;
763  break;
764  }
765  }
766 
767  return returnval;
768  }
769 
770 
771  /**
772  *
773  */
775  {
776  return kNverts;
777  }
778 
779 
780  /**
781  *
782  */
784  {
785  return kNedges;
786  }
787 
788 
789  /**
790  * @brief Determines if a point specified in global coordinates is
791  * located within this tetrahedral geometry.
792  */
794  const Array<OneD, const NekDouble> &gloCoord, NekDouble tol)
795  {
796  Array<OneD,NekDouble> locCoord(GetCoordim(),0.0);
797  return v_ContainsPoint(gloCoord,locCoord,tol);
798 
799  }
800 
801  /**
802  *
803  */
805  Array<OneD, NekDouble> &stdCoord,
806  NekDouble tol)
807  {
808  NekDouble resid;
809  return v_ContainsPoint(gloCoord,stdCoord,tol,resid);
810  }
811 
813  Array<OneD, NekDouble> &stdCoord,
814  NekDouble tol,
815  NekDouble &resid)
816  {
817  ASSERTL1(gloCoord.num_elements() >= 2,
818  "Two dimensional geometry expects at least two coordinates.");
819 
820  resid = GetLocCoords(gloCoord, stdCoord);
821  if (stdCoord[0] >= -(1+tol) && stdCoord[1] >= -(1+tol)
822  && stdCoord[0] + stdCoord[1] <= tol)
823  {
824  return true;
825  }
826  return false;
827  }
828 
830  CurveMap &curvedEdges,
831  CurveMap &curvedFaces)
832  {
833  Geometry::v_Reset(curvedEdges, curvedFaces);
834  CurveMap::iterator it = curvedFaces.find(m_globalID);
835 
836  if (it != curvedFaces.end())
837  {
838  m_curve = it->second;
839  }
840 
841  for (int i = 0; i < 3; ++i)
842  {
843  m_edges[i]->Reset(curvedEdges, curvedFaces);
844  }
845 
846  SetUpXmap();
847  SetUpCoeffs(m_xmap->GetNcoeffs());
848  }
849 
851  {
852  int order0 = m_edges[0]->GetBasis(0)->GetNumModes();
853  int order1 = max(order0,
854  max(m_edges[1]->GetBasis(0)->GetNumModes(),
855  m_edges[2]->GetBasis(0)->GetNumModes()));
856 
857  const LibUtilities::BasisKey B0(
860  const LibUtilities::BasisKey B1(
863 
865  ::AllocateSharedPtr(B0,B1);
866  }
867  }; //end of namespace
868 }; //end of namespace
StdRegions::StdExpansionSharedPtr m_xmap
Definition: Geometry.h:172
virtual int v_NumElmtConnected() const
Definition: TriGeom.cpp:335
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
virtual int v_GetFid() const
Definition: TriGeom.cpp:359
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:833
virtual void v_GenGeomFactors()
Definition: TriGeom.cpp:413
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
STL namespace.
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:452
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:697
virtual StdRegions::Orientation v_GetEorient(const int i) const
Definition: TriGeom.cpp:717
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:386
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:103
virtual NekDouble v_GetLocCoords(const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
Definition: TriGeom.cpp:606
static std::string npts
Definition: InputFld.cpp:43
virtual void v_AddElmtConnected(int gvo_id, int locid)
Definition: TriGeom.cpp:325
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:687
double NekDouble
virtual int v_GetCoordim() const
Definition: TriGeom.cpp:368
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:404
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:377
void Mult(PointGeom &a, PointGeom &b)
Definition: PointGeom.cpp:177
NekDouble GetCoord(const int i, const Array< OneD, const NekDouble > &Lcoord)
Definition: TriGeom.cpp:225
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:727
virtual int v_WhichEdge(SegGeomSharedPtr edge)
Return the edge number of the given edge.
Definition: TriGeom.cpp:751
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:240
Geometry is straight-sided with constant geometric factors.
virtual int v_GetEid(int i) const
Definition: TriGeom.cpp:677
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:315
static StdRegions::Orientation GetFaceOrientation(const TriGeom &face1, const TriGeom &face2)
Definition: TriGeom.cpp:237
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:793
boost::unordered_map< int, CurveSharedPtr > CurveMap
Definition: Curve.hpp:63
virtual int v_GetNumEdges() const
Definition: TriGeom.cpp:783
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:774
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:344
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:218
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:707
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:829