Nektar++
StdExpansion.h
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File Stdexpansion.h
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: Class definition StdExpansion which is the base class
32 // to all expansion shapes
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #ifndef NEKTAR_LIB_STDREGIONS_STANDARDEXPANSION_H
37 #define NEKTAR_LIB_STDREGIONS_STANDARDEXPANSION_H
38 
39 #include <fstream>
40 #include <vector>
41 #include <memory>
42 
43 #include <boost/core/ignore_unused.hpp>
44 
48 #include <StdRegions/IndexMapKey.h>
50 namespace Nektar { namespace LocalRegions { class MatrixKey; class Expansion; } }
51 
52 namespace Nektar
53 {
54  namespace StdRegions
55  {
56 
57  class StdExpansion1D;
58  class StdExpansion2D;
59 
61 
62  /** \brief The base class for all shapes
63  *
64  * This is the lowest level basic class for all shapes and so
65  * contains the definition of common data and common routine to all
66  * elements
67  */
68  class StdExpansion : public std::enable_shared_from_this<StdExpansion>
69  {
70  public:
71 
72  /** \brief Default Constructor */
74 
75  /** \brief Constructor */
76  STD_REGIONS_EXPORT StdExpansion(const int numcoeffs, const int numbases,
80 
81 
82  /** \brief Copy Constructor */
84 
85  /** \brief Destructor */
87 
88 
89  // Standard Expansion Routines Applicable Regardless of Region
90 
91  /** \brief This function returns the number of 1D bases used in
92  * the expansion
93  *
94  * \return returns the number of 1D bases used in the expansion,
95  * which is equal to number dimension of the expansion
96  */
97  inline int GetNumBases() const
98  {
99  return m_base.num_elements();
100  }
101 
102  /** \brief This function gets the shared point to basis
103  *
104  * \return returns the shared pointer to the bases
105  */
107  {
108  return(m_base);
109  }
110 
111  /** \brief This function gets the shared point to basis in
112  * the \a dir direction
113  *
114  * \return returns the shared pointer to the basis in
115  * directin \a dir
116  */
117  inline const LibUtilities::BasisSharedPtr& GetBasis(int dir) const
118  {
119  ASSERTL1(dir < m_base.num_elements(),
120  "dir is larger than number of bases");
121  return(m_base[dir]);
122  }
123 
124  /** \brief This function returns the total number of coefficients
125  * used in the expansion
126  *
127  * \return returns the total number of coefficients (which is
128  * equivalent to the total number of modes) used in the expansion
129  */
130  inline int GetNcoeffs(void) const
131  {
132  return(m_ncoeffs);
133  }
134 
135  /** \brief This function returns the total number of quadrature
136  * points used in the element
137  *
138  * \return returns the total number of quadrature points
139  */
140  inline int GetTotPoints() const
141  {
142  int i;
143  int nqtot = 1;
144 
145  for(i=0; i < m_base.num_elements(); ++i)
146  {
147  nqtot *= m_base[i]->GetNumPoints();
148  }
149 
150  return nqtot;
151  }
152 
153 
154  /** \brief This function returns the type of basis used in the \a dir
155  * direction
156  *
157  * The different types of bases implemented in the code are defined
158  * in the LibUtilities::BasisType enumeration list. As a result, the
159  * function will return one of the types of this enumeration list.
160  *
161  * \param dir the direction
162  * \return returns the type of basis used in the \a dir direction
163  */
164  inline LibUtilities::BasisType GetBasisType(const int dir) const
165  {
166  ASSERTL1(dir < m_base.num_elements(), "dir is larger than m_numbases");
167  return(m_base[dir]->GetBasisType());
168  }
169 
170  /** \brief This function returns the number of expansion modes
171  * in the \a dir direction
172  *
173  * \param dir the direction
174  * \return returns the number of expansion modes in the \a dir
175  * direction
176  */
177  inline int GetBasisNumModes(const int dir) const
178  {
179  ASSERTL1(dir < m_base.num_elements(),"dir is larger than m_numbases");
180  return(m_base[dir]->GetNumModes());
181  }
182 
183 
184  /** \brief This function returns the maximum number of
185  * expansion modes over all local directions
186  *
187  * \return returns the maximum number of expansion modes
188  * over all local directions
189  */
190  inline int EvalBasisNumModesMax(void) const
191  {
192  int i;
193  int returnval = 0;
194 
195  for(i = 0; i < m_base.num_elements(); ++i)
196  {
197  returnval = std::max(returnval, m_base[i]->GetNumModes());
198  }
199 
200  return returnval;
201  }
202 
203  /** \brief This function returns the type of quadrature points used
204  * in the \a dir direction
205  *
206  * The different types of quadrature points implemented in the code
207  * are defined in the LibUtilities::PointsType enumeration list.
208  * As a result, the function will return one of the types of this
209  * enumeration list.
210  *
211  * \param dir the direction
212  * \return returns the type of quadrature points used in the \a dir
213  * direction
214  */
215  inline LibUtilities::PointsType GetPointsType(const int dir) const
216  {
217  ASSERTL1(dir < m_base.num_elements(), "dir is larger than m_numbases");
218  return(m_base[dir]->GetPointsType());
219  }
220 
221  /** \brief This function returns the number of quadrature points
222  * in the \a dir direction
223  *
224  * \param dir the direction
225  * \return returns the number of quadrature points in the \a dir
226  * direction
227  */
228  inline int GetNumPoints(const int dir) const
229  {
230  ASSERTL1(dir < m_base.num_elements() || dir == 0,
231  "dir is larger than m_numbases");
232  return(m_base.num_elements() > 0 ? m_base[dir]->GetNumPoints() : 1);
233  }
234 
235  /** \brief This function returns a pointer to the array containing
236  * the quadrature points in \a dir direction
237  *
238  * \param dir the direction
239  * \return returns a pointer to the array containing
240  * the quadrature points in \a dir direction
241  */
242  inline const Array<OneD, const NekDouble>& GetPoints(const int dir) const
243  {
244  return m_base[dir]->GetZ();
245  }
246 
247 
248  // Wrappers around virtual Functions
249 
250  /** \brief This function returns the number of vertices of the
251  * expansion domain
252  *
253  * This function is a wrapper around the virtual function
254  * \a v_GetNverts()
255  *
256  * \return returns the number of vertices of the expansion domain
257  */
258  int GetNverts() const
259  {
260  return v_GetNverts();
261  }
262 
263  /** \brief This function returns the number of edges of the
264  * expansion domain
265  *
266  * This function is a wrapper around the virtual function
267  * \a v_GetNedges()
268  *
269  * \return returns the number of edges of the expansion domain
270  */
271  int GetNedges() const
272  {
273  return v_GetNedges();
274  }
275 
276  /** \brief This function returns the number of expansion coefficients
277  * belonging to the \a i-th edge
278  *
279  * This function is a wrapper around the virtual function
280  * \a v_GetEdgeNcoeffs()
281  *
282  * \param i specifies which edge
283  * \return returns the number of expansion coefficients belonging to
284  * the \a i-th edge
285  */
286  int GetEdgeNcoeffs(const int i) const
287  {
288  return v_GetEdgeNcoeffs(i);
289  }
290 
291 
293  {
294  return v_GetTotalEdgeIntNcoeffs();
295  }
296 
297  /** \brief This function returns the number of quadrature points
298  * belonging to the \a i-th edge
299  *
300  * This function is a wrapper around the virtual function
301  * \a v_GetEdgeNumPoints()
302  *
303  * \param i specifies which edge
304  * \return returns the number of expansion coefficients belonging to
305  * the \a i-th edge
306  */
307  int GetEdgeNumPoints(const int i) const
308  {
309  return v_GetEdgeNumPoints(i);
310  }
311 
312 
313  int DetCartesianDirOfEdge(const int edge)
314  {
315  return v_DetCartesianDirOfEdge(edge);
316  }
317 
318  const LibUtilities::BasisKey DetEdgeBasisKey(const int i) const
319  {
320  return v_DetEdgeBasisKey(i);
321  }
322 
323  const LibUtilities::BasisKey DetFaceBasisKey(const int i, const int k) const
324  {
325  return v_DetFaceBasisKey(i, k);
326  }
327  /**
328  * \brief This function returns the number of quadrature points
329  * belonging to the \a i-th face.
330  *
331  * This function is a wrapper around the virtual function \a
332  * v_GetFaceNumPoints()
333  *
334  * \param i specifies which face
335  * \return returns the number of expansion coefficients belonging to
336  * the \a i-th face
337  */
338  int GetFaceNumPoints(const int i) const
339  {
340  return v_GetFaceNumPoints(i);
341  }
342 
343  /** \brief This function returns the number of expansion coefficients
344  * belonging to the \a i-th face
345  *
346  * This function is a wrapper around the virtual function
347  * \a v_GetFaceNcoeffs()
348  *
349  * \param i specifies which face
350  * \return returns the number of expansion coefficients belonging to
351  * the \a i-th face
352  */
353  int GetFaceNcoeffs(const int i) const
354  {
355  return v_GetFaceNcoeffs(i);
356  }
357 
358  int GetFaceIntNcoeffs(const int i) const
359  {
360  return v_GetFaceIntNcoeffs(i);
361  }
362 
364  {
365  return v_GetTotalFaceIntNcoeffs();
366  }
367 
368  /** \brief This function returns the number of expansion coefficients
369  * belonging to the \a i-th edge/face
370  *
371  * This function is a wrapper around the virtual function
372  * \a v_GetTraceNcoeffs()
373  *
374  * \param i specifies which edge/face
375  * \return returns the number of expansion coefficients belonging to
376  * the \a i-th edge/face
377  */
378  int GetTraceNcoeffs(const int i) const
379  {
380  return v_GetTraceNcoeffs(i);
381  }
382 
383 
384  LibUtilities::PointsKey GetFacePointsKey(const int i, const int j) const
385  {
386  return v_GetFacePointsKey(i, j);
387  }
388 
389  int NumBndryCoeffs(void) const
390  {
391  return v_NumBndryCoeffs();
392  }
393 
394  int NumDGBndryCoeffs(void) const
395  {
396  return v_NumDGBndryCoeffs();
397  }
398 
399  /** \brief This function returns the type of expansion basis on the
400  * \a i-th edge
401  *
402  * This function is a wrapper around the virtual function
403  * \a v_GetEdgeBasisType()
404  *
405  * The different types of bases implemented in the code are defined
406  * in the LibUtilities::BasisType enumeration list. As a result, the
407  * function will return one of the types of this enumeration list.
408  *
409  * \param i specifies which edge
410  * \return returns the expansion basis on the \a i-th edge
411  */
413  {
414  return v_GetEdgeBasisType(i);
415  }
416 
417  /** \brief This function returns the type of expansion
418  * Nodal point type if defined
419  *
420  * This function is a wrapper around the virtual function
421  * \a v_GetNodalPointsKey()
422  *
423  */
425  {
426  return v_GetNodalPointsKey();
427  };
428 
429  /** \brief This function returns the number of faces of the
430  * expansion domain
431  *
432  * This function is a wrapper around the virtual function
433  * \a v_GetNFaces()
434  *
435  * \return returns the number of faces of the expansion domain
436  */
437  int GetNfaces() const
438  {
439  return v_GetNfaces();
440  }
441 
442  /**
443  * @brief Returns the number of trace elements connected to this
444  * element.
445  *
446  * For example, a quadrilateral has four edges, so this function
447  * would return 4.
448  */
449  int GetNtrace() const
450  {
451  const size_t nBase = m_base.num_elements();
452  return
453  nBase == 1 ? 2 :
454  nBase == 2 ? GetNedges() :
455  nBase == 3 ? GetNfaces() : 0;
456  }
457 
458  /** \brief This function returns the shape of the expansion domain
459  *
460  * This function is a wrapper around the virtual function
461  * \a v_DetShapeType()
462  *
463  * The different shape types implemented in the code are defined
464  * in the ::ShapeType enumeration list. As a result, the
465  * function will return one of the types of this enumeration list.
466  *
467  * \return returns the shape of the expansion domain
468  */
470  {
471  return v_DetShapeType();
472  }
473 
474  std::shared_ptr<StdExpansion> GetStdExp(void) const
475  {
476  return v_GetStdExp();
477  }
478 
479 
480  std::shared_ptr<StdExpansion> GetLinStdExp(void) const
481  {
482  return v_GetLinStdExp();
483  }
484 
485 
486  int GetShapeDimension() const
487  {
488  return v_GetShapeDimension();
489  }
490 
492  {
494  }
495 
496 
498  {
499  return v_IsNodalNonTensorialExp();
500  }
501 
502  /** \brief This function performs the Backward transformation from
503  * coefficient space to physical space
504  *
505  * This function is a wrapper around the virtual function
506  * \a v_BwdTrans()
507  *
508  * Based on the expansion coefficients, this function evaluates the
509  * expansion at the quadrature points. This is equivalent to the
510  * operation \f[ u(\xi_{1i}) =
511  * \sum_{p=0}^{P-1} \hat{u}_p \phi_p(\xi_{1i}) \f] which can be
512  * evaluated as \f$ {\bf u} = {\bf B}^T {\bf \hat{u}} \f$ with
513  * \f${\bf B}[i][j] = \phi_i(\xi_{j})\f$
514  *
515  * This function requires that the coefficient array
516  * \f$\mathbf{\hat{u}}\f$ provided as \a inarray.
517  *
518  * The resulting array
519  * \f$\mathbf{u}[m]=u(\mathbf{\xi}_m)\f$ containing the
520  * expansion evaluated at the quadrature points, is stored
521  * in the \a outarray.
522  *
523  * \param inarray contains the values of the expansion
524  * coefficients (input of the function)
525  *
526  * \param outarray contains the values of the expansion evaluated
527  * at the quadrature points (output of the function)
528  */
529 
531  Array<OneD, NekDouble> &outarray)
532  {
533  v_BwdTrans (inarray, outarray);
534  }
535 
536  /**
537  * @brief This function performs the Forward transformation from
538  * physical space to coefficient space.
539  */
540  inline void FwdTrans (const Array<OneD, const NekDouble>& inarray,
541  Array<OneD, NekDouble> &outarray);
542 
544  Array<OneD, NekDouble> &outarray)
545  {
546  v_FwdTrans_BndConstrained(inarray,outarray);
547  }
548 
549  /** \brief This function integrates the specified function over the
550  * domain
551  *
552  * This function is a wrapper around the virtual function
553  * \a v_Integral()
554  *
555  * Based on the values of the function evaluated at the quadrature
556  * points (which are stored in \a inarray), this function calculates
557  * the integral of this function over the domain. This is
558  * equivalent to the numerical evaluation of the operation
559  * \f[ I=\int u(\mathbf{\xi})d \mathbf{\xi}\f]
560  *
561  * \param inarray values of the function to be integrated evaluated
562  * at the quadrature points (i.e.
563  * \a inarray[m]=\f$u(\mathbf{\xi}_m)\f$)
564  * \return returns the value of the calculated integral
565  *
566  * Inputs:\n
567 
568  - \a inarray: definition of function to be returned at quadrature point
569  of expansion.
570 
571  Outputs:\n
572 
573  - returns \f$\int^1_{-1}\int^1_{-1} u(\xi_1, \xi_2) J[i,j] d
574  \xi_1 d \xi_2 \f$ where \f$inarray[i,j] =
575  u(\xi_{1i},\xi_{2j}) \f$ and \f$ J[i,j] \f$ is the
576  Jacobian evaluated at the quadrature point.
577  *
578  */
580  {
581  return v_Integral(inarray);
582  }
583 
584  /** \brief This function fills the array \a outarray with the
585  * \a mode-th mode of the expansion
586  *
587  * This function is a wrapper around the virtual function
588  * \a v_FillMode()
589  *
590  * The requested mode is evaluated at the quadrature points
591  *
592  * \param mode the mode that should be filled
593  * \param outarray contains the values of the \a mode-th mode of the
594  * expansion evaluated at the quadrature points (output of the
595  * function)
596  */
597  void FillMode(const int mode, Array<OneD, NekDouble> &outarray)
598  {
599  v_FillMode(mode, outarray);
600  }
601 
602  /** \brief this function calculates the inner product of a given
603  * function \a f with the different modes of the expansion
604  *
605  * This function is a wrapper around the virtual function
606  * \a v_IProductWRTBase()
607  *
608  * This is equivalent to the numerical evaluation of
609  * \f[ I[p] = \int \phi_p(\mathbf{x}) f(\mathbf{x}) d\mathbf{x}\f]
610  * \f$ \begin{array}{rcl} I_{pq} = (\phi_q \phi_q, u) & = &
611  \sum_{i=0}^{nq_0} \sum_{j=0}^{nq_1} \phi_p(\xi_{0,i})
612  \phi_q(\xi_{1,j}) w^0_i w^1_j u(\xi_{0,i} \xi_{1,j})
613  J_{i,j}\\ & = & \sum_{i=0}^{nq_0} \phi_p(\xi_{0,i})
614  \sum_{j=0}^{nq_1} \phi_q(\xi_{1,j}) \tilde{u}_{i,j}
615  J_{i,j} \end{array} \f$
616 
617  where
618 
619  \f$ \tilde{u}_{i,j} = w^0_i w^1_j u(\xi_{0,i},\xi_{1,j}) \f$
620 
621  which can be implemented as
622 
623  \f$ f_{qi} = \sum_{j=0}^{nq_1} \phi_q(\xi_{1,j}) \tilde{u}_{i,j} =
624  {\bf B_1 U} \f$
625  \f$ I_{pq} = \sum_{i=0}^{nq_0} \phi_p(\xi_{0,i}) f_{qi} =
626  {\bf B_0 F} \f$
627  *
628  * \param inarray contains the values of the function \a f
629  * evaluated at the quadrature points
630  * \param outarray contains the values of the inner product of \a f
631  * with the different modes, i.e. \f$ outarray[p] = I[p]\f$
632  * (output of the function)
633  */
635  Array<OneD, NekDouble> &outarray)
636  {
637  v_IProductWRTBase(inarray, outarray);
638  }
639 
641  const Array<OneD, const NekDouble>& base,
642  const Array<OneD, const NekDouble>& inarray,
643  Array<OneD, NekDouble> &outarray,
644  int coll_check)
645  {
646  v_IProductWRTBase(base, inarray, outarray, coll_check);
647  }
648 
649 
651  const int dir,
652  const Array<OneD, const NekDouble>& inarray,
653  Array<OneD, NekDouble> &outarray)
654  {
655  v_IProductWRTDerivBase(dir,inarray, outarray);
656  }
657 
659  const Array<OneD, const NekDouble>& direction,
660  const Array<OneD, const NekDouble>& inarray,
661  Array<OneD, NekDouble> &outarray)
662  {
663  v_IProductWRTDirectionalDerivBase(direction, inarray, outarray);
664  }
665 
666  /// \brief Get the element id of this expansion when used
667  /// in a list by returning value of #m_elmt_id
668  inline int GetElmtId()
669  {
670  return m_elmt_id;
671  }
672 
673 
674  /// \brief Set the element id of this expansion when used
675  /// in a list by returning value of #m_elmt_id
676  inline void SetElmtId(const int id)
677  {
678  m_elmt_id = id;
679  }
680 
681  /** \brief this function returns the physical coordinates of the
682  * quadrature points of the expansion
683  *
684  * This function is a wrapper around the virtual function
685  * \a v_GetCoords()
686  *
687  * \param coords an array containing the coordinates of the
688  * quadrature points (output of the function)
689  */
693  {
694  v_GetCoords(coords_1,coords_2,coords_3);
695  }
696 
697  /** \brief given the coordinates of a point of the element in the
698  * local collapsed coordinate system, this function calculates the
699  * physical coordinates of the point
700  *
701  * This function is a wrapper around the virtual function
702  * \a v_GetCoord()
703  *
704  * \param Lcoords the coordinates in the local collapsed
705  * coordinate system
706  * \param coords the physical coordinates (output of the function)
707  */
709  Array<OneD, NekDouble> &coord)
710  {
711  v_GetCoord(Lcoord, coord);
712  }
713 
715  {
716  return m_stdMatrixManager[mkey];
717  }
718 
720  {
721  return m_stdStaticCondMatrixManager[mkey];
722  }
723 
725  {
726  return m_IndexMapManager[ikey];
727  }
728 
730  {
731  return v_GetPhysNormals();
732  }
733 
735  {
736  v_SetPhysNormals(normal);
737  }
738 
739  STD_REGIONS_EXPORT virtual void SetUpPhysNormals(const int edge);
740 
742  {
743  v_NormVectorIProductWRTBase(Fx,outarray);
744  }
745 
747  {
748  v_NormVectorIProductWRTBase(Fx,Fy,outarray);
749  }
750 
752  {
753  v_NormVectorIProductWRTBase(Fx,Fy,Fz,outarray);
754  }
755 
757  {
758  v_NormVectorIProductWRTBase(Fvec, outarray);
759  }
760 
762  {
763  return v_GetLocStaticCondMatrix(mkey);
764  }
765 
767  {
768  return v_DropLocStaticCondMatrix(mkey);
769  }
770 
772  {
773  return v_GetForient(face);
774  }
775 
777  {
778  return v_GetEorient(edge);
779  }
780 
782  Array<OneD, NekDouble> &coeffs,
784  {
785  v_SetCoeffsToOrientation(coeffs, dir);
786  }
787 
791  Array<OneD, NekDouble> &outarray)
792  {
793  v_SetCoeffsToOrientation(dir,inarray,outarray);
794  }
795 
796  int CalcNumberOfCoefficients(const std::vector<unsigned int> &nummodes, int &modes_offset)
797  {
798  return v_CalcNumberOfCoefficients(nummodes,modes_offset);
799  }
800 
801  // virtual functions related to LocalRegions
803  const Array<OneD, const NekDouble> &Lcoord,
804  const Array<OneD, const NekDouble> &physvals);
805 
806 
808  {
809  return v_GetCoordim();
810  }
811 
813  {
814  v_GetBoundaryMap(outarray);
815  }
816 
818  {
819  v_GetInteriorMap(outarray);
820  }
821 
822  int GetVertexMap(const int localVertexId,
823  bool useCoeffPacking = false)
824  {
825  return v_GetVertexMap(localVertexId,useCoeffPacking);
826  }
827 
828  void GetEdgeInteriorMap(const int eid, const Orientation edgeOrient,
829  Array<OneD, unsigned int> &maparray,
830  Array<OneD, int> &signarray)
831  {
832  v_GetEdgeInteriorMap(eid,edgeOrient,maparray,signarray);
833  }
834 
835  void GetFaceNumModes(const int fid, const Orientation faceOrient,
836  int &numModes0,
837  int &numModes1)
838  {
839  v_GetFaceNumModes(fid,faceOrient,numModes0,numModes1);
840  }
841 
842  void GetFaceInteriorMap(const int fid, const Orientation faceOrient,
843  Array<OneD, unsigned int> &maparray,
844  Array<OneD, int> &signarray)
845  {
846  v_GetFaceInteriorMap(fid,faceOrient,maparray,signarray);
847  }
848 
849  void GetEdgeToElementMap(const int eid,
850  const Orientation edgeOrient,
851  Array<OneD, unsigned int> &maparray,
852  Array<OneD, int> &signarray,
853  int P = -1)
854  {
855  v_GetEdgeToElementMap(eid, edgeOrient, maparray, signarray, P);
856  }
857 
858  void GetFaceToElementMap(const int fid, const Orientation faceOrient,
859  Array<OneD, unsigned int> &maparray,
860  Array<OneD, int> &signarray,
861  int nummodesA = -1, int nummodesB = -1)
862  {
863  v_GetFaceToElementMap(fid,faceOrient,maparray,signarray,
864  nummodesA,nummodesB);
865  }
866 
867 
868  /**
869  * @brief Extract the physical values along edge \a edge from \a
870  * inarray into \a outarray following the local edge orientation
871  * and point distribution defined by defined in \a EdgeExp.
872  */
873 
874  void GetEdgePhysVals(const int edge, const Array<OneD,
875  const NekDouble> &inarray,
876  Array<OneD,NekDouble> &outarray)
877  {
878  v_GetEdgePhysVals(edge,inarray,outarray);
879  }
880 
881  void GetEdgePhysVals(const int edge,
882  const std::shared_ptr<StdExpansion> &EdgeExp,
883  const Array<OneD, const NekDouble> &inarray,
884  Array<OneD,NekDouble> &outarray)
885  {
886  v_GetEdgePhysVals(edge,EdgeExp,inarray,outarray);
887  }
888 
889  void GetTracePhysVals(const int edge, const std::shared_ptr<StdExpansion> &EdgeExp, const Array<OneD, const NekDouble> &inarray, Array<OneD,NekDouble> &outarray)
890  {
891  v_GetTracePhysVals(edge,EdgeExp,inarray,outarray);
892  }
893 
894  void GetVertexPhysVals(const int vertex,
895  const Array<OneD, const NekDouble> &inarray,
896  NekDouble &outarray)
897  {
898  v_GetVertexPhysVals(vertex, inarray, outarray);
899  }
900 
901  void GetEdgeInterpVals(const int edge,const Array<OneD,
902  const NekDouble> &inarray,
903  Array<OneD,NekDouble> &outarray)
904  {
905  v_GetEdgeInterpVals(edge, inarray, outarray);
906  }
907 
908  /**
909  * @brief Extract the metric factors to compute the contravariant
910  * fluxes along edge \a edge and stores them into \a outarray
911  * following the local edge orientation (i.e. anticlockwise
912  * convention).
913  */
915  const int edge,
916  Array<OneD, NekDouble> &outarray)
917  {
918  v_GetEdgeQFactors(edge, outarray);
919  }
920 
922  const int face,
923  const std::shared_ptr<StdExpansion> &FaceExp,
924  const Array<OneD, const NekDouble> &inarray,
925  Array<OneD, NekDouble> &outarray,
927  {
928  v_GetFacePhysVals(face, FaceExp, inarray, outarray, orient);
929  }
930 
932  const int edge,
933  Array<OneD, int> &outarray)
934  {
935  v_GetEdgePhysMap(edge, outarray);
936  }
937 
939  const int face,
940  Array<OneD, int> &outarray)
941  {
942  v_GetFacePhysMap(face, outarray);
943  }
944 
946  const Array<OneD, const NekDouble> &inarray,
947  Array<OneD, NekDouble> &outarray)
948  {
949  v_MultiplyByQuadratureMetric(inarray, outarray);
950  }
951 
953  const Array<OneD, const NekDouble> &inarray,
954  Array<OneD, NekDouble> & outarray)
955  {
956  v_MultiplyByStdQuadratureMetric(inarray, outarray);
957  }
958 
959  // Matrix Routines
960 
961  /** \brief this function generates the mass matrix
962  * \f$\mathbf{M}[i][j] =
963  * \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}\f$
964  *
965  * \return returns the mass matrix
966  */
967 
969 
971  Array<OneD,NekDouble> &outarray,
972  const StdMatrixKey &mkey);
973 
975  Array<OneD,NekDouble> &outarray,
976  const StdMatrixKey &mkey)
977  {
978  v_MassMatrixOp(inarray,outarray,mkey);
979  }
980 
982  Array<OneD,NekDouble> &outarray,
983  const StdMatrixKey &mkey)
984  {
985  v_LaplacianMatrixOp(inarray,outarray,mkey);
986  }
987 
988  void ReduceOrderCoeffs(int numMin,
989  const Array<OneD, const NekDouble> &inarray,
990  Array<OneD,NekDouble> &outarray)
991  {
992  v_ReduceOrderCoeffs(numMin,inarray,outarray);
993  }
994 
996  const StdMatrixKey &mkey)
997  {
998  v_SVVLaplacianFilter(array,mkey);
999  }
1000 
1002  const NekDouble alpha,
1003  const NekDouble exponent,
1004  const NekDouble cutoff)
1005  {
1006  v_ExponentialFilter(array, alpha, exponent, cutoff);
1007  }
1008 
1009  void LaplacianMatrixOp(const int k1, const int k2,
1010  const Array<OneD, const NekDouble> &inarray,
1011  Array<OneD,NekDouble> &outarray,
1012  const StdMatrixKey &mkey)
1013  {
1014  v_LaplacianMatrixOp(k1,k2,inarray,outarray,mkey);
1015  }
1016 
1017  void WeakDerivMatrixOp(const int i,
1018  const Array<OneD, const NekDouble> &inarray,
1019  Array<OneD,NekDouble> &outarray,
1020  const StdMatrixKey &mkey)
1021  {
1022  v_WeakDerivMatrixOp(i,inarray,outarray,mkey);
1023  }
1024 
1026  Array<OneD,NekDouble> &outarray,
1027  const StdMatrixKey &mkey)
1028  {
1029  v_WeakDirectionalDerivMatrixOp(inarray,outarray,mkey);
1030  }
1031 
1033  Array<OneD,NekDouble> &outarray,
1034  const StdMatrixKey &mkey)
1035  {
1036  v_MassLevelCurvatureMatrixOp(inarray,outarray,mkey);
1037  }
1038 
1040  Array<OneD,NekDouble> &outarray,
1041  const StdMatrixKey &mkey,
1042  bool addDiffusionTerm = true)
1043  {
1044  v_LinearAdvectionDiffusionReactionMatrixOp(inarray,outarray,mkey,addDiffusionTerm);
1045  }
1046 
1047  /**
1048  * @param inarray Input array @f$ \mathbf{u} @f$.
1049  * @param outarray Output array @f$ \boldsymbol{\nabla^2u}
1050  * + \lambda \boldsymbol{u} @f$.
1051  * @param mkey
1052  */
1054  Array<OneD,NekDouble> &outarray,
1055  const StdMatrixKey &mkey)
1056  {
1057  v_HelmholtzMatrixOp(inarray,outarray,mkey);
1058  }
1059 
1061  {
1062  return v_GenMatrix(mkey);
1063  }
1064 
1066  Array<OneD, NekDouble> &out_d0,
1069  {
1070  v_PhysDeriv (inarray, out_d0, out_d1, out_d2);
1071  }
1072 
1073  void PhysDeriv(const int dir,
1074  const Array<OneD, const NekDouble>& inarray,
1075  Array<OneD, NekDouble> &outarray)
1076  {
1077  v_PhysDeriv (dir, inarray, outarray);
1078  }
1079 
1081  Array<OneD, NekDouble> &out_ds)
1082  {
1083  v_PhysDeriv_s(inarray,out_ds);
1084  }
1085 
1087  Array<OneD, NekDouble>& out_dn)
1088  {
1089  v_PhysDeriv_n(inarray,out_dn);
1090  }
1091 
1093  const Array<OneD, const NekDouble>& direction,
1094  Array<OneD, NekDouble> &outarray)
1095  {
1096  v_PhysDirectionalDeriv (inarray, direction, outarray);
1097  }
1098 
1100  Array<OneD, NekDouble> &out_d0,
1103  {
1104  v_StdPhysDeriv(inarray, out_d0, out_d1, out_d2);
1105  }
1106 
1107  void StdPhysDeriv (const int dir,
1108  const Array<OneD, const NekDouble>& inarray,
1109  Array<OneD, NekDouble> &outarray)
1110  {
1111  v_StdPhysDeriv(dir,inarray,outarray);
1112  }
1113 
1114  void AddRobinMassMatrix(const int edgeid, const Array<OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
1115  {
1116  v_AddRobinMassMatrix(edgeid,primCoeffs,inoutmat);
1117  }
1118 
1119  void AddRobinEdgeContribution(const int edgeid, const Array<OneD, const NekDouble> &primCoeffs, Array<OneD, NekDouble> &coeffs)
1120  {
1121  v_AddRobinEdgeContribution(edgeid, primCoeffs, coeffs);
1122  }
1123 
1124  /** \brief This function evaluates the expansion at a single
1125  * (arbitrary) point of the domain
1126  *
1127  * This function is a wrapper around the virtual function
1128  * \a v_PhysEvaluate()
1129  *
1130  * Based on the value of the expansion at the quadrature
1131  * points provided in \a physvals, this function
1132  * calculates the value of the expansion at an arbitrary
1133  * single points (with coordinates \f$ \mathbf{x_c}\f$
1134  * given by the pointer \a coords). This operation,
1135  * equivalent to \f[ u(\mathbf{x_c}) = \sum_p
1136  * \phi_p(\mathbf{x_c}) \hat{u}_p \f] is evaluated using
1137  * Lagrangian interpolants through the quadrature points:
1138  * \f[ u(\mathbf{x_c}) = \sum_p h_p(\mathbf{x_c}) u_p\f]
1139  *
1140  * \param coords the coordinates of the single point
1141  * \param physvals the interpolated field at the quadrature points
1142  *
1143  * \return returns the value of the expansion at the
1144  * single point
1145  */
1147  const Array<OneD, const NekDouble>& physvals)
1148  {
1149  return v_PhysEvaluate(coords,physvals);
1150  }
1151 
1152 
1153  /** \brief This function evaluates the expansion at a single
1154  * (arbitrary) point of the domain
1155  *
1156  * This function is a wrapper around the virtual function
1157  * \a v_PhysEvaluate()
1158  *
1159  * Based on the value of the expansion at the quadrature
1160  * points provided in \a physvals, this function
1161  * calculates the value of the expansion at an arbitrary
1162  * single points associated with the interpolation
1163  * matrices provided in \f$ I \f$.
1164  *
1165  * \param I an Array of lagrange interpolantes evaluated
1166  * at the coordinate and going through the local physical
1167  * quadrature
1168  * \param physvals the interpolated field at the quadrature points
1169  *
1170  * \return returns the value of the expansion at the
1171  * single point
1172  */
1174  const Array<OneD, const NekDouble >& physvals)
1175  {
1176  return v_PhysEvaluate(I,physvals);
1177  }
1178 
1179 
1180  /**
1181  * \brief Convert local cartesian coordinate \a xi into local
1182  * collapsed coordinates \a eta
1183  **/
1186  {
1187  v_LocCoordToLocCollapsed(xi,eta);
1188  }
1189 
1190  /// \brief Get the element id of this expansion when used
1191  /// in a list by returning value of #m_elmt_id
1192  STD_REGIONS_EXPORT virtual int v_GetElmtId();
1193 
1195 
1197 
1198  STD_REGIONS_EXPORT virtual void v_SetUpPhysNormals(const int edge);
1199 
1200  STD_REGIONS_EXPORT virtual int v_CalcNumberOfCoefficients(const std::vector<unsigned int> &nummodes, int &modes_offset);
1201 
1203 
1205  const Array<OneD, const NekDouble> &Fx,
1206  const Array<OneD, const NekDouble> &Fy,
1207  Array< OneD, NekDouble> &outarray);
1208 
1210 
1212 
1214 
1216 
1217 
1219 
1221 
1222  /** \brief Function to evaluate the discrete \f$ L_\infty\f$
1223  * error \f$ |\epsilon|_\infty = \max |u - u_{exact}|\f$ where \f$
1224  * u_{exact}\f$ is given by the array \a sol.
1225  *
1226  * This function takes the physical value space array \a m_phys as
1227  * approximate solution
1228  *
1229  * \param sol array of solution function at physical quadrature
1230  * points
1231  * \return returns the \f$ L_\infty \f$ error as a NekDouble.
1232  */
1234 
1235  /** \brief Function to evaluate the discrete \f$ L_2\f$ error,
1236  * \f$ | \epsilon |_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2
1237  * dx \right]^{1/2} d\xi_1 \f$ where \f$ u_{exact}\f$ is given by
1238  * the array \a sol.
1239  *
1240  * This function takes the physical value space array \a m_phys as
1241  * approximate solution
1242  *
1243  * \param sol array of solution function at physical quadrature
1244  * points
1245  * \return returns the \f$ L_2 \f$ error as a double.
1246  */
1248 
1249  /** \brief Function to evaluate the discrete \f$ H^1\f$
1250  * error, \f$ | \epsilon |^1_{2} = \left [ \int^1_{-1} [u -
1251  * u_{exact}]^2 + \nabla(u - u_{exact})\cdot\nabla(u -
1252  * u_{exact})\cdot dx \right]^{1/2} d\xi_1 \f$ where \f$
1253  * u_{exact}\f$ is given by the array \a sol.
1254  *
1255  * This function takes the physical value space array
1256  * \a m_phys as approximate solution
1257  *
1258  * \param sol array of solution function at physical quadrature
1259  * points
1260  * \return returns the \f$ H_1 \f$ error as a double.
1261  */
1263 
1264  // I/O routines
1265  const NormalVector & GetEdgeNormal(const int edge) const
1266  {
1267  return v_GetEdgeNormal(edge);
1268  }
1269 
1270  void ComputeEdgeNormal(const int edge)
1271  {
1272  v_ComputeEdgeNormal(edge);
1273  }
1274 
1275  void NegateEdgeNormal(const int edge)
1276  {
1277  v_NegateEdgeNormal(edge);
1278  }
1279 
1280  bool EdgeNormalNegated(const int edge)
1281  {
1282  return v_EdgeNormalNegated(edge);
1283  }
1284 
1285  void ComputeFaceNormal(const int face)
1286  {
1287  v_ComputeFaceNormal(face);
1288  }
1289 
1290  void NegateFaceNormal(const int face)
1291  {
1292  v_NegateFaceNormal(face);
1293  }
1294 
1295  bool FaceNormalNegated(const int face)
1296  {
1297  return v_FaceNormalNegated(face);
1298  }
1299 
1300  void ComputeVertexNormal(const int vertex)
1301  {
1302  v_ComputeVertexNormal(vertex);
1303  }
1304 
1305  void NegateVertexNormal(const int vertex)
1306  {
1307  v_NegateVertexNormal(vertex);
1308  }
1309 
1310  bool VertexNormalNegated(const int vertex)
1311  {
1312  return v_VertexNormalNegated(vertex);
1313  }
1314 
1315  const NormalVector & GetFaceNormal(const int face) const
1316  {
1317  return v_GetFaceNormal(face);
1318  }
1319 
1320  const NormalVector & GetVertexNormal(const int vertex) const
1321  {
1322  return v_GetVertexNormal(vertex);
1323  }
1324 
1325  const NormalVector & GetSurfaceNormal(const int id) const
1326  {
1327  return v_GetSurfaceNormal(id);
1328  }
1329 
1331  {
1333  p.reserve(m_base.num_elements());
1334  for (int i = 0; i < m_base.num_elements(); ++i)
1335  {
1336  p.push_back(m_base[i]->GetPointsKey());
1337  }
1338  return p;
1339  }
1340 
1343  {
1344  return v_GetEdgeInverseBoundaryMap(eid);
1345  }
1346 
1347 
1349  GetFaceInverseBoundaryMap(int fid, StdRegions::Orientation faceOrient = eNoOrientation, int P1=-1, int P2=-1)
1350  {
1351  return v_GetFaceInverseBoundaryMap(fid,faceOrient,P1,P2);
1352  }
1353 
1354 
1359 
1360  {
1361  v_GetInverseBoundaryMaps(vmap,emap,fmap);
1362  }
1363 
1365  const DNekScalMatSharedPtr & m_transformationmatrix)
1366  {
1368  m_transformationmatrix);
1369  }
1370 
1371 
1372  /** \brief This function performs an interpolation from
1373  * the physical space points provided at input into an
1374  * array of equispaced points which are not the collapsed
1375  * coordinate. So for a tetrahedron you will only get a
1376  * tetrahedral number of values.
1377  *
1378  * This is primarily used for output purposes to get a
1379  * better distribution of points more suitable for most
1380  * postprocessing
1381  */
1383  const Array<OneD, const NekDouble> &inarray,
1384  Array<OneD, NekDouble> &outarray,
1385  int npset = -1);
1386 
1387  /** \brief This function provides the connectivity of
1388  * local simplices (triangles or tets) to connect the
1389  * equispaced data points provided by
1390  * PhysInterpToSimplexEquiSpaced
1391  *
1392  * This is a virtual call to the function
1393  * \a v_GetSimplexEquiSpaceConnectivity
1394  */
1396  Array<OneD, int> &conn,
1397  bool standard = true)
1398  {
1399  v_GetSimplexEquiSpacedConnectivity(conn,standard);
1400  }
1401 
1402  /** \brief This function performs a
1403  * projection/interpolation from the equispaced points
1404  * sometimes used in post-processing onto the coefficient
1405  * space
1406  *
1407  * This is primarily used for output purposes to use a
1408  * more even distribution of points more suitable for alot of
1409  * postprocessing
1410  */
1412  const Array<OneD, const NekDouble> &inarray,
1413  Array<OneD, NekDouble> &outarray);
1414 
1415  template<class T>
1416  std::shared_ptr<T> as()
1417  {
1418  return std::dynamic_pointer_cast<T>( shared_from_this() );
1419  }
1420 
1422  Array<OneD, NekDouble> &outarray,
1423  bool multiplybyweights = true)
1424  {
1425  v_IProductWRTBase_SumFac(inarray,outarray,multiplybyweights);
1426  }
1427 
1428  protected:
1429  Array<OneD, LibUtilities::BasisSharedPtr> m_base; /**< Bases needed for the expansion */
1431  int m_ncoeffs; /**< Total number of coefficients used in the expansion */
1435 
1437  {
1438  return v_CreateStdMatrix(mkey);
1439  }
1440 
1441  /** \brief Create the static condensation of a matrix when
1442  using a boundary interior decomposition
1443 
1444  If a matrix system can be represented by
1445  \f$ Mat = \left [ \begin{array}{cc}
1446  A & B \\
1447  C & D \end{array} \right ] \f$
1448  This routine creates a matrix containing the statically
1449  condense system of the form
1450  \f$ Mat = \left [ \begin{array}{cc}
1451  A - B D^{-1} C & B D^{-1} \\
1452  D^{-1} C & D^{-1} \end{array} \right ] \f$
1453  **/
1455 
1456  /** \brief Create an IndexMap which contains mapping information linking any specific
1457  element shape with either its boundaries, edges, faces, verteces, etc.
1458 
1459  The index member of the IndexMapValue struct gives back an integer associated with an entity index
1460  The sign member of the same struct gives back a sign to algebrically apply entities orientation
1461  **/
1463 
1465  Array<OneD, NekDouble> &outarray);
1466 
1468  Array<OneD, NekDouble> &outarray)
1469  {
1470  v_BwdTrans_SumFac(inarray,outarray);
1471  }
1472 
1473 
1475  const int dir,
1476  const Array<OneD, const NekDouble>& inarray,
1477  Array<OneD, NekDouble> &outarray)
1478  {
1479  v_IProductWRTDerivBase_SumFac(dir,inarray,outarray);
1480  }
1481 
1482 
1484  const Array<OneD, const NekDouble>& direction,
1485  const Array<OneD, const NekDouble>& inarray,
1486  Array<OneD, NekDouble> &outarray)
1487  {
1489  inarray,outarray);
1490  }
1491 
1492 
1493  // The term _MatFree denotes that the action of the MatrixOperation
1494  // is done withouth actually using the matrix (which then needs to be stored/calculated).
1495  // Although this does not strictly mean that no matrix operations are involved in the
1496  // evaluation of the operation, we use this term in the same context used as in the following
1497  // paper:
1498  // R. C. Kirby, M. G. Knepley, A. Logg, and L. R. Scott,
1499  // "Optimizing the evaluation of finite element matrices," SISC 27:741-758 (2005)
1501  Array<OneD,NekDouble> &outarray,
1502  const StdMatrixKey &mkey);
1503 
1505  Array<OneD,NekDouble> &outarray,
1506  const StdMatrixKey &mkey);
1507 
1509  Array<OneD,NekDouble> &outarray,
1510  const StdMatrixKey &mkey)
1511  {
1512  v_LaplacianMatrixOp_MatFree(inarray,outarray,mkey);
1513  }
1514 
1516  const Array<OneD, const NekDouble> &inarray,
1517  Array<OneD, NekDouble> &outarray,
1519  {
1520  v_LaplacianMatrixOp_MatFree_Kernel(inarray, outarray, wsp);
1521  }
1522 
1524  Array<OneD,NekDouble> &outarray,
1525  const StdMatrixKey &mkey);
1526 
1527  STD_REGIONS_EXPORT void LaplacianMatrixOp_MatFree(const int k1, const int k2,
1528  const Array<OneD, const NekDouble> &inarray,
1529  Array<OneD,NekDouble> &outarray,
1530  const StdMatrixKey &mkey);
1531 
1533  const Array<OneD, const NekDouble> &inarray,
1534  Array<OneD,NekDouble> &outarray,
1535  const StdMatrixKey &mkey);
1536 
1538  Array<OneD,NekDouble> &outarray,
1539  const StdMatrixKey &mkey);
1540 
1542  Array<OneD,NekDouble> &outarray,
1543  const StdMatrixKey &mkey);
1544 
1546  Array<OneD,NekDouble> &outarray,
1547  const StdMatrixKey &mkey,
1548  bool addDiffusionTerm = true);
1549 
1551  Array<OneD,NekDouble> &outarray,
1552  const StdMatrixKey &mkey)
1553  {
1554  v_HelmholtzMatrixOp_MatFree(inarray,outarray,mkey);
1555  }
1556 
1558  Array<OneD,NekDouble> &outarray,
1559  const StdMatrixKey &mkey);
1560 
1563  Array<OneD, NekDouble> &outarray);
1564 
1566  Array<OneD, NekDouble> &coeffs,
1568 
1570  const Array<OneD, const NekDouble> &Lcoord,
1571  const Array<OneD, const NekDouble> &physvals);
1572 
1573  private:
1574  // Virtual functions
1575  STD_REGIONS_EXPORT virtual int v_GetNverts() const = 0;
1576  STD_REGIONS_EXPORT virtual int v_GetNedges() const;
1577 
1578  STD_REGIONS_EXPORT virtual int v_GetNfaces() const;
1579 
1580  STD_REGIONS_EXPORT virtual int v_NumBndryCoeffs() const;
1581 
1582  STD_REGIONS_EXPORT virtual int v_NumDGBndryCoeffs() const;
1583 
1584  STD_REGIONS_EXPORT virtual int v_GetEdgeNcoeffs(const int i) const;
1585 
1586  STD_REGIONS_EXPORT virtual int v_GetTotalEdgeIntNcoeffs() const;
1587 
1588  STD_REGIONS_EXPORT virtual int v_GetEdgeNumPoints(const int i) const;
1589 
1590  STD_REGIONS_EXPORT virtual int v_DetCartesianDirOfEdge(const int edge);
1591 
1592  STD_REGIONS_EXPORT virtual const LibUtilities::BasisKey v_DetEdgeBasisKey(const int i) const;
1593 
1594  STD_REGIONS_EXPORT virtual const LibUtilities::BasisKey v_DetFaceBasisKey(const int i, const int k) const;
1595 
1596  STD_REGIONS_EXPORT virtual int v_GetFaceNumPoints(const int i) const;
1597 
1598  STD_REGIONS_EXPORT virtual int v_GetFaceNcoeffs(const int i) const;
1599 
1600  STD_REGIONS_EXPORT virtual int v_GetFaceIntNcoeffs(const int i) const;
1601 
1602  STD_REGIONS_EXPORT virtual int v_GetTotalFaceIntNcoeffs() const;
1603 
1604 
1605  STD_REGIONS_EXPORT virtual int v_GetTraceNcoeffs(const int i) const;
1606 
1607  STD_REGIONS_EXPORT virtual LibUtilities::PointsKey v_GetFacePointsKey(const int i, const int j) const;
1608 
1610 
1612 
1614 
1615  STD_REGIONS_EXPORT virtual std::shared_ptr<StdExpansion>
1616  v_GetStdExp(void) const;
1617 
1618  STD_REGIONS_EXPORT virtual std::shared_ptr<StdExpansion>
1619  v_GetLinStdExp(void) const;
1620 
1621  STD_REGIONS_EXPORT virtual int v_GetShapeDimension() const;
1622 
1624 
1626 
1628  Array<OneD, NekDouble> &outarray) = 0;
1629 
1630  /**
1631  * @brief Transform a given function from physical quadrature space
1632  * to coefficient space.
1633  * @see StdExpansion::FwdTrans
1634  */
1636  const Array<OneD, const NekDouble>& inarray,
1637  Array<OneD, NekDouble> &outarray) = 0;
1638 
1639  /**
1640  * @brief Calculates the inner product of a given function \a f
1641  * with the different modes of the expansion
1642  */
1644  const Array<OneD, const NekDouble>& inarray,
1645  Array<OneD, NekDouble> &outarray) = 0;
1646 
1648  const Array<OneD, const NekDouble>& base,
1649  const Array<OneD, const NekDouble>& inarray,
1650  Array<OneD, NekDouble>& outarray,
1651  int coll_check)
1652  {
1653  boost::ignore_unused(base, inarray, outarray, coll_check);
1655  "StdExpansion::v_IProductWRTBase has no (and should have no) implementation");
1656  }
1657 
1658  STD_REGIONS_EXPORT virtual void v_IProductWRTDerivBase (const int dir,
1659  const Array<OneD, const NekDouble>& inarray,
1660  Array<OneD, NekDouble> &outarray);
1662  const Array<OneD, const NekDouble>& inarray,
1663  Array<OneD, NekDouble> &outarray);
1664 
1666  Array<OneD, NekDouble> &outarray);
1667 
1669 
1670  STD_REGIONS_EXPORT virtual void v_PhysDeriv (const Array<OneD, const NekDouble>& inarray,
1671  Array<OneD, NekDouble> &out_d1,
1672  Array<OneD, NekDouble> &out_d2,
1673  Array<OneD, NekDouble> &out_d3);
1674 
1675  STD_REGIONS_EXPORT virtual void v_PhysDeriv_s (const Array<OneD, const NekDouble>& inarray,
1676  Array<OneD, NekDouble> &out_ds);
1677 
1679  Array<OneD, NekDouble>& out_dn);
1680  STD_REGIONS_EXPORT virtual void v_PhysDeriv(const int dir,
1681  const Array<OneD, const NekDouble>& inarray,
1682  Array<OneD, NekDouble> &out_d0);
1683 
1685  const Array<OneD, const NekDouble>& direction,
1686  Array<OneD, NekDouble> &outarray);
1687 
1688  STD_REGIONS_EXPORT virtual void v_StdPhysDeriv (const Array<OneD, const NekDouble>& inarray,
1689  Array<OneD, NekDouble> &out_d1,
1690  Array<OneD, NekDouble> &out_d2,
1691  Array<OneD, NekDouble> &out_d3);
1692 
1693  STD_REGIONS_EXPORT virtual void v_StdPhysDeriv (const int dir,
1694  const Array<OneD, const NekDouble>& inarray,
1695  Array<OneD, NekDouble> &outarray);
1696 
1697  STD_REGIONS_EXPORT virtual void v_AddRobinMassMatrix(const int edgeid, const Array<OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat);
1698 
1699  STD_REGIONS_EXPORT virtual void v_AddRobinEdgeContribution(const int edgeid, const Array<OneD, const NekDouble> &primCoeffs, Array<OneD, NekDouble> &coeffs);
1700 
1702 
1704 
1706  const Array<OneD, const NekDouble>& xi,
1707  Array<OneD, NekDouble>& eta);
1708 
1709 
1710  STD_REGIONS_EXPORT virtual void v_FillMode(const int mode, Array<OneD, NekDouble> &outarray);
1711 
1713 
1715 
1716  STD_REGIONS_EXPORT virtual void v_GetCoords(Array<OneD, NekDouble> &coords_0,
1717  Array<OneD, NekDouble> &coords_1,
1718  Array<OneD, NekDouble> &coords_2);
1719 
1720  STD_REGIONS_EXPORT virtual void v_GetCoord(const Array<OneD, const NekDouble>& Lcoord,
1721  Array<OneD, NekDouble> &coord);
1722 
1723  STD_REGIONS_EXPORT virtual int v_GetCoordim(void);
1724 
1726 
1728 
1729  STD_REGIONS_EXPORT virtual int v_GetVertexMap(int localVertexId,
1730  bool useCoeffPacking = false);
1731 
1732  STD_REGIONS_EXPORT virtual void v_GetEdgeInteriorMap(const int eid, const Orientation edgeOrient,
1733  Array<OneD, unsigned int> &maparray,
1734  Array<OneD, int> &signarray);
1735 
1737  const int fid,
1738  const Orientation faceOrient,
1739  int &numModes0,
1740  int &numModes1);
1741 
1742  STD_REGIONS_EXPORT virtual void v_GetFaceInteriorMap(const int fid, const Orientation faceOrient,
1743  Array<OneD, unsigned int> &maparray,
1744  Array<OneD, int> &signarray);
1745 
1747  const int eid,
1748  const Orientation edgeOrient,
1749  Array<OneD, unsigned int>& maparray,
1750  Array<OneD, int>& signarray,
1751  int P = -1);
1752 
1753  STD_REGIONS_EXPORT virtual void v_GetFaceToElementMap(const int fid, const Orientation faceOrient,
1754  Array<OneD, unsigned int> &maparray,
1755  Array<OneD, int> &signarray,
1756  int nummodesA = -1, int nummodesB = -1);
1757 
1758  /**
1759  * @brief Extract the physical values along edge \a edge from \a
1760  * inarray into \a outarray following the local edge orientation
1761  * and point distribution defined by defined in \a EdgeExp.
1762  */
1763  STD_REGIONS_EXPORT virtual void v_GetEdgePhysVals(const int edge, const Array<OneD, const NekDouble> &inarray, Array<OneD,NekDouble> &outarray);
1764 
1765  STD_REGIONS_EXPORT virtual void v_GetEdgePhysVals(const int edge, const std::shared_ptr<StdExpansion> &EdgeExp, const Array<OneD, const NekDouble> &inarray, Array<OneD,NekDouble> &outarray);
1766 
1767  STD_REGIONS_EXPORT virtual void v_GetTracePhysVals(const int edge, const std::shared_ptr<StdExpansion> &EdgeExp, const Array<OneD, const NekDouble> &inarray, Array<OneD,NekDouble> &outarray, StdRegions::Orientation orient = eNoOrientation);
1768 
1769  STD_REGIONS_EXPORT virtual void v_GetVertexPhysVals(const int vertex, const Array<OneD, const NekDouble> &inarray, NekDouble &outarray);
1770 
1771  STD_REGIONS_EXPORT virtual void v_GetEdgeInterpVals(const int edge,
1772  const Array<OneD, const NekDouble> &inarray,Array<OneD,NekDouble> &outarray);
1773 
1775  const int edge,
1776  Array<OneD, NekDouble> &outarray);
1777 
1779  const int face,
1780  const std::shared_ptr<StdExpansion> &FaceExp,
1781  const Array<OneD, const NekDouble> &inarray,
1782  Array<OneD, NekDouble> &outarray,
1783  StdRegions::Orientation orient);
1784 
1785  STD_REGIONS_EXPORT virtual void v_GetEdgePhysMap(
1786  const int edge,
1787  Array<OneD,int> &outarray);
1788 
1789  STD_REGIONS_EXPORT virtual void v_GetFacePhysMap(
1790  const int face,
1791  Array<OneD,int> &outarray);
1792 
1794  const Array<OneD, const NekDouble> &inarray,
1795  Array<OneD, NekDouble> &outarray);
1796 
1798  const Array<OneD, const NekDouble> &inarray,
1799  Array<OneD, NekDouble> &outarray);
1800 
1802  Array<OneD, NekDouble> &outarray);
1803 
1805  Array<OneD, NekDouble> &outarray, bool multiplybyweights = true);
1806 
1808  const int dir,
1809  const Array<OneD, const NekDouble>& inarray,
1810  Array<OneD, NekDouble> &outarray);
1811 
1812  STD_REGIONS_EXPORT virtual void
1814  const Array<OneD, const NekDouble>& direction,
1815  const Array<OneD, const NekDouble>& inarray,
1816  Array<OneD, NekDouble> &outarray);
1817 
1819  Array<OneD,NekDouble> &outarray,
1820  const StdMatrixKey &mkey);
1821 
1823  Array<OneD,NekDouble> &outarray,
1824  const StdMatrixKey &mkey);
1825 
1827  const StdMatrixKey &mkey);
1828 
1830  Array<OneD, NekDouble> &array,
1831  const NekDouble alpha,
1832  const NekDouble exponent,
1833  const NekDouble cutoff);
1834 
1836  int numMin,
1837  const Array<OneD, const NekDouble> &inarray,
1838  Array<OneD,NekDouble> &outarray);
1839 
1840  STD_REGIONS_EXPORT virtual void v_LaplacianMatrixOp(const int k1, const int k2,
1841  const Array<OneD, const NekDouble> &inarray,
1842  Array<OneD,NekDouble> &outarray,
1843  const StdMatrixKey &mkey);
1844 
1845  STD_REGIONS_EXPORT virtual void v_WeakDerivMatrixOp(const int i,
1846  const Array<OneD, const NekDouble> &inarray,
1847  Array<OneD,NekDouble> &outarray,
1848  const StdMatrixKey &mkey);
1849 
1851  Array<OneD,NekDouble> &outarray,
1852  const StdMatrixKey &mkey);
1853 
1855  Array<OneD,NekDouble> &outarray,
1856  const StdMatrixKey &mkey);
1857 
1859  const NekDouble> &inarray,
1860  Array<OneD,NekDouble> &outarray,
1861  const StdMatrixKey &mkey,
1862  bool addDiffusionTerm=true);
1863 
1865  Array<OneD,NekDouble> &outarray,
1866  const StdMatrixKey &mkey);
1867 
1869  Array<OneD,NekDouble> &outarray,
1870  const StdMatrixKey &mkey);
1871 
1873  const Array<OneD, const NekDouble> &inarray,
1874  Array<OneD, NekDouble> &outarray,
1875  Array<OneD, NekDouble> &wsp);
1876 
1878  Array<OneD,NekDouble> &outarray,
1879  const StdMatrixKey &mkey);
1880 
1881  STD_REGIONS_EXPORT virtual const NormalVector & v_GetEdgeNormal(const int edge) const;
1882 
1883  STD_REGIONS_EXPORT virtual void v_ComputeEdgeNormal(const int edge);
1884 
1885  STD_REGIONS_EXPORT virtual void v_NegateEdgeNormal(const int edge);
1886 
1887  STD_REGIONS_EXPORT virtual bool v_EdgeNormalNegated(const int edge);
1888 
1889  STD_REGIONS_EXPORT virtual void v_ComputeFaceNormal(const int face);
1890 
1891  STD_REGIONS_EXPORT virtual void v_NegateFaceNormal(const int face);
1892 
1893  STD_REGIONS_EXPORT virtual bool v_FaceNormalNegated(const int face);
1894 
1895  STD_REGIONS_EXPORT virtual const NormalVector & v_GetVertexNormal(const int vertex) const;
1896 
1897  STD_REGIONS_EXPORT virtual void v_ComputeVertexNormal(const int vertex);
1898 
1899  STD_REGIONS_EXPORT virtual void v_NegateVertexNormal(const int vertex);
1900 
1901  STD_REGIONS_EXPORT virtual bool v_VertexNormalNegated(const int vertex);
1902 
1903  STD_REGIONS_EXPORT virtual const NormalVector & v_GetFaceNormal(const int face) const;
1904  STD_REGIONS_EXPORT virtual const NormalVector &
1905  v_GetSurfaceNormal(const int id) const;
1906 
1908  v_GetEdgeInverseBoundaryMap(int eid);
1909 
1911  v_GetFaceInverseBoundaryMap(int fid, StdRegions::Orientation faceOrient = eNoOrientation, int P1=-1, int P2=-1);
1912 
1917 
1919 
1921  Array<OneD, int> &conn,
1922  bool standard = true);
1923  };
1924 
1925 
1926  typedef std::shared_ptr<StdExpansion> StdExpansionSharedPtr;
1927  typedef std::vector< StdExpansionSharedPtr > StdExpansionVector;
1928 
1929  /**
1930  * This function is a wrapper around the virtual function
1931  * \a v_FwdTrans()
1932  *
1933  * Given a function evaluated at the quadrature points, this
1934  * function calculates the expansion coefficients such that the
1935  * resulting expansion approximates the original function.
1936  *
1937  * The calculation of the expansion coefficients is done using a
1938  * Galerkin projection. This is equivalent to the operation:
1939  * \f[ \mathbf{\hat{u}} = \mathbf{M}^{-1} \mathbf{I}\f]
1940  * where
1941  * - \f$\mathbf{M}[p][q]= \int\phi_p(\mathbf{\xi})\phi_q(
1942  * \mathbf{\xi}) d\mathbf{\xi}\f$ is the Mass matrix
1943  * - \f$\mathbf{I}[p] = \int\phi_p(\mathbf{\xi}) u(\mathbf{\xi})
1944  * d\mathbf{\xi}\f$
1945  *
1946  * This function takes the array \a inarray as the values of the
1947  * function evaluated at the quadrature points
1948  * (i.e. \f$\mathbf{u}\f$),
1949  * and stores the resulting coefficients \f$\mathbf{\hat{u}}\f$
1950  * in the \a outarray
1951  *
1952  * @param inarray array of the function discretely evaluated at the
1953  * quadrature points
1954  *
1955  * @param outarray array of the function coefficieints
1956  */
1958  Array<OneD, NekDouble> &outarray)
1959  {
1960  v_FwdTrans(inarray,outarray);
1961  }
1962 
1963  } //end of namespace
1964 } //end of namespace
1965 
1966 #endif //STANDARDDEXPANSION_H
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
#define STD_REGIONS_EXPORT
Describes the specification for a Basis.
Definition: Basis.h:50
Defines a specification for a set of points.
Definition: Points.h:60
The base class for all shapes.
Definition: StdExpansion.h:69
void PhysDeriv(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:812
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false)
virtual void v_GetFaceToElementMap(const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int nummodesA=-1, int nummodesB=-1)
virtual ~StdExpansion()
Destructor.
virtual const NormalVector & v_GetEdgeNormal(const int edge) const
virtual void v_PhysDeriv_n(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
int GetNtrace() const
Returns the number of trace elements connected to this element.
Definition: StdExpansion.h:449
virtual void v_PhysDirectionalDeriv(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
Physical derivative along a direction vector.
StdExpansion()
Default Constructor.
virtual StdRegions::Orientation v_GetForient(int face)
virtual void v_GetTracePhysVals(const int edge, const std::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=eNoOrientation)
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:130
void SetPhysNormals(Array< OneD, const NekDouble > &normal)
Definition: StdExpansion.h:734
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
virtual void v_GetCoord(const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
const LibUtilities::BasisSharedPtr & GetBasis(int dir) const
This function gets the shared point to basis in the dir direction.
Definition: StdExpansion.h:117
void GeneralMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_ComputeEdgeNormal(const int edge)
void FillMode(const int mode, Array< OneD, NekDouble > &outarray)
This function fills the array outarray with the mode-th mode of the expansion.
Definition: StdExpansion.h:597
void FwdTrans_BndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:543
virtual NekDouble v_StdPhysEvaluate(const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
const LibUtilities::PointsKey GetNodalPointsKey() const
This function returns the type of expansion Nodal point type if defined.
Definition: StdExpansion.h:424
void PhysInterpToSimplexEquiSpaced(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int npset=-1)
This function performs an interpolation from the physical space points provided at input into an arra...
virtual void v_AddRobinEdgeContribution(const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, Array< OneD, NekDouble > &coeffs)
int GetNedges() const
This function returns the number of edges of the expansion domain.
Definition: StdExpansion.h:271
void GetEdgePhysMap(const int edge, Array< OneD, int > &outarray)
Definition: StdExpansion.h:931
virtual void v_GetFaceInteriorMap(const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
int GetFaceNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th face.
Definition: StdExpansion.h:353
void PhysDeriv_s(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
virtual int v_GetTotalFaceIntNcoeffs() const
void WeakDirectionalDerivMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void WeakDerivMatrixOp_MatFree(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
DNekBlkMatSharedPtr CreateStdStaticCondMatrix(const StdMatrixKey &mkey)
Create the static condensation of a matrix when using a boundary interior decomposition.
void GetInverseBoundaryMaps(Array< OneD, unsigned int > &vmap, Array< OneD, Array< OneD, unsigned int > > &emap, Array< OneD, Array< OneD, unsigned int > > &fmap)
void DropLocStaticCondMatrix(const LocalRegions::MatrixKey &mkey)
Definition: StdExpansion.h:766
virtual void v_SetPhysNormals(Array< OneD, const NekDouble > &normal)
void NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:751
void LaplacianMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_GetEdgeInterpVals(const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:164
void AddRobinEdgeContribution(const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, Array< OneD, NekDouble > &coeffs)
virtual void v_IProductWRTDirectionalDerivBase_SumFac(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual int v_DetCartesianDirOfEdge(const int edge)
void HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual int v_GetEdgeNcoeffs(const int i) const
virtual int v_GetFaceNumPoints(const int i) const
virtual int v_GetNverts() const =0
void GetEdgeInterpVals(const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:901
void BwdTrans_MatOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_GetFacePhysMap(const int face, Array< OneD, int > &outarray)
void MassLevelCurvatureMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2)
void LinearAdvectionDiffusionReactionMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
virtual std::shared_ptr< StdExpansion > v_GetStdExp(void) const
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:714
void LaplacianMatrixOp(const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
Array< OneD, unsigned int > GetFaceInverseBoundaryMap(int fid, StdRegions::Orientation faceOrient=eNoOrientation, int P1=-1, int P2=-1)
const LibUtilities::PointsKeyVector GetPointsKeys() const
virtual void v_AddRobinMassMatrix(const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
void GetFaceInteriorMap(const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
Definition: StdExpansion.h:842
virtual bool v_IsBoundaryInteriorExpansion()
virtual void v_NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
virtual int v_GetTraceNcoeffs(const int i) const
virtual void v_GetEdgePhysMap(const int edge, Array< OneD, int > &outarray)
void NegateVertexNormal(const int vertex)
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
LibUtilities::NekManager< IndexMapKey, IndexMapValues, IndexMapKey::opLess > m_IndexMapManager
virtual void v_GetSimplexEquiSpacedConnectivity(Array< OneD, int > &conn, bool standard=true)
const LibUtilities::BasisKey DetEdgeBasisKey(const int i) const
Definition: StdExpansion.h:318
void GetSimplexEquiSpacedConnectivity(Array< OneD, int > &conn, bool standard=true)
This function provides the connectivity of local simplices (triangles or tets) to connect the equispa...
DNekBlkMatSharedPtr GetStdStaticCondMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:719
void AddRobinMassMatrix(const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
virtual void v_ComputeFaceNormal(const int face)
virtual void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
virtual std::shared_ptr< StdExpansion > v_GetLinStdExp(void) const
void WeakDerivMatrixOp(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
int EvalBasisNumModesMax(void) const
This function returns the maximum number of expansion modes over all local directions.
Definition: StdExpansion.h:190
virtual void v_GetInteriorMap(Array< OneD, unsigned int > &outarray)
virtual void v_SetCoeffsToOrientation(StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_GetEdgeQFactors(const int edge, Array< OneD, NekDouble > &outarray)
NekDouble Linf(const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
Function to evaluate the discrete error where is given by the array sol.
void PhysDirectionalDeriv(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
void EquiSpacedToCoeffs(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs a projection/interpolation from the equispaced points sometimes used in post-p...
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:974
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void GetTracePhysVals(const int edge, const std::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:889
int CalcNumberOfCoefficients(const std::vector< unsigned int > &nummodes, int &modes_offset)
Definition: StdExpansion.h:796
void GetCoord(const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
given the coordinates of a point of the element in the local collapsed coordinate system,...
Definition: StdExpansion.h:708
virtual void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void ComputeEdgeNormal(const int edge)
void StdPhysDeriv(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void PhysDeriv_n(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
virtual void SetUpPhysNormals(const int edge)
virtual void v_WeakDirectionalDerivMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase() const
This function gets the shared point to basis.
Definition: StdExpansion.h:106
void BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix(const LocalRegions::MatrixKey &mkey)
Definition: StdExpansion.h:761
virtual NekDouble v_PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
virtual DNekMatSharedPtr v_BuildInverseTransformationMatrix(const DNekScalMatSharedPtr &m_transformationmatrix)
NekDouble L2(const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
Function to evaluate the discrete error, where is given by the array sol.
virtual void v_GetVertexPhysVals(const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray)
LibUtilities::BasisType GetEdgeBasisType(const int i) const
This function returns the type of expansion basis on the i-th edge.
Definition: StdExpansion.h:412
void GetFaceToElementMap(const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int nummodesA=-1, int nummodesB=-1)
Definition: StdExpansion.h:858
virtual int v_GetShapeDimension() const
void SetElmtId(const int id)
Set the element id of this expansion when used in a list by returning value of m_elmt_id.
Definition: StdExpansion.h:676
bool EdgeNormalNegated(const int edge)
void IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
Definition: StdExpansion.h:822
void SetCoeffsToOrientation(Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
Definition: StdExpansion.h:781
StdRegions::Orientation GetEorient(int edge)
Definition: StdExpansion.h:776
void NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:741
virtual bool v_VertexNormalNegated(const int vertex)
virtual bool v_FaceNormalNegated(const int face)
void GetFacePhysVals(const int face, const std::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=eNoOrientation)
Definition: StdExpansion.h:921
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
IndexMapValuesSharedPtr GetIndexMap(const IndexMapKey &ikey)
Definition: StdExpansion.h:724
bool VertexNormalNegated(const int vertex)
virtual int v_NumBndryCoeffs() const
virtual void v_MassLevelCurvatureMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual int v_GetElmtId()
Get the element id of this expansion when used in a list by returning value of m_elmt_id.
virtual int v_GetEdgeNumPoints(const int i) const
int GetFaceIntNcoeffs(const int i) const
Definition: StdExpansion.h:358
void ComputeFaceNormal(const int face)
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
this function calculates the inner product of a given function f with the different modes of the expa...
Definition: StdExpansion.h:634
void GetEdgePhysVals(const int edge, const std::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:881
virtual void v_MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int GetNverts() const
This function returns the number of vertices of the expansion domain.
Definition: StdExpansion.h:258
std::shared_ptr< T > as()
virtual LibUtilities::ShapeType v_DetShapeType() const
const NormalVector & GetEdgeNormal(const int edge) const
virtual void v_NegateEdgeNormal(const int edge)
virtual bool v_EdgeNormalNegated(const int edge)
void LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:981
void SetCoeffsToOrientation(StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:788
virtual void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray)
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_GetEdgeToElementMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-1)
Array< OneD, unsigned int > GetEdgeInverseBoundaryMap(int eid)
void NormVectorIProductWRTBase(const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:756
virtual const NormalVector & v_GetSurfaceNormal(const int id) const
const NormalVector & GetVertexNormal(const int vertex) const
void GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2=NullNekDouble1DArray, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
this function returns the physical coordinates of the quadrature points of the expansion
Definition: StdExpansion.h:690
int DetCartesianDirOfEdge(const int edge)
Definition: StdExpansion.h:313
void IProductWRTDirectionalDerivBase(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:658
virtual void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int GetElmtId()
Get the element id of this expansion when used in a list by returning value of m_elmt_id.
Definition: StdExpansion.h:668
virtual void v_ExponentialFilter(Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
virtual const NormalVector & v_GetVertexNormal(const int vertex) const
virtual StdRegions::Orientation v_GetEorient(int edge)
virtual void v_NegateFaceNormal(const int face)
NekDouble PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
This function evaluates the expansion at a single (arbitrary) point of the domain.
virtual const Array< OneD, const NekDouble > & v_GetPhysNormals(void)
void IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:650
void MassLevelCurvatureMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:469
virtual void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_SetUpPhysNormals(const int edge)
virtual int v_GetTotalEdgeIntNcoeffs() const
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
void GetEdgeToElementMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-1)
Definition: StdExpansion.h:849
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:817
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space.
Definition: StdExpansion.h:530
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey)
virtual void v_GetInverseBoundaryMaps(Array< OneD, unsigned int > &vmap, Array< OneD, Array< OneD, unsigned int > > &emap, Array< OneD, Array< OneD, unsigned int > > &fmap)
void GetFaceNumModes(const int fid, const Orientation faceOrient, int &numModes0, int &numModes1)
Definition: StdExpansion.h:835
int GetTraceNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th edge/face.
Definition: StdExpansion.h:378
virtual Array< OneD, unsigned int > v_GetEdgeInverseBoundaryMap(int eid)
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)
virtual const LibUtilities::BasisKey v_DetFaceBasisKey(const int i, const int k) const
int GetEdgeNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th edge.
Definition: StdExpansion.h:286
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual int v_NumDGBndryCoeffs() const
virtual void v_FwdTrans_BndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_GetFacePhysVals(const int face, const std::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient)
int GetEdgeNumPoints(const int i) const
This function returns the number of quadrature points belonging to the i-th edge.
Definition: StdExpansion.h:307
virtual void v_PhysDeriv_s(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
virtual const LibUtilities::BasisKey v_DetEdgeBasisKey(const int i) const
void ExponentialFilter(Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
int GetNumBases() const
This function returns the number of 1D bases used in the expansion.
Definition: StdExpansion.h:97
DNekMatSharedPtr BuildInverseTransformationMatrix(const DNekScalMatSharedPtr &m_transformationmatrix)
IndexMapValuesSharedPtr CreateIndexMap(const IndexMapKey &ikey)
Create an IndexMap which contains mapping information linking any specific element shape with either ...
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:215
void NegateFaceNormal(const int face)
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space.
virtual void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray)
virtual LibUtilities::PointsKey v_GetFacePointsKey(const int i, const int j) const
void ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:988
virtual void v_GetEdgeInteriorMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
virtual Array< OneD, unsigned int > v_GetFaceInverseBoundaryMap(int fid, StdRegions::Orientation faceOrient=eNoOrientation, int P1=-1, int P2=-1)
virtual void v_MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLess > m_stdStaticCondMatrixManager
virtual const LibUtilities::PointsKey v_GetNodalPointsKey() const
virtual void v_WeakDerivMatrixOp(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LinearAdvectionDiffusionReactionMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:228
void IProductWRTDirectionalDerivBase_SumFac(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:945
const NormalVector & GetSurfaceNormal(const int id) const
virtual void v_LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
virtual void v_NegateVertexNormal(const int vertex)
LibUtilities::PointsKey GetFacePointsKey(const int i, const int j) const
Definition: StdExpansion.h:384
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)=0
Calculates the inner product of a given function f with the different modes of the expansion.
void StdPhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
void GetEdgeInteriorMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
Definition: StdExpansion.h:828
void MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:952
int GetFaceNumPoints(const int i) const
This function returns the number of quadrature points belonging to the i-th face.
Definition: StdExpansion.h:338
void NegateEdgeNormal(const int edge)
virtual void v_LinearAdvectionDiffusionReactionMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
void HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
NekDouble PhysEvaluate(const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
This function evaluates the expansion at a single (arbitrary) point of the domain.
virtual void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)=0
Transform a given function from physical quadrature space to coefficient space.
virtual void v_DropLocStaticCondMatrix(const LocalRegions::MatrixKey &mkey)
const LibUtilities::BasisKey DetFaceBasisKey(const int i, const int k) const
Definition: StdExpansion.h:323
int GetNfaces() const
This function returns the number of faces of the expansion domain.
Definition: StdExpansion.h:437
StdRegions::Orientation GetForient(int face)
Definition: StdExpansion.h:771
void ComputeVertexNormal(const int vertex)
void WeakDirectionalDerivMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_LaplacianMatrixOp_MatFree_Kernel(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
NekDouble H1(const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
Function to evaluate the discrete error, where is given by the array sol.
virtual int v_CalcNumberOfCoefficients(const std::vector< unsigned int > &nummodes, int &modes_offset)
void GetVertexPhysVals(const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray)
Definition: StdExpansion.h:894
void IProductWRTBase(const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
Definition: StdExpansion.h:640
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
Integrates the specified function over the domain.
void GetFacePhysMap(const int face, Array< OneD, int > &outarray)
Definition: StdExpansion.h:938
void GetEdgePhysVals(const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Extract the physical values along edge edge from inarray into outarray following the local edge orien...
Definition: StdExpansion.h:874
virtual int v_GetFaceNcoeffs(const int i) const
virtual const NormalVector & v_GetFaceNormal(const int face) const
virtual LibUtilities::BasisType v_GetEdgeBasisType(const int i) const
void SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
Definition: StdExpansion.h:995
virtual void v_GetEdgePhysVals(const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Extract the physical values along edge edge from inarray into outarray following the local edge orien...
const Array< OneD, const NekDouble > & GetPhysNormals(void)
Definition: StdExpansion.h:729
void HelmholtzMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)=0
LibUtilities::NekManager< StdMatrixKey, DNekMat, StdMatrixKey::opLess > m_stdMatrixManager
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:177
void PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void LaplacianMatrixOp_MatFree_Kernel(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
std::shared_ptr< StdExpansion > GetStdExp(void) const
Definition: StdExpansion.h:474
const Array< OneD, const NekDouble > & GetPoints(const int dir) const
This function returns a pointer to the array containing the quadrature points in dir direction.
Definition: StdExpansion.h:242
void GetEdgeQFactors(const int edge, Array< OneD, NekDouble > &outarray)
Extract the metric factors to compute the contravariant fluxes along edge edge and stores them into o...
Definition: StdExpansion.h:914
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
DNekMatSharedPtr CreateStdMatrix(const StdMatrixKey &mkey)
virtual void v_ComputeVertexNormal(const int vertex)
virtual void v_IProductWRTDirectionalDerivBase(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, const Array< OneD, NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:746
virtual int v_GetFaceIntNcoeffs(const int i) const
virtual void v_GetFaceNumModes(const int fid, const Orientation faceOrient, int &numModes0, int &numModes1)
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2, Array< OneD, NekDouble > &out_d3)
Calculate the derivative of the physical points.
virtual void v_StdPhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2, Array< OneD, NekDouble > &out_d3)
const NormalVector & GetFaceNormal(const int face) const
NekDouble StdPhysEvaluate(const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
void MassMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void GeneralMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix(const LocalRegions::MatrixKey &mkey)
std::shared_ptr< StdExpansion > GetLinStdExp(void) const
Definition: StdExpansion.h:480
NekDouble Integral(const Array< OneD, const NekDouble > &inarray)
This function integrates the specified function over the domain.
Definition: StdExpansion.h:579
bool FaceNormalNegated(const int face)
std::shared_ptr< Basis > BasisSharedPtr
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:246
std::vector< StdExpansionSharedPtr > StdExpansionVector
std::shared_ptr< IndexMapValues > IndexMapValuesSharedPtr
Definition: IndexMapKey.h:124
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
Array< OneD, Array< OneD, NekDouble > > NormalVector
Definition: StdExpansion.h:58
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:71
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:73
static Array< OneD, NekDouble > NullNekDouble1DArray
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
double NekDouble
P
Definition: main.py:133