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