Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator 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 = 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  int GetShapeDimension() const
482  {
483  return v_GetShapeDimension();
484  }
485 
487  {
489  }
490 
491 
493  {
494  return v_IsNodalNonTensorialExp();
495  }
496 
497  /** \brief This function performs the Backward transformation from
498  * coefficient space to physical space
499  *
500  * This function is a wrapper around the virtual function
501  * \a v_BwdTrans()
502  *
503  * Based on the expansion coefficients, this function evaluates the
504  * expansion at the quadrature points. This is equivalent to the
505  * operation \f[ u(\xi_{1i}) =
506  * \sum_{p=0}^{P-1} \hat{u}_p \phi_p(\xi_{1i}) \f] which can be
507  * evaluated as \f$ {\bf u} = {\bf B}^T {\bf \hat{u}} \f$ with
508  * \f${\bf B}[i][j] = \phi_i(\xi_{j})\f$
509  *
510  * This function requires that the coefficient array
511  * \f$\mathbf{\hat{u}}\f$ provided as \a inarray.
512  *
513  * The resulting array
514  * \f$\mathbf{u}[m]=u(\mathbf{\xi}_m)\f$ containing the
515  * expansion evaluated at the quadrature points, is stored
516  * in the \a outarray.
517  *
518  * \param inarray contains the values of the expansion
519  * coefficients (input of the function)
520  *
521  * \param outarray contains the values of the expansion evaluated
522  * at the quadrature points (output of the function)
523  */
524 
526  Array<OneD, NekDouble> &outarray)
527  {
528  v_BwdTrans (inarray, outarray);
529  }
530 
531  /**
532  * @brief This function performs the Forward transformation from
533  * physical space to coefficient space.
534  */
535  inline void FwdTrans (const Array<OneD, const NekDouble>& inarray,
536  Array<OneD, NekDouble> &outarray);
537 
539  Array<OneD, NekDouble> &outarray)
540  {
541  v_FwdTrans_BndConstrained(inarray,outarray);
542  }
543 
544  /** \brief This function integrates the specified function over the
545  * domain
546  *
547  * This function is a wrapper around the virtual function
548  * \a v_Integral()
549  *
550  * Based on the values of the function evaluated at the quadrature
551  * points (which are stored in \a inarray), this function calculates
552  * the integral of this function over the domain. This is
553  * equivalent to the numerical evaluation of the operation
554  * \f[ I=\int u(\mathbf{\xi})d \mathbf{\xi}\f]
555  *
556  * \param inarray values of the function to be integrated evaluated
557  * at the quadrature points (i.e.
558  * \a inarray[m]=\f$u(\mathbf{\xi}_m)\f$)
559  * \return returns the value of the calculated integral
560  *
561  * Inputs:\n
562 
563  - \a inarray: definition of function to be returned at quadrature point
564  of expansion.
565 
566  Outputs:\n
567 
568  - returns \f$\int^1_{-1}\int^1_{-1} u(\xi_1, \xi_2) J[i,j] d
569  \xi_1 d \xi_2 \f$ where \f$inarray[i,j] =
570  u(\xi_{1i},\xi_{2j}) \f$ and \f$ J[i,j] \f$ is the
571  Jacobian evaluated at the quadrature point.
572  *
573  */
575  {
576  return v_Integral(inarray);
577  }
578 
579  /** \brief This function fills the array \a outarray with the
580  * \a mode-th mode of the expansion
581  *
582  * This function is a wrapper around the virtual function
583  * \a v_FillMode()
584  *
585  * The requested mode is evaluated at the quadrature points
586  *
587  * \param mode the mode that should be filled
588  * \param outarray contains the values of the \a mode-th mode of the
589  * expansion evaluated at the quadrature points (output of the
590  * function)
591  */
592  void FillMode(const int mode, Array<OneD, NekDouble> &outarray)
593  {
594  v_FillMode(mode, outarray);
595  }
596 
597  /** \brief this function calculates the inner product of a given
598  * function \a f with the different modes of the expansion
599  *
600  * This function is a wrapper around the virtual function
601  * \a v_IProductWRTBase()
602  *
603  * This is equivalent to the numerical evaluation of
604  * \f[ I[p] = \int \phi_p(\mathbf{x}) f(\mathbf{x}) d\mathbf{x}\f]
605  * \f$ \begin{array}{rcl} I_{pq} = (\phi_q \phi_q, u) & = &
606  \sum_{i=0}^{nq_0} \sum_{j=0}^{nq_1} \phi_p(\xi_{0,i})
607  \phi_q(\xi_{1,j}) w^0_i w^1_j u(\xi_{0,i} \xi_{1,j})
608  J_{i,j}\\ & = & \sum_{i=0}^{nq_0} \phi_p(\xi_{0,i})
609  \sum_{j=0}^{nq_1} \phi_q(\xi_{1,j}) \tilde{u}_{i,j}
610  J_{i,j} \end{array} \f$
611 
612  where
613 
614  \f$ \tilde{u}_{i,j} = w^0_i w^1_j u(\xi_{0,i},\xi_{1,j}) \f$
615 
616  which can be implemented as
617 
618  \f$ f_{qi} = \sum_{j=0}^{nq_1} \phi_q(\xi_{1,j}) \tilde{u}_{i,j} =
619  {\bf B_1 U} \f$
620  \f$ I_{pq} = \sum_{i=0}^{nq_0} \phi_p(\xi_{0,i}) f_{qi} =
621  {\bf B_0 F} \f$
622  *
623  * \param inarray contains the values of the function \a f
624  * evaluated at the quadrature points
625  * \param outarray contains the values of the inner product of \a f
626  * with the different modes, i.e. \f$ outarray[p] = I[p]\f$
627  * (output of the function)
628  */
630  Array<OneD, NekDouble> &outarray)
631  {
632  v_IProductWRTBase(inarray, outarray);
633  }
634 
636  const Array<OneD, const NekDouble>& base,
637  const Array<OneD, const NekDouble>& inarray,
638  Array<OneD, NekDouble> &outarray,
639  int coll_check)
640  {
641  v_IProductWRTBase(base, inarray, outarray, coll_check);
642  }
643 
644 
645  void IProductWRTDerivBase(const int dir,
646  const Array<OneD, const NekDouble>& inarray,
647  Array<OneD, NekDouble> &outarray)
648  {
649  v_IProductWRTDerivBase(dir,inarray, outarray);
650  }
651 
652  /// \brief Get the element id of this expansion when used
653  /// in a list by returning value of #m_elmt_id
654  inline int GetElmtId()
655  {
656  return m_elmt_id;
657  }
658 
659 
660  /// \brief Set the element id of this expansion when used
661  /// in a list by returning value of #m_elmt_id
662  inline void SetElmtId(const int id)
663  {
664  m_elmt_id = id;
665  }
666 
667  /** \brief this function returns the physical coordinates of the
668  * quadrature points of the expansion
669  *
670  * This function is a wrapper around the virtual function
671  * \a v_GetCoords()
672  *
673  * \param coords an array containing the coordinates of the
674  * quadrature points (output of the function)
675  */
679  {
680  v_GetCoords(coords_1,coords_2,coords_3);
681  }
682 
683  /** \brief given the coordinates of a point of the element in the
684  * local collapsed coordinate system, this function calculates the
685  * physical coordinates of the point
686  *
687  * This function is a wrapper around the virtual function
688  * \a v_GetCoord()
689  *
690  * \param Lcoords the coordinates in the local collapsed
691  * coordinate system
692  * \param coords the physical coordinates (output of the function)
693  */
695  Array<OneD, NekDouble> &coord)
696  {
697  v_GetCoord(Lcoord, coord);
698  }
699 
701  {
702  return m_stdMatrixManager[mkey];
703  }
704 
706  {
707  return m_stdStaticCondMatrixManager[mkey];
708  }
709 
711  {
712  return m_IndexMapManager[ikey];
713  }
714 
716  {
717  return v_GetPhysNormals();
718  }
719 
721  {
722  v_SetPhysNormals(normal);
723  }
724 
725  STD_REGIONS_EXPORT virtual void SetUpPhysNormals(const int edge);
726 
728  {
729  v_NormVectorIProductWRTBase(Fx,outarray);
730  }
731 
733  {
734  v_NormVectorIProductWRTBase(Fx,Fy,outarray);
735  }
736 
738  {
739  v_NormVectorIProductWRTBase(Fx,Fy,Fz,outarray);
740  }
741 
743  {
744  v_NormVectorIProductWRTBase(Fvec, outarray);
745  }
746 
748  {
749  return v_GetLocStaticCondMatrix(mkey);
750  }
751 
753  {
754  return v_DropLocStaticCondMatrix(mkey);
755  }
756 
758  {
759  return v_GetForient(face);
760  }
761 
763  {
764  return v_GetEorient(edge);
765  }
766 
768  {
769  return v_GetPorient(point);
770  }
771 
773  {
774  return v_GetCartesianEorient(edge);
775  }
776 
778  Array<OneD, NekDouble> &coeffs,
780  {
781  v_SetCoeffsToOrientation(coeffs, dir);
782  }
783 
787  Array<OneD, NekDouble> &outarray)
788  {
789  v_SetCoeffsToOrientation(dir,inarray,outarray);
790  }
791 
792  int CalcNumberOfCoefficients(const std::vector<unsigned int> &nummodes, int &modes_offset)
793  {
794  return v_CalcNumberOfCoefficients(nummodes,modes_offset);
795  }
796 
797  void ExtractDataToCoeffs(const NekDouble *data,
798  const std::vector<unsigned int > &nummodes,
799  const int nmodes_offset,
800  NekDouble *coeffs)
801  {
802  v_ExtractDataToCoeffs(data,nummodes,nmodes_offset,coeffs);
803  }
804 
805  // virtual functions related to LocalRegions
807  const Array<OneD, const NekDouble> &Lcoord,
808  const Array<OneD, const NekDouble> &physvals);
809 
810 
812  {
813  return v_GetCoordim();
814  }
815 
817  {
818  v_GetBoundaryMap(outarray);
819  }
820 
822  {
823  v_GetInteriorMap(outarray);
824  }
825 
826  int GetVertexMap(const int localVertexId,
827  bool useCoeffPacking = false)
828  {
829  return v_GetVertexMap(localVertexId,useCoeffPacking);
830  }
831 
832  void GetEdgeInteriorMap(const int eid, const Orientation edgeOrient,
833  Array<OneD, unsigned int> &maparray,
834  Array<OneD, int> &signarray)
835  {
836  v_GetEdgeInteriorMap(eid,edgeOrient,maparray,signarray);
837  }
838 
839  void GetFaceInteriorMap(const int fid, const Orientation faceOrient,
840  Array<OneD, unsigned int> &maparray,
841  Array<OneD, int> &signarray)
842  {
843  v_GetFaceInteriorMap(fid,faceOrient,maparray,signarray);
844  }
845 
846  void GetEdgeToElementMap(const int eid,
847  const Orientation edgeOrient,
848  Array<OneD, unsigned int> &maparray,
849  Array<OneD, int> &signarray,
850  int P = -1)
851  {
852  v_GetEdgeToElementMap(eid, edgeOrient, maparray, signarray, P);
853  }
854 
855  void GetFaceToElementMap(const int fid, const Orientation faceOrient,
856  Array<OneD, unsigned int> &maparray,
857  Array<OneD, int> &signarray,
858  int nummodesA = -1, int nummodesB = -1)
859  {
860  v_GetFaceToElementMap(fid,faceOrient,maparray,signarray,
861  nummodesA,nummodesB);
862  }
863 
864 
865  /**
866  * @brief Extract the physical values along edge \a edge from \a
867  * inarray into \a outarray following the local edge orientation
868  * and point distribution defined by defined in \a EdgeExp.
869  */
870 
871  void GetEdgePhysVals(const int edge, const Array<OneD,
872  const NekDouble> &inarray,
873  Array<OneD,NekDouble> &outarray)
874  {
875  v_GetEdgePhysVals(edge,inarray,outarray);
876  }
877 
878  void GetEdgePhysVals(const int edge,
879  const boost::shared_ptr<StdExpansion> &EdgeExp,
880  const Array<OneD, const NekDouble> &inarray,
881  Array<OneD,NekDouble> &outarray)
882  {
883  v_GetEdgePhysVals(edge,EdgeExp,inarray,outarray);
884  }
885 
886  void GetTracePhysVals(const int edge, const boost::shared_ptr<StdExpansion> &EdgeExp, const Array<OneD, const NekDouble> &inarray, Array<OneD,NekDouble> &outarray)
887  {
888  v_GetTracePhysVals(edge,EdgeExp,inarray,outarray);
889  }
890 
891  void GetVertexPhysVals(const int vertex,
892  const Array<OneD, const NekDouble> &inarray,
893  NekDouble &outarray)
894  {
895  v_GetVertexPhysVals(vertex, inarray, outarray);
896  }
897 
898  void GetEdgeInterpVals(const int edge,const Array<OneD,
899  const NekDouble> &inarray,
900  Array<OneD,NekDouble> &outarray)
901  {
902  v_GetEdgeInterpVals(edge, inarray, outarray);
903  }
904 
905  /**
906  * @brief Extract the metric factors to compute the contravariant
907  * fluxes along edge \a edge and stores them into \a outarray
908  * following the local edge orientation (i.e. anticlockwise
909  * convention).
910  */
912  const int edge,
913  Array<OneD, NekDouble> &outarray)
914  {
915  v_GetEdgeQFactors(edge, outarray);
916  }
917 
919  const int face,
920  const boost::shared_ptr<StdExpansion> &FaceExp,
921  const Array<OneD, const NekDouble> &inarray,
922  Array<OneD, NekDouble> &outarray,
924  {
925  v_GetFacePhysVals(face, FaceExp, inarray, outarray, orient);
926  }
927 
929  const int edge,
930  Array<OneD, int> &outarray)
931  {
932  v_GetEdgePhysMap(edge, outarray);
933  }
934 
936  const int face,
937  Array<OneD, int> &outarray)
938  {
939  v_GetFacePhysMap(face, outarray);
940  }
941 
943  const Array<OneD, const NekDouble> &inarray,
944  Array<OneD, NekDouble> &outarray)
945  {
946  v_MultiplyByQuadratureMetric(inarray, outarray);
947  }
948 
950  const Array<OneD, const NekDouble> &inarray,
951  Array<OneD, NekDouble> & outarray)
952  {
953  v_MultiplyByStdQuadratureMetric(inarray, outarray);
954  }
955 
956  // Matrix Routines
957 
958  /** \brief this function generates the mass matrix
959  * \f$\mathbf{M}[i][j] =
960  * \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}\f$
961  *
962  * \return returns the mass matrix
963  */
964 
966 
968  Array<OneD,NekDouble> &outarray,
969  const StdMatrixKey &mkey);
970 
972  Array<OneD,NekDouble> &outarray,
973  const StdMatrixKey &mkey)
974  {
975  v_MassMatrixOp(inarray,outarray,mkey);
976  }
977 
979  Array<OneD,NekDouble> &outarray,
980  const StdMatrixKey &mkey)
981  {
982  v_LaplacianMatrixOp(inarray,outarray,mkey);
983  }
984 
985  void ReduceOrderCoeffs(int numMin,
986  const Array<OneD, const NekDouble> &inarray,
987  Array<OneD,NekDouble> &outarray)
988  {
989  v_ReduceOrderCoeffs(numMin,inarray,outarray);
990  }
991 
993  const StdMatrixKey &mkey)
994  {
995  v_SVVLaplacianFilter(array,mkey);
996  }
997 
998  void LaplacianMatrixOp(const int k1, const int k2,
999  const Array<OneD, const NekDouble> &inarray,
1000  Array<OneD,NekDouble> &outarray,
1001  const StdMatrixKey &mkey)
1002  {
1003  v_LaplacianMatrixOp(k1,k2,inarray,outarray,mkey);
1004  }
1005 
1006  void WeakDerivMatrixOp(const int i,
1007  const Array<OneD, const NekDouble> &inarray,
1008  Array<OneD,NekDouble> &outarray,
1009  const StdMatrixKey &mkey)
1010  {
1011  v_WeakDerivMatrixOp(i,inarray,outarray,mkey);
1012  }
1013 
1015  Array<OneD,NekDouble> &outarray,
1016  const StdMatrixKey &mkey)
1017  {
1018  v_WeakDirectionalDerivMatrixOp(inarray,outarray,mkey);
1019  }
1020 
1022  Array<OneD,NekDouble> &outarray,
1023  const StdMatrixKey &mkey)
1024  {
1025  v_MassLevelCurvatureMatrixOp(inarray,outarray,mkey);
1026  }
1027 
1029  Array<OneD,NekDouble> &outarray,
1030  const StdMatrixKey &mkey,
1031  bool addDiffusionTerm = true)
1032  {
1033  v_LinearAdvectionDiffusionReactionMatrixOp(inarray,outarray,mkey,addDiffusionTerm);
1034  }
1035 
1036  /**
1037  * @param inarray Input array @f$ \mathbf{u} @f$.
1038  * @param outarray Output array @f$ \boldsymbol{\nabla^2u}
1039  * + \lambda \boldsymbol{u} @f$.
1040  * @param mkey
1041  */
1043  Array<OneD,NekDouble> &outarray,
1044  const StdMatrixKey &mkey)
1045  {
1046  v_HelmholtzMatrixOp(inarray,outarray,mkey);
1047  }
1048 
1050  {
1051  return v_GenMatrix(mkey);
1052  }
1053 
1055  Array<OneD, NekDouble> &out_d0,
1058  {
1059  v_PhysDeriv (inarray, out_d0, out_d1, out_d2);
1060  }
1061 
1062  void PhysDeriv(const int dir,
1063  const Array<OneD, const NekDouble>& inarray,
1064  Array<OneD, NekDouble> &outarray)
1065  {
1066  v_PhysDeriv (dir, inarray, outarray);
1067  }
1068 
1070  Array<OneD, NekDouble> &out_ds)
1071  {
1072  v_PhysDeriv_s(inarray,out_ds);
1073  }
1074 
1076  Array<OneD, NekDouble>& out_dn)
1077  {
1078  v_PhysDeriv_n(inarray,out_dn);
1079  }
1080 
1082  const Array<OneD, const NekDouble>& direction,
1083  Array<OneD, NekDouble> &outarray)
1084  {
1085  v_PhysDirectionalDeriv (inarray, direction, outarray);
1086  }
1087 
1089  Array<OneD, NekDouble> &out_d0,
1092  {
1093  v_StdPhysDeriv(inarray, out_d0, out_d1, out_d2);
1094  }
1095 
1096  void StdPhysDeriv (const int dir,
1097  const Array<OneD, const NekDouble>& inarray,
1098  Array<OneD, NekDouble> &outarray)
1099  {
1100  v_StdPhysDeriv(dir,inarray,outarray);
1101  }
1102 
1103  void AddRobinMassMatrix(const int edgeid, const Array<OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
1104  {
1105  v_AddRobinMassMatrix(edgeid,primCoeffs,inoutmat);
1106  }
1107 
1108  void AddRobinEdgeContribution(const int edgeid, const Array<OneD, const NekDouble> &primCoeffs, Array<OneD, NekDouble> &coeffs)
1109  {
1110  v_AddRobinEdgeContribution(edgeid, primCoeffs, coeffs);
1111  }
1112 
1113  /** \brief This function evaluates the expansion at a single
1114  * (arbitrary) point of the domain
1115  *
1116  * This function is a wrapper around the virtual function
1117  * \a v_PhysEvaluate()
1118  *
1119  * Based on the value of the expansion at the quadrature
1120  * points provided in \a physvals, this function
1121  * calculates the value of the expansion at an arbitrary
1122  * single points (with coordinates \f$ \mathbf{x_c}\f$
1123  * given by the pointer \a coords). This operation,
1124  * equivalent to \f[ u(\mathbf{x_c}) = \sum_p
1125  * \phi_p(\mathbf{x_c}) \hat{u}_p \f] is evaluated using
1126  * Lagrangian interpolants through the quadrature points:
1127  * \f[ u(\mathbf{x_c}) = \sum_p h_p(\mathbf{x_c}) u_p\f]
1128  *
1129  * \param coords the coordinates of the single point
1130  * \param physvals the interpolated field at the quadrature points
1131  *
1132  * \return returns the value of the expansion at the
1133  * single point
1134  */
1136  const Array<OneD, const NekDouble>& physvals)
1137  {
1138  return v_PhysEvaluate(coords,physvals);
1139  }
1140 
1141 
1142  /** \brief This function evaluates the expansion at a single
1143  * (arbitrary) point of the domain
1144  *
1145  * This function is a wrapper around the virtual function
1146  * \a v_PhysEvaluate()
1147  *
1148  * Based on the value of the expansion at the quadrature
1149  * points provided in \a physvals, this function
1150  * calculates the value of the expansion at an arbitrary
1151  * single points associated with the interpolation
1152  * matrices provided in \f$ I \f$.
1153  *
1154  * \param I an Array of lagrange interpolantes evaluated
1155  * at the coordinate and going through the local physical
1156  * quadrature
1157  * \param physvals the interpolated field at the quadrature points
1158  *
1159  * \return returns the value of the expansion at the
1160  * single point
1161  */
1163  const Array<OneD, const NekDouble >& physvals)
1164  {
1165  return v_PhysEvaluate(I,physvals);
1166  }
1167 
1168 
1169  /**
1170  * \brief Convert local cartesian coordinate \a xi into local
1171  * collapsed coordinates \a eta
1172  **/
1175  {
1176  v_LocCoordToLocCollapsed(xi,eta);
1177  }
1178 
1179  const boost::shared_ptr<SpatialDomains::GeomFactors>& GetMetricInfo(void) const
1180  {
1181  return v_GetMetricInfo();
1182  }
1183 
1184  /// \brief Get the element id of this expansion when used
1185  /// in a list by returning value of #m_elmt_id
1186  STD_REGIONS_EXPORT virtual int v_GetElmtId();
1187 
1189 
1191 
1192  STD_REGIONS_EXPORT virtual void v_SetUpPhysNormals(const int edge);
1193 
1194  STD_REGIONS_EXPORT virtual int v_CalcNumberOfCoefficients(const std::vector<unsigned int> &nummodes, int &modes_offset);
1195 
1196  /**
1197  * @brief Unpack data from input file assuming it comes from the
1198  * same expansion type.
1199  * @see StdExpansion::ExtractDataToCoeffs
1200  */
1201  STD_REGIONS_EXPORT virtual void v_ExtractDataToCoeffs(const NekDouble *data,
1202  const std::vector<unsigned int > &nummodes,
1203  const int nmode_offset,
1204  NekDouble *coeffs);
1205 
1207 
1209  const Array<OneD, const NekDouble> &Fx,
1210  const Array<OneD, const NekDouble> &Fy,
1211  Array< OneD, NekDouble> &outarray);
1212 
1214 
1216 
1218 
1220 
1221 
1223 
1225 
1227 
1229 
1230  /** \brief Function to evaluate the discrete \f$ L_\infty\f$
1231  * error \f$ |\epsilon|_\infty = \max |u - u_{exact}|\f$ where \f$
1232  * u_{exact}\f$ is given by the array \a sol.
1233  *
1234  * This function takes the physical value space array \a m_phys as
1235  * approximate solution
1236  *
1237  * \param sol array of solution function at physical quadrature
1238  * points
1239  * \return returns the \f$ L_\infty \f$ error as a NekDouble.
1240  */
1242 
1243  /** \brief Function to evaluate the discrete \f$ L_2\f$ error,
1244  * \f$ | \epsilon |_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2
1245  * dx \right]^{1/2} d\xi_1 \f$ where \f$ u_{exact}\f$ is given by
1246  * the array \a sol.
1247  *
1248  * This function takes the physical value space array \a m_phys as
1249  * approximate solution
1250  *
1251  * \param sol array of solution function at physical quadrature
1252  * points
1253  * \return returns the \f$ L_2 \f$ error as a double.
1254  */
1256 
1257  /** \brief Function to evaluate the discrete \f$ H^1\f$
1258  * error, \f$ | \epsilon |^1_{2} = \left [ \int^1_{-1} [u -
1259  * u_{exact}]^2 + \nabla(u - u_{exact})\cdot\nabla(u -
1260  * u_{exact})\cdot dx \right]^{1/2} d\xi_1 \f$ where \f$
1261  * u_{exact}\f$ is given by the array \a sol.
1262  *
1263  * This function takes the physical value space array
1264  * \a m_phys as approximate solution
1265  *
1266  * \param sol array of solution function at physical quadrature
1267  * points
1268  * \return returns the \f$ H_1 \f$ error as a double.
1269  */
1271 
1272  // I/O routines
1273  const NormalVector & GetEdgeNormal(const int edge) const
1274  {
1275  return v_GetEdgeNormal(edge);
1276  }
1277 
1278  void ComputeEdgeNormal(const int edge)
1279  {
1280  v_ComputeEdgeNormal(edge);
1281  }
1282 
1283  void NegateEdgeNormal(const int edge)
1284  {
1285  v_NegateEdgeNormal(edge);
1286  }
1287 
1288  bool EdgeNormalNegated(const int edge)
1289  {
1290  return v_EdgeNormalNegated(edge);
1291  }
1292 
1293  void ComputeFaceNormal(const int face)
1294  {
1295  v_ComputeFaceNormal(face);
1296  }
1297 
1298  void NegateFaceNormal(const int face)
1299  {
1300  v_NegateFaceNormal(face);
1301  }
1302 
1303  bool FaceNormalNegated(const int face)
1304  {
1305  return v_FaceNormalNegated(face);
1306  }
1307 
1308  void ComputeVertexNormal(const int vertex)
1309  {
1310  v_ComputeVertexNormal(vertex);
1311  }
1312 
1313  const NormalVector & GetFaceNormal(const int face) const
1314  {
1315  return v_GetFaceNormal(face);
1316  }
1317 
1318  const NormalVector & GetVertexNormal(const int vertex) const
1319  {
1320  return v_GetVertexNormal(vertex);
1321  }
1322 
1323  const NormalVector & GetSurfaceNormal(const int id) const
1324  {
1325  return v_GetSurfaceNormal(id);
1326  }
1327 
1329  {
1331  for (int i = 0; i < m_base.num_elements(); ++i)
1332  {
1333  p.push_back(m_base[i]->GetPointsKey());
1334  }
1335  return p;
1336  }
1337 
1340  {
1341  return v_GetEdgeInverseBoundaryMap(eid);
1342  }
1343 
1346  {
1347  return v_GetFaceInverseBoundaryMap(fid,faceOrient);
1348  }
1349 
1351  const DNekScalMatSharedPtr & m_transformationmatrix)
1352  {
1354  m_transformationmatrix);
1355  }
1356 
1357 
1358  /** \brief This function performs an interpolation from
1359  * the physical space points provided at input into an
1360  * array of equispaced points which are not the collapsed
1361  * coordinate. So for a tetrahedron you will only get a
1362  * tetrahedral number of values.
1363  *
1364  * This is primarily used for output purposes to get a
1365  * better distribution of points more suitable for most
1366  * postprocessing
1367  */
1369  const Array<OneD, const NekDouble> &inarray,
1370  Array<OneD, NekDouble> &outarray,
1371  int npset = -1);
1372 
1373  /** \brief This function provides the connectivity of
1374  * local simplices (triangles or tets) to connect the
1375  * equispaced data points provided by
1376  * PhysInterpToSimplexEquiSpaced
1377  *
1378  * This is a virtual call to the function
1379  * \a v_GetSimplexEquiSpaceConnectivity
1380  */
1382  Array<OneD, int> &conn,
1383  bool standard = true)
1384  {
1385  v_GetSimplexEquiSpacedConnectivity(conn,standard);
1386  }
1387 
1388  /** \brief This function performs a
1389  * projection/interpolation from the equispaced points
1390  * sometimes used in post-processing onto the coefficient
1391  * space
1392  *
1393  * This is primarily used for output purposes to use a
1394  * more even distribution of points more suitable for alot of
1395  * postprocessing
1396  */
1398  const Array<OneD, const NekDouble> &inarray,
1399  Array<OneD, NekDouble> &outarray);
1400 
1401  template<class T>
1402  boost::shared_ptr<T> as()
1403  {
1404 #if defined __INTEL_COMPILER && BOOST_VERSION > 105200
1405  typedef typename boost::shared_ptr<T>::element_type E;
1406  E * p = dynamic_cast< E* >( shared_from_this().get() );
1407  ASSERTL1(p, "Cannot perform cast");
1408  return boost::shared_ptr<T>( shared_from_this(), p );
1409 #else
1410  return boost::dynamic_pointer_cast<T>( shared_from_this() );
1411 #endif
1412  }
1413 
1415  Array<OneD, NekDouble> &outarray,
1416  bool multiplybyweights = true)
1417  {
1418  v_IProductWRTBase_SumFac(inarray,outarray,multiplybyweights);
1419  }
1420 
1421  protected:
1422  Array<OneD, LibUtilities::BasisSharedPtr> m_base; /**< Bases needed for the expansion */
1424  int m_ncoeffs; /**< Total number of coefficients used in the expansion */
1428 
1430  {
1431  return v_CreateStdMatrix(mkey);
1432  }
1433 
1434  /** \brief Create the static condensation of a matrix when
1435  using a boundary interior decomposition
1436 
1437  If a matrix system can be represented by
1438  \f$ Mat = \left [ \begin{array}{cc}
1439  A & B \\
1440  C & D \end{array} \right ] \f$
1441  This routine creates a matrix containing the statically
1442  condense system of the form
1443  \f$ Mat = \left [ \begin{array}{cc}
1444  A - B D^{-1} C & B D^{-1} \\
1445  D^{-1} C & D^{-1} \end{array} \right ] \f$
1446  **/
1448 
1449  /** \brief Create an IndexMap which contains mapping information linking any specific
1450  element shape with either its boundaries, edges, faces, verteces, etc.
1451 
1452  The index member of the IndexMapValue struct gives back an integer associated with an entity index
1453  The sign member of the same struct gives back a sign to algebrically apply entities orientation
1454  **/
1456 
1458  Array<OneD, NekDouble> &outarray);
1459 
1461  Array<OneD, NekDouble> &outarray)
1462  {
1463  v_BwdTrans_SumFac(inarray,outarray);
1464  }
1465 
1466 
1467  void IProductWRTDerivBase_SumFac(const int dir,
1468  const Array<OneD, const NekDouble>& inarray,
1469  Array<OneD, NekDouble> &outarray)
1470  {
1471  v_IProductWRTDerivBase_SumFac(dir,inarray,outarray);
1472  }
1473 
1474  // The term _MatFree denotes that the action of the MatrixOperation
1475  // is done withouth actually using the matrix (which then needs to be stored/calculated).
1476  // Although this does not strictly mean that no matrix operations are involved in the
1477  // evaluation of the operation, we use this term in the same context used as in the following
1478  // paper:
1479  // R. C. Kirby, M. G. Knepley, A. Logg, and L. R. Scott,
1480  // "Optimizing the evaluation of finite element matrices," SISC 27:741-758 (2005)
1482  Array<OneD,NekDouble> &outarray,
1483  const StdMatrixKey &mkey);
1484 
1486  Array<OneD,NekDouble> &outarray,
1487  const StdMatrixKey &mkey);
1488 
1490  Array<OneD,NekDouble> &outarray,
1491  const StdMatrixKey &mkey)
1492  {
1493  v_LaplacianMatrixOp_MatFree(inarray,outarray,mkey);
1494  }
1495 
1497  const Array<OneD, const NekDouble> &inarray,
1498  Array<OneD, NekDouble> &outarray,
1500  {
1501  v_LaplacianMatrixOp_MatFree_Kernel(inarray, outarray, wsp);
1502  }
1503 
1505  Array<OneD,NekDouble> &outarray,
1506  const StdMatrixKey &mkey);
1507 
1508  STD_REGIONS_EXPORT void LaplacianMatrixOp_MatFree(const int k1, const int k2,
1509  const Array<OneD, const NekDouble> &inarray,
1510  Array<OneD,NekDouble> &outarray,
1511  const StdMatrixKey &mkey);
1512 
1514  const Array<OneD, const NekDouble> &inarray,
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 
1527  Array<OneD,NekDouble> &outarray,
1528  const StdMatrixKey &mkey,
1529  bool addDiffusionTerm = true);
1530 
1532  Array<OneD,NekDouble> &outarray,
1533  const StdMatrixKey &mkey)
1534  {
1535  v_HelmholtzMatrixOp_MatFree(inarray,outarray,mkey);
1536  }
1537 
1539  Array<OneD,NekDouble> &outarray,
1540  const StdMatrixKey &mkey);
1541 
1544  Array<OneD, NekDouble> &outarray);
1545 
1547  Array<OneD, NekDouble> &coeffs,
1549 
1551  const Array<OneD, const NekDouble> &Lcoord,
1552  const Array<OneD, const NekDouble> &physvals);
1553 
1554  private:
1555  // Virtual functions
1556  STD_REGIONS_EXPORT virtual int v_GetNverts() const = 0;
1557  STD_REGIONS_EXPORT virtual int v_GetNedges() const;
1558 
1559  STD_REGIONS_EXPORT virtual int v_GetNfaces() const;
1560 
1561  STD_REGIONS_EXPORT virtual int v_NumBndryCoeffs() const;
1562 
1563  STD_REGIONS_EXPORT virtual int v_NumDGBndryCoeffs() const;
1564 
1565  STD_REGIONS_EXPORT virtual int v_GetEdgeNcoeffs(const int i) const;
1566 
1567  STD_REGIONS_EXPORT virtual int v_GetTotalEdgeIntNcoeffs() const;
1568 
1569  STD_REGIONS_EXPORT virtual int v_GetEdgeNumPoints(const int i) const;
1570 
1571  STD_REGIONS_EXPORT virtual int v_DetCartesianDirOfEdge(const int edge);
1572 
1573  STD_REGIONS_EXPORT virtual const LibUtilities::BasisKey v_DetEdgeBasisKey(const int i) const;
1574 
1575  STD_REGIONS_EXPORT virtual const LibUtilities::BasisKey v_DetFaceBasisKey(const int i, const int k) const;
1576 
1577  STD_REGIONS_EXPORT virtual int v_GetFaceNumPoints(const int i) const;
1578 
1579  STD_REGIONS_EXPORT virtual int v_GetFaceNcoeffs(const int i) const;
1580 
1581  STD_REGIONS_EXPORT virtual int v_GetFaceIntNcoeffs(const int i) const;
1582 
1583  STD_REGIONS_EXPORT virtual int v_GetTotalFaceIntNcoeffs() const;
1584 
1585 
1586  STD_REGIONS_EXPORT virtual int v_GetTraceNcoeffs(const int i) const;
1587 
1588  STD_REGIONS_EXPORT virtual LibUtilities::PointsKey v_GetFacePointsKey(const int i, const int j) const;
1589 
1591 
1593 
1595 
1596  STD_REGIONS_EXPORT virtual boost::shared_ptr<StdExpansion>
1597  v_GetStdExp(void) const;
1598 
1599  STD_REGIONS_EXPORT virtual int v_GetShapeDimension() const;
1600 
1602 
1604 
1605  STD_REGIONS_EXPORT virtual void v_BwdTrans (const Array<OneD, const NekDouble>& inarray,
1606  Array<OneD, NekDouble> &outarray) = 0;
1607 
1608  /**
1609  * @brief Transform a given function from physical quadrature space
1610  * to coefficient space.
1611  * @see StdExpansion::FwdTrans
1612  */
1613  STD_REGIONS_EXPORT virtual void v_FwdTrans (
1614  const Array<OneD, const NekDouble>& inarray,
1615  Array<OneD, NekDouble> &outarray) = 0;
1616 
1617  /**
1618  * @brief Calculates the inner product of a given function \a f
1619  * with the different modes of the expansion
1620  */
1622  const Array<OneD, const NekDouble>& inarray,
1623  Array<OneD, NekDouble> &outarray) = 0;
1624 
1626  const Array<OneD, const NekDouble>& base,
1627  const Array<OneD, const NekDouble>& inarray,
1628  Array<OneD, NekDouble>& outarray,
1629  int coll_check)
1630  {
1631  ASSERTL0(false, "StdExpansion::v_IProductWRTBase has no (and should have no) implementation");
1632  }
1633 
1634  STD_REGIONS_EXPORT virtual void v_IProductWRTDerivBase (const int dir,
1635  const Array<OneD, const NekDouble>& inarray,
1636  Array<OneD, NekDouble> &outarray);
1637 
1639  Array<OneD, NekDouble> &outarray);
1640 
1642 
1643  STD_REGIONS_EXPORT virtual void v_PhysDeriv (const Array<OneD, const NekDouble>& inarray,
1644  Array<OneD, NekDouble> &out_d1,
1645  Array<OneD, NekDouble> &out_d2,
1646  Array<OneD, NekDouble> &out_d3);
1647 
1648  STD_REGIONS_EXPORT virtual void v_PhysDeriv_s (const Array<OneD, const NekDouble>& inarray,
1649  Array<OneD, NekDouble> &out_ds);
1650 
1652  Array<OneD, NekDouble>& out_dn);
1653  STD_REGIONS_EXPORT virtual void v_PhysDeriv(const int dir,
1654  const Array<OneD, const NekDouble>& inarray,
1655  Array<OneD, NekDouble> &out_d0);
1656 
1658  const Array<OneD, const NekDouble>& direction,
1659  Array<OneD, NekDouble> &outarray);
1660 
1661  STD_REGIONS_EXPORT virtual void v_StdPhysDeriv (const Array<OneD, const NekDouble>& inarray,
1662  Array<OneD, NekDouble> &out_d1,
1663  Array<OneD, NekDouble> &out_d2,
1664  Array<OneD, NekDouble> &out_d3);
1665 
1666  STD_REGIONS_EXPORT virtual void v_StdPhysDeriv (const int dir,
1667  const Array<OneD, const NekDouble>& inarray,
1668  Array<OneD, NekDouble> &outarray);
1669 
1670  STD_REGIONS_EXPORT virtual void v_AddRobinMassMatrix(const int edgeid, const Array<OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat);
1671 
1672  STD_REGIONS_EXPORT virtual void v_AddRobinEdgeContribution(const int edgeid, const Array<OneD, const NekDouble> &primCoeffs, Array<OneD, NekDouble> &coeffs);
1673 
1675 
1677 
1679  const Array<OneD, const NekDouble>& xi,
1680  Array<OneD, NekDouble>& eta);
1681 
1682 
1683  STD_REGIONS_EXPORT virtual void v_FillMode(const int mode, Array<OneD, NekDouble> &outarray);
1684 
1686 
1688 
1689  STD_REGIONS_EXPORT virtual void v_GetCoords(Array<OneD, NekDouble> &coords_0,
1690  Array<OneD, NekDouble> &coords_1,
1691  Array<OneD, NekDouble> &coords_2);
1692 
1693  STD_REGIONS_EXPORT virtual void v_GetCoord(const Array<OneD, const NekDouble>& Lcoord,
1694  Array<OneD, NekDouble> &coord);
1695 
1696  STD_REGIONS_EXPORT virtual int v_GetCoordim(void);
1697 
1699 
1701 
1702  STD_REGIONS_EXPORT virtual int v_GetVertexMap(int localVertexId,
1703  bool useCoeffPacking = false);
1704 
1705  STD_REGIONS_EXPORT virtual void v_GetEdgeInteriorMap(const int eid, const Orientation edgeOrient,
1706  Array<OneD, unsigned int> &maparray,
1707  Array<OneD, int> &signarray);
1708 
1709  STD_REGIONS_EXPORT virtual void v_GetFaceInteriorMap(const int fid, const Orientation faceOrient,
1710  Array<OneD, unsigned int> &maparray,
1711  Array<OneD, int> &signarray);
1712 
1714  const int eid,
1715  const Orientation edgeOrient,
1716  Array<OneD, unsigned int>& maparray,
1717  Array<OneD, int>& signarray,
1718  int P = -1);
1719 
1720  STD_REGIONS_EXPORT virtual void v_GetFaceToElementMap(const int fid, const Orientation faceOrient,
1721  Array<OneD, unsigned int> &maparray,
1722  Array<OneD, int> &signarray,
1723  int nummodesA = -1, int nummodesB = -1);
1724 
1725  /**
1726  * @brief Extract the physical values along edge \a edge from \a
1727  * inarray into \a outarray following the local edge orientation
1728  * and point distribution defined by defined in \a EdgeExp.
1729  */
1730  STD_REGIONS_EXPORT virtual void v_GetEdgePhysVals(const int edge, const Array<OneD, const NekDouble> &inarray, Array<OneD,NekDouble> &outarray);
1731 
1732  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);
1733 
1734  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);
1735 
1736  STD_REGIONS_EXPORT virtual void v_GetVertexPhysVals(const int vertex, const Array<OneD, const NekDouble> &inarray, NekDouble &outarray);
1737 
1738  STD_REGIONS_EXPORT virtual void v_GetEdgeInterpVals(const int edge,
1739  const Array<OneD, const NekDouble> &inarray,Array<OneD,NekDouble> &outarray);
1740 
1742  const int edge,
1743  Array<OneD, NekDouble> &outarray);
1744 
1746  const int face,
1747  const boost::shared_ptr<StdExpansion> &FaceExp,
1748  const Array<OneD, const NekDouble> &inarray,
1749  Array<OneD, NekDouble> &outarray,
1750  StdRegions::Orientation orient);
1751 
1752  STD_REGIONS_EXPORT virtual void v_GetEdgePhysMap(
1753  const int edge,
1754  Array<OneD,int> &outarray);
1755 
1756  STD_REGIONS_EXPORT virtual void v_GetFacePhysMap(
1757  const int face,
1758  Array<OneD,int> &outarray);
1759 
1761  const Array<OneD, const NekDouble> &inarray,
1762  Array<OneD, NekDouble> &outarray);
1763 
1765  const Array<OneD, const NekDouble> &inarray,
1766  Array<OneD, NekDouble> &outarray);
1767 
1768  STD_REGIONS_EXPORT virtual const boost::shared_ptr<SpatialDomains::GeomFactors>& v_GetMetricInfo() const;
1769 
1771  Array<OneD, NekDouble> &outarray);
1772 
1774  Array<OneD, NekDouble> &outarray, bool multiplybyweights = true);
1775 
1776  STD_REGIONS_EXPORT virtual void v_IProductWRTDerivBase_SumFac(const int dir,
1777  const Array<OneD, const NekDouble>& inarray,
1778  Array<OneD, NekDouble> &outarray);
1779 
1781  Array<OneD,NekDouble> &outarray,
1782  const StdMatrixKey &mkey);
1783 
1785  Array<OneD,NekDouble> &outarray,
1786  const StdMatrixKey &mkey);
1787 
1789  const StdMatrixKey &mkey);
1790 
1792  int numMin,
1793  const Array<OneD, const NekDouble> &inarray,
1794  Array<OneD,NekDouble> &outarray);
1795 
1796  STD_REGIONS_EXPORT virtual void v_LaplacianMatrixOp(const int k1, const int k2,
1797  const Array<OneD, const NekDouble> &inarray,
1798  Array<OneD,NekDouble> &outarray,
1799  const StdMatrixKey &mkey);
1800 
1801  STD_REGIONS_EXPORT virtual void v_WeakDerivMatrixOp(const int i,
1802  const Array<OneD, const NekDouble> &inarray,
1803  Array<OneD,NekDouble> &outarray,
1804  const StdMatrixKey &mkey);
1805 
1807  Array<OneD,NekDouble> &outarray,
1808  const StdMatrixKey &mkey);
1809 
1811  Array<OneD,NekDouble> &outarray,
1812  const StdMatrixKey &mkey);
1813 
1815  const NekDouble> &inarray,
1816  Array<OneD,NekDouble> &outarray,
1817  const StdMatrixKey &mkey,
1818  bool addDiffusionTerm=true);
1819 
1821  Array<OneD,NekDouble> &outarray,
1822  const StdMatrixKey &mkey);
1823 
1825  Array<OneD,NekDouble> &outarray,
1826  const StdMatrixKey &mkey);
1827 
1829  const Array<OneD, const NekDouble> &inarray,
1830  Array<OneD, NekDouble> &outarray,
1831  Array<OneD, NekDouble> &wsp);
1832 
1834  Array<OneD,NekDouble> &outarray,
1835  const StdMatrixKey &mkey);
1836 
1837  STD_REGIONS_EXPORT virtual const NormalVector & v_GetEdgeNormal(const int edge) const;
1838 
1839  STD_REGIONS_EXPORT virtual void v_ComputeEdgeNormal(const int edge);
1840 
1841  STD_REGIONS_EXPORT virtual void v_NegateEdgeNormal(const int edge);
1842 
1843  STD_REGIONS_EXPORT virtual bool v_EdgeNormalNegated(const int edge);
1844 
1845  STD_REGIONS_EXPORT virtual void v_ComputeFaceNormal(const int face);
1846 
1847  STD_REGIONS_EXPORT virtual void v_NegateFaceNormal(const int face);
1848 
1849  STD_REGIONS_EXPORT virtual bool v_FaceNormalNegated(const int face);
1850 
1851  STD_REGIONS_EXPORT virtual const NormalVector & v_GetVertexNormal(const int vertex) const;
1852 
1853  STD_REGIONS_EXPORT virtual void v_ComputeVertexNormal(const int vertex);
1854 
1855  STD_REGIONS_EXPORT virtual const NormalVector & v_GetFaceNormal(const int face) const;
1856  STD_REGIONS_EXPORT virtual const NormalVector &
1857  v_GetSurfaceNormal(const int id) const;
1858 
1860  v_GetEdgeInverseBoundaryMap(int eid);
1861 
1864 
1866 
1868  Array<OneD, int> &conn,
1869  bool standard = true);
1870  };
1871 
1872 
1873  typedef boost::shared_ptr<StdExpansion> StdExpansionSharedPtr;
1874  typedef std::vector< StdExpansionSharedPtr > StdExpansionVector;
1876 
1877  /**
1878  * This function is a wrapper around the virtual function
1879  * \a v_FwdTrans()
1880  *
1881  * Given a function evaluated at the quadrature points, this
1882  * function calculates the expansion coefficients such that the
1883  * resulting expansion approximates the original function.
1884  *
1885  * The calculation of the expansion coefficients is done using a
1886  * Galerkin projection. This is equivalent to the operation:
1887  * \f[ \mathbf{\hat{u}} = \mathbf{M}^{-1} \mathbf{I}\f]
1888  * where
1889  * - \f$\mathbf{M}[p][q]= \int\phi_p(\mathbf{\xi})\phi_q(
1890  * \mathbf{\xi}) d\mathbf{\xi}\f$ is the Mass matrix
1891  * - \f$\mathbf{I}[p] = \int\phi_p(\mathbf{\xi}) u(\mathbf{\xi})
1892  * d\mathbf{\xi}\f$
1893  *
1894  * This function takes the array \a inarray as the values of the
1895  * function evaluated at the quadrature points
1896  * (i.e. \f$\mathbf{u}\f$),
1897  * and stores the resulting coefficients \f$\mathbf{\hat{u}}\f$
1898  * in the \a outarray
1899  *
1900  * @param inarray array of the function discretely evaluated at the
1901  * quadrature points
1902  *
1903  * @param outarray array of the function coefficieints
1904  */
1906  Array<OneD, NekDouble> &outarray)
1907  {
1908  v_FwdTrans(inarray,outarray);
1909  }
1910 
1911  } //end of namespace
1912 } //end of namespace
1913 
1914 #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)
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:161
virtual void v_GetInteriorMap(Array< OneD, unsigned int > &outarray)
const NormalVector & GetEdgeNormal(const int edge) const
StdRegions::Orientation GetCartesianEorient(int edge)
Definition: StdExpansion.h:772
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)
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:220
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:654
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:592
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:971
int GetNfaces() const
This function returns the number of faces of the expansion domain.
Definition: StdExpansion.h:438
static Array< OneD, NekDouble > NullNekDouble1DArray
virtual void v_ExtractDataToCoeffs(const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmode_offset, NekDouble *coeffs)
Unpack data from input file assuming it comes from the same expansion type.
StdRegions::Orientation GetPorient(int point)
Definition: StdExpansion.h:767
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:715
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:710
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:942
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:832
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:935
void MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:949
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:629
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:574
void SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
Definition: StdExpansion.h:992
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:752
void GetEdgePhysVals(const int edge, const boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:878
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:727
StdRegions::Orientation GetEorient(int edge)
Definition: StdExpansion.h:762
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:757
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:928
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:911
virtual void v_NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix(const LocalRegions::MatrixKey &mkey)
Definition: StdExpansion.h:747
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:985
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:700
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:676
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
void NormVectorIProductWRTBase(const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:742
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:720
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:732
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:705
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:826
void FwdTrans_BndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:538
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:855
void GetVertexPhysVals(const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray)
Definition: StdExpansion.h:891
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)
Definition: StdExpansion.h:998
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:694
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:821
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:662
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:918
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:886
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:645
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...
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
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:871
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:978
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:846
#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:737
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:792
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:525
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:777
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:635
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:898
virtual int v_GetNfaces() const
void SetCoeffsToOrientation(StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:784
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
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:191
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)
void ExtractDataToCoeffs(const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs)
Definition: StdExpansion.h:797
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:816
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:839
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)