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 
46 namespace Nektar
47 {
48  namespace SpatialDomains
49  {
50  /**
51  *
52  */
54  {
56  }
57 
58 
59  /**
60  *
61  */
62  QuadGeom::QuadGeom(int id, const int coordim):
63  Geometry2D(coordim)
64  {
65 
68 
69  m_globalID = m_fid = id;
70 
72  SetUpCoeffs(m_xmap->GetNcoeffs());
73  }
74 
75  /**
76  *
77  */
78  QuadGeom::QuadGeom(const int id,
79  const PointGeomSharedPtr verts[],
80  const SegGeomSharedPtr edges[],
81  const StdRegions::Orientation eorient[]):
82  Geometry2D(verts[0]->GetCoordim()),
83  m_fid(id)
84  {
85  m_globalID = id;
87 
88  /// Copy the vert shared pointers.
89  m_verts.insert(m_verts.begin(), verts, verts+QuadGeom::kNverts);
90 
91  /// Copy the edge shared pointers.
92  m_edges.insert(m_edges.begin(), edges, edges+QuadGeom::kNedges);
93 
94 
95  for (int j=0; j<kNedges; ++j)
96  {
97  m_eorient[j] = eorient[j];
98  }
99 
100  m_coordim = verts[0]->GetCoordim();
101  ASSERTL0(m_coordim > 1,
102  "Cannot call function with dim == 1");
103 
104  int order0 = max(edges[0]->GetBasis(0)->GetNumModes(),
105  edges[2]->GetBasis(0)->GetNumModes());
106  int points0 = max(edges[0]->GetBasis(0)->GetNumPoints(),
107  edges[2]->GetBasis(0)->GetNumPoints());
108  int order1 = max(edges[1]->GetBasis(0)->GetNumModes(),
109  edges[3]->GetBasis(0)->GetNumModes());
110  int points1 = max(edges[1]->GetBasis(0)->GetNumPoints(),
111  edges[3]->GetBasis(0)->GetNumPoints());
112 
117 
119  SetUpCoeffs(m_xmap->GetNcoeffs());
120  }
121 
122 
123  /**
124  *
125  */
126  QuadGeom::QuadGeom(const int id,
127  const SegGeomSharedPtr edges[],
128  const StdRegions::Orientation eorient[],
129  const CurveSharedPtr &curve) :
130  Geometry2D(edges[0]->GetVertex(0)->GetCoordim()),
131  m_fid(id)
132  {
133  int j;
134 
135  m_globalID = m_fid;
136 
138 
139  /// Copy the edge shared pointers.
140  m_edges.insert(m_edges.begin(), edges, edges+QuadGeom::kNedges);
141 
142  for(j=0; j <kNedges; ++j)
143  {
144  if(eorient[j] == StdRegions::eForwards)
145  {
146  m_verts.push_back(edges[j]->GetVertex(0));
147  }
148  else
149  {
150  m_verts.push_back(edges[j]->GetVertex(1));
151  }
152  }
153 
154  for (j=0; j<kNedges; ++j)
155  {
156  m_eorient[j] = eorient[j];
157  }
158 
159  m_coordim = edges[0]->GetVertex(0)->GetCoordim();
160  ASSERTL0(m_coordim > 1,
161  "Cannot call function with dim == 1");
162 
163  int order0 = max(edges[0]->GetBasis(0)->GetNumModes(),
164  edges[2]->GetBasis(0)->GetNumModes());
165  int points0 = max(edges[0]->GetBasis(0)->GetNumPoints(),
166  edges[2]->GetBasis(0)->GetNumPoints());
167  int order1 = max(edges[1]->GetBasis(0)->GetNumModes(),
168  edges[3]->GetBasis(0)->GetNumModes());
169  int points1 = max(edges[1]->GetBasis(0)->GetNumPoints(),
170  edges[3]->GetBasis(0)->GetNumPoints());
171 
176 
178  SetUpCoeffs(m_xmap->GetNcoeffs());
179 
180  for(int i = 0; i < m_coordim; ++i)
181  {
182  int npts = curve->m_points.size();
183  int nEdgePts = (int)sqrt(static_cast<NekDouble>(npts));
184  Array<OneD,NekDouble> tmp(npts);
185  LibUtilities::PointsKey curveKey(nEdgePts, curve->m_ptype);
186 
187  // Sanity checks:
188  // - Curved faces should have square number of points;
189  // - Each edge should have sqrt(npts) points.
190  ASSERTL0(nEdgePts*nEdgePts == npts,
191  "NUMPOINTS should be a square number");
192 
193  for (j = 0; j < kNedges; ++j)
194  {
195  ASSERTL0(edges[j]->GetXmap()->GetNcoeffs() == nEdgePts,
196  "Number of edge points does not correspond "
197  "to number of face points.");
198  }
199 
200  for (j = 0; j < npts; ++j)
201  {
202  tmp[j] = (curve->m_points[j]->GetPtr())[i];
203  }
204 
205  // Interpolate curve points to GLL points
206  Array<OneD, NekDouble> tmp2(points0*points1);
207  LibUtilities::Interp2D(curveKey,curveKey,tmp,
208  B0.GetPointsKey(),B1.GetPointsKey(),
209  tmp2);
210 
211  // Forwards transform to get coefficient space.
212  m_xmap->FwdTrans(tmp2, m_coeffs[i]);
213  }
214  }
215 
216 
217  /**
218  *
219  */
220  QuadGeom::QuadGeom(const int id,
221  const SegGeomSharedPtr edges[],
222  const StdRegions::Orientation eorient[]):
223  Geometry2D(edges[0]->GetVertex(0)->GetCoordim()),
224  m_fid(id)
225  {
226  int j;
227 
228  m_globalID = m_fid;
230 
231  /// Copy the edge shared pointers.
232  m_edges.insert(m_edges.begin(), edges, edges+QuadGeom::kNedges);
233 
234  for(j=0; j <kNedges; ++j)
235  {
236  if(eorient[j] == StdRegions::eForwards)
237  {
238  m_verts.push_back(edges[j]->GetVertex(0));
239  }
240  else
241  {
242  m_verts.push_back(edges[j]->GetVertex(1));
243  }
244  }
245 
246  for (j=0; j<kNedges; ++j)
247  {
248  m_eorient[j] = eorient[j];
249  }
250 
251  m_coordim = edges[0]->GetVertex(0)->GetCoordim();
252  ASSERTL0(m_coordim > 1,
253  "Cannot call function with dim == 1");
254 
255  int order0 = max(edges[0]->GetBasis(0)->GetNumModes(),
256  edges[2]->GetBasis(0)->GetNumModes());
257  int points0 = max(edges[0]->GetBasis(0)->GetNumPoints(),
258  edges[2]->GetBasis(0)->GetNumPoints());
259  int order1 = max(edges[1]->GetBasis(0)->GetNumModes(),
260  edges[3]->GetBasis(0)->GetNumModes());
261  int points1 = max(edges[1]->GetBasis(0)->GetNumPoints(),
262  edges[3]->GetBasis(0)->GetNumPoints());
263 
268 
270  SetUpCoeffs(m_xmap->GetNcoeffs());
271  }
272 
273 
274  /**
275  *
276  */
278  {
279  // From Geometry
281 
282  // From QuadFaceComponent
283  m_fid = in.m_fid;
284  m_globalID = m_fid;
285  m_ownVerts = in.m_ownVerts;
286  std::list<CompToElmt>::const_iterator def;
287  for(def = in.m_elmtMap.begin(); def != in.m_elmtMap.end(); def++)
288  {
289  m_elmtMap.push_back(*def);
290  }
291 
292  // From QuadGeom
293  m_verts = in.m_verts;
294  m_edges = in.m_edges;
295  for (int i = 0; i < kNedges; i++)
296  {
297  m_eorient[i] = in.m_eorient[i];
298  }
299  m_ownData = in.m_ownData;
300  }
301 
302 
303  /**
304  *
305  */
307  {
308  }
309 
310 
311  /**
312  *
313  */
315  const Array<OneD, const NekDouble> &Lcoord)
316  {
318  "Geometry is not in physical space");
319 
320  Array<OneD, NekDouble> tmp(m_xmap->GetTotPoints());
321  m_xmap->BwdTrans(m_coeffs[i], tmp);
322 
323  return m_xmap->PhysEvaluate(Lcoord, tmp);
324  }
325 
327  const QuadGeom &face1,
328  const QuadGeom &face2)
329  {
330  return GetFaceOrientation(face1.m_verts, face2.m_verts);
331  }
332 
333  /**
334  * Calculate the orientation of face2 to face1 (note this is
335  * not face1 to face2!).
336  */
338  const PointGeomVector &face1,
339  const PointGeomVector &face2)
340  {
341  int i, j, vmap[4] = {-1,-1,-1,-1};
342  NekDouble x, y, z, x1, y1, z1, cx = 0.0, cy = 0.0, cz = 0.0;
343 
344  // For periodic faces, we calculate the vector between the centre
345  // points of the two faces. (For connected faces this will be
346  // zero). We can then use this to determine alignment later in the
347  // algorithm.
348  for (i = 0; i < 4; ++i)
349  {
350  cx += (*face2[i])(0) - (*face1[i])(0);
351  cy += (*face2[i])(1) - (*face1[i])(1);
352  cz += (*face2[i])(2) - (*face1[i])(2);
353  }
354  cx /= 4;
355  cy /= 4;
356  cz /= 4;
357 
358  // Now construct a mapping which takes us from the vertices of one
359  // face to the other. That is, vertex j of face2 corresponds to
360  // vertex vmap[j] of face1.
361  for (i = 0; i < 4; ++i)
362  {
363  x = (*face1[i])(0);
364  y = (*face1[i])(1);
365  z = (*face1[i])(2);
366  for (j = 0; j < 4; ++j)
367  {
368  x1 = (*face2[j])(0)-cx;
369  y1 = (*face2[j])(1)-cy;
370  z1 = (*face2[j])(2)-cz;
371  if (sqrt((x1-x)*(x1-x)+(y1-y)*(y1-y)+(z1-z)*(z1-z)) < 1e-8)
372  {
373  vmap[j] = i;
374  break;
375  }
376  }
377  }
378 
379  // Use the mapping to determine the eight alignment options between
380  // faces.
381  if (vmap[1] == (vmap[0]+1) % 4)
382  {
383  switch (vmap[0])
384  {
385  case 0:
387  break;
388  case 1:
390  break;
391  case 2:
393  break;
394  case 3:
396  break;
397  }
398  }
399  else
400  {
401  switch (vmap[0])
402  {
403  case 0:
405  break;
406  case 1:
408  break;
409  case 2:
411  break;
412  case 3:
414  break;
415  }
416  }
417 
418  ASSERTL0(false, "unable to determine face orientation");
420  }
421 
422 
423  /**
424  *
425  */
426  void QuadGeom::v_AddElmtConnected(int gvo_id, int locid)
427  {
428  CompToElmt ee(gvo_id,locid);
429  m_elmtMap.push_back(ee);
430  }
431 
432 
433  /**
434  *
435  */
437  {
438  return int(m_elmtMap.size());
439  }
440 
441 
442  /**
443  *
444  */
445  bool QuadGeom::v_IsElmtConnected(int gvo_id, int locid) const
446  {
447  std::list<CompToElmt>::const_iterator def;
448  CompToElmt ee(gvo_id,locid);
449 
450  def = find(m_elmtMap.begin(),m_elmtMap.end(),ee);
451 
452  // Found the element connectivity object in the list
453  if(def != m_elmtMap.end())
454  {
455  return(true);
456  }
457 
458  return(false);
459  }
460 
461 
462  /**
463  *
464  */
465  int QuadGeom::v_GetFid() const
466  {
467  return m_fid;
468  }
469 
470 
471  /**
472  *
473  */
475  {
476  return m_coordim;
477  }
478 
479 
480  /**
481  *
482  */
484  {
485  return m_xmap->GetBasis(i);
486  }
487 
488 
489  /**
490  *
491  */
493  {
494  ASSERTL1(i <= 3,"edge is out of range");
495 
496  if (i == 0 || i == 2)
497  {
498  return m_xmap->GetBasis(0);
499  }
500  else
501  {
502  return m_xmap->GetBasis(1);
503  }
504  }
505 
506  /**
507  *
508  */
509  NekDouble QuadGeom::v_GetCoord(const int i, const Array<OneD, const NekDouble> &Lcoord)
510  {
511  return GetCoord(i,Lcoord);
512  }
513 
514 
515  /**
516  * Set up GeoFac for this geometry using Coord quadrature distribution
517  */
519  {
521  {
522  int i;
523  GeomType Gtype = eRegular;
524 
526 
527  // We will first check whether we have a regular or deformed
528  // geometry. We will define regular as those cases where the
529  // Jacobian and the metric terms of the derivative are constants
530  // (i.e. not coordinate dependent)
531 
532  // Check to see if expansions are linear
533  // If not linear => deformed geometry
534  for(i = 0; i < m_coordim; ++i)
535  {
536  if((m_xmap->GetBasisNumModes(0) != 2)||
537  (m_xmap->GetBasisNumModes(1) != 2))
538  {
539  Gtype = eDeformed;
540  }
541  }
542 
543  // For linear expansions, the mapping from standard to local
544  // element is given by the relation:
545  // x_i = 0.25 * [ ( x_i^A + x_i^B + x_i^C + x_i^D) +
546  // (-x_i^A + x_i^B + x_i^C - x_i^D)*xi_1 +
547  // (-x_i^A - x_i^B + x_i^C + x_i^D)*xi_2 +
548  // ( x_i^A - x_i^B + x_i^C - x_i^D)*xi_1*xi_2 ]
549  //
550  // The jacobian of the transformation and the metric terms
551  // dxi_i/dx_j, involve only terms of the form dx_i/dxi_j (both
552  // for coordim == 2 or 3). Inspecting the formula above, it can
553  // be appreciated that the derivatives dx_i/dxi_j will be
554  // constant, if the coefficient of the non-linear term is zero.
555  //
556  // That is why for regular geometry, we require
557  //
558  // x_i^A - x_i^B + x_i^C - x_i^D = 0
559  //
560  // or equivalently
561  //
562  // x_i^A - x_i^B = x_i^D - x_i^C
563  //
564  // This corresponds to quadrilaterals which are paralellograms.
565  if(Gtype == eRegular)
566  {
567  for(i = 0; i < m_coordim; i++)
568  {
569  if( fabs( (*m_verts[0])(i) - (*m_verts[1])(i) +
570  (*m_verts[2])(i) - (*m_verts[3])(i) )
572  {
573  Gtype = eDeformed;
574  break;
575  }
576  }
577  }
578 
580  Gtype, m_coordim, m_xmap, m_coeffs);
582  }
583  }
584 
585 
586  /**
587  *
588  */
590  {
591  m_ownData = true;
592  }
593 
594 
595  /**
596  * Note verts and edges are listed according to anticlockwise
597  * convention but points in _coeffs have to be in array format from
598  * left to right.
599  */
601  {
602  // check to see if geometry structure is already filled
603  if(m_state != ePtsFilled)
604  {
605  int i,j,k;
606  int nEdgeCoeffs;
607 
608  Array<OneD, unsigned int> mapArray;
609  Array<OneD, int> signArray;
610 
611  for(i = 0; i < kNedges; i++)
612  {
613  m_edges[i]->FillGeom();
614  m_xmap->GetEdgeToElementMap(i,m_eorient[i],
615  mapArray,signArray);
616 
617  nEdgeCoeffs = m_edges[i]->GetXmap()->GetNcoeffs();
618 
619  for(j = 0; j < m_coordim; j++)
620  {
621  for(k = 0; k < nEdgeCoeffs; k++)
622  {
623  m_coeffs[j][mapArray[k]]
624  = signArray[k]*(m_edges[i]->GetCoeffs(j))[k];
625  }
626  }
627  }
628 
630  }
631  }
632 
633  /**
634  *
635  */
636  NekDouble QuadGeom::v_GetLocCoords(const Array<OneD, const NekDouble> &coords,
637  Array<OneD,NekDouble> &Lcoords)
638  {
639  NekDouble resid = 0.0;
640  if(GetMetricInfo()->GetGtype() == eRegular)
641  {
642  NekDouble coords2 = (m_coordim == 3)? coords[2]: 0.0;
643  PointGeom dv1, dv2, norm, orth1, orth2;
644  PointGeom xin(m_coordim,0,coords[0],coords[1],coords2);
645 
646  // Calculate edge vectors from 0-1 and 0-3 edges.
647  dv1.Sub(*m_verts[1],*m_verts[0]);
648  dv2.Sub(*m_verts[3],*m_verts[0]);
649 
650  // Obtain normal to plane in which dv1 and dv2 lie
651  norm.Mult(dv1,dv2);
652 
653  // Obtain vector which are normal to dv1 and dv2.
654  orth1.Mult(norm,dv1);
655  orth2.Mult(norm,dv2);
656 
657  // Start with vector of desired points minus vertex_0
658  xin -= *m_verts[0];
659 
660  // Calculate length using L/|dv1| = (x-v0).n1/(dv1.n1) for coordiante 1
661  // Then rescale to [-1,1].
662  Lcoords[0] = xin.dot(orth2)/dv1.dot(orth2);
663  Lcoords[0] = 2*Lcoords[0]-1;
664  Lcoords[1] = xin.dot(orth1)/dv2.dot(orth1);
665  Lcoords[1] = 2*Lcoords[1]-1;
666  }
667  else
668  {
670 
671  // Determine nearest point of coords to values in m_xmap
672  int npts = m_xmap->GetTotPoints();
673  Array<OneD, NekDouble> ptsx(npts), ptsy(npts);
674  Array<OneD, NekDouble> tmpx(npts), tmpy(npts);
675 
676  m_xmap->BwdTrans(m_coeffs[0], ptsx);
677  m_xmap->BwdTrans(m_coeffs[1], ptsy);
678 
679  const Array<OneD, const NekDouble> za = m_xmap->GetPoints(0);
680  const Array<OneD, const NekDouble> zb = m_xmap->GetPoints(1);
681 
682  //guess the first local coords based on nearest point
683  Vmath::Sadd(npts, -coords[0], ptsx,1,tmpx,1);
684  Vmath::Sadd(npts, -coords[1], ptsy,1,tmpy,1);
685  Vmath::Vmul (npts, tmpx,1,tmpx,1,tmpx,1);
686  Vmath::Vvtvp(npts, tmpy,1,tmpy,1,tmpx,1,tmpx,1);
687 
688  int min_i = Vmath::Imin(npts,tmpx,1);
689 
690  Lcoords[0] = za[min_i%za.num_elements()];
691  Lcoords[1] = zb[min_i/za.num_elements()];
692 
693  // Perform newton iteration to find local coordinates
694  NewtonIterationForLocCoord(coords, ptsx, ptsy, Lcoords,resid);
695  }
696  return resid;
697  }
698 
699 
700  /**
701  *
702  */
703  int QuadGeom::v_GetEid(int i) const
704  {
705  ASSERTL2((i >=0) && (i <= 3),"Edge id must be between 0 and 3");
706  return m_edges[i]->GetEid();
707  }
708 
709 
710  /**
711  *
712  */
713  int QuadGeom::v_GetVid(int i) const
714  {
715  ASSERTL2((i >=0) && (i <= 3),"Verted id must be between 0 and 3");
716  return m_verts[i]->GetVid();
717  }
718 
719 
720  /**
721  *
722  */
724  {
725  ASSERTL2((i >=0) && (i <= 3),"Vertex id must be between 0 and 3");
726  return m_verts[i];
727  }
728 
729 
730  /**
731  *
732  */
734  {
735  ASSERTL2((i >=0) && (i <= 3),"Edge id must be between 0 and 3");
736  return m_edges[i];
737  }
738 
739 
740  /**
741  *
742  */
744  {
745  ASSERTL2((i >=0) && (i <= 3),"Edge id must be between 0 and 3");
746  return m_eorient[i];
747  }
748 
749 
750  /**
751  *
752  */
754  {
755  ASSERTL2((i >=0) && (i <= 3),"Edge id must be between 0 and 3");
756  if(i < 2)
757  {
758  return m_eorient[i];
759  }
760  else
761  {
763  {
764  return StdRegions::eBackwards;
765  }
766  else
767  {
768  return StdRegions::eForwards;
769  }
770  }
771  }
772 
773 
774  /**
775  *
776  */
778  {
779  int returnval = -1;
780 
781  SegGeomVector::iterator edgeIter;
782  int i;
783 
784  for (i=0,edgeIter = m_edges.begin(); edgeIter != m_edges.end(); ++edgeIter,++i)
785  {
786  if (*edgeIter == edge)
787  {
788  returnval = i;
789  break;
790  }
791  }
792 
793  return returnval;
794  }
795 
796 
797  /**
798  *
799  */
801  {
802  return kNverts;
803  }
804 
805 
806  /**
807  *
808  */
810  {
811  return kNedges;
812  }
813 
814 
816  const Array<OneD, const NekDouble> &gloCoord, NekDouble tol)
817  {
818  Array<OneD,NekDouble> locCoord(GetCoordim(),0.0);
819  return v_ContainsPoint(gloCoord,locCoord,tol);
820 
821  }
822  /**
823  *
824  */
825  bool QuadGeom::v_ContainsPoint(const Array<OneD, const NekDouble> &gloCoord,
826  Array<OneD, NekDouble> &stdCoord,
827  NekDouble tol)
828  {
829  NekDouble resid;
830  return v_ContainsPoint(gloCoord,stdCoord,tol,resid);
831  }
832 
833  bool QuadGeom::v_ContainsPoint(const Array<OneD, const NekDouble> &gloCoord,
834  Array<OneD, NekDouble> &stdCoord,
835  NekDouble tol,
836  NekDouble &resid)
837  {
838  ASSERTL1(gloCoord.num_elements() >= 2,
839  "Two dimensional geometry expects at least two coordinates.");
840 
841  resid = GetLocCoords(gloCoord, stdCoord);
842  if (stdCoord[0] >= -(1+tol) && stdCoord[1] >= -(1+tol)
843  && stdCoord[0] <= (1+tol) && stdCoord[1] <= (1+tol))
844  {
845  return true;
846  }
847  return false;
848  }
849  }; //end of namespace
850 }; //end of namespace