Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 // 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: Expansion list 3D definition
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #include <iomanip>
37 #include <MultiRegions/ExpList3D.h>
38 
39 #include <LocalRegions/HexExp.h>
40 #include <LocalRegions/PrismExp.h>
41 #include <LocalRegions/PyrExp.h>
42 #include <LocalRegions/TetExp.h>
43 
46 
47 using namespace std;
48 
49 namespace Nektar
50 {
51  namespace MultiRegions
52  {
53 
54  ExpList3D::ExpList3D(): ExpList()
55  {
56  SetExpType(e3D);
57  }
58 
60  {
61  SetExpType(e3D);
62  }
63 
65  const std::vector<unsigned int> &eIDs,
66  const bool DeclareCoeffPhysArrays,
67  const Collections::ImplementationType ImpType):
68  ExpList(In, eIDs, DeclareCoeffPhysArrays)
69  {
70  SetExpType(e3D);
71 
72  // Setup Default optimisation information.
73  int nel = GetExpSize();
76 
78 
79  if (DeclareCoeffPhysArrays)
80  {
81  // Set up m_coeffs, m_phys.
84  }
85 
87  CreateCollections(ImpType);
88  }
89 
91  {
92  }
93 
94 
96  const LibUtilities::BasisKey &TBa,
97  const LibUtilities::BasisKey &TBb,
98  const LibUtilities::BasisKey &TBc,
99  const LibUtilities::BasisKey &HBa,
100  const LibUtilities::BasisKey &HBb,
101  const LibUtilities::BasisKey &HBc,
102  const SpatialDomains::MeshGraphSharedPtr &graph3D,
103  const LibUtilities::PointsType TetNb,
104  const Collections::ImplementationType ImpType):
105  ExpList(pSession,graph3D)
106  {
107  SetExpType(e3D);
108 
113 
114  const SpatialDomains::ExpansionMap &expansions = graph3D->GetExpansions();
115 
116  SpatialDomains::ExpansionMap::const_iterator expIt;
117  for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
118  {
123 
124  if((TetGeom = boost::dynamic_pointer_cast<SpatialDomains::TetGeom>(expIt->second->m_geomShPtr)))
125  {
127  {
128 // Ntet = MemoryManager<LocalRegions::NodalTetExp>::AllocateSharedPtr(TetBa,TetBb,TetBc,TetNb,TetGeom);
129 // (*m_exp).push_back(Ntet);
130  }
131  else
132  {
134  (*m_exp).push_back(tet);
135  }
136 
138 
139  m_npoints += TBa.GetNumPoints()*TBb.GetNumPoints()*TBc.GetNumPoints();
140  }
141 /*
142  else if((PrismGeom = boost::dynamic_pointer_cast<SpatialDomains::PrismGeom>(expansions[i]->m_geomShPtr)))
143  {
144  prism = MemoryManager<LocalRegions::PrismExp>::AllocateSharedPtr(Ba,Bb,Bc,PrismGeom);
145  (*m_exp).push_back(prism);
146 
147  m_ncoeffs += StdRegions::StdPrismData::getNumberOfCoefficients(Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes());
148  m_npoints += Ba.GetNumPoints()*Bb.GetNumPoints()*Bc.GetNumPoints();
149 
150  }
151  else if((PyrGeom = boost::dynamic_pointer_cast<SpatialDomains::PyrGeom>(expansions[i]->m_geomShPtr)))
152  {
153  pyramid = MemoryManager<LocalRegions::PyrExp>::AllocateSharedPtr(Ba,Bb,Bc,PyrGeom);
154  (*m_exp).push_back(pyramid);
155 
156  m_ncoeffs += StdRegions::StdPyrData::getNumberOfCoefficients(Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes());
157  m_npoints += Ba.GetNumPoints()*Bb.GetNumPoints()*Bc.GetNumPoints();
158 
159  }
160 */
161  else if((HexGeom = boost::dynamic_pointer_cast<SpatialDomains::HexGeom>(expIt->second->m_geomShPtr)))
162  {
163  hex = MemoryManager<LocalRegions::HexExp>::AllocateSharedPtr(HBa,HBb,HBc, HexGeom);
164  (*m_exp).push_back(hex);
165 
166  m_ncoeffs += HBa.GetNumModes()*HBb.GetNumModes()*HBc.GetNumModes();
167  m_npoints += HBa.GetNumPoints()*HBb.GetNumPoints()*HBc.GetNumPoints();
168  }
169  else
170  {
171  ASSERTL0(false,"dynamic cast to a proper Geometry3D failed");
172  }
173 
174  }
175 
176  // Setup Default optimisation information.
177  int nel = GetExpSize();
180 
182 
183  // Set up m_coeffs, m_phys.
186 
188  CreateCollections(ImpType);
189  }
190 
191  /**
192  * Given a mesh \a graph3D, containing information about the domain and
193  * the spectral/hp element expansion, this constructor fills the list
194  * of local expansions \texttt{m_exp} with the proper expansions,
195  * calculates the total number of quadrature points
196  * \f$\boldsymbol{x}_i\f$ and the local expansion coefficients
197  * \f$\hat{u}^e_n\f$ and allocates memory for the arrays #m_coeffs and
198  * #m_phys.
199  *
200  * @param graph3D A mesh, containing information about the domain
201  * and the spectral/hp element expansion.
202  */
204  const SpatialDomains::MeshGraphSharedPtr &graph3D,
205  const std::string &variable,
206  const Collections::ImplementationType ImpType):
207  ExpList(pSession,graph3D)
208  {
209  SetExpType(e3D);
210 
211  int elmtid = 0;
216 
217  const SpatialDomains::ExpansionMap &expansions
218  = graph3D->GetExpansions(variable);
219 
220  SpatialDomains::ExpansionMap::const_iterator expIt;
221  for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
222  {
227 
228  if((TetGeom = boost::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 = boost::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 = boost::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 = boost::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  SpatialDomains::ExpansionMap::const_iterator expmap = expansions.find(i);
359  ASSERTL1(expmap != expansions.end(), "Unable to find expansion.");
360  const SpatialDomains::ExpansionShPtr exp = expmap->second;
361 
362  if((TetGeom = boost::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 = boost::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 = boost::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 = boost::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  int i,j,k;
483  int nquad0 = (*m_exp)[expansion]->GetNumPoints(0);
484  int nquad1 = (*m_exp)[expansion]->GetNumPoints(1);
485  int nquad2 = (*m_exp)[expansion]->GetNumPoints(2);
486  int ntot = nquad0*nquad1*nquad2;
487  int ntotminus = (nquad0-1)*(nquad1-1)*(nquad2-1);
488 
489  Array<OneD,NekDouble> coords[3];
490  coords[0] = Array<OneD,NekDouble>(ntot);
491  coords[1] = Array<OneD,NekDouble>(ntot);
492  coords[2] = Array<OneD,NekDouble>(ntot);
493  (*m_exp)[expansion]->GetCoords(coords[0],coords[1],coords[2]);
494 
495  outfile << " <Piece NumberOfPoints=\""
496  << ntot << "\" NumberOfCells=\""
497  << ntotminus << "\">" << endl;
498  outfile << " <Points>" << endl;
499  outfile << " <DataArray type=\"Float64\" "
500  << "NumberOfComponents=\"3\" format=\"ascii\">" << endl;
501  outfile << " ";
502  for (i = 0; i < ntot; ++i)
503  {
504  for (j = 0; j < 3; ++j)
505  {
506  outfile << setprecision(8) << scientific
507  << (float)coords[j][i] << " ";
508  }
509  outfile << endl;
510  }
511 
512  outfile << endl;
513  outfile << " </DataArray>" << endl;
514  outfile << " </Points>" << endl;
515  outfile << " <Cells>" << endl;
516  outfile << " <DataArray type=\"Int32\" "
517  << "Name=\"connectivity\" format=\"ascii\">" << endl;
518  for (i = 0; i < nquad0-1; ++i)
519  {
520  for (j = 0; j < nquad1-1; ++j)
521  {
522  for (k = 0; k < nquad2-1; ++k)
523  {
524  outfile << k*nquad0*nquad1 + j*nquad0 + i << " "
525  << k*nquad0*nquad1 + j*nquad0 + i + 1 << " "
526  << k*nquad0*nquad1 + (j+1)*nquad0 + i + 1 << " "
527  << k*nquad0*nquad1 + (j+1)*nquad0 + i << " "
528  << (k+1)*nquad0*nquad1 + j*nquad0 + i << " "
529  << (k+1)*nquad0*nquad1 + j*nquad0 + i + 1 << " "
530  << (k+1)*nquad0*nquad1 + (j+1)*nquad0 + i + 1 << " "
531  << (k+1)*nquad0*nquad1 + (j+1)*nquad0 + i << " " << endl;
532  }
533  }
534  }
535  outfile << endl;
536  outfile << " </DataArray>" << endl;
537  outfile << " <DataArray type=\"Int32\" "
538  << "Name=\"offsets\" format=\"ascii\">" << endl;
539  for (i = 0; i < ntotminus; ++i)
540  {
541  outfile << i*8+8 << " ";
542  }
543  outfile << endl;
544  outfile << " </DataArray>" << endl;
545  outfile << " <DataArray type=\"UInt8\" "
546  << "Name=\"types\" format=\"ascii\">" << endl;
547  for (i = 0; i < ntotminus; ++i)
548  {
549  outfile << "12 ";
550  }
551  outfile << endl;
552  outfile << " </DataArray>" << endl;
553  outfile << " </Cells>" << endl;
554  outfile << " <PointData>" << endl;
555  }
556 
558  {
559  int i, j;
560  for (i = 0; i < m_exp->size(); ++i)
561  {
562  for (j = 0; j < (*m_exp)[i]->GetNfaces(); ++j)
563  {
564  (*m_exp)[i]->ComputeFaceNormal(j);
565  }
566  }
567  }
568 
570  const Array<OneD, NekDouble> &inarray,
571  Array<OneD, NekDouble> &outarray)
572  {
573  int cnt,cnt1;
574 
575  cnt = cnt1 = 0;
576  for(int i = 0; i < GetExpSize(); ++i)
577  {
578  // get new points key
579  int pt0 = (*m_exp)[i]->GetNumPoints(0);
580  int pt1 = (*m_exp)[i]->GetNumPoints(1);
581  int pt2 = (*m_exp)[i]->GetNumPoints(2);
582  int npt0 = (int) pt0*scale;
583  int npt1 = (int) pt1*scale;
584  int npt2 = (int) pt2*scale;
585 
586  LibUtilities::PointsKey newPointsKey0(npt0,(*m_exp)[i]->GetPointsType(0));
587  LibUtilities::PointsKey newPointsKey1(npt1,(*m_exp)[i]->GetPointsType(1));
588  LibUtilities::PointsKey newPointsKey2(npt2,(*m_exp)[i]->GetPointsType(2));
589 
590  // Interpolate points;
591  LibUtilities::Interp3D((*m_exp)[i]->GetBasis(0)->GetPointsKey(),
592  (*m_exp)[i]->GetBasis(1)->GetPointsKey(),
593  (*m_exp)[i]->GetBasis(2)->GetPointsKey(),
594  &inarray[cnt], newPointsKey0,
595  newPointsKey1, newPointsKey2,
596  &outarray[cnt1]);
597 
598  cnt += pt0*pt1*pt2;
599  cnt1 += npt0*npt1*npt2;
600  }
601  }
602 
604  const Array<OneD, NekDouble> &inarray,
605  Array<OneD, NekDouble> &outarray)
606  {
607  int cnt,cnt1;
608 
609  cnt = cnt1 = 0;
610  for(int i = 0; i < GetExpSize(); ++i)
611  {
612  // get new points key
613  int pt0 = (*m_exp)[i]->GetNumPoints(0);
614  int pt1 = (*m_exp)[i]->GetNumPoints(1);
615  int pt2 = (*m_exp)[i]->GetNumPoints(2);
616  int npt0 = (int) pt0*scale;
617  int npt1 = (int) pt1*scale;
618  int npt2 = (int) pt2*scale;
619 
620  LibUtilities::PointsKey newPointsKey0(npt0,(*m_exp)[i]->GetPointsType(0));
621  LibUtilities::PointsKey newPointsKey1(npt1,(*m_exp)[i]->GetPointsType(1));
622  LibUtilities::PointsKey newPointsKey2(npt2,(*m_exp)[i]->GetPointsType(2));
623 
624  // Project points;
626  newPointsKey1,
627  newPointsKey2,
628  &inarray[cnt],
629  (*m_exp)[i]->GetBasis(0)->GetPointsKey(),
630  (*m_exp)[i]->GetBasis(1)->GetPointsKey(),
631  (*m_exp)[i]->GetBasis(2)->GetPointsKey(),
632  &outarray[cnt1]);
633 
634  cnt += npt0*npt1*npt2;
635  cnt1 += pt0*pt1*pt2;
636  }
637 
638  }
639  } //end of namespace
640 } //end of namespace
641 
boost::shared_ptr< PyrGeom > PyrGeomSharedPtr
Definition: PyrGeom.h:84
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
virtual void v_SetUpPhysNormals()
Set up the normals on each expansion.
Definition: ExpList3D.cpp:557
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:1052
BasisType GetBasisType() const
Return type of expansion basis.
Definition: Basis.h:139
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:54
STL namespace.
boost::shared_ptr< HexExp > HexExpSharedPtr
Definition: HexExp.h:57
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:1015
boost::shared_ptr< HexGeom > HexGeomSharedPtr
Definition: HexGeom.h:110
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:186
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:998
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2046
boost::shared_ptr< TetExp > TetExpSharedPtr
Definition: TetExp.h:228
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
virtual void v_PhysGalerkinProjection1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList3D.cpp:603
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:101
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1036
boost::shared_ptr< PyrExp > PyrExpSharedPtr
Definition: PyrExp.h:164
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:976
int GetNumPoints() const
Return points order at which basis is defined.
Definition: Basis.h:128
void SetCoeffPhysOffsets()
Definition of the total number of degrees of freedom and quadrature points and offsets to access data...
Definition: ExpList.cpp:247
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:969
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:58
double NekDouble
boost::shared_ptr< PrismExp > PrismExpSharedPtr
Definition: PrismExp.h:222
virtual void v_ReadGlobalOptimizationParameters()
Definition: ExpList3D.cpp:457
virtual ~ExpList3D()
Destructor.
Definition: ExpList3D.cpp:90
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:186
boost::shared_ptr< Expansion > ExpansionShPtr
Definition: MeshGraph.h:173
boost::shared_ptr< PrismGeom > PrismGeomSharedPtr
Definition: PrismGeom.h:109
boost::shared_ptr< TetGeom > TetGeomSharedPtr
Definition: TetGeom.h:106
Abstraction of a three-dimensional multi-elemental expansion which is merely a collection of local ex...
Definition: ExpList3D.h:49
ExpList3D()
Default constructor.
Definition: ExpList3D.cpp:54
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:3151
int GetNumModes() const
Returns the order of the basis.
Definition: Basis.h:84
Lagrange for SEM basis .
Definition: BasisType.h:53
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:277
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:442
Describes the specification for a Basis.
Definition: Basis.h:50
virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip)
Definition: ExpList3D.cpp:480
std::map< int, ExpansionShPtr > ExpansionMap
Definition: MeshGraph.h:174
virtual void v_PhysInterp1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList3D.cpp:569