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