Nektar++
StdExpansion.h
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File Stdexpansion.h
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // 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_GetFaceNcoeffs()
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  LibUtilities::PointsKey GetFacePointsKey(const int i, const int j) const
370  {
371  return v_GetFacePointsKey(i, j);
372  }
373 
374  int NumBndryCoeffs(void) const
375  {
376  return v_NumBndryCoeffs();
377  }
378 
379  int NumDGBndryCoeffs(void) const
380  {
381  return v_NumDGBndryCoeffs();
382  }
383 
384  /** \brief This function returns the type of expansion basis on the
385  * \a i-th edge
386  *
387  * This function is a wrapper around the virtual function
388  * \a v_GetEdgeBasisType()
389  *
390  * The different types of bases implemented in the code are defined
391  * in the LibUtilities::BasisType enumeration list. As a result, the
392  * function will return one of the types of this enumeration list.
393  *
394  * \param i specifies which edge
395  * \return returns the expansion basis on the \a i-th edge
396  */
398  {
399  return v_GetEdgeBasisType(i);
400  }
401 
402  /** \brief This function returns the type of expansion
403  * Nodal point type if defined
404  *
405  * This function is a wrapper around the virtual function
406  * \a v_GetNodalPointsKey()
407  *
408  */
410  {
411  return v_GetNodalPointsKey();
412  };
413 
414  /** \brief This function returns the number of faces of the
415  * expansion domain
416  *
417  * This function is a wrapper around the virtual function
418  * \a v_GetNFaces()
419  *
420  * \return returns the number of faces of the expansion domain
421  */
422  int GetNfaces() const
423  {
424  return v_GetNfaces();
425  }
426 
427  /**
428  * @brief Returns the number of trace elements connected to this
429  * element.
430  *
431  * For example, a quadrilateral has four edges, so this function
432  * would return 4.
433  */
434  int GetNtrace() const
435  {
436  const int nBase = m_base.num_elements();
437  return
438  nBase == 1 ? 2 :
439  nBase == 2 ? GetNedges() :
440  nBase == 3 ? GetNfaces() : 0;
441  }
442 
443  /** \brief This function returns the shape of the expansion domain
444  *
445  * This function is a wrapper around the virtual function
446  * \a v_DetShapeType()
447  *
448  * The different shape types implemented in the code are defined
449  * in the ::ShapeType enumeration list. As a result, the
450  * function will return one of the types of this enumeration list.
451  *
452  * \return returns the shape of the expansion domain
453  */
455  {
456  return v_DetShapeType();
457  }
458 
459  boost::shared_ptr<StdExpansion> GetStdExp(void) const
460  {
461  return v_GetStdExp();
462  }
463 
464 
465  int GetShapeDimension() const
466  {
467  return v_GetShapeDimension();
468  }
469 
471  {
473  }
474 
475 
477  {
478  return v_IsNodalNonTensorialExp();
479  }
480 
481  /** \brief This function performs the Backward transformation from
482  * coefficient space to physical space
483  *
484  * This function is a wrapper around the virtual function
485  * \a v_BwdTrans()
486  *
487  * Based on the expansion coefficients, this function evaluates the
488  * expansion at the quadrature points. This is equivalent to the
489  * operation \f[ u(\xi_{1i}) =
490  * \sum_{p=0}^{P-1} \hat{u}_p \phi_p(\xi_{1i}) \f] which can be
491  * evaluated as \f$ {\bf u} = {\bf B}^T {\bf \hat{u}} \f$ with
492  * \f${\bf B}[i][j] = \phi_i(\xi_{j})\f$
493  *
494  * This function requires that the coefficient array
495  * \f$\mathbf{\hat{u}}\f$ provided as \a inarray.
496  *
497  * The resulting array
498  * \f$\mathbf{u}[m]=u(\mathbf{\xi}_m)\f$ containing the
499  * expansion evaluated at the quadrature points, is stored
500  * in the \a outarray.
501  *
502  * \param inarray contains the values of the expansion
503  * coefficients (input of the function)
504  *
505  * \param outarray contains the values of the expansion evaluated
506  * at the quadrature points (output of the function)
507  */
508 
510  Array<OneD, NekDouble> &outarray)
511  {
512  v_BwdTrans (inarray, outarray);
513  }
514 
515  /**
516  * @brief This function performs the Forward transformation from
517  * physical space to coefficient space.
518  */
519  inline void FwdTrans (const Array<OneD, const NekDouble>& inarray,
520  Array<OneD, NekDouble> &outarray);
521 
523  Array<OneD, NekDouble> &outarray)
524  {
525  v_FwdTrans_BndConstrained(inarray,outarray);
526  }
527 
528  /** \brief This function integrates the specified function over the
529  * domain
530  *
531  * This function is a wrapper around the virtual function
532  * \a v_Integral()
533  *
534  * Based on the values of the function evaluated at the quadrature
535  * points (which are stored in \a inarray), this function calculates
536  * the integral of this function over the domain. This is
537  * equivalent to the numerical evaluation of the operation
538  * \f[ I=\int u(\mathbf{\xi})d \mathbf{\xi}\f]
539  *
540  * \param inarray values of the function to be integrated evaluated
541  * at the quadrature points (i.e.
542  * \a inarray[m]=\f$u(\mathbf{\xi}_m)\f$)
543  * \return returns the value of the calculated integral
544  *
545  * Inputs:\n
546 
547  - \a inarray: definition of function to be returned at quadrature point
548  of expansion.
549 
550  Outputs:\n
551 
552  - returns \f$\int^1_{-1}\int^1_{-1} u(\xi_1, \xi_2) J[i,j] d
553  \xi_1 d \xi_2 \f$ where \f$inarray[i,j] =
554  u(\xi_{1i},\xi_{2j}) \f$ and \f$ J[i,j] \f$ is the
555  Jacobian evaluated at the quadrature point.
556  *
557  */
559  {
560  return v_Integral(inarray);
561  }
562 
563  /** \brief This function fills the array \a outarray with the
564  * \a mode-th mode of the expansion
565  *
566  * This function is a wrapper around the virtual function
567  * \a v_FillMode()
568  *
569  * The requested mode is evaluated at the quadrature points
570  *
571  * \param mode the mode that should be filled
572  * \param outarray contains the values of the \a mode-th mode of the
573  * expansion evaluated at the quadrature points (output of the
574  * function)
575  */
576  void FillMode(const int mode, Array<OneD, NekDouble> &outarray)
577  {
578  v_FillMode(mode, outarray);
579  }
580 
581  /** \brief this function calculates the inner product of a given
582  * function \a f with the different modes of the expansion
583  *
584  * This function is a wrapper around the virtual function
585  * \a v_IProductWRTBase()
586  *
587  * This is equivalent to the numerical evaluation of
588  * \f[ I[p] = \int \phi_p(\mathbf{x}) f(\mathbf{x}) d\mathbf{x}\f]
589  * \f$ \begin{array}{rcl} I_{pq} = (\phi_q \phi_q, u) & = &
590  \sum_{i=0}^{nq_0} \sum_{j=0}^{nq_1} \phi_p(\xi_{0,i})
591  \phi_q(\xi_{1,j}) w^0_i w^1_j u(\xi_{0,i} \xi_{1,j})
592  J_{i,j}\\ & = & \sum_{i=0}^{nq_0} \phi_p(\xi_{0,i})
593  \sum_{j=0}^{nq_1} \phi_q(\xi_{1,j}) \tilde{u}_{i,j}
594  J_{i,j} \end{array} \f$
595 
596  where
597 
598  \f$ \tilde{u}_{i,j} = w^0_i w^1_j u(\xi_{0,i},\xi_{1,j}) \f$
599 
600  which can be implemented as
601 
602  \f$ f_{qi} = \sum_{j=0}^{nq_1} \phi_q(\xi_{1,j}) \tilde{u}_{i,j} =
603  {\bf B_1 U} \f$
604  \f$ I_{pq} = \sum_{i=0}^{nq_0} \phi_p(\xi_{0,i}) f_{qi} =
605  {\bf B_0 F} \f$
606  *
607  * \param inarray contains the values of the function \a f
608  * evaluated at the quadrature points
609  * \param outarray contains the values of the inner product of \a f
610  * with the different modes, i.e. \f$ outarray[p] = I[p]\f$
611  * (output of the function)
612  */
614  Array<OneD, NekDouble> &outarray)
615  {
616  v_IProductWRTBase(inarray, outarray);
617  }
618 
620  const Array<OneD, const NekDouble>& base,
621  const Array<OneD, const NekDouble>& inarray,
622  Array<OneD, NekDouble> &outarray,
623  int coll_check)
624  {
625  v_IProductWRTBase(base, inarray, outarray, coll_check);
626  }
627 
628 
629  void IProductWRTDerivBase(const int dir,
630  const Array<OneD, const NekDouble>& inarray,
631  Array<OneD, NekDouble> &outarray)
632  {
633  v_IProductWRTDerivBase(dir,inarray, outarray);
634  }
635 
636  /// \brief Get the element id of this expansion when used
637  /// in a list by returning value of #m_elmt_id
638  inline int GetElmtId()
639  {
640  return m_elmt_id;
641  }
642 
643 
644  /// \brief Set the element id of this expansion when used
645  /// in a list by returning value of #m_elmt_id
646  inline void SetElmtId(const int id)
647  {
648  m_elmt_id = id;
649  }
650 
651  /** \brief this function returns the physical coordinates of the
652  * quadrature points of the expansion
653  *
654  * This function is a wrapper around the virtual function
655  * \a v_GetCoords()
656  *
657  * \param coords an array containing the coordinates of the
658  * quadrature points (output of the function)
659  */
663  {
664  v_GetCoords(coords_1,coords_2,coords_3);
665  }
666 
667  /** \brief given the coordinates of a point of the element in the
668  * local collapsed coordinate system, this function calculates the
669  * physical coordinates of the point
670  *
671  * This function is a wrapper around the virtual function
672  * \a v_GetCoord()
673  *
674  * \param Lcoords the coordinates in the local collapsed
675  * coordinate system
676  * \param coords the physical coordinates (output of the function)
677  */
679  Array<OneD, NekDouble> &coord)
680  {
681  v_GetCoord(Lcoord, coord);
682  }
683 
685  {
686  return m_stdMatrixManager[mkey];
687  }
688 
690  {
691  return m_stdStaticCondMatrixManager[mkey];
692  }
693 
695  {
696  return m_IndexMapManager[ikey];
697  }
698 
700  {
701  return v_GetPhysNormals();
702  }
703 
705  {
706  v_SetPhysNormals(normal);
707  }
708 
709  STD_REGIONS_EXPORT virtual void SetUpPhysNormals(const int edge);
710 
712  {
713  v_NormVectorIProductWRTBase(Fx,Fy,outarray);
714  }
715 
717  {
718  v_NormVectorIProductWRTBase(Fx,Fy,Fz,outarray);
719  }
720 
722  {
723  return v_GetLocStaticCondMatrix(mkey);
724  }
725 
727  {
728  return v_DropLocStaticCondMatrix(mkey);
729  }
730 
732  {
733  return v_GetForient(face);
734  }
735 
737  {
738  return v_GetEorient(edge);
739  }
740 
742  {
743  return v_GetPorient(point);
744  }
745 
747  {
748  return v_GetCartesianEorient(edge);
749  }
750 
752  Array<OneD, NekDouble> &coeffs,
754  {
755  v_SetCoeffsToOrientation(coeffs, dir);
756  }
757 
761  Array<OneD, NekDouble> &outarray)
762  {
763  v_SetCoeffsToOrientation(dir,inarray,outarray);
764  }
765 
766  int CalcNumberOfCoefficients(const std::vector<unsigned int> &nummodes, int &modes_offset)
767  {
768  return v_CalcNumberOfCoefficients(nummodes,modes_offset);
769  }
770 
771  void ExtractDataToCoeffs(const NekDouble *data,
772  const std::vector<unsigned int > &nummodes,
773  const int nmodes_offset,
774  NekDouble *coeffs)
775  {
776  v_ExtractDataToCoeffs(data,nummodes,nmodes_offset,coeffs);
777  }
778 
779  // virtual functions related to LocalRegions
781  const Array<OneD, const NekDouble> &Lcoord,
782  const Array<OneD, const NekDouble> &physvals);
783 
784 
786  {
787  return v_GetCoordim();
788  }
789 
791  {
792  v_GetBoundaryMap(outarray);
793  }
794 
796  {
797  v_GetInteriorMap(outarray);
798  }
799 
800  int GetVertexMap(const int localVertexId,
801  bool useCoeffPacking = false)
802  {
803  return v_GetVertexMap(localVertexId,useCoeffPacking);
804  }
805 
806  void GetEdgeInteriorMap(const int eid, const Orientation edgeOrient,
807  Array<OneD, unsigned int> &maparray,
808  Array<OneD, int> &signarray)
809  {
810  v_GetEdgeInteriorMap(eid,edgeOrient,maparray,signarray);
811  }
812 
813  void GetFaceInteriorMap(const int fid, const Orientation faceOrient,
814  Array<OneD, unsigned int> &maparray,
815  Array<OneD, int> &signarray)
816  {
817  v_GetFaceInteriorMap(fid,faceOrient,maparray,signarray);
818  }
819 
820  void GetEdgeToElementMap(const int eid, const Orientation edgeOrient,
821  Array<OneD, unsigned int> &maparray,
822  Array<OneD, int> &signarray)
823  {
824  v_GetEdgeToElementMap(eid,edgeOrient,maparray,signarray);
825  }
826 
827  void GetFaceToElementMap(const int fid, const Orientation faceOrient,
828  Array<OneD, unsigned int> &maparray,
829  Array<OneD, int> &signarray,
830  int nummodesA = -1, int nummodesB = -1)
831  {
832  v_GetFaceToElementMap(fid,faceOrient,maparray,signarray,
833  nummodesA,nummodesB);
834  }
835 
836 
837  /**
838  * @brief Extract the physical values along edge \a edge from \a
839  * inarray into \a outarray following the local edge orientation
840  * and point distribution defined by defined in \a EdgeExp.
841  */
842 
843  void GetEdgePhysVals(const int edge, const Array<OneD,
844  const NekDouble> &inarray,
845  Array<OneD,NekDouble> &outarray)
846  {
847  v_GetEdgePhysVals(edge,inarray,outarray);
848  }
849 
850  void GetEdgePhysVals(const int edge,
851  const boost::shared_ptr<StdExpansion> &EdgeExp,
852  const Array<OneD, const NekDouble> &inarray,
853  Array<OneD,NekDouble> &outarray)
854  {
855  v_GetEdgePhysVals(edge,EdgeExp,inarray,outarray);
856  }
857 
858  void GetTracePhysVals(const int edge, const boost::shared_ptr<StdExpansion> &EdgeExp, const Array<OneD, const NekDouble> &inarray, Array<OneD,NekDouble> &outarray)
859  {
860  v_GetTracePhysVals(edge,EdgeExp,inarray,outarray);
861  }
862 
863  void GetVertexPhysVals(const int vertex,
864  const Array<OneD, const NekDouble> &inarray,
865  NekDouble &outarray)
866  {
867  v_GetVertexPhysVals(vertex, inarray, outarray);
868  }
869 
870  void GetEdgeInterpVals(const int edge,const Array<OneD,
871  const NekDouble> &inarray,
872  Array<OneD,NekDouble> &outarray)
873  {
874  v_GetEdgeInterpVals(edge, inarray, outarray);
875  }
876 
877  /**
878  * @brief Extract the metric factors to compute the contravariant
879  * fluxes along edge \a edge and stores them into \a outarray
880  * following the local edge orientation (i.e. anticlockwise
881  * convention).
882  */
884  const int edge,
885  Array<OneD, NekDouble> &outarray)
886  {
887  v_GetEdgeQFactors(edge, outarray);
888  }
889 
890 
891 
893  const int face,
894  const boost::shared_ptr<StdExpansion> &FaceExp,
895  const Array<OneD, const NekDouble> &inarray,
896  Array<OneD, NekDouble> &outarray,
898  {
899  v_GetFacePhysVals(face, FaceExp, inarray, outarray, orient);
900  }
901 
903  const Array<OneD, const NekDouble> &inarray,
904  Array<OneD, NekDouble> &outarray)
905  {
906  v_MultiplyByQuadratureMetric(inarray, outarray);
907  }
908 
910  const Array<OneD, const NekDouble> &inarray,
911  Array<OneD, NekDouble> & outarray)
912  {
913  v_MultiplyByStdQuadratureMetric(inarray, outarray);
914  }
915 
916  // Matrix Routines
917 
918  /** \brief this function generates the mass matrix
919  * \f$\mathbf{M}[i][j] =
920  * \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}\f$
921  *
922  * \return returns the mass matrix
923  */
924 
926 
928  Array<OneD,NekDouble> &outarray,
929  const StdMatrixKey &mkey);
930 
932  Array<OneD,NekDouble> &outarray,
933  const StdMatrixKey &mkey)
934  {
935  v_MassMatrixOp(inarray,outarray,mkey);
936  }
937 
939  Array<OneD,NekDouble> &outarray,
940  const StdMatrixKey &mkey)
941  {
942  v_LaplacianMatrixOp(inarray,outarray,mkey);
943  }
944 
945  void ReduceOrderCoeffs(int numMin,
946  const Array<OneD, const NekDouble> &inarray,
947  Array<OneD,NekDouble> &outarray)
948  {
949  v_ReduceOrderCoeffs(numMin,inarray,outarray);
950  }
951 
953  const StdMatrixKey &mkey)
954  {
955  v_SVVLaplacianFilter(array,mkey);
956  }
957 
958  void LaplacianMatrixOp(const int k1, const int k2,
959  const Array<OneD, const NekDouble> &inarray,
960  Array<OneD,NekDouble> &outarray,
961  const StdMatrixKey &mkey)
962  {
963  v_LaplacianMatrixOp(k1,k2,inarray,outarray,mkey);
964  }
965 
966  void WeakDerivMatrixOp(const int i,
967  const Array<OneD, const NekDouble> &inarray,
968  Array<OneD,NekDouble> &outarray,
969  const StdMatrixKey &mkey)
970  {
971  v_WeakDerivMatrixOp(i,inarray,outarray,mkey);
972  }
973 
975  Array<OneD,NekDouble> &outarray,
976  const StdMatrixKey &mkey)
977  {
978  v_WeakDirectionalDerivMatrixOp(inarray,outarray,mkey);
979  }
980 
982  Array<OneD,NekDouble> &outarray,
983  const StdMatrixKey &mkey)
984  {
985  v_MassLevelCurvatureMatrixOp(inarray,outarray,mkey);
986  }
987 
989  Array<OneD,NekDouble> &outarray,
990  const StdMatrixKey &mkey,
991  bool addDiffusionTerm = true)
992  {
993  v_LinearAdvectionDiffusionReactionMatrixOp(inarray,outarray,mkey,addDiffusionTerm);
994  }
995 
996  /**
997  * @param inarray Input array @f$ \mathbf{u} @f$.
998  * @param outarray Output array @f$ \boldsymbol{\nabla^2u}
999  * + \lambda \boldsymbol{u} @f$.
1000  * @param mkey
1001  */
1003  Array<OneD,NekDouble> &outarray,
1004  const StdMatrixKey &mkey)
1005  {
1006  v_HelmholtzMatrixOp(inarray,outarray,mkey);
1007  }
1008 
1010  {
1011  return v_GenMatrix(mkey);
1012  }
1013 
1015  Array<OneD, NekDouble> &out_d0,
1018  {
1019  v_PhysDeriv (inarray, out_d0, out_d1, out_d2);
1020  }
1021 
1022  void PhysDeriv(const int dir,
1023  const Array<OneD, const NekDouble>& inarray,
1024  Array<OneD, NekDouble> &outarray)
1025  {
1026  v_PhysDeriv (dir, inarray, outarray);
1027  }
1028 
1030  Array<OneD, NekDouble> &out_ds)
1031  {
1032  v_PhysDeriv_s(inarray,out_ds);
1033  }
1034 
1036  Array<OneD, NekDouble>& out_dn)
1037  {
1038  v_PhysDeriv_n(inarray,out_dn);
1039  }
1040 
1042  const Array<OneD, const NekDouble>& direction,
1043  Array<OneD, NekDouble> &outarray)
1044  {
1045  v_PhysDirectionalDeriv (inarray, direction, outarray);
1046  }
1047 
1049  Array<OneD, NekDouble> &out_d0,
1052  {
1053  v_StdPhysDeriv(inarray, out_d0, out_d1, out_d2);
1054  }
1055 
1056  void StdPhysDeriv (const int dir,
1057  const Array<OneD, const NekDouble>& inarray,
1058  Array<OneD, NekDouble> &outarray)
1059  {
1060  v_StdPhysDeriv(dir,inarray,outarray);
1061  }
1062 
1063  void AddRobinMassMatrix(const int edgeid, const Array<OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
1064  {
1065  v_AddRobinMassMatrix(edgeid,primCoeffs,inoutmat);
1066  }
1067 
1068  void AddRobinEdgeContribution(const int edgeid, const Array<OneD, const NekDouble> &primCoeffs, Array<OneD, NekDouble> &coeffs)
1069  {
1070  v_AddRobinEdgeContribution(edgeid, primCoeffs, coeffs);
1071  }
1072 
1073  /** \brief This function evaluates the expansion at a single
1074  * (arbitrary) point of the domain
1075  *
1076  * This function is a wrapper around the virtual function
1077  * \a v_PhysEvaluate()
1078  *
1079  * Based on the value of the expansion at the quadrature
1080  * points provided in \a physvals, this function
1081  * calculates the value of the expansion at an arbitrary
1082  * single points (with coordinates \f$ \mathbf{x_c}\f$
1083  * given by the pointer \a coords). This operation,
1084  * equivalent to \f[ u(\mathbf{x_c}) = \sum_p
1085  * \phi_p(\mathbf{x_c}) \hat{u}_p \f] is evaluated using
1086  * Lagrangian interpolants through the quadrature points:
1087  * \f[ u(\mathbf{x_c}) = \sum_p h_p(\mathbf{x_c}) u_p\f]
1088  *
1089  * \param coords the coordinates of the single point
1090  * \param physvals the interpolated field at the quadrature points
1091  *
1092  * \return returns the value of the expansion at the
1093  * single point
1094  */
1096  const Array<OneD, const NekDouble>& physvals)
1097  {
1098  return v_PhysEvaluate(coords,physvals);
1099  }
1100 
1101 
1102  /** \brief This function evaluates the expansion at a single
1103  * (arbitrary) point of the domain
1104  *
1105  * This function is a wrapper around the virtual function
1106  * \a v_PhysEvaluate()
1107  *
1108  * Based on the value of the expansion at the quadrature
1109  * points provided in \a physvals, this function
1110  * calculates the value of the expansion at an arbitrary
1111  * single points associated with the interpolation
1112  * matrices provided in \f$ I \f$.
1113  *
1114  * \param I an Array of lagrange interpolantes evaluated
1115  * at the coordinate and going through the local physical
1116  * quadrature
1117  * \param physvals the interpolated field at the quadrature points
1118  *
1119  * \return returns the value of the expansion at the
1120  * single point
1121  */
1123  const Array<OneD, const NekDouble >& physvals)
1124  {
1125  return v_PhysEvaluate(I,physvals);
1126  }
1127 
1128 
1129  /**
1130  * \brief Convert local cartesian coordinate \a xi into local
1131  * collapsed coordinates \a eta
1132  **/
1135  {
1136  v_LocCoordToLocCollapsed(xi,eta);
1137  }
1138 
1139  const boost::shared_ptr<SpatialDomains::GeomFactors>& GetMetricInfo(void) const
1140  {
1141  return v_GetMetricInfo();
1142  }
1143 
1144  /// \brief Get the element id of this expansion when used
1145  /// in a list by returning value of #m_elmt_id
1146  STD_REGIONS_EXPORT virtual int v_GetElmtId();
1147 
1149 
1151 
1152  STD_REGIONS_EXPORT virtual void v_SetUpPhysNormals(const int edge);
1153 
1154  STD_REGIONS_EXPORT virtual int v_CalcNumberOfCoefficients(const std::vector<unsigned int> &nummodes, int &modes_offset);
1155 
1156  /**
1157  * @brief Unpack data from input file assuming it comes from the
1158  * same expansion type.
1159  * @see StdExpansion::ExtractDataToCoeffs
1160  */
1161  STD_REGIONS_EXPORT virtual void v_ExtractDataToCoeffs(const NekDouble *data,
1162  const std::vector<unsigned int > &nummodes,
1163  const int nmode_offset,
1164  NekDouble *coeffs);
1165 
1167 
1169 
1171 
1173 
1174 
1176 
1178 
1180 
1182 
1183  /** \brief Function to evaluate the discrete \f$ L_\infty\f$
1184  * error \f$ |\epsilon|_\infty = \max |u - u_{exact}|\f$ where \f$
1185  * u_{exact}\f$ is given by the array \a sol.
1186  *
1187  * This function takes the physical value space array \a m_phys as
1188  * approximate solution
1189  *
1190  * \param sol array of solution function at physical quadrature
1191  * points
1192  * \return returns the \f$ L_\infty \f$ error as a NekDouble.
1193  */
1195 
1196  /** \brief Function to evaluate the discrete \f$ L_2\f$ error,
1197  * \f$ | \epsilon |_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2
1198  * dx \right]^{1/2} d\xi_1 \f$ where \f$ u_{exact}\f$ is given by
1199  * the array \a sol.
1200  *
1201  * This function takes the physical value space array \a m_phys as
1202  * approximate solution
1203  *
1204  * \param sol array of solution function at physical quadrature
1205  * points
1206  * \return returns the \f$ L_2 \f$ error as a double.
1207  */
1209 
1210  /** \brief Function to evaluate the discrete \f$ H^1\f$
1211  * error, \f$ | \epsilon |^1_{2} = \left [ \int^1_{-1} [u -
1212  * u_{exact}]^2 + \nabla(u - u_{exact})\cdot\nabla(u -
1213  * u_{exact})\cdot dx \right]^{1/2} d\xi_1 \f$ where \f$
1214  * u_{exact}\f$ is given by the array \a sol.
1215  *
1216  * This function takes the physical value space array
1217  * \a m_phys as approximate solution
1218  *
1219  * \param sol array of solution function at physical quadrature
1220  * points
1221  * \return returns the \f$ H_1 \f$ error as a double.
1222  */
1224 
1225  // I/O routines
1226  const NormalVector & GetEdgeNormal(const int edge) const
1227  {
1228  return v_GetEdgeNormal(edge);
1229  }
1230 
1231  void ComputeEdgeNormal(const int edge)
1232  {
1233  v_ComputeEdgeNormal(edge);
1234  }
1235 
1236  void NegateEdgeNormal(const int edge)
1237  {
1238  v_NegateEdgeNormal(edge);
1239  }
1240 
1241  bool EdgeNormalNegated(const int edge)
1242  {
1243  return v_EdgeNormalNegated(edge);
1244  }
1245 
1246  void ComputeFaceNormal(const int face)
1247  {
1248  v_ComputeFaceNormal(face);
1249  }
1250 
1251  void NegateFaceNormal(const int face)
1252  {
1253  v_NegateFaceNormal(face);
1254  }
1255 
1256  void ComputeVertexNormal(const int vertex)
1257  {
1258  v_ComputeVertexNormal(vertex);
1259  }
1260 
1261  const NormalVector & GetFaceNormal(const int face) const
1262  {
1263  return v_GetFaceNormal(face);
1264  }
1265 
1266  const NormalVector & GetVertexNormal(const int vertex) const
1267  {
1268  return v_GetVertexNormal(vertex);
1269  }
1270 
1271  const NormalVector & GetSurfaceNormal(const int id) const
1272  {
1273  return v_GetSurfaceNormal(id);
1274  }
1275 
1277  {
1279  for (int i = 0; i < m_base.num_elements(); ++i)
1280  {
1281  p.push_back(m_base[i]->GetPointsKey());
1282  }
1283  return p;
1284  }
1285 
1288  {
1289  return v_GetEdgeInverseBoundaryMap(eid);
1290  }
1291 
1294  {
1295  return v_GetFaceInverseBoundaryMap(fid,faceOrient);
1296  }
1297 
1299  const DNekScalMatSharedPtr & m_transformationmatrix)
1300  {
1302  m_transformationmatrix);
1303  }
1304 
1305 
1306  /** \brief This function performs an interpolation from
1307  * the physical space points provided at input into an
1308  * array of equispaced points which are not the collapsed
1309  * coordinate. So for a tetrahedron you will only get a
1310  * tetrahedral number of values.
1311  *
1312  * This is primarily used for output purposes to get a
1313  * better distribution of points more suitable for most
1314  * postprocessing
1315  */
1317  const Array<OneD, const NekDouble> &inarray,
1318  Array<OneD, NekDouble> &outarray);
1319 
1320  /** \brief This function provides the connectivity of
1321  * local simplices (triangles or tets) to connect the
1322  * equispaced data points provided by
1323  * PhysInterpToSimplexEquiSpaced
1324  *
1325  * This is a virtual call to the function
1326  * \a v_GetSimplexEquiSpaceConnectivity
1327  */
1329  Array<OneD, int> &conn,
1330  bool standard = true)
1331  {
1332  v_GetSimplexEquiSpacedConnectivity(conn,standard);
1333  }
1334 
1335  template<class T>
1336  boost::shared_ptr<T> as()
1337  {
1338 #if defined __INTEL_COMPILER && BOOST_VERSION > 105200
1339  typedef typename boost::shared_ptr<T>::element_type E;
1340  E * p = dynamic_cast< E* >( shared_from_this().get() );
1341  ASSERTL1(p, "Cannot perform cast");
1342  return boost::shared_ptr<T>( shared_from_this(), p );
1343 #else
1344  return boost::dynamic_pointer_cast<T>( shared_from_this() );
1345 #endif
1346  }
1347 
1349  Array<OneD, NekDouble> &outarray,
1350  bool multiplybyweights = true)
1351  {
1352  v_IProductWRTBase_SumFac(inarray,outarray,multiplybyweights);
1353  }
1354 
1355  protected:
1356  Array<OneD, LibUtilities::BasisSharedPtr> m_base; /**< Bases needed for the expansion */
1358  int m_ncoeffs; /**< Total number of coefficients used in the expansion */
1362 
1364  {
1365  return v_CreateStdMatrix(mkey);
1366  }
1367 
1368  /** \brief Create the static condensation of a matrix when
1369  using a boundary interior decomposition
1370 
1371  If a matrix system can be represented by
1372  \f$ Mat = \left [ \begin{array}{cc}
1373  A & B \\
1374  C & D \end{array} \right ] \f$
1375  This routine creates a matrix containing the statically
1376  condense system of the form
1377  \f$ Mat = \left [ \begin{array}{cc}
1378  A - B D^{-1} C & B D^{-1} \\
1379  D^{-1} C & D^{-1} \end{array} \right ] \f$
1380  **/
1382 
1383  /** \brief Create an IndexMap which contains mapping information linking any specific
1384  element shape with either its boundaries, edges, faces, verteces, etc.
1385 
1386  The index member of the IndexMapValue struct gives back an integer associated with an entity index
1387  The sign member of the same struct gives back a sign to algebrically apply entities orientation
1388  **/
1390 
1392  Array<OneD, NekDouble> &outarray);
1393 
1395  Array<OneD, NekDouble> &outarray)
1396  {
1397  v_BwdTrans_SumFac(inarray,outarray);
1398  }
1399 
1400 
1401  void IProductWRTDerivBase_SumFac(const int dir,
1402  const Array<OneD, const NekDouble>& inarray,
1403  Array<OneD, NekDouble> &outarray)
1404  {
1405  v_IProductWRTDerivBase_SumFac(dir,inarray,outarray);
1406  }
1407 
1408  // The term _MatFree denotes that the action of the MatrixOperation
1409  // is done withouth actually using the matrix (which then needs to be stored/calculated).
1410  // Although this does not strictly mean that no matrix operations are involved in the
1411  // evaluation of the operation, we use this term in the same context used as in the following
1412  // paper:
1413  // R. C. Kirby, M. G. Knepley, A. Logg, and L. R. Scott,
1414  // "Optimizing the evaluation of finite element matrices," SISC 27:741-758 (2005)
1416  Array<OneD,NekDouble> &outarray,
1417  const StdMatrixKey &mkey);
1418 
1420  Array<OneD,NekDouble> &outarray,
1421  const StdMatrixKey &mkey);
1422 
1424  Array<OneD,NekDouble> &outarray,
1425  const StdMatrixKey &mkey)
1426  {
1427  v_LaplacianMatrixOp_MatFree(inarray,outarray,mkey);
1428  }
1429 
1431  const Array<OneD, const NekDouble> &inarray,
1432  Array<OneD, NekDouble> &outarray,
1434  {
1435  v_LaplacianMatrixOp_MatFree_Kernel(inarray, outarray, wsp);
1436  }
1437 
1439  Array<OneD,NekDouble> &outarray,
1440  const StdMatrixKey &mkey);
1441 
1442  STD_REGIONS_EXPORT void LaplacianMatrixOp_MatFree(const int k1, const int k2,
1443  const Array<OneD, const NekDouble> &inarray,
1444  Array<OneD,NekDouble> &outarray,
1445  const StdMatrixKey &mkey);
1446 
1448  const Array<OneD, const NekDouble> &inarray,
1449  Array<OneD,NekDouble> &outarray,
1450  const StdMatrixKey &mkey);
1451 
1453  Array<OneD,NekDouble> &outarray,
1454  const StdMatrixKey &mkey);
1455 
1457  Array<OneD,NekDouble> &outarray,
1458  const StdMatrixKey &mkey);
1459 
1461  Array<OneD,NekDouble> &outarray,
1462  const StdMatrixKey &mkey,
1463  bool addDiffusionTerm = true);
1464 
1466  Array<OneD,NekDouble> &outarray,
1467  const StdMatrixKey &mkey)
1468  {
1469  v_HelmholtzMatrixOp_MatFree(inarray,outarray,mkey);
1470  }
1471 
1473  Array<OneD,NekDouble> &outarray,
1474  const StdMatrixKey &mkey);
1475 
1478  Array<OneD, NekDouble> &outarray);
1479 
1481  Array<OneD, NekDouble> &coeffs,
1483 
1485  const Array<OneD, const NekDouble> &Lcoord,
1486  const Array<OneD, const NekDouble> &physvals);
1487 
1488  private:
1489  // Virtual functions
1490  STD_REGIONS_EXPORT virtual int v_GetNverts() const = 0;
1491  STD_REGIONS_EXPORT virtual int v_GetNedges() const;
1492 
1493  STD_REGIONS_EXPORT virtual int v_GetNfaces() const;
1494 
1495  STD_REGIONS_EXPORT virtual int v_NumBndryCoeffs() const;
1496 
1497  STD_REGIONS_EXPORT virtual int v_NumDGBndryCoeffs() const;
1498 
1499  STD_REGIONS_EXPORT virtual int v_GetEdgeNcoeffs(const int i) const;
1500 
1501  STD_REGIONS_EXPORT virtual int v_GetTotalEdgeIntNcoeffs() const;
1502 
1503  STD_REGIONS_EXPORT virtual int v_GetEdgeNumPoints(const int i) const;
1504 
1505  STD_REGIONS_EXPORT virtual int v_DetCartesianDirOfEdge(const int edge);
1506 
1507  STD_REGIONS_EXPORT virtual const LibUtilities::BasisKey v_DetEdgeBasisKey(const int i) const;
1508 
1509  STD_REGIONS_EXPORT virtual const LibUtilities::BasisKey v_DetFaceBasisKey(const int i, const int k) const;
1510 
1511  STD_REGIONS_EXPORT virtual int v_GetFaceNumPoints(const int i) const;
1512 
1513  STD_REGIONS_EXPORT virtual int v_GetFaceNcoeffs(const int i) const;
1514 
1515  STD_REGIONS_EXPORT virtual int v_GetFaceIntNcoeffs(const int i) const;
1516 
1517  STD_REGIONS_EXPORT virtual int v_GetTotalFaceIntNcoeffs() const;
1518 
1519  STD_REGIONS_EXPORT virtual LibUtilities::PointsKey v_GetFacePointsKey(const int i, const int j) const;
1520 
1522 
1524 
1526 
1527  STD_REGIONS_EXPORT virtual boost::shared_ptr<StdExpansion>
1528  v_GetStdExp(void) const;
1529 
1530  STD_REGIONS_EXPORT virtual int v_GetShapeDimension() const;
1531 
1533 
1535 
1536  STD_REGIONS_EXPORT virtual void v_BwdTrans (const Array<OneD, const NekDouble>& inarray,
1537  Array<OneD, NekDouble> &outarray) = 0;
1538 
1539  /**
1540  * @brief Transform a given function from physical quadrature space
1541  * to coefficient space.
1542  * @see StdExpansion::FwdTrans
1543  */
1544  STD_REGIONS_EXPORT virtual void v_FwdTrans (
1545  const Array<OneD, const NekDouble>& inarray,
1546  Array<OneD, NekDouble> &outarray) = 0;
1547 
1548  /**
1549  * @brief Calculates the inner product of a given function \a f
1550  * with the different modes of the expansion
1551  */
1553  const Array<OneD, const NekDouble>& inarray,
1554  Array<OneD, NekDouble> &outarray) = 0;
1555 
1557  const Array<OneD, const NekDouble>& base,
1558  const Array<OneD, const NekDouble>& inarray,
1559  Array<OneD, NekDouble>& outarray,
1560  int coll_check)
1561  {
1562  ASSERTL0(false, "StdExpansion::v_IProductWRTBase has no (and should have no) implementation");
1563  }
1564 
1565  STD_REGIONS_EXPORT virtual void v_IProductWRTDerivBase (const int dir,
1566  const Array<OneD, const NekDouble>& inarray,
1567  Array<OneD, NekDouble> &outarray);
1568 
1570  Array<OneD, NekDouble> &outarray);
1571 
1573 
1574  STD_REGIONS_EXPORT virtual void v_PhysDeriv (const Array<OneD, const NekDouble>& inarray,
1575  Array<OneD, NekDouble> &out_d1,
1576  Array<OneD, NekDouble> &out_d2,
1577  Array<OneD, NekDouble> &out_d3);
1578 
1579  STD_REGIONS_EXPORT virtual void v_PhysDeriv_s (const Array<OneD, const NekDouble>& inarray,
1580  Array<OneD, NekDouble> &out_ds);
1581 
1583  Array<OneD, NekDouble>& out_dn);
1584  STD_REGIONS_EXPORT virtual void v_PhysDeriv(const int dir,
1585  const Array<OneD, const NekDouble>& inarray,
1586  Array<OneD, NekDouble> &out_d0);
1587 
1589  const Array<OneD, const NekDouble>& direction,
1590  Array<OneD, NekDouble> &outarray);
1591 
1592  STD_REGIONS_EXPORT virtual void v_StdPhysDeriv (const Array<OneD, const NekDouble>& inarray,
1593  Array<OneD, NekDouble> &out_d1,
1594  Array<OneD, NekDouble> &out_d2,
1595  Array<OneD, NekDouble> &out_d3);
1596 
1597  STD_REGIONS_EXPORT virtual void v_StdPhysDeriv (const int dir,
1598  const Array<OneD, const NekDouble>& inarray,
1599  Array<OneD, NekDouble> &outarray);
1600 
1601  STD_REGIONS_EXPORT virtual void v_AddRobinMassMatrix(const int edgeid, const Array<OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat);
1602 
1603  STD_REGIONS_EXPORT virtual void v_AddRobinEdgeContribution(const int edgeid, const Array<OneD, const NekDouble> &primCoeffs, Array<OneD, NekDouble> &coeffs);
1604 
1606 
1608 
1610  const Array<OneD, const NekDouble>& xi,
1611  Array<OneD, NekDouble>& eta);
1612 
1613 
1614  STD_REGIONS_EXPORT virtual void v_FillMode(const int mode, Array<OneD, NekDouble> &outarray);
1615 
1617 
1619 
1620  STD_REGIONS_EXPORT virtual void v_GetCoords(Array<OneD, NekDouble> &coords_0,
1621  Array<OneD, NekDouble> &coords_1,
1622  Array<OneD, NekDouble> &coords_2);
1623 
1624  STD_REGIONS_EXPORT virtual void v_GetCoord(const Array<OneD, const NekDouble>& Lcoord,
1625  Array<OneD, NekDouble> &coord);
1626 
1627  STD_REGIONS_EXPORT virtual int v_GetCoordim(void);
1628 
1630 
1632 
1633  STD_REGIONS_EXPORT virtual int v_GetVertexMap(int localVertexId,
1634  bool useCoeffPacking = false);
1635 
1636  STD_REGIONS_EXPORT virtual void v_GetEdgeInteriorMap(const int eid, const Orientation edgeOrient,
1637  Array<OneD, unsigned int> &maparray,
1638  Array<OneD, int> &signarray);
1639 
1640  STD_REGIONS_EXPORT virtual void v_GetFaceInteriorMap(const int fid, const Orientation faceOrient,
1641  Array<OneD, unsigned int> &maparray,
1642  Array<OneD, int> &signarray);
1643 
1644  STD_REGIONS_EXPORT virtual void v_GetEdgeToElementMap(const int eid, const Orientation edgeOrient,
1645  Array<OneD, unsigned int> &maparray,
1646  Array<OneD, int> &signarray);
1647 
1648  STD_REGIONS_EXPORT virtual void v_GetFaceToElementMap(const int fid, const Orientation faceOrient,
1649  Array<OneD, unsigned int> &maparray,
1650  Array<OneD, int> &signarray,
1651  int nummodesA = -1, int nummodesB = -1);
1652 
1653  /**
1654  * @brief Extract the physical values along edge \a edge from \a
1655  * inarray into \a outarray following the local edge orientation
1656  * and point distribution defined by defined in \a EdgeExp.
1657  */
1658  STD_REGIONS_EXPORT virtual void v_GetEdgePhysVals(const int edge, const Array<OneD, const NekDouble> &inarray, Array<OneD,NekDouble> &outarray);
1659 
1660  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);
1661 
1662  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);
1663 
1664  STD_REGIONS_EXPORT virtual void v_GetVertexPhysVals(const int vertex, const Array<OneD, const NekDouble> &inarray, NekDouble &outarray);
1665 
1666  STD_REGIONS_EXPORT virtual void v_GetEdgeInterpVals(const int edge,
1667  const Array<OneD, const NekDouble> &inarray,Array<OneD,NekDouble> &outarray);
1668 
1670  const int edge,
1671  Array<OneD, NekDouble> &outarray);
1672 
1674  const int face,
1675  const boost::shared_ptr<StdExpansion> &FaceExp,
1676  const Array<OneD, const NekDouble> &inarray,
1677  Array<OneD, NekDouble> &outarray,
1678  StdRegions::Orientation orient);
1679 
1681  const Array<OneD, const NekDouble> &inarray,
1682  Array<OneD, NekDouble> &outarray);
1683 
1685  const Array<OneD, const NekDouble> &inarray,
1686  Array<OneD, NekDouble> &outarray);
1687 
1688  STD_REGIONS_EXPORT virtual const boost::shared_ptr<SpatialDomains::GeomFactors>& v_GetMetricInfo() const;
1689 
1691  Array<OneD, NekDouble> &outarray);
1692 
1694  Array<OneD, NekDouble> &outarray, bool multiplybyweights = true);
1695 
1696  STD_REGIONS_EXPORT virtual void v_IProductWRTDerivBase_SumFac(const int dir,
1697  const Array<OneD, const NekDouble>& inarray,
1698  Array<OneD, NekDouble> &outarray);
1699 
1701  Array<OneD,NekDouble> &outarray,
1702  const StdMatrixKey &mkey);
1703 
1705  Array<OneD,NekDouble> &outarray,
1706  const StdMatrixKey &mkey);
1707 
1709  const StdMatrixKey &mkey);
1710 
1712  int numMin,
1713  const Array<OneD, const NekDouble> &inarray,
1714  Array<OneD,NekDouble> &outarray);
1715 
1716  STD_REGIONS_EXPORT virtual void v_LaplacianMatrixOp(const int k1, const int k2,
1717  const Array<OneD, const NekDouble> &inarray,
1718  Array<OneD,NekDouble> &outarray,
1719  const StdMatrixKey &mkey);
1720 
1721  STD_REGIONS_EXPORT virtual void v_WeakDerivMatrixOp(const int i,
1722  const Array<OneD, const NekDouble> &inarray,
1723  Array<OneD,NekDouble> &outarray,
1724  const StdMatrixKey &mkey);
1725 
1727  Array<OneD,NekDouble> &outarray,
1728  const StdMatrixKey &mkey);
1729 
1731  Array<OneD,NekDouble> &outarray,
1732  const StdMatrixKey &mkey);
1733 
1735  const NekDouble> &inarray,
1736  Array<OneD,NekDouble> &outarray,
1737  const StdMatrixKey &mkey,
1738  bool addDiffusionTerm=true);
1739 
1741  Array<OneD,NekDouble> &outarray,
1742  const StdMatrixKey &mkey);
1743 
1745  Array<OneD,NekDouble> &outarray,
1746  const StdMatrixKey &mkey);
1747 
1749  const Array<OneD, const NekDouble> &inarray,
1750  Array<OneD, NekDouble> &outarray,
1751  Array<OneD, NekDouble> &wsp);
1752 
1754  Array<OneD,NekDouble> &outarray,
1755  const StdMatrixKey &mkey);
1756 
1757  STD_REGIONS_EXPORT virtual const NormalVector & v_GetEdgeNormal(const int edge) const;
1758 
1759  STD_REGIONS_EXPORT virtual void v_ComputeEdgeNormal(const int edge);
1760 
1761  STD_REGIONS_EXPORT virtual void v_NegateEdgeNormal(const int edge);
1762 
1763  STD_REGIONS_EXPORT virtual bool v_EdgeNormalNegated(const int edge);
1764 
1765  STD_REGIONS_EXPORT virtual void v_ComputeFaceNormal(const int face);
1766 
1767  STD_REGIONS_EXPORT virtual void v_NegateFaceNormal(const int face);
1768 
1769  STD_REGIONS_EXPORT virtual const NormalVector & v_GetVertexNormal(const int vertex) const;
1770 
1771  STD_REGIONS_EXPORT virtual void v_ComputeVertexNormal(const int vertex);
1772 
1773  STD_REGIONS_EXPORT virtual const NormalVector & v_GetFaceNormal(const int face) const;
1774  STD_REGIONS_EXPORT virtual const NormalVector &
1775  v_GetSurfaceNormal(const int id) const;
1776 
1778  v_GetEdgeInverseBoundaryMap(int eid);
1779 
1782 
1784 
1786  Array<OneD, int> &conn,
1787  bool standard = true);
1788  };
1789 
1790 
1791  typedef boost::shared_ptr<StdExpansion> StdExpansionSharedPtr;
1792  typedef std::vector< StdExpansionSharedPtr > StdExpansionVector;
1794 
1795  /**
1796  * This function is a wrapper around the virtual function
1797  * \a v_FwdTrans()
1798  *
1799  * Given a function evaluated at the quadrature points, this
1800  * function calculates the expansion coefficients such that the
1801  * resulting expansion approximates the original function.
1802  *
1803  * The calculation of the expansion coefficients is done using a
1804  * Galerkin projection. This is equivalent to the operation:
1805  * \f[ \mathbf{\hat{u}} = \mathbf{M}^{-1} \mathbf{I}\f]
1806  * where
1807  * - \f$\mathbf{M}[p][q]= \int\phi_p(\mathbf{\xi})\phi_q(
1808  * \mathbf{\xi}) d\mathbf{\xi}\f$ is the Mass matrix
1809  * - \f$\mathbf{I}[p] = \int\phi_p(\mathbf{\xi}) u(\mathbf{\xi})
1810  * d\mathbf{\xi}\f$
1811  *
1812  * This function takes the array \a inarray as the values of the
1813  * function evaluated at the quadrature points
1814  * (i.e. \f$\mathbf{u}\f$),
1815  * and stores the resulting coefficients \f$\mathbf{\hat{u}}\f$
1816  * in the \a outarray
1817  *
1818  * @param inarray array of the function discretely evaluated at the
1819  * quadrature points
1820  *
1821  * @param outarray array of the function coefficieints
1822  */
1824  Array<OneD, NekDouble> &outarray)
1825  {
1826  v_FwdTrans(inarray,outarray);
1827  }
1828 
1829  } //end of namespace
1830 } //end of namespace
1831 
1832 #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:454
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
virtual void v_NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
virtual void v_GetInteriorMap(Array< OneD, unsigned int > &outarray)
const NormalVector & GetEdgeNormal(const int edge) const
StdRegions::Orientation GetCartesianEorient(int edge)
Definition: StdExpansion.h:746
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:459
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:638
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:576
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:931
int GetNfaces() const
This function returns the number of faces of the expansion domain.
Definition: StdExpansion.h:422
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:741
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:699
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:694
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:902
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:806
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...
void GetEdgeToElementMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
Definition: StdExpansion.h:820
virtual void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray)
void MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:909
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
void NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:711
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:613
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:558
void SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
Definition: StdExpansion.h:952
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:726
void GetEdgePhysVals(const int edge, const boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:850
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
StdRegions::Orientation GetEorient(int edge)
Definition: StdExpansion.h:736
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.
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:409
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:731
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 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:883
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix(const LocalRegions::MatrixKey &mkey)
Definition: StdExpansion.h:721
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:945
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:684
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:660
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
const NormalVector & GetFaceNormal(const int face) const
void IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
void PhysInterpToSimplexEquiSpaced(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs an interpolation from the physical space points provided at input into an arra...
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:704
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)
Definition: StdExpansion.h:966
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)
Definition: StdExpansion.h:981
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:141
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:397
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:689
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:800
void FwdTrans_BndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:522
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:827
void GetVertexPhysVals(const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray)
Definition: StdExpansion.h:863
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:958
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:678
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:795
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:369
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:646
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:892
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:858
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:629
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)
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:843
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:938
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)
#define STD_REGIONS_EXPORT
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:716
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:766
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:509
int GetNtrace() const
Returns the number of trace elements connected to this element.
Definition: StdExpansion.h:434
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:751
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:619
boost::shared_ptr< Basis > BasisSharedPtr
void WeakDirectionalDerivMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:974
virtual void v_GetEdgeToElementMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
void GetEdgeInterpVals(const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:870
virtual int v_GetNfaces() const
void SetCoeffsToOrientation(StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:758
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:771
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:790
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:813
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)
void LinearAdvectionDiffusionReactionMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
Definition: StdExpansion.h:988
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)