Nektar++
ExpList3DHomogeneous1D.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File ExpList3DHomogeneous1D.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: An ExpList which is homogeneous in 1 direction and so
32 // uses much of the functionality from a ExpList2D and its daughters
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
37 #include <MultiRegions/ExpList2D.h>
38 
39 using namespace std;
40 
41 namespace Nektar
42 {
43  namespace MultiRegions
44  {
45  // Forward declaration for typedefs
46  ExpList3DHomogeneous1D::ExpList3DHomogeneous1D():
48  {
50  }
51 
54  const LibUtilities::BasisKey &HomoBasis,
55  const NekDouble lhom,
56  const bool useFFT,
57  const bool dealiasing):
58  ExpListHomogeneous1D(pSession,HomoBasis,lhom,useFFT,dealiasing)
59  {
61  }
62 
63  // Constructor for ExpList3DHomogeneous1D to act as a Explist2D field
66  const LibUtilities::BasisKey &HomoBasis,
67  const NekDouble lhom,
68  const bool useFFT,
69  const bool dealiasing,
71  const std::string &var,
72  const Collections::ImplementationType ImpType):
73  ExpListHomogeneous1D(pSession,HomoBasis,lhom,useFFT,dealiasing)
74  {
76 
77  GenExpList3DHomogeneous1D(graph2D->GetExpansions(var),ImpType);
78  }
79 
80  // Constructor for ExpList3DHomogeneous1D to act as a Explist2D field
83  const LibUtilities::BasisKey &HomoBasis,
84  const NekDouble lhom,
85  const bool useFFT,
86  const bool dealiasing,
87  const SpatialDomains::ExpansionMap &expansions,
88  const Collections::ImplementationType ImpType):
89  ExpListHomogeneous1D(pSession,HomoBasis,lhom,useFFT,dealiasing)
90  {
92 
93  GenExpList3DHomogeneous1D(expansions,ImpType);
94  }
95 
97  const SpatialDomains::ExpansionMap &expansions,
98  const Collections::ImplementationType ImpType)
99  {
100  int n,j,nel;
101  bool False = false;
102  ExpList2DSharedPtr plane_zero;
103 
104  // note that nzplanes can be larger than nzmodes
105  m_planes[0] = plane_zero = MemoryManager<ExpList2D>::AllocateSharedPtr(m_session, expansions, False,ImpType);
106 
108  nel = m_planes[0]->GetExpSize();
109 
110  for(j = 0; j < nel; ++j)
111  {
112  (*m_exp).push_back(m_planes[0]->GetExp(j));
113  }
114 
115  for(n = 1; n < m_planes.num_elements(); ++n)
116  {
118  for(j = 0; j < nel; ++j)
119  {
120  (*m_exp).push_back((*m_exp)[j]);
121  }
122  }
123 
124  // Setup Default optimisation information.
125  nel = GetExpSize();
128 
129  SetCoeffPhys();
130  }
131 
132 
133  /**
134  * @param In ExpList3DHomogeneous1D object to copy.
135  */
136  ExpList3DHomogeneous1D::ExpList3DHomogeneous1D(const ExpList3DHomogeneous1D &In, bool DeclarePlanesSetCoeffPhys):
138  {
140 
141  if(DeclarePlanesSetCoeffPhys)
142  {
143  bool False = false;
144  ExpList2DSharedPtr zero_plane = std::dynamic_pointer_cast<ExpList2D> (In.m_planes[0]);
145 
146  for(int n = 0; n < m_planes.num_elements(); ++n)
147  {
149  }
150 
151  SetCoeffPhys();
152  }
153  }
154 
155  /**
156  * @param In ExpList3DHomogeneous1D object to copy.
157  * @param eIDs Id of elements that should be copied.
158  */
160  const ExpList3DHomogeneous1D &In,
161  const std::vector<unsigned int> &eIDs,
162  bool DeclarePlanesSetCoeffPhys,
163  const Collections::ImplementationType ImpType):
164  ExpListHomogeneous1D(In, eIDs)
165  {
167 
168  if(DeclarePlanesSetCoeffPhys)
169  {
170  bool False = false;
171  std::vector<unsigned int> eIDsPlane;
172  int nel = eIDs.size()/m_planes.num_elements();
173 
174  for (int i = 0; i < nel; ++i)
175  {
176  eIDsPlane.push_back(eIDs[i]);
177  }
178 
179  ExpList2DSharedPtr zero_plane_old =
180  std::dynamic_pointer_cast<ExpList2D> (In.m_planes[0]);
181 
182  ExpList2DSharedPtr zero_plane =
183  MemoryManager<ExpList2D>::AllocateSharedPtr(*(zero_plane_old), eIDsPlane, ImpType);
184 
185  for(int n = 0; n < m_planes.num_elements(); ++n)
186  {
188  }
189 
190  SetCoeffPhys();
191  }
192  }
193 
194  /**
195  * Destructor
196  */
198  {
199  }
200 
202  {
203  int i,n,cnt;
204  int ncoeffs_per_plane = m_planes[0]->GetNcoeffs();
205  int npoints_per_plane = m_planes[0]->GetTotPoints();
206 
207  int nzplanes = m_planes.num_elements();
208 
209  // Set total coefficients and points
210  m_ncoeffs = ncoeffs_per_plane*nzplanes;
211  m_npoints = npoints_per_plane*nzplanes;
212 
215 
216  int nel = m_planes[0]->GetExpSize();
217  m_coeff_offset = Array<OneD,int>(nel*nzplanes);
218  m_phys_offset = Array<OneD,int>(nel*nzplanes);
219  Array<OneD, NekDouble> tmparray;
220 
221  for(cnt = n = 0; n < nzplanes; ++n)
222  {
223  m_planes[n]->SetCoeffsArray(tmparray= m_coeffs + ncoeffs_per_plane*n);
224  m_planes[n]->SetPhysArray(tmparray = m_phys + npoints_per_plane*n);
225 
226  for(i = 0; i < nel; ++i)
227  {
228  m_coeff_offset[cnt] = m_planes[n]->GetCoeff_Offset(i) + n*ncoeffs_per_plane;
229  m_phys_offset[cnt++] = m_planes[n]->GetPhys_Offset(i) + n*npoints_per_plane;
230  }
231  }
232  }
233 
238  {
239  int n;
240  Array<OneD, NekDouble> tmp_xc;
241  int nzplanes = m_planes.num_elements();
242  int npoints = GetTotPoints(eid);
243 
244  (*m_exp)[eid]->GetCoords(xc0,xc1);
245 
246  // Fill z-direction
248  Array<OneD, NekDouble> local_pts(m_planes.num_elements());
249 
250  for(n = 0; n < m_planes.num_elements(); n++)
251  {
252  local_pts[n] = pts[m_transposition->GetPlaneID(n)];
253  }
254 
255  Array<OneD, NekDouble> z(nzplanes);
256 
257  Vmath::Smul(nzplanes,m_lhom/2.0,local_pts,1,z,1);
258  Vmath::Sadd(nzplanes,m_lhom/2.0,z,1,z,1);
259 
260  for(n = 0; n < nzplanes; ++n)
261  {
262  Vmath::Fill(npoints,z[n],tmp_xc = xc2 + npoints*n,1);
263  if(n)
264  {
265  Vmath::Vcopy(npoints,xc0,1,tmp_xc = xc0+npoints*n,1);
266  Vmath::Vcopy(npoints,xc1,1,tmp_xc = xc1+npoints*n,1);
267  }
268  }
269  }
270 
271  /**
272  * The operation calls the 2D plane coordinates through the
273  * function ExpList#GetCoords and then evaluated the third
274  * coordinate using the member \a m_lhom
275  *
276  * @param coord_0 After calculation, the \f$x_1\f$ coordinate
277  * will be stored in this array.
278  *
279  * @param coord_1 After calculation, the \f$x_2\f$ coordinate
280  * will be stored in this array.
281  *
282  * @param coord_2 After calculation, the \f$x_3\f$ coordinate
283  * will be stored in this array. This
284  * coordinate is evaluated using the
285  * predefined value \a m_lhom
286  */
290  {
291  int n;
292  Array<OneD, NekDouble> tmp_xc;
293  int nzplanes = m_planes.num_elements();
294  int npoints = m_planes[0]->GetTotPoints();
295 
296  m_planes[0]->GetCoords(xc0,xc1);
297 
298  // Fill z-direction
300 
301  Array<OneD, NekDouble> local_pts(m_planes.num_elements());
302 
303  for(n = 0; n < m_planes.num_elements(); n++)
304  {
305  local_pts[n] = pts[m_transposition->GetPlaneID(n)];
306  }
307 
308  Array<OneD, NekDouble> z(nzplanes);
309 
310  Vmath::Smul(nzplanes,m_lhom/2.0,local_pts,1,z,1);
311  Vmath::Sadd(nzplanes,m_lhom/2.0,z,1,z,1);
312 
313  for(n = 0; n < nzplanes; ++n)
314  {
315  Vmath::Fill(npoints,z[n],tmp_xc = xc2 + npoints*n,1);
316  if(n)
317  {
318  Vmath::Vcopy(npoints,xc0,1,tmp_xc = xc0+npoints*n,1);
319  Vmath::Vcopy(npoints,xc1,1,tmp_xc = xc1+npoints*n,1);
320  }
321  }
322  }
323 
324  void ExpList3DHomogeneous1D::v_WriteTecplotConnectivity(std::ostream &outfile, int expansion)
325  {
326  ASSERTL0(expansion == -1, "Multi-zone output not supported for homogeneous expansions.");
327 
328  const int nPtsPlane = m_planes[0]->GetNpoints();
329  const int nElmt = m_planes[0]->GetExpSize();
330  const int nPlanes = m_planes.num_elements();
331 
332  int cnt = 0;
333  int cnt2 = 0;
334  for (int i = 0; i < nElmt; ++i)
335  {
336  const int np0 = (*m_exp)[i]->GetNumPoints(0);
337  const int np1 = (*m_exp)[i]->GetNumPoints(1);
338 
339  for (int n = 1; n < nPlanes; ++n)
340  {
341  const int o1 = (n-1) * nPtsPlane;
342  const int o2 = n * nPtsPlane;
343  for (int j = 1; j < np1; ++j)
344  {
345  for(int k = 1; k < np0; ++k)
346  {
347  outfile << cnt + (j-1)*np0 + (k-1) + o1 + 1 << " ";
348  outfile << cnt + (j-1)*np0 + (k-1) + o2 + 1 << " ";
349  outfile << cnt + (j-1)*np0 + k + o2 + 1 << " ";
350  outfile << cnt + (j-1)*np0 + k + o1 + 1 << " ";
351  outfile << cnt + j *np0 + (k-1) + o1 + 1 << " ";
352  outfile << cnt + j *np0 + (k-1) + o2 + 1 << " ";
353  outfile << cnt + j *np0 + k + o2 + 1 << " ";
354  outfile << cnt + j *np0 + k + o1 + 1 << endl;
355  cnt2++;
356  }
357  }
358  }
359 
360  cnt += np0*np1;
361  }
362  }
363 
364 
366  std::ostream &outfile,
367  int expansion,
368  int istrip)
369  {
370  // If there is only one plane (e.g. HalfMode), we write a 2D plane.
371  if (m_planes.num_elements() == 1)
372  {
373  m_planes[0]->WriteVtkPieceHeader(outfile, expansion);
374  return;
375  }
376 
377  // If we are using Fourier points, output extra plane to fill domain
378  int outputExtraPlane = 0;
379  if ( m_homogeneousBasis->GetBasisType() == LibUtilities::eFourier
380  && m_homogeneousBasis->GetPointsType() ==
382  {
383  outputExtraPlane = 1;
384  }
385 
386  int i,j,k;
387  int nq0 = (*m_exp)[expansion]->GetNumPoints(0);
388  int nq1 = (*m_exp)[expansion]->GetNumPoints(1);
389  int nq2 = m_planes.num_elements() + outputExtraPlane;
390  int ntot = nq0*nq1*nq2;
391  int ntotminus = (nq0-1)*(nq1-1)*(nq2-1);
392 
393  Array<OneD,NekDouble> coords[3];
394  coords[0] = Array<OneD,NekDouble>(ntot);
395  coords[1] = Array<OneD,NekDouble>(ntot);
396  coords[2] = Array<OneD,NekDouble>(ntot);
397  GetCoords(expansion,coords[0],coords[1],coords[2]);
398 
399  if (outputExtraPlane)
400  {
401  // Copy coords[0] and coords[1] to extra plane
403  Vmath::Vcopy (nq0*nq1, coords[0], 1,
404  tmp = coords[0] + (nq2-1)*nq0*nq1, 1);
405  Vmath::Vcopy (nq0*nq1, coords[1], 1,
406  tmp = coords[1] + (nq2-1)*nq0*nq1, 1);
407  // Fill coords[2] for extra plane
408  NekDouble z = coords[2][nq0*nq1*m_planes.num_elements()-1] +
409  (coords[2][nq0*nq1] - coords[2][0]);
410  Vmath::Fill(nq0*nq1, z, tmp = coords[2] + (nq2-1)*nq0*nq1, 1);
411  }
412 
413  NekDouble DistStrip;
414  m_session->LoadParameter("DistStrip", DistStrip, 0);
415  // Reset the z-coords for homostrips
416  for(int i = 0; i < ntot; i++)
417  {
418  coords[2][i] += istrip*DistStrip;
419  }
420 
421  outfile << " <Piece NumberOfPoints=\""
422  << ntot << "\" NumberOfCells=\""
423  << ntotminus << "\">" << endl;
424  outfile << " <Points>" << endl;
425  outfile << " <DataArray type=\"Float64\" "
426  << "NumberOfComponents=\"3\" format=\"ascii\">" << endl;
427  outfile << " ";
428  for (i = 0; i < ntot; ++i)
429  {
430  for (j = 0; j < 3; ++j)
431  {
432  outfile << coords[j][i] << " ";
433  }
434  outfile << endl;
435  }
436  outfile << endl;
437  outfile << " </DataArray>" << endl;
438  outfile << " </Points>" << endl;
439  outfile << " <Cells>" << endl;
440  outfile << " <DataArray type=\"Int32\" "
441  << "Name=\"connectivity\" format=\"ascii\">" << endl;
442  for (i = 0; i < nq0-1; ++i)
443  {
444  for (j = 0; j < nq1-1; ++j)
445  {
446  for (k = 0; k < nq2-1; ++k)
447  {
448  outfile << k*nq0*nq1 + j*nq0 + i << " "
449  << k*nq0*nq1 + j*nq0 + i + 1 << " "
450  << k*nq0*nq1 + (j+1)*nq0 + i + 1 << " "
451  << k*nq0*nq1 + (j+1)*nq0 + i << " "
452  << (k+1)*nq0*nq1 + j*nq0 + i << " "
453  << (k+1)*nq0*nq1 + j*nq0 + i + 1 << " "
454  << (k+1)*nq0*nq1 + (j+1)*nq0 + i + 1 << " "
455  << (k+1)*nq0*nq1 + (j+1)*nq0 + i << endl;
456  }
457  }
458  }
459  outfile << endl;
460  outfile << " </DataArray>" << endl;
461  outfile << " <DataArray type=\"Int32\" "
462  << "Name=\"offsets\" format=\"ascii\">" << endl;
463  for (i = 0; i < ntotminus; ++i)
464  {
465  outfile << i*8+8 << " ";
466  }
467  outfile << endl;
468  outfile << " </DataArray>" << endl;
469  outfile << " <DataArray type=\"UInt8\" "
470  << "Name=\"types\" format=\"ascii\">" << endl;
471  for (i = 0; i < ntotminus; ++i)
472  {
473  outfile << "12 ";
474  }
475  outfile << endl;
476  outfile << " </DataArray>" << endl;
477  outfile << " </Cells>" << endl;
478  outfile << " <PointData>" << endl;
479  }
480 
482  const Array<OneD, const NekDouble> &inarray,
483  const Array<OneD, const NekDouble> &soln)
484  {
485  int cnt = 0;
486  NekDouble errL2, err = 0.0;
488  Array<OneD, NekDouble> local_w(m_planes.num_elements());
489 
490  for(int n = 0; n < m_planes.num_elements(); n++)
491  {
492  local_w[n] = w[m_transposition->GetPlaneID(n)];
493  }
494 
495  if (soln == NullNekDouble1DArray)
496  {
497  for(int n = 0; n < m_planes.num_elements(); ++n)
498  {
499  errL2 = m_planes[n]->L2(inarray + cnt);
500  cnt += m_planes[n]->GetTotPoints();
501  err += errL2*errL2*local_w[n]*m_lhom*0.5;
502  }
503  }
504  else
505  {
506  for(int n = 0; n < m_planes.num_elements(); ++n)
507  {
508  errL2 = m_planes[n]->L2(inarray + cnt, soln + cnt);
509  cnt += m_planes[n]->GetTotPoints();
510  err += errL2*errL2*local_w[n]*m_lhom*0.5;
511  }
512  }
513  m_comm->GetColumnComm()->AllReduce(err, LibUtilities::ReduceSum);
514 
515  return sqrt(err);
516  }
517 
519  {
520  Array<OneD, NekDouble> energy(m_planes.num_elements()/2);
521  NekDouble area = 0.0;
522 
523  // Calculate total area of elements.
524  for (int n = 0; n < m_planes[0]->GetExpSize(); ++n)
525  {
526  Array<OneD, NekDouble> inarray(m_planes[0]->GetExp(n)->GetTotPoints(), 1.0);
527  area += m_planes[0]->GetExp(n)->Integral(inarray);
528  }
529 
530  m_comm->GetRowComm()->AllReduce(area, LibUtilities::ReduceSum);
531 
532  // Calculate L2 norm of real/imaginary planes.
533  for (int n = 0; n < m_planes.num_elements(); n += 2)
534  {
535  NekDouble err;
536 
537  energy[n/2] = 0;
538 
539  for(int i = 0; i < m_planes[n]->GetExpSize(); ++i)
540  {
541  LocalRegions::ExpansionSharedPtr exp = m_planes[n]->GetExp( i );
542  Array<OneD, NekDouble> phys(exp->GetTotPoints());
543  exp->BwdTrans(m_planes[n]->GetCoeffs()+m_planes[n]->GetCoeff_Offset(i),
544  phys);
545  err = exp->L2(phys);
546  energy[n/2] += err*err;
547 
548  exp = m_planes[n+1]->GetExp(i);
549  exp->BwdTrans(m_planes[n+1]->GetCoeffs()+m_planes[n+1]->GetCoeff_Offset(i),
550  phys);
551  err = exp->L2(phys);
552  energy[n/2] += err*err;
553  }
554 
555  m_comm->GetRowComm()->AllReduce(energy[n/2], LibUtilities::ReduceSum);
556  energy[n/2] /= 2.0*area;
557  }
558 
559  return energy;
560  }
561  } //end of namespace
562 } //end of namespace
563 
564 
Abstraction of a two-dimensional multi-elemental expansion which is merely a collection of local expa...
virtual void v_WriteTecplotConnectivity(std::ostream &outfile, int expansion)
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:163
virtual Array< OneD, const NekDouble > v_HomogeneousEnergy(void)
static Array< OneD, NekDouble > NullNekDouble1DArray
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1090
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:1106
LibUtilities::TranspositionSharedPtr m_transposition
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45
const Array< OneD, const NekDouble > & GetCoeffs() const
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
Definition: ExpList.h:2062
STL namespace.
NekDouble m_lhom
Width of homogeneous direction.
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:1069
Fourier Expansion .
Definition: BasisType.h:53
std::shared_ptr< ExpList2D > ExpList2DSharedPtr
Shared pointer to an ExpList2D object.
Definition: ExpList2D.h:48
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
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:2191
virtual void v_GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2)
void SetCoeffPhys(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1101
LibUtilities::BasisSharedPtr m_homogeneousBasis
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:65
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:216
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1104
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:1030
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
Array< OneD, ExpListSharedPtr > m_planes
double NekDouble
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.cpp:318
virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip)
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:65
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:1020
void GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1=NullNekDouble1DArray, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
This function calculates the coordinates of all the elemental quadrature points . ...
Abstraction of a two-dimensional multi-elemental expansion which is merely a collection of local expa...
Definition: ExpList2D.h:57
void GenExpList3DHomogeneous1D(const SpatialDomains::ExpansionMap &expansions, const Collections::ImplementationType ImpType)
virtual NekDouble v_L2(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:279
int GetCoeff_Offset(int n) const
Get the start offset position for a global list of m_coeffs correspoinding to element n...
Definition: ExpList.h:2200
int GetTotPoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1608
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
Abstraction of a two-dimensional multi-elemental expansion which is merely a collection of local expa...
Describes the specification for a Basis.
Definition: Basis.h:49
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::map< int, ExpansionShPtr > ExpansionMap
Definition: MeshGraph.h:152