Nektar++
Expansion.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File Expansion.cpp
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: File for Expansion routines
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <boost/core/ignore_unused.hpp>
36 
38 #include <LocalRegions/Expansion.h>
39 #include <LocalRegions/MatrixKey.h>
40 
41 using namespace std;
42 
43 namespace Nektar
44 {
45  namespace LocalRegions
46  {
47  Expansion::Expansion(SpatialDomains::GeometrySharedPtr pGeom) :
48  m_indexMapManager
49  (std::bind(&Expansion::CreateIndexMap,this, std::placeholders::_1),
50  std::string("ExpansionIndexMap")),
51  m_geom(pGeom),
52  m_metricinfo(m_geom->GetGeomFactors()),
53  m_elementTraceLeft(-1),
54  m_elementTraceRight(-1)
55  {
56  if (!m_metricinfo)
57  {
58  return;
59  }
60 
61  if (!m_metricinfo->IsValid())
62  {
63  int nDim = m_base.size();
64  string type = "regular";
65  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
66  {
67  type = "deformed";
68  }
69 
70  stringstream err;
71  err << nDim << "D " << type << " Jacobian not positive "
72  << "(element ID = " << m_geom->GetGlobalID() << ") "
73  << "(first vertex ID = " << m_geom->GetVid(0) << ")";
74  NEKERROR(ErrorUtil::ewarning, err.str());
75  }
76  }
77 
79  StdExpansion(pSrc),
80  m_indexMapManager(pSrc.m_indexMapManager),
81  m_geom(pSrc.m_geom),
82  m_metricinfo(pSrc.m_metricinfo)
83  {
84  }
85 
87  {
88  }
89 
91  {
92  return v_GetLocMatrix(mkey);
93  }
94 
96  const DNekScalMatSharedPtr &r_bnd,
97  const StdRegions::MatrixType matrixType)
98  {
99  return v_BuildTransformationMatrix(r_bnd,matrixType);
100  }
101 
102 
104  const DNekScalMatSharedPtr &r_bnd)
105  {
106  return v_BuildVertexMatrix(r_bnd);
107  }
108 
109 
111  const NekDouble *data,
112  const std::vector<unsigned int > &nummodes,
113  const int nmodes_offset,
114  NekDouble *coeffs,
115  std::vector<LibUtilities::BasisType> &fromType)
116  {
117  v_ExtractDataToCoeffs(data,nummodes,nmodes_offset,coeffs,fromType);
118  }
119 
121  const int edge,
122  const std::shared_ptr<Expansion> &EdgeExp,
125  Array<OneD, NekDouble> &outarray)
126  {
127  v_AddEdgeNormBoundaryInt(edge, EdgeExp, Fx, Fy, outarray);
128  }
129 
131  const int edge,
132  const std::shared_ptr<Expansion> &EdgeExp,
134  Array<OneD, NekDouble> &outarray)
135  {
136  v_AddEdgeNormBoundaryInt(edge, EdgeExp, Fn, outarray);
137  }
138 
140  const int face,
141  const std::shared_ptr<Expansion> &FaceExp,
143  Array<OneD, NekDouble> &outarray)
144  {
145  v_AddFaceNormBoundaryInt(face, FaceExp, Fn, outarray);
146  }
147 
149  const int dir,
150  const Array<OneD, const NekDouble>& inarray,
152  Array<OneD, Array<OneD, NekDouble> > &coeffs,
153  Array<OneD, NekDouble> &outarray)
154  {
155  v_DGDeriv(dir, inarray, EdgeExp, coeffs, outarray);
156  }
157 
159  const Array<OneD, Array<OneD, NekDouble > > &vec)
160  {
161  return v_VectorFlux(vec);
162  }
163 
165  const StdRegions::ConstFactorMap &factors,
166  const StdRegions::VarCoeffMap &varcoeffs)
167  {
168  MatrixKey mkey(mtype, DetShapeType(), *this, factors, varcoeffs);
169  return GetLocMatrix(mkey);
170  }
171 
173  {
174  return m_geom;
175  }
176 
178  {
179  // Clear metrics
180  m_metrics.clear();
181 
182  // Regenerate geometry factors
183  m_metricinfo = m_geom->GetGeomFactors();
184  }
185 
187  {
188  IndexMapValuesSharedPtr returnval;
189 
190  IndexMapType itype = ikey.GetIndexMapType();
191 
192  int entity = ikey.GetIndexEntity();
193 
195 
198 
199  switch(itype)
200  {
201  case eEdgeToElement:
202  {
203  GetTraceToElementMap(entity,map,sign,orient);
204  }
205  break;
206  case eFaceToElement:
207  {
208  GetTraceToElementMap(entity,map,sign,orient);
209  }
210  break;
211  case eEdgeInterior:
212  {
213  ASSERTL0(false,"Boundary Index Map not implemented yet.");
214  //v_GetEdgeInteriorMap(entity,orient,map,sign);
215  }
216  break;
217  case eFaceInterior:
218  {
219  ASSERTL0(false,"Boundary Index Map not implemented yet.");
220  //v_GetFaceInteriorMap(entity,orient,map,sign);
221  }
222  break;
223  case eBoundary:
224  {
225  ASSERTL0(false,"Boundary Index Map not implemented yet.");
226  }
227  break;
228  case eVertex:
229  {
230  ASSERTL0(false,"Vertex Index Map not implemented yet.");
231  }
232  break;
233  default:
234  {
235  ASSERTL0(false,"The Index Map you are requiring "
236  "is not between the possible options.");
237  }
238  }
239 
240  returnval = MemoryManager<IndexMapValues>::AllocateSharedPtr(map.size());
241 
242  for(int i = 0; i < map.size(); i++)
243  {
244  (*returnval)[i].index = map[i];
245  (*returnval)[i].sign = sign[i];
246  }
247 
248  return returnval;
249  }
250 
252  {
253  return m_metricinfo;
254  }
255 
256 
258  {
259  boost::ignore_unused(mkey);
260  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
262  }
263 
265  Array<OneD, NekDouble> &outarray)
266  {
267  const int nqtot = GetTotPoints();
268 
269  if (m_metrics.count(eMetricQuadrature) == 0)
270  {
272  }
273 
274  Vmath::Vmul(nqtot, m_metrics[eMetricQuadrature], 1, inarray, 1, outarray, 1);
275  }
276 
278  const Array<OneD,
279  const NekDouble>& inarray,
280  Array<OneD, NekDouble> &outarray)
281  {
282  const int nqtot = GetTotPoints();
283 
284  if (m_metrics.count(eMetricQuadrature) == 0)
285  {
287  }
288 
289  Vmath::Vdiv(nqtot, inarray,
291  1, outarray, 1);
292  }
293 
295  {
297  }
298 
300  {
301  unsigned int nqtot = GetTotPoints();
302  SpatialDomains::GeomType type = m_metricinfo->GetGtype();
304  if (type == SpatialDomains::eRegular ||
306  {
308  }
309  else
310  {
312  }
313 
316  }
317 
319  Array<OneD, NekDouble> &coords_0,
320  Array<OneD, NekDouble> &coords_1,
321  Array<OneD, NekDouble> &coords_2)
322  {
323  ASSERTL1(m_geom, "m_geom not defined");
324 
325  // get physical points defined in Geom
326  m_geom->FillGeom();
327 
328  const int expDim = m_base.size();
329  int nqGeom = 1;
330  bool doCopy = true;
331 
334 
335  for (int i = 0; i < expDim; ++i)
336  {
337  CBasis[i] = m_geom->GetXmap()->GetBasis(i);
338  nqGeom *= CBasis[i]->GetNumPoints();
339  doCopy = doCopy && m_base[i]->GetBasisKey().SamePoints(
340  CBasis[i]->GetBasisKey());
341  }
342 
343  tmp[0] = coords_0;
344  tmp[1] = coords_1;
345  tmp[2] = coords_2;
346 
347  if (doCopy)
348  {
349  for (int i = 0; i < m_geom->GetCoordim(); ++i)
350  {
351  m_geom->GetXmap()->BwdTrans(m_geom->GetCoeffs(i), tmp[i]);
352  }
353  }
354  else
355  {
356  for (int i = 0; i < m_geom->GetCoordim(); ++i)
357  {
358  Array<OneD, NekDouble> tmpGeom(nqGeom);
359  m_geom->GetXmap()->BwdTrans(m_geom->GetCoeffs(i), tmpGeom);
360 
361  switch (expDim)
362  {
363  case 1:
364  {
366  CBasis[0]->GetPointsKey(), &tmpGeom[0],
367  m_base[0]->GetPointsKey(), &tmp[i][0]);
368  break;
369  }
370  case 2:
371  {
373  CBasis[0]->GetPointsKey(),
374  CBasis[1]->GetPointsKey(),
375  &tmpGeom[0],
376  m_base[0]->GetPointsKey(),
377  m_base[1]->GetPointsKey(),
378  &tmp[i][0]);
379  break;
380  }
381  case 3:
382  {
384  CBasis[0]->GetPointsKey(),
385  CBasis[1]->GetPointsKey(),
386  CBasis[2]->GetPointsKey(),
387  &tmpGeom[0],
388  m_base[0]->GetPointsKey(),
389  m_base[1]->GetPointsKey(),
390  m_base[2]->GetPointsKey(),
391  &tmp[i][0]);
392  break;
393  }
394  }
395  }
396  }
397  }
398 
401  const Array<OneD, const NekDouble> &direction,
403  {
404  int shapedim = dfdir.size();
405  int coordim = GetCoordim();
406  int nqtot = direction.size()/coordim;
407 
408  for(int j = 0; j < shapedim; j++)
409  {
410  dfdir[j] = Array<OneD, NekDouble>(nqtot,0.0);
411  for (int k = 0; k < coordim; k++)
412  {
413  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
414  {
415  Vmath::Vvtvp(nqtot,
416  &df[shapedim*k+j][0], 1,
417  &direction[k*nqtot], 1,
418  &dfdir[j][0], 1,
419  &dfdir[j][0], 1);
420  }
421  else
422  {
423  Vmath::Svtvp(nqtot,
424  df[shapedim*k+j][0],
425  &direction[k*nqtot], 1,
426  &dfdir[j][0], 1,
427  &dfdir[j][0], 1);
428  }
429  }
430  }
431  }
432 
433  // Get Moving frames
435  const int dir,
436  const int shapedim,
437  const StdRegions::VarCoeffMap &varcoeffs)
438  {
439  int coordim = GetCoordim();
440 
441  int nquad0, nquad1, nquad2;
442  int nqtot=1;
443  switch(shapedim)
444  {
445  case 2:
446  {
447  nquad0 = m_base[0]->GetNumPoints();
448  nquad1 = m_base[1]->GetNumPoints();
449  nqtot = nquad0*nquad1;
450  break;
451  }
452  case 3:
453  {
454  nquad0 = m_base[0]->GetNumPoints();
455  nquad1 = m_base[1]->GetNumPoints();
456  nquad2 = m_base[2]->GetNumPoints();
457  nqtot = nquad0 * nquad1 * nquad2;
458  break;
459  }
460  default:
461  break;
462  }
463 
464  StdRegions::VarCoeffType MMFCoeffs[15] =
465  {
481  };
482 
483  Array<OneD, NekDouble> MF(coordim*nqtot);
484  Array<OneD, NekDouble> tmp(nqtot);
485  for (int k = 0; k < coordim; k++)
486  {
487  StdRegions::VarCoeffMap::const_iterator MFdir =
488  varcoeffs.find(MMFCoeffs[5*dir+k]);
489  tmp = MFdir->second;
490 
491  Vmath::Vcopy(nqtot, &tmp[0], 1, &MF[k*nqtot], 1);
492  }
493 
494  return MF;
495  }
496 
497  // Get magnitude of MF
499  const int dir,
500  const StdRegions::VarCoeffMap &varcoeffs)
501  {
502  int indxDiv = 3;
503 
504  StdRegions::VarCoeffType MMFCoeffs[15] =
505  {
521  };
522 
523  StdRegions::VarCoeffMap::const_iterator MFdir =
524  varcoeffs.find(MMFCoeffs[5*dir+indxDiv]);
525  Array<OneD, NekDouble> MFDiv = MFdir->second;
526 
527  return MFDiv;
528  }
529 
530  // Get magnitude of MF
532  const int dir,
533  const StdRegions::VarCoeffMap &varcoeffs)
534  {
535  int indxMag = 4;
536 
537  StdRegions::VarCoeffType MMFCoeffs[15] =
538  {
554  };
555 
556  StdRegions::VarCoeffMap::const_iterator MFdir = varcoeffs.find(MMFCoeffs[5*dir+indxMag]);
557  Array<OneD, NekDouble> MFmag = MFdir->second;
558 
559  return MFmag;
560  }
561 
562 
564  const DNekScalMatSharedPtr &r_bnd,
565  const StdRegions::MatrixType matrixType)
566  {
567  boost::ignore_unused(r_bnd, matrixType);
568  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
569  return NullDNekMatSharedPtr;
570  }
571 
573  const DNekScalMatSharedPtr &r_bnd)
574  {
575  boost::ignore_unused(r_bnd);
576  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
577  return NullDNekMatSharedPtr;
578  }
579 
581  const NekDouble *data,
582  const std::vector<unsigned int > &nummodes,
583  const int nmodes_offset,
584  NekDouble *coeffs,
585  std::vector<LibUtilities::BasisType> &fromType)
586  {
587  boost::ignore_unused(data, nummodes, nmodes_offset,
588  coeffs, fromType);
589  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
590  }
591 
593  const int edge,
594  const std::shared_ptr<Expansion> &EdgeExp,
597  Array<OneD, NekDouble> &outarray)
598  {
599  boost::ignore_unused(edge, EdgeExp, Fx, Fy, outarray);
600  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
601  }
602 
604  const int edge,
605  const std::shared_ptr<Expansion> &EdgeExp,
607  Array<OneD, NekDouble> &outarray)
608  {
609  boost::ignore_unused(edge, EdgeExp, Fn, outarray);
610  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
611  }
612 
614  const int face,
615  const std::shared_ptr<Expansion> &FaceExp,
617  Array<OneD, NekDouble> &outarray)
618  {
619  boost::ignore_unused(face, FaceExp, Fn, outarray);
620  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
621  }
622 
624  const int dir,
625  const Array<OneD, const NekDouble>& inarray,
627  Array<OneD, Array<OneD, NekDouble> > &coeffs,
628  Array<OneD, NekDouble> &outarray)
629  {
630  boost::ignore_unused(dir, inarray, EdgeExp, coeffs, outarray);
631  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
632  }
633 
635  const Array<OneD, Array<OneD, NekDouble > > &vec)
636  {
637  boost::ignore_unused(vec);
638  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
639  return 0.0;
640  }
641 
643  {
644  boost::ignore_unused(trace);
645  return StdRegions::eForwards;
646  }
647 
650  Array<OneD, NekDouble> &outarray)
651  {
652  boost::ignore_unused(dir, inarray, outarray);
653  NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
654  }
655 
656  void Expansion::v_GetTraceQFactors(const int trace,
657  Array<OneD, NekDouble> &outarray)
658  {
659  boost::ignore_unused(trace, outarray);
661  "Method does not exist for this shape or library");
662  }
663 
665  (const int trace,
666  const StdRegions::StdExpansionSharedPtr &TraceExp,
667  const Array<OneD, const NekDouble> &inarray,
668  Array<OneD, NekDouble> &outarray,
670  {
671  boost::ignore_unused(trace,TraceExp,inarray,outarray,orient);
673  "Method does not exist for this shape or library" );
674  }
675 
676  void Expansion::v_GetTracePhysMap(const int edge,
677  Array<OneD, int> &outarray)
678  {
679  boost::ignore_unused(edge, outarray);
681  "Method does not exist for this shape or library" );
682  }
683 
685  const StdRegions::Orientation orient,
686  Array<OneD, int> &idmap,
687  const int nq0, const int nq1)
688  {
689  boost::ignore_unused(orient,idmap,nq0,nq1);
691  "Method does not exist for this shape or library" );
692  }
693 
694  const NormalVector & Expansion::v_GetTraceNormal(const int id) const
695  {
696  boost::ignore_unused(id);
697  ASSERTL0(false, "Cannot get trace normals for this expansion.");
698  static NormalVector result;
699  return result;
700  }
701 
703  {
704  boost::ignore_unused(id);
705  ASSERTL0(false, "Cannot compute trace normal for this expansion.");
706  }
707 
709  {
710  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
711  return NullNekDouble1DArray;
712  }
713 
714 
716  {
717  boost::ignore_unused(normal);
718  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
719  }
720 
721  void Expansion::v_SetUpPhysNormals(const int edge)
722  {
723  boost::ignore_unused(edge);
724  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
725  }
726 
728  const int trace,
729  const Array<OneD, const NekDouble > &primCoeffs,
730  DNekMatSharedPtr &inoutmat)
731  {
732  boost::ignore_unused(trace,primCoeffs,inoutmat);
733  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
734  }
735 
737  const int traceid,
738  const Array<OneD, const NekDouble> &primCoeffs,
739  const Array<OneD, NekDouble> &incoeffs,
740  Array<OneD, NekDouble> &coeffs)
741  {
742  boost::ignore_unused(traceid,primCoeffs,incoeffs, coeffs);
743  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
744  }
745 
747  GetElmtBndNormDirElmtLen(const int nbnd) const
748  {
749  auto x = m_elmtBndNormDirElmtLen.find(nbnd);
750  ASSERTL0 (x != m_elmtBndNormDirElmtLen.end(),
751  "m_elmtBndNormDirElmtLen normal not computed.");
752  return x->second;
753  }
754 
756  const int dir,
757  const Array<OneD, const NekDouble> &inarray,
758  Array<OneD, Array<OneD, NekDouble> > &outarray)
759  {
760  boost::ignore_unused(dir,inarray,outarray);
761  NEKERROR(ErrorUtil::efatal, "v_AlignVectorToCollapsedDir is not defined");
762  }
763  } //end of namespace
764 } //end of namespace
765 
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:15
std::map< int, Array< OneD, NekDouble > > m_elmtBndNormDirElmtLen
the element length in each element boundary(Vertex, edge or face) normal direction calculated based o...
Definition: Expansion.h:284
virtual void v_SetPhysNormals(Array< OneD, const NekDouble > &normal)
Definition: Expansion.cpp:715
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:172
IndexMapValuesSharedPtr CreateIndexMap(const IndexMapKey &ikey)
Definition: Expansion.cpp:186
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:272
Array< OneD, NekDouble > v_GetMF(const int dir, const int shapedim, const StdRegions::VarCoeffMap &varcoeffs)
Definition: Expansion.cpp:434
const SpatialDomains::GeomFactorsSharedPtr & GetMetricInfo() const
Definition: Expansion.cpp:251
DNekMatSharedPtr BuildVertexMatrix(const DNekScalMatSharedPtr &r_bnd)
Definition: Expansion.cpp:103
virtual void v_DivideByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:277
virtual void v_AddRobinTraceContribution(const int traceid, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs)
Definition: Expansion.cpp:736
void ComputeGmatcdotMF(const Array< TwoD, const NekDouble > &df, const Array< OneD, const NekDouble > &direction, Array< OneD, Array< OneD, NekDouble > > &dfdir)
Definition: Expansion.cpp:399
virtual void v_GetTraceQFactors(const int trace, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:656
virtual void v_ComputeTraceNormal(const int id)
Definition: Expansion.cpp:702
virtual DNekMatSharedPtr v_BuildVertexMatrix(const DNekScalMatSharedPtr &r_bnd)
Definition: Expansion.cpp:572
virtual void v_SetCoeffsToOrientation(StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:648
virtual DNekScalMatSharedPtr v_GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:257
virtual void v_AddFaceNormBoundaryInt(const int face, const std::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:613
virtual void v_DGDeriv(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &coeffs, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:623
Array< OneD, NekDouble > v_GetMFMag(const int dir, const StdRegions::VarCoeffMap &varcoeffs)
Definition: Expansion.cpp:531
virtual DNekMatSharedPtr v_BuildTransformationMatrix(const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
Definition: Expansion.cpp:563
const Array< OneD, const NekDouble > & GetElmtBndNormDirElmtLen(const int nbnd) const
Definition: Expansion.cpp:747
virtual void v_AddRobinMassMatrix(const int face, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
Definition: Expansion.cpp:727
Array< OneD, NekDouble > v_GetMFDiv(const int dir, const StdRegions::VarCoeffMap &varcoeffs)
Definition: Expansion.cpp:498
void AddEdgeNormBoundaryInt(const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:120
void AddFaceNormBoundaryInt(const int face, const std::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:139
virtual NekDouble v_VectorFlux(const Array< OneD, Array< OneD, NekDouble > > &vec)
Definition: Expansion.cpp:634
virtual const Array< OneD, const NekDouble > & v_GetPhysNormals(void)
Definition: Expansion.cpp:708
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:273
virtual StdRegions::Orientation v_GetTraceOrient(int trace)
Definition: Expansion.cpp:642
void DGDeriv(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &coeffs, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:148
void ExtractDataToCoeffs(const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType)
Definition: Expansion.cpp:110
virtual void v_ReOrientTracePhysMap(const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1=-1)
Definition: Expansion.cpp:684
virtual void v_MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:264
virtual void v_ComputeLaplacianMetric()
Definition: Expansion.h:301
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3)
Definition: Expansion.cpp:318
virtual void v_ExtractDataToCoeffs(const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType)
Definition: Expansion.cpp:580
virtual void v_AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: Expansion.cpp:755
virtual void v_AddEdgeNormBoundaryInt(const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:592
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:90
virtual void v_SetUpPhysNormals(const int id)
Definition: Expansion.cpp:721
virtual void v_GetTracePhysMap(const int edge, Array< OneD, int > &outarray)
Definition: Expansion.cpp:676
virtual const NormalVector & v_GetTraceNormal(const int id) const
Definition: Expansion.cpp:694
virtual void v_GetTracePhysVals(const int trace, const StdRegions::StdExpansionSharedPtr &TraceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient)
Definition: Expansion.cpp:665
Expansion(SpatialDomains::GeometrySharedPtr pGeom)
Definition: Expansion.cpp:47
DNekMatSharedPtr BuildTransformationMatrix(const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
Definition: Expansion.cpp:95
NekDouble VectorFlux(const Array< OneD, Array< OneD, NekDouble > > &vec)
Definition: Expansion.cpp:158
StdRegions::Orientation GetIndexOrientation() const
Definition: IndexMapKey.h:93
IndexMapType GetIndexMapType() const
Definition: IndexMapKey.h:88
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:134
const LibUtilities::PointsKeyVector GetPointsKeys() const
virtual void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void GetTraceToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
Definition: StdExpansion.h:703
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:376
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Interp1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
this function interpolates a 1D function evaluated at the quadrature points of the basis fbasis0 to ...
Definition: Interp.cpp:53
void Interp3D(const BasisKey &fbasis0, const BasisKey &fbasis1, const BasisKey &fbasis2, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, const BasisKey &tbasis2, Array< OneD, NekDouble > &to)
this function interpolates a 3D function evaluated at the quadrature points of the 3D basis,...
Definition: Interp.cpp:185
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:115
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:246
std::shared_ptr< IndexMapValues > IndexMapValuesSharedPtr
Definition: IndexMapKey.h:126
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:62
GeomType
Indicates the type of element geometry.
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eMovingRegular
Currently unused.
@ eDeformed
Geometry is curved or has non-constant factors.
std::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:53
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:272
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:314
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
static DNekMatSharedPtr NullDNekMatSharedPtr
Definition: NekTypeDefs.hpp:78
static Array< OneD, NekDouble > NullNekDouble1DArray
static DNekScalMatSharedPtr NullDNekScalMatSharedPtr
Definition: NekTypeDefs.hpp:79
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
double NekDouble
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:192
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:565
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:513
void Vdiv(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:257
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199