Nektar++
ExpList3D.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File ExpList3D.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: Expansion list 3D definition
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <iomanip>
36 
37 #include <boost/core/ignore_unused.hpp>
38 
39 #include <MultiRegions/ExpList3D.h>
40 
41 #include <LocalRegions/HexExp.h>
42 #include <LocalRegions/PrismExp.h>
43 #include <LocalRegions/PyrExp.h>
44 #include <LocalRegions/TetExp.h>
45 
48 
49 using namespace std;
50 
51 namespace Nektar
52 {
53  namespace MultiRegions
54  {
55 
56  ExpList3D::ExpList3D(): ExpList()
57  {
58  SetExpType(e3D);
59  }
60 
62  {
63  SetExpType(e3D);
64  }
65 
67  const std::vector<unsigned int> &eIDs,
68  const bool DeclareCoeffPhysArrays,
69  const Collections::ImplementationType ImpType):
70  ExpList(In, eIDs, DeclareCoeffPhysArrays)
71  {
72  SetExpType(e3D);
73 
74  // Setup Default optimisation information.
75  int nel = GetExpSize();
78 
80 
81  if (DeclareCoeffPhysArrays)
82  {
83  // Set up m_coeffs, m_phys.
86  }
87 
89  CreateCollections(ImpType);
90  }
91 
93  {
94  }
95 
96 
98  const LibUtilities::BasisKey &TBa,
99  const LibUtilities::BasisKey &TBb,
100  const LibUtilities::BasisKey &TBc,
101  const LibUtilities::BasisKey &HBa,
102  const LibUtilities::BasisKey &HBb,
103  const LibUtilities::BasisKey &HBc,
104  const SpatialDomains::MeshGraphSharedPtr &graph3D,
105  const LibUtilities::PointsType TetNb,
106  const Collections::ImplementationType ImpType):
107  ExpList(pSession,graph3D)
108  {
109  SetExpType(e3D);
110 
115 
116  const SpatialDomains::ExpansionMap &expansions = graph3D->GetExpansions();
117 
118  for (auto &expIt : expansions)
119  {
124 
125  if((TetGeom = std::dynamic_pointer_cast<SpatialDomains::TetGeom>(expIt.second->m_geomShPtr)))
126  {
128  {
129 // Ntet = MemoryManager<LocalRegions::NodalTetExp>::AllocateSharedPtr(TetBa,TetBb,TetBc,TetNb,TetGeom);
130 // (*m_exp).push_back(Ntet);
131  }
132  else
133  {
135  (*m_exp).push_back(tet);
136  }
137 
139 
140  m_npoints += TBa.GetNumPoints()*TBb.GetNumPoints()*TBc.GetNumPoints();
141  }
142 /*
143  else if((PrismGeom = std::dynamic_pointer_cast<SpatialDomains::PrismGeom>(expansions[i]->m_geomShPtr)))
144  {
145  prism = MemoryManager<LocalRegions::PrismExp>::AllocateSharedPtr(Ba,Bb,Bc,PrismGeom);
146  (*m_exp).push_back(prism);
147 
148  m_ncoeffs += StdRegions::StdPrismData::getNumberOfCoefficients(Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes());
149  m_npoints += Ba.GetNumPoints()*Bb.GetNumPoints()*Bc.GetNumPoints();
150 
151  }
152  else if((PyrGeom = std::dynamic_pointer_cast<SpatialDomains::PyrGeom>(expansions[i]->m_geomShPtr)))
153  {
154  pyramid = MemoryManager<LocalRegions::PyrExp>::AllocateSharedPtr(Ba,Bb,Bc,PyrGeom);
155  (*m_exp).push_back(pyramid);
156 
157  m_ncoeffs += StdRegions::StdPyrData::getNumberOfCoefficients(Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes());
158  m_npoints += Ba.GetNumPoints()*Bb.GetNumPoints()*Bc.GetNumPoints();
159 
160  }
161 */
162  else if((HexGeom = std::dynamic_pointer_cast<SpatialDomains::HexGeom>(expIt.second->m_geomShPtr)))
163  {
164  hex = MemoryManager<LocalRegions::HexExp>::AllocateSharedPtr(HBa,HBb,HBc, HexGeom);
165  (*m_exp).push_back(hex);
166 
167  m_ncoeffs += HBa.GetNumModes()*HBb.GetNumModes()*HBc.GetNumModes();
168  m_npoints += HBa.GetNumPoints()*HBb.GetNumPoints()*HBc.GetNumPoints();
169  }
170  else
171  {
172  ASSERTL0(false,"dynamic cast to a proper Geometry3D failed");
173  }
174 
175  }
176 
177  // Setup Default optimisation information.
178  int nel = GetExpSize();
181 
183 
184  // Set up m_coeffs, m_phys.
187 
189  CreateCollections(ImpType);
190  }
191 
192  /**
193  * Given a mesh \a graph3D, containing information about the domain and
194  * the spectral/hp element expansion, this constructor fills the list
195  * of local expansions \texttt{m_exp} with the proper expansions,
196  * calculates the total number of quadrature points
197  * \f$\boldsymbol{x}_i\f$ and the local expansion coefficients
198  * \f$\hat{u}^e_n\f$ and allocates memory for the arrays #m_coeffs and
199  * #m_phys.
200  *
201  * @param graph3D A mesh, containing information about the domain
202  * and the spectral/hp element expansion.
203  */
205  const SpatialDomains::MeshGraphSharedPtr &graph3D,
206  const std::string &variable,
207  const Collections::ImplementationType ImpType):
208  ExpList(pSession,graph3D)
209  {
210  SetExpType(e3D);
211 
212  int elmtid = 0;
217 
218  const SpatialDomains::ExpansionMap &expansions
219  = graph3D->GetExpansions(variable);
220 
221  for (auto &expIt : expansions)
222  {
227 
228  if((TetGeom = std::dynamic_pointer_cast<
229  SpatialDomains::TetGeom>(expIt.second->m_geomShPtr)))
230  {
232  = expIt.second->m_basisKeyVector[0];
234  = expIt.second->m_basisKeyVector[1];
236  = expIt.second->m_basisKeyVector[2];
237 
240  {
241  ASSERTL0(false,"LocalRegions::NodalTetExp is not "
242  "implemented yet");
243  }
244  else
245  {
247  ::AllocateSharedPtr(TetBa,TetBb,TetBc,
248  TetGeom);
249  tet->SetElmtId(elmtid++);
250  (*m_exp).push_back(tet);
251  }
252  }
253  else if((PrismGeom = std::dynamic_pointer_cast<SpatialDomains
254  ::PrismGeom>(expIt.second->m_geomShPtr)))
255  {
256  LibUtilities::BasisKey PrismBa
257  = expIt.second->m_basisKeyVector[0];
258  LibUtilities::BasisKey PrismBb
259  = expIt.second->m_basisKeyVector[1];
260  LibUtilities::BasisKey PrismBc
261  = expIt.second->m_basisKeyVector[2];
262 
264  ::AllocateSharedPtr(PrismBa,PrismBb,
265  PrismBc,PrismGeom);
266  prism->SetElmtId(elmtid++);
267  (*m_exp).push_back(prism);
268  }
269  else if((PyrGeom = std::dynamic_pointer_cast<
270  SpatialDomains::PyrGeom>(expIt.second->m_geomShPtr)))
271  {
273  = expIt.second->m_basisKeyVector[0];
275  = expIt.second->m_basisKeyVector[1];
277  = expIt.second->m_basisKeyVector[2];
278 
280  ::AllocateSharedPtr(PyrBa,PyrBb,PyrBc,
281  PyrGeom);
282  pyramid->SetElmtId(elmtid++);
283  (*m_exp).push_back(pyramid);
284  }
285  else if((HexGeom = std::dynamic_pointer_cast<
286  SpatialDomains::HexGeom>(expIt.second->m_geomShPtr)))
287  {
289  = expIt.second->m_basisKeyVector[0];
291  = expIt.second->m_basisKeyVector[1];
293  = expIt.second->m_basisKeyVector[2];
294 
296  ::AllocateSharedPtr(HexBa,HexBb,HexBc,
297  HexGeom);
298  hex->SetElmtId(elmtid++);
299  (*m_exp).push_back(hex);
300  }
301  else
302  {
303  ASSERTL0(false,"dynamic cast to a proper Geometry3D "
304  "failed");
305  }
306 
307  }
308 
309  // Setup Default optimisation information.
310  int nel = GetExpSize();
313 
315 
316  // Set up m_coeffs, m_phys.
319 
321  CreateCollections(ImpType);
322  }
323 
324  /**
325  * Given an expansion vector \a expansions, containing
326  * information about the domain and the spectral/hp element
327  * expansion, this constructor fills the list of local
328  * expansions \texttt{m_exp} with the proper expansions,
329  * calculates the total number of quadrature points
330  * \f$\boldsymbol{x}_i\f$ and the local expansion coefficients
331  * \f$\hat{u}^e_n\f$ and allocates memory for the arrays
332  * #m_coeffs and #m_phys.
333  *
334  * @param expansions An expansion vector, containing
335  * information about the domain and the
336  * spectral/hp element expansion.
337  */
339  const Collections::ImplementationType ImpType):
340  ExpList()
341  {
342  SetExpType(e3D);
343 
344  int elmtid = 0;
349 
350 
351  for(int i = 0; i < expansions.size(); ++i)
352  {
357 
358  auto expmap = expansions.find(i);
359  ASSERTL1(expmap != expansions.end(), "Unable to find expansion.");
360  const SpatialDomains::ExpansionShPtr exp = expmap->second;
361 
362  if((TetGeom = std::dynamic_pointer_cast<
363  SpatialDomains::TetGeom>(exp->m_geomShPtr)))
364  {
366  = exp->m_basisKeyVector[0];
368  = exp->m_basisKeyVector[1];
370  = exp->m_basisKeyVector[2];
371 
374  {
375  ASSERTL0(false,"LocalRegions::NodalTetExp is not "
376  "implemented yet");
377  }
378  else
379  {
381  ::AllocateSharedPtr(TetBa,TetBb,TetBc,
382  TetGeom);
383  tet->SetElmtId(elmtid++);
384  (*m_exp).push_back(tet);
385  }
386  }
387  else if((PrismGeom = std::dynamic_pointer_cast<
388  SpatialDomains::PrismGeom>(exp->m_geomShPtr)))
389  {
390  LibUtilities::BasisKey PrismBa
391  = exp->m_basisKeyVector[0];
392  LibUtilities::BasisKey PrismBb
393  = exp->m_basisKeyVector[1];
394  LibUtilities::BasisKey PrismBc
395  = exp->m_basisKeyVector[2];
396 
398  ::AllocateSharedPtr(PrismBa,PrismBb,
399  PrismBc,PrismGeom);
400  prism->SetElmtId(elmtid++);
401  (*m_exp).push_back(prism);
402  }
403  else if((PyrGeom = std::dynamic_pointer_cast<
404  SpatialDomains::PyrGeom>(exp->m_geomShPtr)))
405  {
407  = exp->m_basisKeyVector[0];
409  = exp->m_basisKeyVector[1];
411  = exp->m_basisKeyVector[2];
412 
414  ::AllocateSharedPtr(PyrBa,PyrBb,PyrBc,
415  PyrGeom);
416  pyramid->SetElmtId(elmtid++);
417  (*m_exp).push_back(pyramid);
418  }
419  else if((HexGeom = std::dynamic_pointer_cast<
420  SpatialDomains::HexGeom>(exp->m_geomShPtr)))
421  {
423  = exp->m_basisKeyVector[0];
425  = exp->m_basisKeyVector[1];
427  = exp->m_basisKeyVector[2];
428 
430  ::AllocateSharedPtr(HexBa,HexBb,HexBc,
431  HexGeom);
432  hex->SetElmtId(elmtid++);
433  (*m_exp).push_back(hex);
434  }
435  else
436  {
437  ASSERTL0(false,"dynamic cast to a proper Geometry3D "
438  "failed");
439  }
440 
441  }
442 
443  // Setup Default optimisation information.
444  int nel = GetExpSize();
447 
449 
450  // Set up m_coeffs, m_phys.
453 
454  CreateCollections(ImpType);
455  }
456 
458  {
459  Array<OneD, int> NumShape(4,0);
460 
461  for(int i = 0; i < GetExpSize(); ++i)
462  {
463  switch ((*m_exp)[i]->DetShapeType())
464  {
465  case LibUtilities::eTetrahedron: NumShape[0]++; break;
466  case LibUtilities::ePyramid: NumShape[1]++; break;
467  case LibUtilities::ePrism: NumShape[2]++; break;
468  case LibUtilities::eHexahedron: NumShape[3]++; break;
469  default:
470  ASSERTL0(false, "Unknown expansion type.");
471  break;
472  }
473  }
474 
475  int three = 3;
477  ::AllocateSharedPtr(m_session,three,NumShape);
478  }
479 
480  void ExpList3D::v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip)
481  {
482  boost::ignore_unused(istrip);
483 
484  int i,j,k;
485  int nquad0 = (*m_exp)[expansion]->GetNumPoints(0);
486  int nquad1 = (*m_exp)[expansion]->GetNumPoints(1);
487  int nquad2 = (*m_exp)[expansion]->GetNumPoints(2);
488  int ntot = nquad0*nquad1*nquad2;
489  int ntotminus = (nquad0-1)*(nquad1-1)*(nquad2-1);
490 
491  Array<OneD,NekDouble> coords[3];
492  coords[0] = Array<OneD,NekDouble>(ntot);
493  coords[1] = Array<OneD,NekDouble>(ntot);
494  coords[2] = Array<OneD,NekDouble>(ntot);
495  (*m_exp)[expansion]->GetCoords(coords[0],coords[1],coords[2]);
496 
497  outfile << " <Piece NumberOfPoints=\""
498  << ntot << "\" NumberOfCells=\""
499  << ntotminus << "\">" << endl;
500  outfile << " <Points>" << endl;
501  outfile << " <DataArray type=\"Float64\" "
502  << "NumberOfComponents=\"3\" format=\"ascii\">" << endl;
503  outfile << " ";
504  for (i = 0; i < ntot; ++i)
505  {
506  for (j = 0; j < 3; ++j)
507  {
508  outfile << setprecision(8) << scientific
509  << (float)coords[j][i] << " ";
510  }
511  outfile << endl;
512  }
513 
514  outfile << endl;
515  outfile << " </DataArray>" << endl;
516  outfile << " </Points>" << endl;
517  outfile << " <Cells>" << endl;
518  outfile << " <DataArray type=\"Int32\" "
519  << "Name=\"connectivity\" format=\"ascii\">" << endl;
520  for (i = 0; i < nquad0-1; ++i)
521  {
522  for (j = 0; j < nquad1-1; ++j)
523  {
524  for (k = 0; k < nquad2-1; ++k)
525  {
526  outfile << k*nquad0*nquad1 + j*nquad0 + i << " "
527  << k*nquad0*nquad1 + j*nquad0 + i + 1 << " "
528  << k*nquad0*nquad1 + (j+1)*nquad0 + i + 1 << " "
529  << k*nquad0*nquad1 + (j+1)*nquad0 + i << " "
530  << (k+1)*nquad0*nquad1 + j*nquad0 + i << " "
531  << (k+1)*nquad0*nquad1 + j*nquad0 + i + 1 << " "
532  << (k+1)*nquad0*nquad1 + (j+1)*nquad0 + i + 1 << " "
533  << (k+1)*nquad0*nquad1 + (j+1)*nquad0 + i << " " << endl;
534  }
535  }
536  }
537  outfile << endl;
538  outfile << " </DataArray>" << endl;
539  outfile << " <DataArray type=\"Int32\" "
540  << "Name=\"offsets\" format=\"ascii\">" << endl;
541  for (i = 0; i < ntotminus; ++i)
542  {
543  outfile << i*8+8 << " ";
544  }
545  outfile << endl;
546  outfile << " </DataArray>" << endl;
547  outfile << " <DataArray type=\"UInt8\" "
548  << "Name=\"types\" format=\"ascii\">" << endl;
549  for (i = 0; i < ntotminus; ++i)
550  {
551  outfile << "12 ";
552  }
553  outfile << endl;
554  outfile << " </DataArray>" << endl;
555  outfile << " </Cells>" << endl;
556  outfile << " <PointData>" << endl;
557  }
558 
560  {
561  int i, j;
562  for (i = 0; i < m_exp->size(); ++i)
563  {
564  for (j = 0; j < (*m_exp)[i]->GetNfaces(); ++j)
565  {
566  (*m_exp)[i]->ComputeFaceNormal(j);
567  }
568  }
569  }
570 
572  const Array<OneD, NekDouble> &inarray,
573  Array<OneD, NekDouble> &outarray)
574  {
575  int cnt,cnt1;
576 
577  cnt = cnt1 = 0;
578  for(int i = 0; i < GetExpSize(); ++i)
579  {
580  // get new points key
581  int pt0 = (*m_exp)[i]->GetNumPoints(0);
582  int pt1 = (*m_exp)[i]->GetNumPoints(1);
583  int pt2 = (*m_exp)[i]->GetNumPoints(2);
584  int npt0 = (int) pt0*scale;
585  int npt1 = (int) pt1*scale;
586  int npt2 = (int) pt2*scale;
587 
588  LibUtilities::PointsKey newPointsKey0(npt0,(*m_exp)[i]->GetPointsType(0));
589  LibUtilities::PointsKey newPointsKey1(npt1,(*m_exp)[i]->GetPointsType(1));
590  LibUtilities::PointsKey newPointsKey2(npt2,(*m_exp)[i]->GetPointsType(2));
591 
592  // Interpolate points;
593  LibUtilities::Interp3D((*m_exp)[i]->GetBasis(0)->GetPointsKey(),
594  (*m_exp)[i]->GetBasis(1)->GetPointsKey(),
595  (*m_exp)[i]->GetBasis(2)->GetPointsKey(),
596  &inarray[cnt], newPointsKey0,
597  newPointsKey1, newPointsKey2,
598  &outarray[cnt1]);
599 
600  cnt += pt0*pt1*pt2;
601  cnt1 += npt0*npt1*npt2;
602  }
603  }
604 
606  const Array<OneD, NekDouble> &inarray,
607  Array<OneD, NekDouble> &outarray)
608  {
609  int cnt,cnt1;
610 
611  cnt = cnt1 = 0;
612  for(int i = 0; i < GetExpSize(); ++i)
613  {
614  // get new points key
615  int pt0 = (*m_exp)[i]->GetNumPoints(0);
616  int pt1 = (*m_exp)[i]->GetNumPoints(1);
617  int pt2 = (*m_exp)[i]->GetNumPoints(2);
618  int npt0 = (int) pt0*scale;
619  int npt1 = (int) pt1*scale;
620  int npt2 = (int) pt2*scale;
621 
622  LibUtilities::PointsKey newPointsKey0(npt0,(*m_exp)[i]->GetPointsType(0));
623  LibUtilities::PointsKey newPointsKey1(npt1,(*m_exp)[i]->GetPointsType(1));
624  LibUtilities::PointsKey newPointsKey2(npt2,(*m_exp)[i]->GetPointsType(2));
625 
626  // Project points;
628  newPointsKey1,
629  newPointsKey2,
630  &inarray[cnt],
631  (*m_exp)[i]->GetBasis(0)->GetPointsKey(),
632  (*m_exp)[i]->GetBasis(1)->GetPointsKey(),
633  (*m_exp)[i]->GetBasis(2)->GetPointsKey(),
634  &outarray[cnt1]);
635 
636  cnt += npt0*npt1*npt2;
637  cnt1 += pt0*pt1*pt2;
638  }
639 
640  }
641  } //end of namespace
642 } //end of namespace
643 
int GetNumPoints() const
Return points order at which basis is defined.
Definition: Basis.h:133
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:163
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1090
virtual void v_SetUpPhysNormals()
Set up the normals on each expansion.
Definition: ExpList3D.cpp:559
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:1106
std::shared_ptr< TetGeom > TetGeomSharedPtr
Definition: TetGeom.h:88
std::shared_ptr< HexExp > HexExpSharedPtr
Definition: HexExp.h:56
BasisType GetBasisType() const
Return type of expansion basis.
Definition: Basis.h:144
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:55
STL namespace.
std::shared_ptr< TetExp > TetExpSharedPtr
Definition: TetExp.h:226
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:1069
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:194
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:1052
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2170
virtual void v_PhysGalerkinProjection1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList3D.cpp:605
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:103
int GetNumModes() const
Returns the order of the basis.
Definition: Basis.h:86
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:1030
void SetCoeffPhysOffsets()
Definition of the total number of degrees of freedom and quadrature points and offsets to access data...
Definition: ExpList.cpp:249
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1023
void PhysGalerkinProject3D(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)
Defines a specification for a set of points.
Definition: Points.h:59
double NekDouble
std::shared_ptr< Expansion > ExpansionShPtr
Definition: MeshGraph.h:151
virtual void v_ReadGlobalOptimizationParameters()
Definition: ExpList3D.cpp:457
std::shared_ptr< PyrExp > PyrExpSharedPtr
Definition: PyrExp.h:191
std::shared_ptr< PyrGeom > PyrGeomSharedPtr
Definition: PyrGeom.h:80
virtual ~ExpList3D()
Destructor.
Definition: ExpList3D.cpp:92
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
std::shared_ptr< PrismExp > PrismExpSharedPtr
Definition: PrismExp.h:220
std::shared_ptr< HexGeom > HexGeomSharedPtr
Definition: HexGeom.h:90
Abstraction of a three-dimensional multi-elemental expansion which is merely a collection of local ex...
Definition: ExpList3D.h:48
ExpList3D()
Default constructor.
Definition: ExpList3D.cpp:56
void CreateCollections(Collections::ImplementationType ImpType=Collections::eNoImpType)
Construct collections of elements containing a single element type and polynomial order from the list...
Definition: ExpList.cpp:3295
std::shared_ptr< PrismGeom > PrismGeomSharedPtr
Definition: PrismGeom.h:88
Lagrange for SEM basis .
Definition: BasisType.h:54
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:279
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
Describes the specification for a Basis.
Definition: Basis.h:49
std::shared_ptr< SessionReader > SessionReaderSharedPtr
virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip)
Definition: ExpList3D.cpp:480
std::map< int, ExpansionShPtr > ExpansionMap
Definition: MeshGraph.h:152
virtual void v_PhysInterp1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList3D.cpp:571