Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 // 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 
39 
40 #include <StdRegions/StdQuadExp.h>
41 #include <SpatialDomains/SegGeom.h>
42 #include <SpatialDomains/Curve.hpp>
44 
45 using namespace std;
46 
47 namespace Nektar
48 {
49  namespace SpatialDomains
50  {
51  /**
52  *
53  */
54  QuadGeom::QuadGeom()
55  {
56  m_shapeType = LibUtilities::eQuadrilateral;
57  }
58 
59 
60  /**
61  *
62  */
63  QuadGeom::QuadGeom(const int id,
64  const PointGeomSharedPtr verts[],
65  const SegGeomSharedPtr edges[],
66  const StdRegions::Orientation eorient[]):
67  Geometry2D(verts[0]->GetCoordim()),
68  m_fid(id)
69  {
70  m_globalID = id;
72 
73  /// Copy the vert shared pointers.
74  m_verts.insert(m_verts.begin(), verts, verts+QuadGeom::kNverts);
75 
76  /// Copy the edge shared pointers.
77  m_edges.insert(m_edges.begin(), edges, edges+QuadGeom::kNedges);
78 
79 
80  for (int j=0; j<kNedges; ++j)
81  {
82  m_eorient[j] = eorient[j];
83  }
84 
85  m_coordim = verts[0]->GetCoordim();
86  ASSERTL0(m_coordim > 1,
87  "Cannot call function with dim == 1");
88 
89  SetUpXmap();
90  SetUpCoeffs(m_xmap->GetNcoeffs());
91  }
92 
93 
94  /**
95  *
96  */
97  QuadGeom::QuadGeom(const int id,
98  const SegGeomSharedPtr edges[],
99  const StdRegions::Orientation eorient[],
100  const CurveSharedPtr &curve) :
101  Geometry2D(edges[0]->GetVertex(0)->GetCoordim()),
102  m_fid(id),
103  m_curve(curve)
104  {
105  int j;
106 
107  m_globalID = m_fid;
108 
110 
111  /// Copy the edge shared pointers.
112  m_edges.insert(m_edges.begin(), edges, edges+QuadGeom::kNedges);
113 
114  for(j=0; j <kNedges; ++j)
115  {
116  if(eorient[j] == StdRegions::eForwards)
117  {
118  m_verts.push_back(edges[j]->GetVertex(0));
119  }
120  else
121  {
122  m_verts.push_back(edges[j]->GetVertex(1));
123  }
124  }
125 
126  for (j=0; j<kNedges; ++j)
127  {
128  m_eorient[j] = eorient[j];
129  }
130 
131  m_coordim = edges[0]->GetVertex(0)->GetCoordim();
132  ASSERTL0(m_coordim > 1,
133  "Cannot call function with dim == 1");
134 
135  SetUpXmap();
136  SetUpCoeffs(m_xmap->GetNcoeffs());
137  }
138 
139 
140  /**
141  *
142  */
143  QuadGeom::QuadGeom(const int id,
144  const SegGeomSharedPtr edges[],
145  const StdRegions::Orientation eorient[]):
146  Geometry2D(edges[0]->GetVertex(0)->GetCoordim()),
147  m_fid(id)
148  {
149  int j;
150 
151  m_globalID = m_fid;
153 
154  /// Copy the edge shared pointers.
155  m_edges.insert(m_edges.begin(), edges, edges+QuadGeom::kNedges);
156 
157  for(j=0; j <kNedges; ++j)
158  {
159  if(eorient[j] == StdRegions::eForwards)
160  {
161  m_verts.push_back(edges[j]->GetVertex(0));
162  }
163  else
164  {
165  m_verts.push_back(edges[j]->GetVertex(1));
166  }
167  }
168 
169  for (j=0; j<kNedges; ++j)
170  {
171  m_eorient[j] = eorient[j];
172  }
173 
174  m_coordim = edges[0]->GetVertex(0)->GetCoordim();
175  ASSERTL0(m_coordim > 1,
176  "Cannot call function with dim == 1");
177 
178  SetUpXmap();
179  SetUpCoeffs(m_xmap->GetNcoeffs());
180  }
181 
182 
183  /**
184  *
185  */
187  {
188  // From Geometry
190 
191  // From QuadFaceComponent
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 QuadGeom
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 
220  {
221  int order0 = max(m_edges[0]->GetBasis(0)->GetNumModes(),
222  m_edges[2]->GetBasis(0)->GetNumModes());
223  int order1 = max(m_edges[1]->GetBasis(0)->GetNumModes(),
224  m_edges[3]->GetBasis(0)->GetNumModes());
225 
226  const LibUtilities::BasisKey B0(
230  const LibUtilities::BasisKey B1(
234 
236  ::AllocateSharedPtr(B0,B1);
237  }
238 
239  /**
240  *
241  */
243  const Array<OneD, const NekDouble> &Lcoord)
244  {
246  "Geometry is not in physical space");
247 
248  Array<OneD, NekDouble> tmp(m_xmap->GetTotPoints());
249  m_xmap->BwdTrans(m_coeffs[i], tmp);
250 
251  return m_xmap->PhysEvaluate(Lcoord, tmp);
252  }
253 
255  const QuadGeom &face1,
256  const QuadGeom &face2)
257  {
258  return GetFaceOrientation(face1.m_verts, face2.m_verts);
259  }
260 
261  /**
262  * Calculate the orientation of face2 to face1 (note this is
263  * not face1 to face2!).
264  */
266  const PointGeomVector &face1,
267  const PointGeomVector &face2)
268  {
269  int i, j, vmap[4] = {-1,-1,-1,-1};
270  NekDouble x, y, z, x1, y1, z1, cx = 0.0, cy = 0.0, cz = 0.0;
271 
272  // For periodic faces, we calculate the vector between the centre
273  // points of the two faces. (For connected faces this will be
274  // zero). We can then use this to determine alignment later in the
275  // algorithm.
276  for (i = 0; i < 4; ++i)
277  {
278  cx += (*face2[i])(0) - (*face1[i])(0);
279  cy += (*face2[i])(1) - (*face1[i])(1);
280  cz += (*face2[i])(2) - (*face1[i])(2);
281  }
282  cx /= 4;
283  cy /= 4;
284  cz /= 4;
285 
286  // Now construct a mapping which takes us from the vertices of one
287  // face to the other. That is, vertex j of face2 corresponds to
288  // vertex vmap[j] of face1.
289  for (i = 0; i < 4; ++i)
290  {
291  x = (*face1[i])(0);
292  y = (*face1[i])(1);
293  z = (*face1[i])(2);
294  for (j = 0; j < 4; ++j)
295  {
296  x1 = (*face2[j])(0)-cx;
297  y1 = (*face2[j])(1)-cy;
298  z1 = (*face2[j])(2)-cz;
299  if (sqrt((x1-x)*(x1-x)+(y1-y)*(y1-y)+(z1-z)*(z1-z)) < 1e-8)
300  {
301  vmap[j] = i;
302  break;
303  }
304  }
305  }
306 
307  // Use the mapping to determine the eight alignment options between
308  // faces.
309  if (vmap[1] == (vmap[0]+1) % 4)
310  {
311  switch (vmap[0])
312  {
313  case 0:
315  break;
316  case 1:
318  break;
319  case 2:
321  break;
322  case 3:
324  break;
325  }
326  }
327  else
328  {
329  switch (vmap[0])
330  {
331  case 0:
333  break;
334  case 1:
336  break;
337  case 2:
339  break;
340  case 3:
342  break;
343  }
344  }
345 
346  ASSERTL0(false, "unable to determine face orientation");
348  }
349 
350 
351  /**
352  *
353  */
354  void QuadGeom::v_AddElmtConnected(int gvo_id, int locid)
355  {
356  CompToElmt ee(gvo_id,locid);
357  m_elmtMap.push_back(ee);
358  }
359 
360 
361  /**
362  *
363  */
365  {
366  return int(m_elmtMap.size());
367  }
368 
369 
370  /**
371  *
372  */
373  bool QuadGeom::v_IsElmtConnected(int gvo_id, int locid) const
374  {
375  std::list<CompToElmt>::const_iterator def;
376  CompToElmt ee(gvo_id,locid);
377 
378  def = find(m_elmtMap.begin(),m_elmtMap.end(),ee);
379 
380  // Found the element connectivity object in the list
381  if(def != m_elmtMap.end())
382  {
383  return(true);
384  }
385 
386  return(false);
387  }
388 
389 
390  /**
391  *
392  */
393  int QuadGeom::v_GetFid() const
394  {
395  return m_fid;
396  }
397 
398 
399  /**
400  *
401  */
403  {
404  return m_coordim;
405  }
406 
407 
408  /**
409  *
410  */
412  {
413  return m_xmap->GetBasis(i);
414  }
415 
416 
417  /**
418  *
419  */
421  {
422  ASSERTL1(i <= 3,"edge is out of range");
423 
424  if (i == 0 || i == 2)
425  {
426  return m_xmap->GetBasis(0);
427  }
428  else
429  {
430  return m_xmap->GetBasis(1);
431  }
432  }
433 
434  /**
435  *
436  */
438  {
439  return GetCoord(i,Lcoord);
440  }
441 
442 
443  /**
444  * Set up GeoFac for this geometry using Coord quadrature distribution
445  */
447  {
449  {
450  int i;
451  GeomType Gtype = eRegular;
452 
454 
455  // We will first check whether we have a regular or deformed
456  // geometry. We will define regular as those cases where the
457  // Jacobian and the metric terms of the derivative are constants
458  // (i.e. not coordinate dependent)
459 
460  // Check to see if expansions are linear
461  // If not linear => deformed geometry
462  for(i = 0; i < m_coordim; ++i)
463  {
464  if((m_xmap->GetBasisNumModes(0) != 2)||
465  (m_xmap->GetBasisNumModes(1) != 2))
466  {
467  Gtype = eDeformed;
468  }
469  }
470 
471  // For linear expansions, the mapping from standard to local
472  // element is given by the relation:
473  // x_i = 0.25 * [ ( x_i^A + x_i^B + x_i^C + x_i^D) +
474  // (-x_i^A + x_i^B + x_i^C - x_i^D)*xi_1 +
475  // (-x_i^A - x_i^B + x_i^C + x_i^D)*xi_2 +
476  // ( x_i^A - x_i^B + x_i^C - x_i^D)*xi_1*xi_2 ]
477  //
478  // The jacobian of the transformation and the metric terms
479  // dxi_i/dx_j, involve only terms of the form dx_i/dxi_j (both
480  // for coordim == 2 or 3). Inspecting the formula above, it can
481  // be appreciated that the derivatives dx_i/dxi_j will be
482  // constant, if the coefficient of the non-linear term is zero.
483  //
484  // That is why for regular geometry, we require
485  //
486  // x_i^A - x_i^B + x_i^C - x_i^D = 0
487  //
488  // or equivalently
489  //
490  // x_i^A - x_i^B = x_i^D - x_i^C
491  //
492  // This corresponds to quadrilaterals which are paralellograms.
493  if(Gtype == eRegular)
494  {
495  for(i = 0; i < m_coordim; i++)
496  {
497  if( fabs( (*m_verts[0])(i) - (*m_verts[1])(i) +
498  (*m_verts[2])(i) - (*m_verts[3])(i) )
500  {
501  Gtype = eDeformed;
502  break;
503  }
504  }
505  }
506 
508  Gtype, m_coordim, m_xmap, m_coeffs);
510  }
511  }
512 
513 
514  /**
515  *
516  */
518  {
519  m_ownData = true;
520  }
521 
522 
523  /**
524  * Note verts and edges are listed according to anticlockwise
525  * convention but points in _coeffs have to be in array format from
526  * left to right.
527  */
529  {
530  // check to see if geometry structure is already filled
531  if(m_state != ePtsFilled)
532  {
533  int i,j,k;
534  int nEdgeCoeffs;
535 
536  if (m_curve)
537  {
538  int npts = m_curve->m_points.size();
539  int nEdgePts = (int)sqrt(static_cast<NekDouble>(npts));
540  Array<OneD, NekDouble> tmp(npts);
541  Array<OneD, NekDouble> tmp2(m_xmap->GetTotPoints());
542  LibUtilities::PointsKey curveKey(
543  nEdgePts, m_curve->m_ptype);
544 
545  // Sanity checks:
546  // - Curved faces should have square number of points;
547  // - Each edge should have sqrt(npts) points.
548  ASSERTL0(nEdgePts * nEdgePts == npts,
549  "NUMPOINTS should be a square number in"
550  " quadrilteral "
551  + boost::lexical_cast<string>(m_globalID));
552 
553  for (i = 0; i < kNedges; ++i)
554  {
555  ASSERTL0(
556  m_edges[i]->GetXmap()->GetNcoeffs() == nEdgePts,
557  "Number of edge points does not correspond to "
558  "number of face points in quadrilateral "
559  + boost::lexical_cast<string>(m_globalID));
560  }
561 
562  for (i = 0; i < m_coordim; ++i)
563  {
564  for (j = 0; j < npts; ++j)
565  {
566  tmp[j] = (m_curve->m_points[j]->GetPtr())[i];
567  }
568 
569  // Interpolate m_curve points to GLL points
571  curveKey, curveKey, tmp,
572  m_xmap->GetBasis(0)->GetPointsKey(),
573  m_xmap->GetBasis(1)->GetPointsKey(),
574  tmp2);
575 
576  // Forwards transform to get coefficient space.
577  m_xmap->FwdTrans(tmp2, m_coeffs[i]);
578  }
579  }
580 
581  // Now fill in edges.
582  Array<OneD, unsigned int> mapArray;
583  Array<OneD, int> signArray;
584 
585  for(i = 0; i < kNedges; i++)
586  {
587  m_edges[i]->FillGeom();
588  m_xmap->GetEdgeToElementMap(i,m_eorient[i],
589  mapArray,signArray);
590 
591  nEdgeCoeffs = m_edges[i]->GetXmap()->GetNcoeffs();
592 
593  for(j = 0; j < m_coordim; j++)
594  {
595  for(k = 0; k < nEdgeCoeffs; k++)
596  {
597  m_coeffs[j][mapArray[k]]
598  = signArray[k]*(m_edges[i]->GetCoeffs(j))[k];
599  }
600  }
601  }
602 
604  }
605  }
606 
607  /**
608  *
609  */
611  Array<OneD,NekDouble> &Lcoords)
612  {
613  NekDouble resid = 0.0;
614  if(GetMetricInfo()->GetGtype() == eRegular)
615  {
616  NekDouble coords2 = (m_coordim == 3)? coords[2]: 0.0;
617  PointGeom dv1, dv2, norm, orth1, orth2;
618  PointGeom xin(m_coordim,0,coords[0],coords[1],coords2);
619 
620  // Calculate edge vectors from 0-1 and 0-3 edges.
621  dv1.Sub(*m_verts[1],*m_verts[0]);
622  dv2.Sub(*m_verts[3],*m_verts[0]);
623 
624  // Obtain normal to plane in which dv1 and dv2 lie
625  norm.Mult(dv1,dv2);
626 
627  // Obtain vector which are normal to dv1 and dv2.
628  orth1.Mult(norm,dv1);
629  orth2.Mult(norm,dv2);
630 
631  // Start with vector of desired points minus vertex_0
632  xin -= *m_verts[0];
633 
634  // Calculate length using L/|dv1| = (x-v0).n1/(dv1.n1) for coordiante 1
635  // Then rescale to [-1,1].
636  Lcoords[0] = xin.dot(orth2)/dv1.dot(orth2);
637  Lcoords[0] = 2*Lcoords[0]-1;
638  Lcoords[1] = xin.dot(orth1)/dv2.dot(orth1);
639  Lcoords[1] = 2*Lcoords[1]-1;
640  }
641  else
642  {
644 
645  // Determine nearest point of coords to values in m_xmap
646  int npts = m_xmap->GetTotPoints();
647  Array<OneD, NekDouble> ptsx(npts), ptsy(npts);
648  Array<OneD, NekDouble> tmpx(npts), tmpy(npts);
649 
650  m_xmap->BwdTrans(m_coeffs[0], ptsx);
651  m_xmap->BwdTrans(m_coeffs[1], ptsy);
652 
653  const Array<OneD, const NekDouble> za = m_xmap->GetPoints(0);
654  const Array<OneD, const NekDouble> zb = m_xmap->GetPoints(1);
655 
656  //guess the first local coords based on nearest point
657  Vmath::Sadd(npts, -coords[0], ptsx,1,tmpx,1);
658  Vmath::Sadd(npts, -coords[1], ptsy,1,tmpy,1);
659  Vmath::Vmul (npts, tmpx,1,tmpx,1,tmpx,1);
660  Vmath::Vvtvp(npts, tmpy,1,tmpy,1,tmpx,1,tmpx,1);
661 
662  int min_i = Vmath::Imin(npts,tmpx,1);
663 
664  Lcoords[0] = za[min_i%za.num_elements()];
665  Lcoords[1] = zb[min_i/za.num_elements()];
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 QuadGeom::v_GetEid(int i) const
678  {
679  ASSERTL2((i >=0) && (i <= 3),"Edge id must be between 0 and 3");
680  return m_edges[i]->GetEid();
681  }
682 
683 
684  /**
685  *
686  */
687  int QuadGeom::v_GetVid(int i) const
688  {
689  ASSERTL2((i >=0) && (i <= 3),"Verted id must be between 0 and 3");
690  return m_verts[i]->GetVid();
691  }
692 
693 
694  /**
695  *
696  */
698  {
699  ASSERTL2((i >=0) && (i <= 3),"Vertex id must be between 0 and 3");
700  return m_verts[i];
701  }
702 
703 
704  /**
705  *
706  */
708  {
709  ASSERTL2((i >=0) && (i <= 3),"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 <= 3),"Edge id must be between 0 and 3");
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 
790  const Array<OneD, const NekDouble> &gloCoord, NekDouble tol)
791  {
792  Array<OneD,NekDouble> locCoord(GetCoordim(),0.0);
793  return v_ContainsPoint(gloCoord,locCoord,tol);
794 
795  }
796  /**
797  *
798  */
800  Array<OneD, NekDouble> &stdCoord,
801  NekDouble tol)
802  {
803  NekDouble resid;
804  return v_ContainsPoint(gloCoord,stdCoord,tol,resid);
805  }
806 
808  Array<OneD, NekDouble> &stdCoord,
809  NekDouble tol,
810  NekDouble &resid)
811  {
812  ASSERTL1(gloCoord.num_elements() >= 2,
813  "Two dimensional geometry expects at least two coordinates.");
814 
815  resid = GetLocCoords(gloCoord, stdCoord);
816  if (stdCoord[0] >= -(1+tol) && stdCoord[1] >= -(1+tol)
817  && stdCoord[0] <= (1+tol) && stdCoord[1] <= (1+tol))
818  {
819  return true;
820  }
821  return false;
822  }
823 
824  void QuadGeom::v_Reset(CurveMap &curvedEdges,
825  CurveMap &curvedFaces)
826  {
827  Geometry::v_Reset(curvedEdges, curvedFaces);
828  CurveMap::iterator it = curvedFaces.find(m_globalID);
829 
830  if (it != curvedFaces.end())
831  {
832  m_curve = it->second;
833  }
834 
835  for (int i = 0; i < 4; ++i)
836  {
837  m_edges[i]->Reset(curvedEdges, curvedFaces);
838  }
839 
840  SetUpXmap();
841  SetUpCoeffs(m_xmap->GetNcoeffs());
842  }
843  }; //end of namespace
844 }; //end of namespace
StdRegions::StdExpansionSharedPtr m_xmap
Definition: Geometry.h:172
virtual const Geometry1DSharedPtr v_GetEdge(int i) const
Definition: QuadGeom.cpp:707
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
virtual int v_GetCoordim() const
Definition: QuadGeom.cpp:402
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
virtual bool v_ContainsPoint(const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0)
Definition: QuadGeom.cpp:789
std::vector< PointGeomSharedPtr > PointGeomVector
Definition: Geometry3D.h:61
GeomFactorsSharedPtr m_geomFactors
Definition: Geometry.h:170
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 int v_NumElmtConnected() const
Definition: QuadGeom.cpp:364
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
virtual PointGeomSharedPtr v_GetVertex(int i) const
Definition: QuadGeom.cpp:697
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< Curve > CurveSharedPtr
Definition: Curve.hpp:62
virtual bool v_IsElmtConnected(int gvo_id, int locid) const
Definition: QuadGeom.cpp:373
virtual int v_GetVid(int i) const
Definition: QuadGeom.cpp:687
const LibUtilities::BasisSharedPtr GetBasis(const int i)
Return the j-th basis of the i-th co-ordinate dimension.
Definition: Geometry.h:475
virtual int v_GetNumEdges() const
Definition: QuadGeom.cpp:783
std::list< CompToElmt > m_elmtMap
Definition: QuadGeom.h:110
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
NekDouble GetCoord(const int i, const Array< OneD, const NekDouble > &Lcoord)
Definition: QuadGeom.cpp:242
static StdRegions::Orientation GetFaceOrientation(const QuadGeom &face1, const QuadGeom &face2)
Get the orientation of face1.
Definition: QuadGeom.cpp:254
static const NekDouble kNekZeroTol
virtual NekDouble v_GetCoord(const int i, const Array< OneD, const NekDouble > &Lcoord)
Definition: QuadGeom.cpp:437
virtual void v_AddElmtConnected(int gvo_id, int locid)
Definition: QuadGeom.cpp:354
boost::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:60
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 void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state and remove allocated GeomFactors.
Definition: QuadGeom.cpp:824
virtual StdRegions::Orientation v_GetCartesianEorient(const int i) const
Definition: QuadGeom.cpp:727
virtual const LibUtilities::BasisSharedPtr v_GetEdgeBasis(const int i)
Definition: QuadGeom.cpp:420
static std::string npts
Definition: InputFld.cpp:43
virtual int v_GetNumVerts() const
Definition: QuadGeom.cpp:774
virtual int v_GetFid() const
Definition: QuadGeom.cpp:393
Defines a specification for a set of points.
Definition: Points.h:58
double NekDouble
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
PointGeomSharedPtr GetVertex(int i) const
Definition: Geometry.h:348
virtual int v_WhichEdge(SegGeomSharedPtr edge)
Return the edge number of the given edge.
Definition: QuadGeom.cpp:751
void Mult(PointGeom &a, PointGeom &b)
Definition: PointGeom.cpp:177
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
virtual void v_FillGeom()
Put all quadrature information into edge structure.
Definition: QuadGeom.cpp:528
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.
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
virtual int v_GetEid(int i) const
Definition: QuadGeom.cpp:677
boost::unordered_map< int, CurveSharedPtr > CurveMap
Definition: Curve.hpp:63
Geometric information has been generated.
void SetUpCoeffs(const int nCoeffs)
Initialise the m_coeffs array.
Definition: Geometry.h:484
StdRegions::Orientation m_eorient[kNedges]
Definition: QuadGeom.h:107
GeomFactorsSharedPtr GetMetricInfo()
Definition: Geometry.h:299
GeomType
Indicates the type of element geometry.
bool m_ownData
Boolean indicating whether object owns the data.
Definition: QuadGeom.h:193
boost::shared_ptr< Basis > BasisSharedPtr
virtual NekDouble v_GetLocCoords(const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
Definition: QuadGeom.cpp:610
GeomState m_state
enum identifier to determine if quad points are filled
Definition: Geometry.h:175
virtual StdRegions::Orientation v_GetEorient(const int i) const
Definition: QuadGeom.cpp:717
virtual const LibUtilities::BasisSharedPtr v_GetBasis(const int i)
Definition: QuadGeom.cpp:411
#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
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