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 // 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 namespace Nektar
48 {
49  namespace MultiRegions
50  {
51 
53  {
54  SetExpType(e3D);
55  }
56 
58  {
59  SetExpType(e3D);
60  }
61 
63  {
64  }
65 
66 
68  const LibUtilities::BasisKey &TBa,
69  const LibUtilities::BasisKey &TBb,
70  const LibUtilities::BasisKey &TBc,
71  const LibUtilities::BasisKey &HBa,
72  const LibUtilities::BasisKey &HBb,
73  const LibUtilities::BasisKey &HBc,
75  const LibUtilities::PointsType TetNb):
76  ExpList(pSession,graph3D)
77  {
78  SetExpType(e3D);
79 
84 
85  const SpatialDomains::ExpansionMap &expansions = graph3D->GetExpansions();
86 
87  SpatialDomains::ExpansionMap::const_iterator expIt;
88  for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
89  {
94 
95  if((TetGeom = boost::dynamic_pointer_cast<SpatialDomains::TetGeom>(expIt->second->m_geomShPtr)))
96  {
98  {
99 // Ntet = MemoryManager<LocalRegions::NodalTetExp>::AllocateSharedPtr(TetBa,TetBb,TetBc,TetNb,TetGeom);
100 // (*m_exp).push_back(Ntet);
101  }
102  else
103  {
105  (*m_exp).push_back(tet);
106  }
107 
109 
110  m_npoints += TBa.GetNumPoints()*TBb.GetNumPoints()*TBc.GetNumPoints();
111  }
112 /*
113  else if((PrismGeom = boost::dynamic_pointer_cast<SpatialDomains::PrismGeom>(expansions[i]->m_geomShPtr)))
114  {
115  prism = MemoryManager<LocalRegions::PrismExp>::AllocateSharedPtr(Ba,Bb,Bc,PrismGeom);
116  (*m_exp).push_back(prism);
117 
118  m_ncoeffs += StdRegions::StdPrismData::getNumberOfCoefficients(Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes());
119  m_npoints += Ba.GetNumPoints()*Bb.GetNumPoints()*Bc.GetNumPoints();
120 
121  }
122  else if((PyrGeom = boost::dynamic_pointer_cast<SpatialDomains::PyrGeom>(expansions[i]->m_geomShPtr)))
123  {
124  pyramid = MemoryManager<LocalRegions::PyrExp>::AllocateSharedPtr(Ba,Bb,Bc,PyrGeom);
125  (*m_exp).push_back(pyramid);
126 
127  m_ncoeffs += StdRegions::StdPyrData::getNumberOfCoefficients(Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes());
128  m_npoints += Ba.GetNumPoints()*Bb.GetNumPoints()*Bc.GetNumPoints();
129 
130  }
131 */
132  else if((HexGeom = boost::dynamic_pointer_cast<SpatialDomains::HexGeom>(expIt->second->m_geomShPtr)))
133  {
134  hex = MemoryManager<LocalRegions::HexExp>::AllocateSharedPtr(HBa,HBb,HBc, HexGeom);
135  (*m_exp).push_back(hex);
136 
137  m_ncoeffs += HBa.GetNumModes()*HBb.GetNumModes()*HBc.GetNumModes();
138  m_npoints += HBa.GetNumPoints()*HBb.GetNumPoints()*HBc.GetNumPoints();
139  }
140  else
141  {
142  ASSERTL0(false,"dynamic cast to a proper Geometry3D failed");
143  }
144 
145  }
146 
147  // Setup Default optimisation information.
148  int nel = GetExpSize();
151 
152  SetCoeffPhys();
153 
156  }
157 
158  /**
159  * Given a mesh \a graph3D, containing information about the domain and
160  * the spectral/hp element expansion, this constructor fills the list
161  * of local expansions \texttt{m_exp} with the proper expansions,
162  * calculates the total number of quadrature points
163  * \f$\boldsymbol{x}_i\f$ and the local expansion coefficients
164  * \f$\hat{u}^e_n\f$ and allocates memory for the arrays #m_coeffs and
165  * #m_phys.
166  *
167  * @param graph3D A mesh, containing information about the domain
168  * and the spectral/hp element expansion.
169  */
171  const SpatialDomains::MeshGraphSharedPtr &graph3D,
172  const std::string &variable) :
173  ExpList(pSession,graph3D)
174  {
175  SetExpType(e3D);
176 
177  int elmtid = 0;
182 
183  const SpatialDomains::ExpansionMap &expansions
184  = graph3D->GetExpansions(variable);
185 
186  SpatialDomains::ExpansionMap::const_iterator expIt;
187  for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
188  {
193 
194  if((TetGeom = boost::dynamic_pointer_cast<
195  SpatialDomains::TetGeom>(expIt->second->m_geomShPtr)))
196  {
198  = expIt->second->m_basisKeyVector[0];
200  = expIt->second->m_basisKeyVector[1];
202  = expIt->second->m_basisKeyVector[2];
203 
206  {
207  ASSERTL0(false,"LocalRegions::NodalTetExp is not "
208  "implemented yet");
209  }
210  else
211  {
213  ::AllocateSharedPtr(TetBa,TetBb,TetBc,
214  TetGeom);
215  tet->SetElmtId(elmtid++);
216  (*m_exp).push_back(tet);
217  }
218  }
219  else if((PrismGeom = boost::dynamic_pointer_cast<SpatialDomains
220  ::PrismGeom>(expIt->second->m_geomShPtr)))
221  {
222  LibUtilities::BasisKey PrismBa
223  = expIt->second->m_basisKeyVector[0];
224  LibUtilities::BasisKey PrismBb
225  = expIt->second->m_basisKeyVector[1];
226  LibUtilities::BasisKey PrismBc
227  = expIt->second->m_basisKeyVector[2];
228 
230  ::AllocateSharedPtr(PrismBa,PrismBb,
231  PrismBc,PrismGeom);
232  prism->SetElmtId(elmtid++);
233  (*m_exp).push_back(prism);
234  }
235  else if((PyrGeom = boost::dynamic_pointer_cast<
236  SpatialDomains::PyrGeom>(expIt->second->m_geomShPtr)))
237  {
239  = expIt->second->m_basisKeyVector[0];
241  = expIt->second->m_basisKeyVector[1];
243  = expIt->second->m_basisKeyVector[2];
244 
246  ::AllocateSharedPtr(PyrBa,PyrBb,PyrBc,
247  PyrGeom);
248  pyramid->SetElmtId(elmtid++);
249  (*m_exp).push_back(pyramid);
250  }
251  else if((HexGeom = boost::dynamic_pointer_cast<
252  SpatialDomains::HexGeom>(expIt->second->m_geomShPtr)))
253  {
255  = expIt->second->m_basisKeyVector[0];
257  = expIt->second->m_basisKeyVector[1];
259  = expIt->second->m_basisKeyVector[2];
260 
262  ::AllocateSharedPtr(HexBa,HexBb,HexBc,
263  HexGeom);
264  hex->SetElmtId(elmtid++);
265  (*m_exp).push_back(hex);
266  }
267  else
268  {
269  ASSERTL0(false,"dynamic cast to a proper Geometry3D "
270  "failed");
271  }
272 
273  }
274 
275  // Setup Default optimisation information.
276  int nel = GetExpSize();
279 
280  SetCoeffPhys();
283  }
284 
285  /**
286  * Given an expansion vector \a expansions, containing
287  * information about the domain and the spectral/hp element
288  * expansion, this constructor fills the list of local
289  * expansions \texttt{m_exp} with the proper expansions,
290  * calculates the total number of quadrature points
291  * \f$\boldsymbol{x}_i\f$ and the local expansion coefficients
292  * \f$\hat{u}^e_n\f$ and allocates memory for the arrays
293  * #m_coeffs and #m_phys.
294  *
295  * @param expansions An expansion vector, containing
296  * information about the domain and the
297  * spectral/hp element expansion.
298  */
300  ExpList()
301  {
302  SetExpType(e3D);
303 
304  int elmtid = 0;
309 
310 
311  for(int i = 0; i < expansions.size(); ++i)
312  {
317 
318  SpatialDomains::ExpansionMap::const_iterator expmap = expansions.find(i);
319  ASSERTL1(expmap != expansions.end(), "Unable to find expansion.");
320  const SpatialDomains::ExpansionShPtr exp = expmap->second;
321 
322  if((TetGeom = boost::dynamic_pointer_cast<
323  SpatialDomains::TetGeom>(exp->m_geomShPtr)))
324  {
326  = exp->m_basisKeyVector[0];
328  = exp->m_basisKeyVector[1];
330  = exp->m_basisKeyVector[2];
331 
334  {
335  ASSERTL0(false,"LocalRegions::NodalTetExp is not "
336  "implemented yet");
337  }
338  else
339  {
341  ::AllocateSharedPtr(TetBa,TetBb,TetBc,
342  TetGeom);
343  tet->SetElmtId(elmtid++);
344  (*m_exp).push_back(tet);
345  }
346  }
347  else if((PrismGeom = boost::dynamic_pointer_cast<
348  SpatialDomains::PrismGeom>(exp->m_geomShPtr)))
349  {
350  LibUtilities::BasisKey PrismBa
351  = exp->m_basisKeyVector[0];
352  LibUtilities::BasisKey PrismBb
353  = exp->m_basisKeyVector[1];
354  LibUtilities::BasisKey PrismBc
355  = exp->m_basisKeyVector[2];
356 
358  ::AllocateSharedPtr(PrismBa,PrismBb,
359  PrismBc,PrismGeom);
360  prism->SetElmtId(elmtid++);
361  (*m_exp).push_back(prism);
362  }
363  else if((PyrGeom = boost::dynamic_pointer_cast<
364  SpatialDomains::PyrGeom>(exp->m_geomShPtr)))
365  {
367  = exp->m_basisKeyVector[0];
369  = exp->m_basisKeyVector[1];
371  = exp->m_basisKeyVector[2];
372 
374  ::AllocateSharedPtr(PyrBa,PyrBb,PyrBc,
375  PyrGeom);
376  pyramid->SetElmtId(elmtid++);
377  (*m_exp).push_back(pyramid);
378  }
379  else if((HexGeom = boost::dynamic_pointer_cast<
380  SpatialDomains::HexGeom>(exp->m_geomShPtr)))
381  {
383  = exp->m_basisKeyVector[0];
385  = exp->m_basisKeyVector[1];
387  = exp->m_basisKeyVector[2];
388 
390  ::AllocateSharedPtr(HexBa,HexBb,HexBc,
391  HexGeom);
392  hex->SetElmtId(elmtid++);
393  (*m_exp).push_back(hex);
394  }
395  else
396  {
397  ASSERTL0(false,"dynamic cast to a proper Geometry3D "
398  "failed");
399  }
400 
401  }
402 
403  // Setup Default optimisation information.
404  int nel = GetExpSize();
407 
408  SetCoeffPhys();
410  }
411 
412  /**
413  * Set up the storage for the concatenated list of
414  * coefficients and physical evaluations at the quadrature
415  * points. Each expansion (local element) is processed in turn
416  * to determine the number of coefficients and physical data
417  * points it contributes to the domain. Three arrays,
418  * #m_coeff_offset, #m_phys_offset and #m_offset_elmt_id, are
419  * also initialised and updated to store the data offsets of
420  * each element in the #m_coeffs and #m_phys arrays, and the
421  * element id that each consecutive block is associated
422  * respectively.
423  */
425  {
426  int i;
427 
428  // Set up offset information and array sizes
432 
433  m_ncoeffs = m_npoints = 0;
434 
435  for(i = 0; i < m_exp->size(); ++i)
436  {
438  m_phys_offset [i] = m_npoints;
439  m_offset_elmt_id[i] = i;
440  m_ncoeffs += (*m_exp)[i]->GetNcoeffs();
441  m_npoints += (*m_exp)[i]->GetTotPoints();
442  }
443 
446  }
447 
449  {
450  Array<OneD, int> NumShape(4,0);
451 
452  for(int i = 0; i < GetExpSize(); ++i)
453  {
454  switch ((*m_exp)[i]->DetShapeType())
455  {
456  case LibUtilities::eTetrahedron: NumShape[0]++; break;
457  case LibUtilities::ePyramid: NumShape[1]++; break;
458  case LibUtilities::ePrism: NumShape[2]++; break;
459  case LibUtilities::eHexahedron: NumShape[3]++; break;
460  default:
461  ASSERTL0(false, "Unknown expansion type.");
462  break;
463  }
464  }
465 
466  int three = 3;
468  ::AllocateSharedPtr(m_session,three,NumShape);
469  }
470 
471  void ExpList3D::v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip)
472  {
473  int i,j,k;
474  int nquad0 = (*m_exp)[expansion]->GetNumPoints(0);
475  int nquad1 = (*m_exp)[expansion]->GetNumPoints(1);
476  int nquad2 = (*m_exp)[expansion]->GetNumPoints(2);
477  int ntot = nquad0*nquad1*nquad2;
478  int ntotminus = (nquad0-1)*(nquad1-1)*(nquad2-1);
479 
480  Array<OneD,NekDouble> coords[3];
481  coords[0] = Array<OneD,NekDouble>(ntot);
482  coords[1] = Array<OneD,NekDouble>(ntot);
483  coords[2] = Array<OneD,NekDouble>(ntot);
484  (*m_exp)[expansion]->GetCoords(coords[0],coords[1],coords[2]);
485 
486  outfile << " <Piece NumberOfPoints=\""
487  << ntot << "\" NumberOfCells=\""
488  << ntotminus << "\">" << endl;
489  outfile << " <Points>" << endl;
490  outfile << " <DataArray type=\"Float64\" "
491  << "NumberOfComponents=\"3\" format=\"ascii\">" << endl;
492  outfile << " ";
493  for (i = 0; i < ntot; ++i)
494  {
495  for (j = 0; j < 3; ++j)
496  {
497  outfile << setprecision(8) << scientific
498  << (float)coords[j][i] << " ";
499  }
500  outfile << endl;
501  }
502 
503  outfile << endl;
504  outfile << " </DataArray>" << endl;
505  outfile << " </Points>" << endl;
506  outfile << " <Cells>" << endl;
507  outfile << " <DataArray type=\"Int32\" "
508  << "Name=\"connectivity\" format=\"ascii\">" << endl;
509  for (i = 0; i < nquad0-1; ++i)
510  {
511  for (j = 0; j < nquad1-1; ++j)
512  {
513  for (k = 0; k < nquad2-1; ++k)
514  {
515  outfile << k*nquad0*nquad1 + j*nquad0 + i << " "
516  << k*nquad0*nquad1 + j*nquad0 + i + 1 << " "
517  << k*nquad0*nquad1 + (j+1)*nquad0 + i + 1 << " "
518  << k*nquad0*nquad1 + (j+1)*nquad0 + i << " "
519  << (k+1)*nquad0*nquad1 + j*nquad0 + i << " "
520  << (k+1)*nquad0*nquad1 + j*nquad0 + i + 1 << " "
521  << (k+1)*nquad0*nquad1 + (j+1)*nquad0 + i + 1 << " "
522  << (k+1)*nquad0*nquad1 + (j+1)*nquad0 + i << " " << endl;
523  }
524  }
525  }
526  outfile << endl;
527  outfile << " </DataArray>" << endl;
528  outfile << " <DataArray type=\"Int32\" "
529  << "Name=\"offsets\" format=\"ascii\">" << endl;
530  for (i = 0; i < ntotminus; ++i)
531  {
532  outfile << i*8+8 << " ";
533  }
534  outfile << endl;
535  outfile << " </DataArray>" << endl;
536  outfile << " <DataArray type=\"UInt8\" "
537  << "Name=\"types\" format=\"ascii\">" << endl;
538  for (i = 0; i < ntotminus; ++i)
539  {
540  outfile << "12 ";
541  }
542  outfile << endl;
543  outfile << " </DataArray>" << endl;
544  outfile << " </Cells>" << endl;
545  outfile << " <PointData>" << endl;
546  }
547 
549  {
550  int i, j;
551  for (i = 0; i < m_exp->size(); ++i)
552  {
553  for (j = 0; j < (*m_exp)[i]->GetNfaces(); ++j)
554  {
555  (*m_exp)[i]->ComputeFaceNormal(j);
556  }
557  }
558  }
559 
561  const Array<OneD, NekDouble> &inarray,
562  Array<OneD, NekDouble> &outarray)
563  {
564  int cnt,cnt1;
565 
566  cnt = cnt1 = 0;
567  for(int i = 0; i < GetExpSize(); ++i)
568  {
569  // get new points key
570  int pt0 = (*m_exp)[i]->GetNumPoints(0);
571  int pt1 = (*m_exp)[i]->GetNumPoints(1);
572  int pt2 = (*m_exp)[i]->GetNumPoints(2);
573  int npt0 = (int) pt0*scale;
574  int npt1 = (int) pt1*scale;
575  int npt2 = (int) pt2*scale;
576 
577  LibUtilities::PointsKey newPointsKey0(npt0,(*m_exp)[i]->GetPointsType(0));
578  LibUtilities::PointsKey newPointsKey1(npt1,(*m_exp)[i]->GetPointsType(1));
579  LibUtilities::PointsKey newPointsKey2(npt2,(*m_exp)[i]->GetPointsType(2));
580 
581  // Interpolate points;
582  LibUtilities::Interp3D((*m_exp)[i]->GetBasis(0)->GetPointsKey(),
583  (*m_exp)[i]->GetBasis(1)->GetPointsKey(),
584  (*m_exp)[i]->GetBasis(2)->GetPointsKey(),
585  &inarray[cnt], newPointsKey0,
586  newPointsKey1, newPointsKey2,
587  &outarray[cnt1]);
588 
589  cnt += pt0*pt1*pt2;
590  cnt1 += npt0*npt1*npt2;
591  }
592  }
593 
595  const Array<OneD, NekDouble> &inarray,
596  Array<OneD, NekDouble> &outarray)
597  {
598  int cnt,cnt1;
599 
600  cnt = cnt1 = 0;
601  for(int i = 0; i < GetExpSize(); ++i)
602  {
603  // get new points key
604  int pt0 = (*m_exp)[i]->GetNumPoints(0);
605  int pt1 = (*m_exp)[i]->GetNumPoints(1);
606  int pt2 = (*m_exp)[i]->GetNumPoints(2);
607  int npt0 = (int) pt0*scale;
608  int npt1 = (int) pt1*scale;
609  int npt2 = (int) pt2*scale;
610 
611  LibUtilities::PointsKey newPointsKey0(npt0,(*m_exp)[i]->GetPointsType(0));
612  LibUtilities::PointsKey newPointsKey1(npt1,(*m_exp)[i]->GetPointsType(1));
613  LibUtilities::PointsKey newPointsKey2(npt2,(*m_exp)[i]->GetPointsType(2));
614 
615  // Project points;
617  newPointsKey1,
618  newPointsKey2,
619  &inarray[cnt],
620  (*m_exp)[i]->GetBasis(0)->GetPointsKey(),
621  (*m_exp)[i]->GetBasis(1)->GetPointsKey(),
622  (*m_exp)[i]->GetBasis(2)->GetPointsKey(),
623  &outarray[cnt1]);
624 
625  cnt += npt0*npt1*npt2;
626  cnt1 += pt0*pt1*pt2;
627  }
628 
629  }
630  } //end of namespace
631 } //end of namespace
632 
boost::shared_ptr< PyrGeom > PyrGeomSharedPtr
Definition: PyrGeom.h:83
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
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:548
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:971
BasisType GetBasisType() const
Return type of expansion basis.
Definition: Basis.h:139
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:54
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:926
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:909
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1837
boost::shared_ptr< TetExp > TetExpSharedPtr
Definition: TetExp.h:234
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:50
virtual void v_PhysGalerkinProjection1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList3D.cpp:594
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:958
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:947
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:961
boost::shared_ptr< PyrExp > PyrExpSharedPtr
Definition: PyrExp.h:163
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:887
Array< OneD, int > m_offset_elmt_id
Array containing the element id m_offset_elmt_id[n] that the n^th consecutive block of data in m_coef...
Definition: ExpList.h:969
int GetNumPoints() const
Return points order at which basis is defined.
Definition: Basis.h:128
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:880
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:226
virtual void v_ReadGlobalOptimizationParameters()
Definition: ExpList3D.cpp:448
virtual ~ExpList3D()
Destructor.
Definition: ExpList3D.cpp:62
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:178
boost::shared_ptr< Expansion > ExpansionShPtr
Definition: MeshGraph.h:170
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:52
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:2723
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:211
void SetCoeffPhys(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Definition: ExpList3D.cpp:424
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:165
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:432
Describes the specification for a Basis.
Definition: Basis.h:50
virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip)
Definition: ExpList3D.cpp:471
std::map< int, ExpansionShPtr > ExpansionMap
Definition: MeshGraph.h:171
virtual void v_PhysInterp1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList3D.cpp:560