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  : ExpList(In, eIDs, DeclareCoeffPhysArrays)
68  {
69  SetExpType(e3D);
70 
71  // Setup Default optimisation information.
72  int nel = GetExpSize();
75 
76  SetCoeffPhys();
77 
80  }
81 
83  {
84  }
85 
86 
88  const LibUtilities::BasisKey &TBa,
89  const LibUtilities::BasisKey &TBb,
90  const LibUtilities::BasisKey &TBc,
91  const LibUtilities::BasisKey &HBa,
92  const LibUtilities::BasisKey &HBb,
93  const LibUtilities::BasisKey &HBc,
95  const LibUtilities::PointsType TetNb):
96  ExpList(pSession,graph3D)
97  {
98  SetExpType(e3D);
99 
104 
105  const SpatialDomains::ExpansionMap &expansions = graph3D->GetExpansions();
106 
107  SpatialDomains::ExpansionMap::const_iterator expIt;
108  for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
109  {
114 
115  if((TetGeom = boost::dynamic_pointer_cast<SpatialDomains::TetGeom>(expIt->second->m_geomShPtr)))
116  {
118  {
119 // Ntet = MemoryManager<LocalRegions::NodalTetExp>::AllocateSharedPtr(TetBa,TetBb,TetBc,TetNb,TetGeom);
120 // (*m_exp).push_back(Ntet);
121  }
122  else
123  {
125  (*m_exp).push_back(tet);
126  }
127 
129 
130  m_npoints += TBa.GetNumPoints()*TBb.GetNumPoints()*TBc.GetNumPoints();
131  }
132 /*
133  else if((PrismGeom = boost::dynamic_pointer_cast<SpatialDomains::PrismGeom>(expansions[i]->m_geomShPtr)))
134  {
135  prism = MemoryManager<LocalRegions::PrismExp>::AllocateSharedPtr(Ba,Bb,Bc,PrismGeom);
136  (*m_exp).push_back(prism);
137 
138  m_ncoeffs += StdRegions::StdPrismData::getNumberOfCoefficients(Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes());
139  m_npoints += Ba.GetNumPoints()*Bb.GetNumPoints()*Bc.GetNumPoints();
140 
141  }
142  else if((PyrGeom = boost::dynamic_pointer_cast<SpatialDomains::PyrGeom>(expansions[i]->m_geomShPtr)))
143  {
144  pyramid = MemoryManager<LocalRegions::PyrExp>::AllocateSharedPtr(Ba,Bb,Bc,PyrGeom);
145  (*m_exp).push_back(pyramid);
146 
147  m_ncoeffs += StdRegions::StdPyrData::getNumberOfCoefficients(Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes());
148  m_npoints += Ba.GetNumPoints()*Bb.GetNumPoints()*Bc.GetNumPoints();
149 
150  }
151 */
152  else if((HexGeom = boost::dynamic_pointer_cast<SpatialDomains::HexGeom>(expIt->second->m_geomShPtr)))
153  {
154  hex = MemoryManager<LocalRegions::HexExp>::AllocateSharedPtr(HBa,HBb,HBc, HexGeom);
155  (*m_exp).push_back(hex);
156 
157  m_ncoeffs += HBa.GetNumModes()*HBb.GetNumModes()*HBc.GetNumModes();
158  m_npoints += HBa.GetNumPoints()*HBb.GetNumPoints()*HBc.GetNumPoints();
159  }
160  else
161  {
162  ASSERTL0(false,"dynamic cast to a proper Geometry3D failed");
163  }
164 
165  }
166 
167  // Setup Default optimisation information.
168  int nel = GetExpSize();
171 
172  SetCoeffPhys();
173 
176  }
177 
178  /**
179  * Given a mesh \a graph3D, containing information about the domain and
180  * the spectral/hp element expansion, this constructor fills the list
181  * of local expansions \texttt{m_exp} with the proper expansions,
182  * calculates the total number of quadrature points
183  * \f$\boldsymbol{x}_i\f$ and the local expansion coefficients
184  * \f$\hat{u}^e_n\f$ and allocates memory for the arrays #m_coeffs and
185  * #m_phys.
186  *
187  * @param graph3D A mesh, containing information about the domain
188  * and the spectral/hp element expansion.
189  */
191  const SpatialDomains::MeshGraphSharedPtr &graph3D,
192  const std::string &variable) :
193  ExpList(pSession,graph3D)
194  {
195  SetExpType(e3D);
196 
197  int elmtid = 0;
202 
203  const SpatialDomains::ExpansionMap &expansions
204  = graph3D->GetExpansions(variable);
205 
206  SpatialDomains::ExpansionMap::const_iterator expIt;
207  for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
208  {
213 
214  if((TetGeom = boost::dynamic_pointer_cast<
215  SpatialDomains::TetGeom>(expIt->second->m_geomShPtr)))
216  {
218  = expIt->second->m_basisKeyVector[0];
220  = expIt->second->m_basisKeyVector[1];
222  = expIt->second->m_basisKeyVector[2];
223 
226  {
227  ASSERTL0(false,"LocalRegions::NodalTetExp is not "
228  "implemented yet");
229  }
230  else
231  {
233  ::AllocateSharedPtr(TetBa,TetBb,TetBc,
234  TetGeom);
235  tet->SetElmtId(elmtid++);
236  (*m_exp).push_back(tet);
237  }
238  }
239  else if((PrismGeom = boost::dynamic_pointer_cast<SpatialDomains
240  ::PrismGeom>(expIt->second->m_geomShPtr)))
241  {
242  LibUtilities::BasisKey PrismBa
243  = expIt->second->m_basisKeyVector[0];
244  LibUtilities::BasisKey PrismBb
245  = expIt->second->m_basisKeyVector[1];
246  LibUtilities::BasisKey PrismBc
247  = expIt->second->m_basisKeyVector[2];
248 
250  ::AllocateSharedPtr(PrismBa,PrismBb,
251  PrismBc,PrismGeom);
252  prism->SetElmtId(elmtid++);
253  (*m_exp).push_back(prism);
254  }
255  else if((PyrGeom = boost::dynamic_pointer_cast<
256  SpatialDomains::PyrGeom>(expIt->second->m_geomShPtr)))
257  {
259  = expIt->second->m_basisKeyVector[0];
261  = expIt->second->m_basisKeyVector[1];
263  = expIt->second->m_basisKeyVector[2];
264 
266  ::AllocateSharedPtr(PyrBa,PyrBb,PyrBc,
267  PyrGeom);
268  pyramid->SetElmtId(elmtid++);
269  (*m_exp).push_back(pyramid);
270  }
271  else if((HexGeom = boost::dynamic_pointer_cast<
272  SpatialDomains::HexGeom>(expIt->second->m_geomShPtr)))
273  {
275  = expIt->second->m_basisKeyVector[0];
277  = expIt->second->m_basisKeyVector[1];
279  = expIt->second->m_basisKeyVector[2];
280 
282  ::AllocateSharedPtr(HexBa,HexBb,HexBc,
283  HexGeom);
284  hex->SetElmtId(elmtid++);
285  (*m_exp).push_back(hex);
286  }
287  else
288  {
289  ASSERTL0(false,"dynamic cast to a proper Geometry3D "
290  "failed");
291  }
292 
293  }
294 
295  // Setup Default optimisation information.
296  int nel = GetExpSize();
299 
300  SetCoeffPhys();
303  }
304 
305  /**
306  * Given an expansion vector \a expansions, containing
307  * information about the domain and the spectral/hp element
308  * expansion, this constructor fills the list of local
309  * expansions \texttt{m_exp} with the proper expansions,
310  * calculates the total number of quadrature points
311  * \f$\boldsymbol{x}_i\f$ and the local expansion coefficients
312  * \f$\hat{u}^e_n\f$ and allocates memory for the arrays
313  * #m_coeffs and #m_phys.
314  *
315  * @param expansions An expansion vector, containing
316  * information about the domain and the
317  * spectral/hp element expansion.
318  */
320  ExpList()
321  {
322  SetExpType(e3D);
323 
324  int elmtid = 0;
329 
330 
331  for(int i = 0; i < expansions.size(); ++i)
332  {
337 
338  SpatialDomains::ExpansionMap::const_iterator expmap = expansions.find(i);
339  ASSERTL1(expmap != expansions.end(), "Unable to find expansion.");
340  const SpatialDomains::ExpansionShPtr exp = expmap->second;
341 
342  if((TetGeom = boost::dynamic_pointer_cast<
343  SpatialDomains::TetGeom>(exp->m_geomShPtr)))
344  {
346  = exp->m_basisKeyVector[0];
348  = exp->m_basisKeyVector[1];
350  = exp->m_basisKeyVector[2];
351 
354  {
355  ASSERTL0(false,"LocalRegions::NodalTetExp is not "
356  "implemented yet");
357  }
358  else
359  {
361  ::AllocateSharedPtr(TetBa,TetBb,TetBc,
362  TetGeom);
363  tet->SetElmtId(elmtid++);
364  (*m_exp).push_back(tet);
365  }
366  }
367  else if((PrismGeom = boost::dynamic_pointer_cast<
368  SpatialDomains::PrismGeom>(exp->m_geomShPtr)))
369  {
370  LibUtilities::BasisKey PrismBa
371  = exp->m_basisKeyVector[0];
372  LibUtilities::BasisKey PrismBb
373  = exp->m_basisKeyVector[1];
374  LibUtilities::BasisKey PrismBc
375  = exp->m_basisKeyVector[2];
376 
378  ::AllocateSharedPtr(PrismBa,PrismBb,
379  PrismBc,PrismGeom);
380  prism->SetElmtId(elmtid++);
381  (*m_exp).push_back(prism);
382  }
383  else if((PyrGeom = boost::dynamic_pointer_cast<
384  SpatialDomains::PyrGeom>(exp->m_geomShPtr)))
385  {
387  = exp->m_basisKeyVector[0];
389  = exp->m_basisKeyVector[1];
391  = exp->m_basisKeyVector[2];
392 
394  ::AllocateSharedPtr(PyrBa,PyrBb,PyrBc,
395  PyrGeom);
396  pyramid->SetElmtId(elmtid++);
397  (*m_exp).push_back(pyramid);
398  }
399  else if((HexGeom = boost::dynamic_pointer_cast<
400  SpatialDomains::HexGeom>(exp->m_geomShPtr)))
401  {
403  = exp->m_basisKeyVector[0];
405  = exp->m_basisKeyVector[1];
407  = exp->m_basisKeyVector[2];
408 
410  ::AllocateSharedPtr(HexBa,HexBb,HexBc,
411  HexGeom);
412  hex->SetElmtId(elmtid++);
413  (*m_exp).push_back(hex);
414  }
415  else
416  {
417  ASSERTL0(false,"dynamic cast to a proper Geometry3D "
418  "failed");
419  }
420 
421  }
422 
423  // Setup Default optimisation information.
424  int nel = GetExpSize();
427 
428  SetCoeffPhys();
430  }
431 
432  /**
433  * Set up the storage for the concatenated list of
434  * coefficients and physical evaluations at the quadrature
435  * points. Each expansion (local element) is processed in turn
436  * to determine the number of coefficients and physical data
437  * points it contributes to the domain. Three arrays,
438  * #m_coeff_offset, #m_phys_offset and #m_offset_elmt_id, are
439  * also initialised and updated to store the data offsets of
440  * each element in the #m_coeffs and #m_phys arrays, and the
441  * element id that each consecutive block is associated
442  * respectively.
443  */
445  {
446  int i;
447 
448  // Set up offset information and array sizes
452 
453  m_ncoeffs = m_npoints = 0;
454 
455  for(i = 0; i < m_exp->size(); ++i)
456  {
458  m_phys_offset [i] = m_npoints;
459  m_offset_elmt_id[i] = i;
460  m_ncoeffs += (*m_exp)[i]->GetNcoeffs();
461  m_npoints += (*m_exp)[i]->GetTotPoints();
462  }
463 
466  }
467 
469  {
470  Array<OneD, int> NumShape(4,0);
471 
472  for(int i = 0; i < GetExpSize(); ++i)
473  {
474  switch ((*m_exp)[i]->DetShapeType())
475  {
476  case LibUtilities::eTetrahedron: NumShape[0]++; break;
477  case LibUtilities::ePyramid: NumShape[1]++; break;
478  case LibUtilities::ePrism: NumShape[2]++; break;
479  case LibUtilities::eHexahedron: NumShape[3]++; break;
480  default:
481  ASSERTL0(false, "Unknown expansion type.");
482  break;
483  }
484  }
485 
486  int three = 3;
488  ::AllocateSharedPtr(m_session,three,NumShape);
489  }
490 
491  void ExpList3D::v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip)
492  {
493  int i,j,k;
494  int nquad0 = (*m_exp)[expansion]->GetNumPoints(0);
495  int nquad1 = (*m_exp)[expansion]->GetNumPoints(1);
496  int nquad2 = (*m_exp)[expansion]->GetNumPoints(2);
497  int ntot = nquad0*nquad1*nquad2;
498  int ntotminus = (nquad0-1)*(nquad1-1)*(nquad2-1);
499 
500  Array<OneD,NekDouble> coords[3];
501  coords[0] = Array<OneD,NekDouble>(ntot);
502  coords[1] = Array<OneD,NekDouble>(ntot);
503  coords[2] = Array<OneD,NekDouble>(ntot);
504  (*m_exp)[expansion]->GetCoords(coords[0],coords[1],coords[2]);
505 
506  outfile << " <Piece NumberOfPoints=\""
507  << ntot << "\" NumberOfCells=\""
508  << ntotminus << "\">" << endl;
509  outfile << " <Points>" << endl;
510  outfile << " <DataArray type=\"Float64\" "
511  << "NumberOfComponents=\"3\" format=\"ascii\">" << endl;
512  outfile << " ";
513  for (i = 0; i < ntot; ++i)
514  {
515  for (j = 0; j < 3; ++j)
516  {
517  outfile << setprecision(8) << scientific
518  << (float)coords[j][i] << " ";
519  }
520  outfile << endl;
521  }
522 
523  outfile << endl;
524  outfile << " </DataArray>" << endl;
525  outfile << " </Points>" << endl;
526  outfile << " <Cells>" << endl;
527  outfile << " <DataArray type=\"Int32\" "
528  << "Name=\"connectivity\" format=\"ascii\">" << endl;
529  for (i = 0; i < nquad0-1; ++i)
530  {
531  for (j = 0; j < nquad1-1; ++j)
532  {
533  for (k = 0; k < nquad2-1; ++k)
534  {
535  outfile << k*nquad0*nquad1 + j*nquad0 + i << " "
536  << k*nquad0*nquad1 + j*nquad0 + i + 1 << " "
537  << k*nquad0*nquad1 + (j+1)*nquad0 + i + 1 << " "
538  << k*nquad0*nquad1 + (j+1)*nquad0 + i << " "
539  << (k+1)*nquad0*nquad1 + j*nquad0 + i << " "
540  << (k+1)*nquad0*nquad1 + j*nquad0 + i + 1 << " "
541  << (k+1)*nquad0*nquad1 + (j+1)*nquad0 + i + 1 << " "
542  << (k+1)*nquad0*nquad1 + (j+1)*nquad0 + i << " " << endl;
543  }
544  }
545  }
546  outfile << endl;
547  outfile << " </DataArray>" << endl;
548  outfile << " <DataArray type=\"Int32\" "
549  << "Name=\"offsets\" format=\"ascii\">" << endl;
550  for (i = 0; i < ntotminus; ++i)
551  {
552  outfile << i*8+8 << " ";
553  }
554  outfile << endl;
555  outfile << " </DataArray>" << endl;
556  outfile << " <DataArray type=\"UInt8\" "
557  << "Name=\"types\" format=\"ascii\">" << endl;
558  for (i = 0; i < ntotminus; ++i)
559  {
560  outfile << "12 ";
561  }
562  outfile << endl;
563  outfile << " </DataArray>" << endl;
564  outfile << " </Cells>" << endl;
565  outfile << " <PointData>" << endl;
566  }
567 
569  {
570  int i, j;
571  for (i = 0; i < m_exp->size(); ++i)
572  {
573  for (j = 0; j < (*m_exp)[i]->GetNfaces(); ++j)
574  {
575  (*m_exp)[i]->ComputeFaceNormal(j);
576  }
577  }
578  }
579 
581  const Array<OneD, NekDouble> &inarray,
582  Array<OneD, NekDouble> &outarray)
583  {
584  int cnt,cnt1;
585 
586  cnt = cnt1 = 0;
587  for(int i = 0; i < GetExpSize(); ++i)
588  {
589  // get new points key
590  int pt0 = (*m_exp)[i]->GetNumPoints(0);
591  int pt1 = (*m_exp)[i]->GetNumPoints(1);
592  int pt2 = (*m_exp)[i]->GetNumPoints(2);
593  int npt0 = (int) pt0*scale;
594  int npt1 = (int) pt1*scale;
595  int npt2 = (int) pt2*scale;
596 
597  LibUtilities::PointsKey newPointsKey0(npt0,(*m_exp)[i]->GetPointsType(0));
598  LibUtilities::PointsKey newPointsKey1(npt1,(*m_exp)[i]->GetPointsType(1));
599  LibUtilities::PointsKey newPointsKey2(npt2,(*m_exp)[i]->GetPointsType(2));
600 
601  // Interpolate points;
602  LibUtilities::Interp3D((*m_exp)[i]->GetBasis(0)->GetPointsKey(),
603  (*m_exp)[i]->GetBasis(1)->GetPointsKey(),
604  (*m_exp)[i]->GetBasis(2)->GetPointsKey(),
605  &inarray[cnt], newPointsKey0,
606  newPointsKey1, newPointsKey2,
607  &outarray[cnt1]);
608 
609  cnt += pt0*pt1*pt2;
610  cnt1 += npt0*npt1*npt2;
611  }
612  }
613 
615  const Array<OneD, NekDouble> &inarray,
616  Array<OneD, NekDouble> &outarray)
617  {
618  int cnt,cnt1;
619 
620  cnt = cnt1 = 0;
621  for(int i = 0; i < GetExpSize(); ++i)
622  {
623  // get new points key
624  int pt0 = (*m_exp)[i]->GetNumPoints(0);
625  int pt1 = (*m_exp)[i]->GetNumPoints(1);
626  int pt2 = (*m_exp)[i]->GetNumPoints(2);
627  int npt0 = (int) pt0*scale;
628  int npt1 = (int) pt1*scale;
629  int npt2 = (int) pt2*scale;
630 
631  LibUtilities::PointsKey newPointsKey0(npt0,(*m_exp)[i]->GetPointsType(0));
632  LibUtilities::PointsKey newPointsKey1(npt1,(*m_exp)[i]->GetPointsType(1));
633  LibUtilities::PointsKey newPointsKey2(npt2,(*m_exp)[i]->GetPointsType(2));
634 
635  // Project points;
637  newPointsKey1,
638  newPointsKey2,
639  &inarray[cnt],
640  (*m_exp)[i]->GetBasis(0)->GetPointsKey(),
641  (*m_exp)[i]->GetBasis(1)->GetPointsKey(),
642  (*m_exp)[i]->GetBasis(2)->GetPointsKey(),
643  &outarray[cnt1]);
644 
645  cnt += npt0*npt1*npt2;
646  cnt1 += pt0*pt1*pt2;
647  }
648 
649  }
650  } //end of namespace
651 } //end of namespace
652 
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:568
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:1060
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:2054
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:614
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1047
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
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1050
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
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:1058
int GetNumPoints() const
Return points order at which basis is defined.
Definition: Basis.h:128
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:468
virtual ~ExpList3D()
Destructor.
Definition: ExpList3D.cpp:82
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:3129
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:253
void SetCoeffPhys(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Definition: ExpList3D.cpp:444
#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:491
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:580