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), m_metricinfo(m_geom->GetGeomFactors()),
52  m_elementTraceLeft(-1), m_elementTraceRight(-1)
53 {
54  if (!m_metricinfo)
55  {
56  return;
57  }
58 
59  if (!m_metricinfo->IsValid())
60  {
61  int nDim = m_base.size();
62  string type = "regular";
63  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
64  {
65  type = "deformed";
66  }
67 
68  stringstream err;
69  err << nDim << "D " << type << " Jacobian not positive "
70  << "(element ID = " << m_geom->GetGlobalID() << ") "
71  << "(first vertex ID = " << m_geom->GetVid(0) << ")";
72  NEKERROR(ErrorUtil::ewarning, err.str());
73  }
74 
75  m_traceExp.empty();
76 }
77 
79  : StdExpansion(pSrc), m_indexMapManager(pSrc.m_indexMapManager),
80  m_geom(pSrc.m_geom), m_metricinfo(pSrc.m_metricinfo)
81 {
82 }
83 
85 {
86 }
87 
89  const LocalRegions::MatrixKey &mkey)
90 {
91  return v_GetLocMatrix(mkey);
92 }
93 
95 {
96  return v_DropLocMatrix(mkey);
97 }
98 
100  const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
101 {
102  return v_BuildTransformationMatrix(r_bnd, matrixType);
103 }
104 
106 {
107  return v_BuildVertexMatrix(r_bnd);
108 }
109 
111  const NekDouble *data, const std::vector<unsigned int> &nummodes,
112  const int nmodes_offset, NekDouble *coeffs,
113  std::vector<LibUtilities::BasisType> &fromType)
114 {
115  v_ExtractDataToCoeffs(data, nummodes, nmodes_offset, coeffs, fromType);
116 }
117 
119  const int edge, const std::shared_ptr<Expansion> &EdgeExp,
122 {
123  v_AddEdgeNormBoundaryInt(edge, EdgeExp, Fx, Fy, outarray);
124 }
125 
127  const int edge, const std::shared_ptr<Expansion> &EdgeExp,
129 {
130  v_AddEdgeNormBoundaryInt(edge, EdgeExp, Fn, outarray);
131 }
132 
134  const int face, const std::shared_ptr<Expansion> &FaceExp,
136 {
137  v_AddFaceNormBoundaryInt(face, FaceExp, Fn, outarray);
138 }
139 
140 void Expansion::DGDeriv(const int dir,
141  const Array<OneD, const NekDouble> &inarray,
144  Array<OneD, NekDouble> &outarray)
145 {
146  v_DGDeriv(dir, inarray, EdgeExp, coeffs, outarray);
147 }
148 
150 {
151  return v_VectorFlux(vec);
152 }
153 
155  Array<OneD, Array<OneD, NekDouble>> &factors,
156  Array<OneD, Array<OneD, NekDouble>> &d0factors,
157  Array<OneD, Array<OneD, NekDouble>> &d1factors)
158 {
159  return v_NormalTraceDerivFactors(factors, d0factors, d1factors);
160 }
161 
163  const StdRegions::MatrixType mtype,
164  const StdRegions::ConstFactorMap &factors,
165  const StdRegions::VarCoeffMap &varcoeffs)
166 {
167  MatrixKey mkey(mtype, DetShapeType(), *this, factors, varcoeffs);
168  return GetLocMatrix(mkey);
169 }
170 
172 {
173  return m_geom;
174 }
175 
177 {
178  // Clear metrics
179  m_metrics.clear();
180 
181  // Regenerate geometry factors
182  m_metricinfo = m_geom->GetGeomFactors();
183 }
184 
186 {
187  IndexMapValuesSharedPtr returnval;
188 
189  IndexMapType itype = ikey.GetIndexMapType();
190 
191  int entity = ikey.GetIndexEntity();
192 
194 
197 
198  switch (itype)
199  {
200  case eEdgeToElement:
201  {
202  GetTraceToElementMap(entity, map, sign, orient);
203  }
204  break;
205  case eFaceToElement:
206  {
207  GetTraceToElementMap(entity, map, sign, orient);
208  }
209  break;
210  case eEdgeInterior:
211  {
212  ASSERTL0(false, "Boundary Index Map not implemented yet.");
213  // v_GetEdgeInteriorMap(entity,orient,map,sign);
214  }
215  break;
216  case eFaceInterior:
217  {
218  ASSERTL0(false, "Boundary Index Map not implemented yet.");
219  // v_GetFaceInteriorMap(entity,orient,map,sign);
220  }
221  break;
222  case eBoundary:
223  {
224  ASSERTL0(false, "Boundary Index Map not implemented yet.");
225  }
226  break;
227  case eVertex:
228  {
229  ASSERTL0(false, "Vertex Index Map not implemented yet.");
230  }
231  break;
232  default:
233  {
234  ASSERTL0(false, "The Index Map you are requiring "
235  "is not between the possible options.");
236  }
237  }
238 
239  returnval = MemoryManager<IndexMapValues>::AllocateSharedPtr(map.size());
240 
241  for (int i = 0; i < map.size(); i++)
242  {
243  (*returnval)[i].index = map[i];
244  (*returnval)[i].sign = sign[i];
245  }
246 
247  return returnval;
248 }
249 
251 {
252  return m_metricinfo;
253 }
254 
256 {
257  std::map<int, NormalVector>::const_iterator x;
258  x = m_traceNormals.find(id);
259 
260  // if edge normal not defined compute it
261  if (x == m_traceNormals.end())
262  {
264  x = m_traceNormals.find(id);
265  }
266  return x->second;
267 }
268 
270  const LocalRegions::MatrixKey &mkey)
271 {
272  boost::ignore_unused(mkey);
273  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
275 }
276 
278 {
279  DNekScalBlkMatSharedPtr returnval;
280 
282  "Geometric information is not set up");
283 
284  // set up block matrix system
285  unsigned int nbdry = NumBndryCoeffs();
286  unsigned int nint = (unsigned int)(m_ncoeffs - nbdry);
287  unsigned int exp_size[] = {nbdry, nint};
288  unsigned int nblks = 2;
290  nblks, nblks, exp_size, exp_size);
291  // Really need a constructor which takes Arrays
292  NekDouble factor = 1.0;
293 
294  switch (mkey.GetMatrixType())
295  {
296  // this can only use stdregions statically condensed system
297  // for mass matrix
298  case StdRegions::eMass:
299  if ((m_metricinfo->GetGtype() == SpatialDomains::eDeformed) ||
300  (mkey.GetNVarCoeff()))
301  {
302  factor = 1.0;
303  goto UseLocRegionsMatrix;
304  }
305  else
306  {
307  factor = (m_metricinfo->GetJac(GetPointsKeys()))[0];
308  goto UseStdRegionsMatrix;
309  }
310  break;
311  default: // use Deformed case for both
312  // regular and deformed geometries
313  factor = 1.0;
314  goto UseLocRegionsMatrix;
315  break;
316  UseStdRegionsMatrix:
317  {
318  NekDouble invfactor = 1.0 / factor;
319  NekDouble one = 1.0;
322  DNekMatSharedPtr Asubmat;
323 
324  returnval->SetBlock(
325  0, 0,
327  factor, Asubmat = mat->GetBlock(0, 0)));
328  returnval->SetBlock(
329  0, 1,
331  one, Asubmat = mat->GetBlock(0, 1)));
332  returnval->SetBlock(
333  1, 0,
335  factor, Asubmat = mat->GetBlock(1, 0)));
336  returnval->SetBlock(
337  1, 1,
339  invfactor, Asubmat = mat->GetBlock(1, 1)));
340  }
341  break;
342  UseLocRegionsMatrix:
343  {
344  int i, j;
345  NekDouble invfactor = 1.0 / factor;
346  NekDouble one = 1.0;
347  DNekScalMat &mat = *GetLocMatrix(mkey);
350  DNekMatSharedPtr B =
352  DNekMatSharedPtr C =
354  DNekMatSharedPtr D =
356 
357  Array<OneD, unsigned int> bmap(nbdry);
358  Array<OneD, unsigned int> imap(nint);
359  GetBoundaryMap(bmap);
360  GetInteriorMap(imap);
361 
362  for (i = 0; i < nbdry; ++i)
363  {
364  for (j = 0; j < nbdry; ++j)
365  {
366  (*A)(i, j) = mat(bmap[i], bmap[j]);
367  }
368 
369  for (j = 0; j < nint; ++j)
370  {
371  (*B)(i, j) = mat(bmap[i], imap[j]);
372  }
373  }
374 
375  for (i = 0; i < nint; ++i)
376  {
377  for (j = 0; j < nbdry; ++j)
378  {
379  (*C)(i, j) = mat(imap[i], bmap[j]);
380  }
381 
382  for (j = 0; j < nint; ++j)
383  {
384  (*D)(i, j) = mat(imap[i], imap[j]);
385  }
386  }
387 
388  // Calculate static condensed system
389  if (nint)
390  {
391  D->Invert();
392  (*B) = (*B) * (*D);
393  (*A) = (*A) - (*B) * (*C);
394  }
395 
397 
398  returnval->SetBlock(
399  0, 0,
400  Atmp =
402  returnval->SetBlock(
403  0, 1,
405  returnval->SetBlock(
406  1, 0,
407  Atmp =
409  returnval->SetBlock(
410  1, 1,
412  D));
413  }
414  }
415  return returnval;
416 }
417 
419 {
420  boost::ignore_unused(mkey);
421  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
422 }
423 
425  const Array<OneD, const NekDouble> &inarray,
426  Array<OneD, NekDouble> &outarray)
427 {
428  const int nqtot = GetTotPoints();
429 
430  if (m_metrics.count(eMetricQuadrature) == 0)
431  {
433  }
434 
435  Vmath::Vmul(nqtot, m_metrics[eMetricQuadrature], 1, inarray, 1, outarray,
436  1);
437 }
438 
440  const Array<OneD, const NekDouble> &inarray,
441  Array<OneD, NekDouble> &outarray)
442 {
443  const int nqtot = GetTotPoints();
444 
445  if (m_metrics.count(eMetricQuadrature) == 0)
446  {
448  }
449 
450  Vmath::Vdiv(nqtot, inarray, 1, m_metrics[eMetricQuadrature], 1, outarray,
451  1);
452 }
453 
455 {
457 }
458 
460 {
461  unsigned int nqtot = GetTotPoints();
462  SpatialDomains::GeomType type = m_metricinfo->GetGtype();
464  if (type == SpatialDomains::eRegular ||
466  {
468  Array<OneD, NekDouble>(nqtot, m_metricinfo->GetJac(p)[0]);
469  }
470  else
471  {
473  }
474 
477 }
478 
480 {
481  int nquad = GetTotPoints();
482  int ntraces = GetNtraces();
483  int ndir = m_base.size();
484 
486  Array<OneD, NekDouble> phys(nquad);
487 
488  Array<OneD, Array<OneD, int>> traceids(ntraces);
489 
490  int tottracepts = 0;
491  for (int i = 0; i < ntraces; ++i)
492  {
493  GetTracePhysMap(i, traceids[i]);
494  tottracepts += GetTraceNumPoints(i);
495  }
496 
497  // initialise array to null so can call for
498  // differnt dimensions
500 
501  DerivMat = Array<OneD, DNekMatSharedPtr>(ndir);
502 
503  for (int i = 0; i < ndir; ++i)
504  {
505  Deriv[i] = Array<OneD, NekDouble>(nquad);
506  DerivMat[i] =
508  }
509 
510  for (int i = 0; i < m_ncoeffs; ++i)
511  {
512  Vmath::Zero(m_ncoeffs, coeffs, 1);
513  coeffs[i] = 1.0;
514  BwdTrans(coeffs, phys);
515 
516  // dphi_i/d\xi_1, dphi_i/d\xi_2 dphi_i/d\xi_3
517  StdPhysDeriv(phys, Deriv[0], Deriv[1], Deriv[2]);
518 
519  int cnt = 0;
520  for (int j = 0; j < ntraces; ++j)
521  {
522  int nTracePts = GetTraceNumPoints(j);
523  for (int k = 0; k < nTracePts; ++k)
524  {
525  for (int d = 0; d < ndir; ++d)
526  {
527  (*DerivMat[d])(i, cnt + k) = Deriv[d][traceids[j][k]];
528  }
529  }
530  cnt += nTracePts;
531  }
532  }
533 }
534 
536  Array<OneD, NekDouble> &coords_1,
537  Array<OneD, NekDouble> &coords_2)
538 {
539  ASSERTL1(m_geom, "m_geom not defined");
540 
541  // get physical points defined in Geom
542  m_geom->FillGeom();
543 
544  const int expDim = m_base.size();
545  int nqGeom = 1;
546  bool doCopy = true;
547 
550 
551  for (int i = 0; i < expDim; ++i)
552  {
553  CBasis[i] = m_geom->GetXmap()->GetBasis(i);
554  nqGeom *= CBasis[i]->GetNumPoints();
555  doCopy = doCopy &&
556  m_base[i]->GetBasisKey().SamePoints(CBasis[i]->GetBasisKey());
557  }
558 
559  tmp[0] = coords_0;
560  tmp[1] = coords_1;
561  tmp[2] = coords_2;
562 
563  if (doCopy)
564  {
565  for (int i = 0; i < m_geom->GetCoordim(); ++i)
566  {
567  m_geom->GetXmap()->BwdTrans(m_geom->GetCoeffs(i), tmp[i]);
568  }
569  }
570  else
571  {
572  for (int i = 0; i < m_geom->GetCoordim(); ++i)
573  {
574  Array<OneD, NekDouble> tmpGeom(nqGeom);
575  m_geom->GetXmap()->BwdTrans(m_geom->GetCoeffs(i), tmpGeom);
576 
577  switch (expDim)
578  {
579  case 1:
580  {
582  CBasis[0]->GetPointsKey(), &tmpGeom[0],
583  m_base[0]->GetPointsKey(), &tmp[i][0]);
584  break;
585  }
586  case 2:
587  {
589  CBasis[0]->GetPointsKey(), CBasis[1]->GetPointsKey(),
590  &tmpGeom[0], m_base[0]->GetPointsKey(),
591  m_base[1]->GetPointsKey(), &tmp[i][0]);
592  break;
593  }
594  case 3:
595  {
597  CBasis[0]->GetPointsKey(), CBasis[1]->GetPointsKey(),
598  CBasis[2]->GetPointsKey(), &tmpGeom[0],
599  m_base[0]->GetPointsKey(), m_base[1]->GetPointsKey(),
600  m_base[2]->GetPointsKey(), &tmp[i][0]);
601  break;
602  }
603  }
604  }
605  }
606 }
607 
609  const Array<OneD, const NekDouble> &direction,
611 {
612  int shapedim = dfdir.size();
613  int coordim = GetCoordim();
614  int nqtot = direction.size() / coordim;
615 
616  for (int j = 0; j < shapedim; j++)
617  {
618  dfdir[j] = Array<OneD, NekDouble>(nqtot, 0.0);
619  for (int k = 0; k < coordim; k++)
620  {
621  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
622  {
623  Vmath::Vvtvp(nqtot, &df[shapedim * k + j][0], 1,
624  &direction[k * nqtot], 1, &dfdir[j][0], 1,
625  &dfdir[j][0], 1);
626  }
627  else
628  {
629  Vmath::Svtvp(nqtot, df[shapedim * k + j][0],
630  &direction[k * nqtot], 1, &dfdir[j][0], 1,
631  &dfdir[j][0], 1);
632  }
633  }
634  }
635 }
636 
637 // Get Moving frames
639  const int dir, const int shapedim, const StdRegions::VarCoeffMap &varcoeffs)
640 {
641  int coordim = GetCoordim();
642 
643  int nquad0, nquad1, nquad2;
644  int nqtot = 1;
645  switch (shapedim)
646  {
647  case 2:
648  {
649  nquad0 = m_base[0]->GetNumPoints();
650  nquad1 = m_base[1]->GetNumPoints();
651  nqtot = nquad0 * nquad1;
652  break;
653  }
654  case 3:
655  {
656  nquad0 = m_base[0]->GetNumPoints();
657  nquad1 = m_base[1]->GetNumPoints();
658  nquad2 = m_base[2]->GetNumPoints();
659  nqtot = nquad0 * nquad1 * nquad2;
660  break;
661  }
662  default:
663  break;
664  }
665 
666  StdRegions::VarCoeffType MMFCoeffs[15] = {
675 
676  Array<OneD, NekDouble> MF(coordim * nqtot);
677  Array<OneD, NekDouble> tmp(nqtot);
678  for (int k = 0; k < coordim; k++)
679  {
680  StdRegions::VarCoeffMap::const_iterator MFdir =
681  varcoeffs.find(MMFCoeffs[5 * dir + k]);
682  tmp = MFdir->second.GetValue();
683 
684  Vmath::Vcopy(nqtot, &tmp[0], 1, &MF[k * nqtot], 1);
685  }
686 
687  return MF;
688 }
689 
690 // Get magnitude of MF
692  const int dir, const StdRegions::VarCoeffMap &varcoeffs)
693 {
694  int indxDiv = 3;
695 
696  StdRegions::VarCoeffType MMFCoeffs[15] = {
705 
706  StdRegions::VarCoeffMap::const_iterator MFdir =
707  varcoeffs.find(MMFCoeffs[5 * dir + indxDiv]);
708  Array<OneD, NekDouble> MFDiv = MFdir->second.GetValue();
709 
710  return MFDiv;
711 }
712 
713 // Get magnitude of MF
715  const int dir, const StdRegions::VarCoeffMap &varcoeffs)
716 {
717  int indxMag = 4;
718 
719  StdRegions::VarCoeffType MMFCoeffs[15] = {
728 
729  StdRegions::VarCoeffMap::const_iterator MFdir =
730  varcoeffs.find(MMFCoeffs[5 * dir + indxMag]);
731  Array<OneD, NekDouble> MFmag = MFdir->second.GetValue();
732 
733  return MFmag;
734 }
735 
737  const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
738 {
739  boost::ignore_unused(r_bnd, matrixType);
740  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
741  return NullDNekMatSharedPtr;
742 }
743 
745  const DNekScalMatSharedPtr &r_bnd)
746 {
747  boost::ignore_unused(r_bnd);
748  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
749  return NullDNekMatSharedPtr;
750 }
751 
753  const NekDouble *data, const std::vector<unsigned int> &nummodes,
754  const int nmodes_offset, NekDouble *coeffs,
755  std::vector<LibUtilities::BasisType> &fromType)
756 {
757  boost::ignore_unused(data, nummodes, nmodes_offset, coeffs, fromType);
758  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
759 }
760 
762  const int edge, const std::shared_ptr<Expansion> &EdgeExp,
765 {
766  boost::ignore_unused(edge, EdgeExp, Fx, Fy, outarray);
767  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
768 }
769 
771  const int edge, const std::shared_ptr<Expansion> &EdgeExp,
773 {
774  boost::ignore_unused(edge, EdgeExp, Fn, outarray);
775  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
776 }
777 
779  const int face, const std::shared_ptr<Expansion> &FaceExp,
781 {
782  boost::ignore_unused(face, FaceExp, Fn, outarray);
783  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
784 }
785 
786 void Expansion::v_DGDeriv(const int dir,
787  const Array<OneD, const NekDouble> &inarray,
790  Array<OneD, NekDouble> &outarray)
791 {
792  boost::ignore_unused(dir, inarray, EdgeExp, coeffs, outarray);
793  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
794 }
795 
797  const Array<OneD, Array<OneD, NekDouble>> &vec)
798 {
799  boost::ignore_unused(vec);
801  "This function is only valid for "
802  "shape expansion in LocalRegions, not parant class");
803  return 0.0;
804 }
805 
807  Array<OneD, Array<OneD, NekDouble>> &factors,
808  Array<OneD, Array<OneD, NekDouble>> &d0factors,
809  Array<OneD, Array<OneD, NekDouble>> &d1factors)
810 {
811  boost::ignore_unused(factors, d0factors, d1factors);
813  "This function is only valid for "
814  "shape expansion in LocalRegions, not parant class");
815 }
816 
818 {
819  boost::ignore_unused(trace);
820  return StdRegions::eForwards;
821 }
822 
825  Array<OneD, NekDouble> &outarray)
826 {
827  boost::ignore_unused(dir, inarray, outarray);
828  NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
829 }
830 
831 void Expansion::v_GetTraceQFactors(const int trace,
832  Array<OneD, NekDouble> &outarray)
833 {
834  boost::ignore_unused(trace, outarray);
836  "Method does not exist for this shape or library");
837 }
838 
840  const int trace, const StdRegions::StdExpansionSharedPtr &TraceExp,
841  const Array<OneD, const NekDouble> &inarray,
843 {
844  boost::ignore_unused(trace, TraceExp, inarray, outarray, orient);
846  "Method does not exist for this shape or library");
847 }
848 
849 void Expansion::v_GetTracePhysMap(const int edge, Array<OneD, int> &outarray)
850 {
851  boost::ignore_unused(edge, outarray);
853  "Method does not exist for this shape or library");
854 }
855 
857  Array<OneD, int> &idmap, const int nq0,
858  const int nq1)
859 {
860  boost::ignore_unused(orient, idmap, nq0, nq1);
862  "Method does not exist for this shape or library");
863 }
864 
865 void Expansion::v_GenTraceExp(const int traceid, ExpansionSharedPtr &exp)
866 {
867  boost::ignore_unused(traceid, exp);
869  "Method does not exist for this shape or library");
870 }
871 
873 {
874  boost::ignore_unused(id);
875  ASSERTL0(false, "Cannot compute trace normal for this expansion.");
876 }
877 
879 {
880  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
881  return NullNekDouble1DArray;
882 }
883 
885 {
886  boost::ignore_unused(normal);
887  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
888 }
889 
890 void Expansion::v_SetUpPhysNormals(const int edge)
891 {
892  boost::ignore_unused(edge);
893  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
894 }
895 
897  const int trace, const Array<OneD, const NekDouble> &primCoeffs,
898  DNekMatSharedPtr &inoutmat)
899 {
900  boost::ignore_unused(trace, primCoeffs, inoutmat);
901  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
902 }
903 
905  const int traceid, const Array<OneD, const NekDouble> &primCoeffs,
906  const Array<OneD, NekDouble> &incoeffs, Array<OneD, NekDouble> &coeffs)
907 {
908  boost::ignore_unused(traceid, primCoeffs, incoeffs, coeffs);
909  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
910 }
911 
912 void Expansion::v_TraceNormLen(const int traceid, NekDouble &h, NekDouble &p)
913 {
914  boost::ignore_unused(traceid, h, p);
915  NEKERROR(ErrorUtil::efatal, "This method has not been defined");
916 }
917 
919  const int nbnd) const
920 {
921  auto x = m_elmtBndNormDirElmtLen.find(nbnd);
923  "m_elmtBndNormDirElmtLen normal not computed.");
924  return x->second;
925 }
926 
928  const int dir, const Array<OneD, const NekDouble> &inarray,
929  Array<OneD, Array<OneD, NekDouble>> &outarray)
930 {
931  boost::ignore_unused(dir, inarray, outarray);
932  NEKERROR(ErrorUtil::efatal, "v_AlignVectorToCollapsedDir is not defined");
933 }
934 } // namespace LocalRegions
935 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#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:249
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:272
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:49
std::map< int, NormalVector > m_traceNormals
Definition: Expansion.h:278
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:288
virtual void v_SetPhysNormals(Array< OneD, const NekDouble > &normal)
Definition: Expansion.cpp:884
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:171
virtual void v_TraceNormLen(const int traceid, NekDouble &h, NekDouble &p)
Definition: Expansion.cpp:912
IndexMapValuesSharedPtr CreateIndexMap(const IndexMapKey &ikey)
Definition: Expansion.cpp:185
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:275
void GetTracePhysMap(const int edge, Array< OneD, int > &outarray)
Definition: Expansion.h:211
void DropLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:94
const SpatialDomains::GeomFactorsSharedPtr & GetMetricInfo() const
Definition: Expansion.cpp:250
DNekMatSharedPtr BuildVertexMatrix(const DNekScalMatSharedPtr &r_bnd)
Definition: Expansion.cpp:105
void ComputeGmatcdotMF(const Array< TwoD, const NekDouble > &df, const Array< OneD, const NekDouble > &direction, Array< OneD, Array< OneD, NekDouble >> &dfdir)
Definition: Expansion.cpp:608
virtual void v_DivideByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:439
virtual void v_GenTraceExp(const int traceid, ExpansionSharedPtr &exp)
Definition: Expansion.cpp:865
NekDouble VectorFlux(const Array< OneD, Array< OneD, NekDouble >> &vec)
Definition: Expansion.cpp:149
virtual void v_SetCoeffsToOrientation(StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: Expansion.cpp:823
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:904
virtual void v_AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
Definition: Expansion.cpp:927
virtual void v_GetTraceQFactors(const int trace, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:831
virtual void v_ComputeTraceNormal(const int id)
Definition: Expansion.cpp:872
DNekScalBlkMatSharedPtr CreateStaticCondMatrix(const MatrixKey &mkey)
Definition: Expansion.cpp:277
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
Definition: Expansion.cpp:535
Array< OneD, NekDouble > GetMFMag(const int dir, const StdRegions::VarCoeffMap &varcoeffs)
Definition: Expansion.cpp:714
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:786
virtual DNekMatSharedPtr v_BuildVertexMatrix(const DNekScalMatSharedPtr &r_bnd)
Definition: Expansion.cpp:744
virtual DNekScalMatSharedPtr v_GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:269
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:778
void NormalTraceDerivFactors(Array< OneD, Array< OneD, NekDouble >> &factors, Array< OneD, Array< OneD, NekDouble >> &d0factors, Array< OneD, Array< OneD, NekDouble >> &d1factors)
Definition: Expansion.cpp:154
virtual DNekMatSharedPtr v_BuildTransformationMatrix(const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
Definition: Expansion.cpp:736
const Array< OneD, const NekDouble > & GetElmtBndNormDirElmtLen(const int nbnd) const
Definition: Expansion.cpp:918
virtual void v_AddRobinMassMatrix(const int face, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
Definition: Expansion.cpp:896
std::map< int, ExpansionWeakPtr > m_traceExp
Definition: Expansion.h:274
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:118
virtual NekDouble v_VectorFlux(const Array< OneD, Array< OneD, NekDouble >> &vec)
Definition: Expansion.cpp:796
void AddFaceNormBoundaryInt(const int face, const std::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:133
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:140
virtual void v_DropLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:418
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:276
virtual StdRegions::Orientation v_GetTraceOrient(int trace)
Definition: Expansion.cpp:817
void StdDerivBaseOnTraceMat(Array< OneD, DNekMatSharedPtr > &DerivMat)
Definition: Expansion.cpp:479
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_MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: Expansion.cpp:424
virtual void v_ReOrientTracePhysMap(const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1=-1)
Definition: Expansion.cpp:856
virtual void v_ComputeLaplacianMetric()
Definition: Expansion.h:313
virtual const Array< OneD, const NekDouble > & v_GetPhysNormals()
Definition: Expansion.cpp:878
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:752
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:761
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:88
virtual void v_SetUpPhysNormals(const int id)
Definition: Expansion.cpp:890
virtual void v_GetTracePhysMap(const int edge, Array< OneD, int > &outarray)
Definition: Expansion.cpp:849
Array< OneD, NekDouble > GetMFDiv(const int dir, const StdRegions::VarCoeffMap &varcoeffs)
Definition: Expansion.cpp:691
const NormalVector & GetTraceNormal(const int id)
Definition: Expansion.cpp:255
Array< OneD, NekDouble > GetMF(const int dir, const int shapedim, const StdRegions::VarCoeffMap &varcoeffs)
Definition: Expansion.cpp:638
virtual void v_NormalTraceDerivFactors(Array< OneD, Array< OneD, NekDouble >> &factors, Array< OneD, Array< OneD, NekDouble >> &d0factors, Array< OneD, Array< OneD, NekDouble >> &d1factors)
Definition: Expansion.cpp:806
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:839
Expansion(SpatialDomains::GeometrySharedPtr pGeom)
Definition: Expansion.cpp:47
DNekMatSharedPtr BuildTransformationMatrix(const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
Definition: Expansion.cpp:99
StdRegions::Orientation GetIndexOrientation() const
Definition: IndexMapKey.h:96
IndexMapType GetIndexMapType() const
Definition: IndexMapKey.h:91
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:675
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
int GetTraceNumPoints(const int i) const
This function returns the number of quadrature points belonging to the i-th trace.
Definition: StdExpansion.h:287
const LibUtilities::PointsKeyVector GetPointsKeys() const
DNekBlkMatSharedPtr GetStdStaticCondMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:614
int GetNtraces() const
Returns the number of trace elements connected to this element.
Definition: StdExpansion.h:357
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:690
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:373
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:680
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space.
Definition: StdExpansion.h:430
void StdPhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Definition: StdExpansion.h:881
Array< OneD, LibUtilities::BasisSharedPtr > m_base
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:85
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:52
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:167
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:106
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:250
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
std::shared_ptr< IndexMapValues > IndexMapValuesSharedPtr
Definition: IndexMapKey.h:128
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.
@ eNoGeomType
No type defined.
@ 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< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:399
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:343
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:77
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:79
static DNekMatSharedPtr NullDNekMatSharedPtr
Definition: NekTypeDefs.hpp:83
static Array< OneD, NekDouble > NullNekDouble1DArray
static DNekScalMatSharedPtr NullDNekScalMatSharedPtr
Definition: NekTypeDefs.hpp:84
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
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:209
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:622
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:574
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:284
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255