Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
StdSegExp.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File StdSegExp.cpp
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // 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: Routines within Standard Segment Expansions
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #include <StdRegions/StdSegExp.h>
38 
39 
40 using namespace std;
41 
42 namespace Nektar
43 {
44  namespace StdRegions
45  {
46 
47  /** \brief Default constructor */
48 
49  StdSegExp::StdSegExp()
50  {
51  }
52 
53 
54  /** \brief Constructor using BasisKey class for quadrature points and
55  * order definition
56  *
57  * \param Ba BasisKey class definition containing order and quadrature
58  * points.
59  */
60 
61  StdSegExp::StdSegExp(const LibUtilities::BasisKey &Ba):
62  StdExpansion(Ba.GetNumModes(), 1, Ba),
63  StdExpansion1D(Ba.GetNumModes(),Ba)
64  {
65  }
66 
67 
68  /** \brief Copy Constructor */
69 
71  StdExpansion(T),
73  {
74  }
75 
76 
78  {
79  }
80 
81  /** \brief Return Shape of region, using ShapeType enum list.
82  * i.e. Segment
83  */
85  {
87  }
88 
90  {
91 
92  bool returnval = false;
93 
95  {
96  returnval = true;
97  }
98 
100  {
101  returnval = true;
102  }
103 
104  return returnval;
105  }
106 
107 
108 
109 
110  //---------------------------------------------------------------------
111  // Integration Methods
112  //---------------------------------------------------------------------
113 
114  /** \brief Integrate the physical point list \a inarray over region
115  * and return the value
116  *
117  * \param inarray definition of function to be integrated evauluated at
118  * quadrature point of expansion.
119  * \return returns \f$\int^1_{-1} u(\xi_1)d \xi_1 \f$ where \f$inarray[i]
120  * = u(\xi_{1i}) \f$
121  */
123  {
124  NekDouble Int = 0.0;
125  int nquad0 = m_base[0]->GetNumPoints();
126  Array<OneD, NekDouble> tmp(nquad0);
127  Array<OneD, const NekDouble> z = m_base[0]->GetZ();
128  Array<OneD, const NekDouble> w0 = m_base[0]->GetW();
129 
130  // multiply by integration constants
131  Vmath::Vmul(nquad0, inarray, 1, w0, 1, tmp, 1);
132 
133  Int = Vmath::Vsum(nquad0, tmp, 1);
134 
135  return Int;
136  }
137 
138 
139 
140 
141  //---------------------------------------------------------------------
142  // Differentiation Methods
143  //---------------------------------------------------------------------
144 
145 
146  /** \brief Evaluate the derivative \f$ d/d{\xi_1} \f$ at the physical
147  * quadrature points given by \a inarray and return in \a outarray.
148  *
149  * This is a wrapper around StdExpansion1D::Tensor_Deriv
150  * \param inarray array of a function evaluated at the quadrature points
151  * \param outarray the resulting array of the derivative \f$
152  * du/d_{\xi_1}|_{\xi_{1i}} \f$ will be stored in the array \a outarra
153  */
154 
156  Array<OneD, NekDouble> &out_d0,
157  Array<OneD, NekDouble> &out_d1,
158  Array<OneD, NekDouble> &out_d2)
159  {
160  PhysTensorDeriv(inarray,out_d0);
161  }
162 
163 
164  void StdSegExp::v_PhysDeriv(const int dir,
165  const Array<OneD, const NekDouble>& inarray,
166  Array<OneD, NekDouble> &outarray)
167  {
168  ASSERTL1(dir==0,"input dir is out of range");
169  PhysTensorDeriv(inarray,outarray);
170  // PhysDeriv(inarray, outarray);
171  }
172 
174  const Array<OneD, const NekDouble>& inarray,
175  Array<OneD, NekDouble> &out_d0,
176  Array<OneD, NekDouble> &out_d1,
177  Array<OneD, NekDouble> &out_d2)
178  {
179  PhysTensorDeriv(inarray,out_d0);
180  // PhysDeriv(inarray, out_d0);
181  }
182 
184  const int dir,
185  const Array<OneD, const NekDouble>& inarray,
186  Array<OneD, NekDouble> &outarray)
187  {
188  ASSERTL1(dir==0,"input dir is out of range");
189  PhysTensorDeriv(inarray,outarray);
190  // PhysDeriv(inarray, outarray);
191  }
192 
193 
194 
195  //---------------------------------------------------------------------
196  // Transforms
197  //---------------------------------------------------------------------
198 
199  /** \brief Backward transform from coefficient space given
200  * in \a inarray and evaluate at the physical quadrature
201  * points \a outarray
202  *
203  * Operation can be evaluated as \f$ u(\xi_{1i}) =
204  * \sum_{p=0}^{order-1} \hat{u}_p \phi_p(\xi_{1i}) \f$ or equivalently
205  * \f$ {\bf u} = {\bf B}^T {\bf \hat{u}} \f$ where
206  * \f${\bf B}[i][j] = \phi_i(\xi_{1j}), \mbox{\_coeffs}[p] = {\bf
207  * \hat{u}}[p] \f$
208  *
209  * The function takes the coefficient array \a inarray as
210  * input for the transformation
211  *
212  * \param inarray: the coeffficients of the expansion
213  *
214  * \param outarray: the resulting array of the values of the function at
215  * the physical quadrature points will be stored in the array \a outarray
216  */
217 
219  const Array<OneD, const NekDouble>& inarray,
220  Array<OneD, NekDouble> &outarray)
221  {
222  int nquad = m_base[0]->GetNumPoints();
223 
224  if(m_base[0]->Collocation())
225  {
226  Vmath::Vcopy(nquad, inarray, 1, outarray, 1);
227  }
228  else
229  {
230 
231 #ifdef NEKTAR_USING_DIRECT_BLAS_CALLS
232 
233  Blas::Dgemv('N',nquad,m_base[0]->GetNumModes(),1.0, (m_base[0]->GetBdata()).get(),
234  nquad,&inarray[0],1,0.0,&outarray[0],1);
235 
236 #else //NEKTAR_USING_DIRECT_BLAS_CALLS
237 
239  NekVector<NekDouble> out(nquad,outarray,eWrapper);
240  NekMatrix<NekDouble> B(nquad,m_ncoeffs,m_base[0]->GetBdata(),eWrapper);
241  out = B * in;
242 
243 #endif //NEKTAR_USING_DIRECT_BLAS_CALLS
244 
245  }
246  }
247 
249  int numMin,
250  const Array<OneD, const NekDouble> &inarray,
251  Array<OneD, NekDouble> &outarray)
252  {
253  int n_coeffs = inarray.num_elements();
254 
255  Array<OneD, NekDouble> coeff(n_coeffs);
256  Array<OneD, NekDouble> coeff_tmp(n_coeffs,0.0);
259 
260  int nmodes0 = m_base[0]->GetNumModes();
261  int numMax = nmodes0;
262 
263  Vmath::Vcopy(n_coeffs,inarray,1,coeff_tmp,1);
264 
265  const LibUtilities::PointsKey Pkey0(
267 
268  LibUtilities::BasisKey b0(m_base[0]->GetBasisType(),nmodes0,Pkey0);
269 
270  LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A,nmodes0,Pkey0);
271 
273  b0, coeff_tmp, bortho0, coeff);
274 
275  Vmath::Zero(n_coeffs,coeff_tmp,1);
276 
277  Vmath::Vcopy(numMin,
278  tmp = coeff,1,
279  tmp2 = coeff_tmp,1);
280 
282  bortho0, coeff_tmp, b0, outarray);
283  }
284 
285  /**
286  * \brief Forward transform from physical quadrature space stored in \a
287  * inarray and evaluate the expansion coefficients and store in \a
288  * outarray
289  *
290  * Perform a forward transform using a Galerkin projection by taking the
291  * inner product of the physical points and multiplying by the inverse
292  * of the mass matrix using the Solve method of the standard matrix
293  * container holding the local mass matrix, i.e. \f$ {\bf \hat{u}} =
294  * {\bf M}^{-1} {\bf I} \f$ where \f$ {\bf I}[p] = \int^1_{-1}
295  * \phi_p(\xi_1) u(\xi_1) d\xi_1 \f$
296  *
297  * This function stores the expansion coefficients calculated by the
298  * transformation in the coefficient space array \a outarray
299  *
300  * \param inarray: array of physical quadrature points to be transformed
301  * \param outarray: the coeffficients of the expansion
302  */
304  Array<OneD, NekDouble> &outarray)
305  {
306  if(m_base[0]->Collocation())
307  {
308  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
309  }
310  else
311  {
312  v_IProductWRTBase(inarray,outarray);
313 
314  // get Mass matrix inverse
315  StdMatrixKey masskey(eInvMass,v_DetShapeType(),*this);
316  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
317 
318  NekVector<NekDouble> in(m_ncoeffs,outarray,eCopy);
320 
321  out = (*matsys)*in;
322  }
323  }
324 
326  const Array<OneD, const NekDouble>& inarray,
327  Array<OneD, NekDouble> &outarray)
328  {
329  if(m_base[0]->Collocation())
330  {
331  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
332  }
333  else
334  {
335  int nInteriorDofs = m_ncoeffs-2;
336  int offset;
337 
338  switch(m_base[0]->GetBasisType())
339  {
341  {
342  offset = 1;
343  }
344  break;
346  {
347  nInteriorDofs = m_ncoeffs;
348  offset = 0;
349  }
350  break;
353  {
354  offset = 2;
355  }
356  break;
357  default:
358  ASSERTL0(false,"This type of FwdTrans is not defined for this expansion type");
359  }
360 
361  fill(outarray.get(), outarray.get()+m_ncoeffs, 0.0 );
362 
364  {
365  outarray[GetVertexMap(0)] = inarray[0];
366  outarray[GetVertexMap(1)] = inarray[m_base[0]->GetNumPoints()-1];
367 
368  if(m_ncoeffs>2)
369  {
370  // ideally, we would like to have tmp0 to be replaced by
371  // outarray (currently MassMatrixOp does not allow aliasing)
374 
375  StdMatrixKey masskey(eMass,v_DetShapeType(),*this);
376  MassMatrixOp(outarray,tmp0,masskey);
377  v_IProductWRTBase(inarray,tmp1);
378 
379  Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
380 
381  // get Mass matrix inverse (only of interior DOF)
382  DNekMatSharedPtr matsys =
383  (m_stdStaticCondMatrixManager[masskey])-> GetBlock(1,1);
384 
385  Blas::Dgemv('N',nInteriorDofs,nInteriorDofs,1.0,
386  &(matsys->GetPtr())[0],nInteriorDofs,tmp1.get()
387  +offset,1,0.0,outarray.get()+offset,1);
388  }
389  }
390  else
391  {
392  StdSegExp::v_FwdTrans(inarray, outarray);
393  }
394  }
395 
396  }
397 
398 
400  Array<OneD, NekDouble> &outarray)
401  {
402  v_BwdTrans(inarray, outarray);
403  }
404 
405 
406 
407  //---------------------------------------------------------------------
408  // Inner product functions
409  //---------------------------------------------------------------------
410 
411 
412 
413  /** \brief Inner product of \a inarray over region with respect to
414  * expansion basis \a base and return in \a outarray
415  *
416  * Calculate \f$ I[p] = \int^{1}_{-1} \phi_p(\xi_1) u(\xi_1) d\xi_1
417  * = \sum_{i=0}^{nq-1} \phi_p(\xi_{1i}) u(\xi_{1i}) w_i \f$ where
418  * \f$ outarray[p] = I[p], inarray[i] = u(\xi_{1i}), base[p*nq+i] =
419  * \phi_p(\xi_{1i}) \f$.
420  *
421  * \param base an array defining the local basis for the inner product
422  * usually passed from Basis->GetBdata() or Basis->GetDbdata()
423  * \param inarray: physical point array of function to be integrated
424  * \f$ u(\xi_1) \f$
425  * \param coll_check flag to identify when a Basis->Collocation() call
426  * should be performed to see if this is a GLL_Lagrange basis with a
427  * collocation property. (should be set to 0 if taking the inner
428  * product with respect to the derivative of basis)
429  * \param outarray the values of the inner product with respect to
430  * each basis over region will be stored in the array \a outarray as
431  * output of the function
432  */
433 
435  const Array<OneD, const NekDouble> &base,
436  const Array<OneD, const NekDouble> &inarray,
437  Array<OneD, NekDouble> &outarray,
438  int coll_check)
439  {
440  int nquad = m_base[0]->GetNumPoints();
441  Array<OneD, NekDouble> tmp(nquad);
442  Array<OneD, const NekDouble> w = m_base[0]->GetW();
443 
444  Vmath::Vmul(nquad, inarray, 1, w, 1, tmp, 1);
445 
446  /* Comment below was a bug for collocated basis
447  if(coll_check&&m_base[0]->Collocation())
448  {
449  Vmath::Vcopy(nquad, tmp, 1, outarray, 1);
450  }
451  else
452  {
453  Blas::Dgemv('T',nquad,m_ncoeffs,1.0,base.get(),nquad,
454  &tmp[0],1,0.0,outarray.get(),1);
455  }*/
456 
457  // Correct implementation
458  Blas::Dgemv('T',nquad,m_ncoeffs,1.0,base.get(),nquad,
459  &tmp[0],1,0.0,outarray.get(),1);
460  }
461 
462  /** \brief Inner product of \a inarray over region with respect to the
463  * expansion basis (this)->m_base[0] and return in \a outarray
464  *
465  * Wrapper call to StdSegExp::IProductWRTBase
466  * \param inarray array of function values evaluated at the physical
467  * collocation points
468  * \param outarray the values of the inner product with respect to
469  * each basis over region will be stored in the array \a outarray as
470  * output of the function
471  */
473  Array<OneD, NekDouble> &outarray)
474  {
475  v_IProductWRTBase(m_base[0]->GetBdata(),inarray,outarray,1);
476  }
477 
479  const int dir,
480  const Array<OneD, const NekDouble>& inarray,
481  Array<OneD, NekDouble> & outarray)
482  {
483  ASSERTL1(dir >= 0 && dir < 1,"input dir is out of range");
484  v_IProductWRTBase(m_base[0]->GetDbdata(),inarray,outarray,1);
485  }
486 
488  const Array<OneD, const NekDouble>& inarray,
489  Array<OneD, NekDouble> &outarray,
490  bool multiplybyweights)
491  {
492  int nquad = m_base[0]->GetNumPoints();
493  Array<OneD, NekDouble> tmp(nquad);
494  Array<OneD, const NekDouble> w = m_base[0]->GetW();
495  Array<OneD, const NekDouble> base = m_base[0]->GetBdata();
496 
497  if(multiplybyweights)
498  {
499  Vmath::Vmul(nquad, inarray, 1, w, 1, tmp, 1);
500 
501  Blas::Dgemv('T',nquad,m_ncoeffs,1.0,base.get(),nquad,
502  &tmp[0],1,0.0,outarray.get(),1);
503  }
504  else
505  {
506  Blas::Dgemv('T',nquad,m_ncoeffs,1.0,base.get(),nquad,
507  &inarray[0],1,0.0,outarray.get(),1);
508  }
509  }
510 
511 
512 
513 
514 
515  //----------------------------
516  // Evaluation
517  //----------------------------
518 
519  void StdSegExp::v_FillMode(const int mode, Array<OneD, NekDouble> &outarray)
520  {
521  int nquad = m_base[0]->GetNumPoints();
522  const NekDouble * base = m_base[0]->GetBdata().get();
523 
524  ASSERTL2(mode <= m_ncoeffs,
525  "calling argument mode is larger than total expansion order");
526 
527  Vmath::Vcopy(nquad,(NekDouble *)base+mode*nquad,1, &outarray[0],1);
528  }
529 
530 
532  const Array<OneD, const NekDouble>& coords,
533  const Array<OneD, const NekDouble>& physvals)
534  {
535  return StdExpansion1D::v_PhysEvaluate(coords, physvals);
536  }
537 
539  const Array<OneD, const NekDouble> &inarray,
540  Array<OneD, NekDouble> &outarray,
541  const StdMatrixKey &mkey)
542  {
543  int nquad = m_base[0]->GetNumPoints();
544 
545  Array<OneD,NekDouble> physValues(nquad);
546  Array<OneD,NekDouble> dPhysValuesdx(nquad);
547 
548  v_BwdTrans(inarray,physValues);
549 
550  // Laplacian matrix operation
551  v_PhysDeriv(physValues,dPhysValuesdx);
552  v_IProductWRTBase(m_base[0]->GetDbdata(),dPhysValuesdx,outarray,1);
553  }
554 
555 
557  const Array<OneD, const NekDouble> &inarray,
558  Array<OneD, NekDouble> &outarray,
559  const StdMatrixKey &mkey)
560  {
561  int nquad = m_base[0]->GetNumPoints();
562 
563  Array<OneD,NekDouble> physValues(nquad);
564  Array<OneD,NekDouble> dPhysValuesdx(nquad);
566 
567  v_BwdTrans(inarray,physValues);
568 
569  // mass matrix operation
570  v_IProductWRTBase((m_base[0]->GetBdata()),physValues,wsp,1);
571 
572  // Laplacian matrix operation
573  v_PhysDeriv(physValues,dPhysValuesdx);
574  v_IProductWRTBase(m_base[0]->GetDbdata(),dPhysValuesdx,outarray,1);
575  Blas::Daxpy(m_ncoeffs, mkey.GetConstFactor(eFactorLambda), wsp.get(), 1, outarray.get(), 1);
576  }
577 
578  //up to here
580  const Array<OneD, const NekDouble> &inarray,
581  Array<OneD, NekDouble> &outarray)
582  {
583  int nquad0 = m_base[0]->GetNumPoints();
584 
585  const Array<OneD, const NekDouble>& w0 = m_base[0]->GetW();
586 
587  Vmath::Vmul(nquad0, inarray.get(),1,
588  w0.get(),1,outarray.get(),1);
589  }
590 
592  Array<OneD, NekDouble> &coords_0,
593  Array<OneD, NekDouble> &coords_1,
594  Array<OneD, NekDouble> &coords_2)
595  {
596  Blas::Dcopy(GetNumPoints(0),(m_base[0]->GetZ()).get(),
597  1,&coords_0[0],1);
598  }
599 
600 
601 
602 
603 
604  //---------------------------------------------------------------------
605  // Helper functions
606  //---------------------------------------------------------------------
607 
609  {
610  return 2;
611  }
612 
614  {
615  return 2;
616  }
617 
619  {
620  return 2;
621  }
622 
624  const std::vector<unsigned int> &nummodes,
625  int &modes_offset)
626  {
627  int nmodes = nummodes[modes_offset];
628  modes_offset += 1;
629 
630  return nmodes;
631  }
632 
633  //---------------------------------------------------------------------
634  // Wrapper functions
635  //---------------------------------------------------------------------
636 
638  {
639  DNekMatSharedPtr Mat;
640  MatrixType mattype;
641 
642  switch(mattype = mkey.GetMatrixType())
643  {
644  case eFwdTrans:
645  {
647  StdMatrixKey iprodkey(eIProductWRTBase,v_DetShapeType(),*this);
648  DNekMat &Iprod = *GetStdMatrix(iprodkey);
649  StdMatrixKey imasskey(eInvMass,v_DetShapeType(),*this);
650  DNekMat &Imass = *GetStdMatrix(imasskey);
651 
652  (*Mat) = Imass*Iprod;
653  }
654  break;
655  default:
656  {
658 
659  if(mattype == eMass)
660  {
661  // For Fourier basis set the imaginary component
662  // of mean mode to have a unit diagonal component
663  // in mass matrix
665  {
666  (*Mat)(1,1) = 1.0;
667  }
668  }
669  }
670  break;
671  }
672 
673  return Mat;
674  }
675 
676 
678  {
679  return v_GenMatrix(mkey);
680  }
681 
682  //---------------------------------------------------------------------
683  // Mappings
684  //---------------------------------------------------------------------
685 
686 
688  {
689  if(outarray.num_elements() != NumBndryCoeffs())
690  {
692  }
693  const LibUtilities::BasisType Btype = GetBasisType(0);
694  int nummodes = m_base[0]->GetNumModes();
695 
696  outarray[0] = 0;
697 
698  switch(Btype)
699  {
704  outarray[1]= nummodes-1;
705  break;
708  outarray[1] = 1;
709  break;
710  default:
711  ASSERTL0(0,"Mapping array is not defined for this expansion");
712  break;
713  }
714  }
715 
717  {
718  int i;
719  if(outarray.num_elements()!=GetNcoeffs()-NumBndryCoeffs())
720  {
722  }
723  const LibUtilities::BasisType Btype = GetBasisType(0);
724 
725  switch(Btype)
726  {
731  for(i = 0 ; i < GetNcoeffs()-2;i++)
732  {
733  outarray[i] = i+1;
734  }
735  break;
738  for(i = 0 ; i < GetNcoeffs()-2;i++)
739  {
740  outarray[i] = i+2;
741  }
742  break;
743  default:
744  ASSERTL0(0,"Mapping array is not defined for this expansion");
745  break;
746  }
747  }
748 
749  int StdSegExp::v_GetVertexMap(int localVertexId,bool useCoeffPacking)
750  {
751  ASSERTL0((localVertexId==0)||(localVertexId==1),"local vertex id"
752  "must be between 0 or 1");
753 
754  int localDOF = localVertexId;
755 
757  (localVertexId==1) )
758  {
759  localDOF = m_base[0]->GetNumModes()-1;
760  }
761  return localDOF;
762  }
763 
764 
765  }//end namespace
766 }//end namespace
767 
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
Integrate the physical point list inarray over region and return the value.
Definition: StdSegExp.cpp:122
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdSegExp.cpp:478
NekDouble GetConstFactor(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:122
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
virtual void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdSegExp.cpp:687
virtual int v_NumDGBndryCoeffs() const
Definition: StdSegExp.cpp:618
virtual void v_StdPhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Definition: StdSegExp.cpp:173
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2)
Definition: StdSegExp.cpp:591
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:82
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:971
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Backward transform from coefficient space given in inarray and evaluate at the physical quadrature po...
Definition: StdSegExp.cpp:218
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdSegExp.cpp:399
void InterpCoeff1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
Definition: InterpCoeff.cpp:47
virtual int v_NumBndryCoeffs() const
Definition: StdSegExp.cpp:613
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:229
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Forward transform from physical quadrature space stored in inarray and evaluate the expansion coeffic...
Definition: StdSegExp.cpp:303
Principle Modified Functions .
Definition: BasisType.h:49
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:54
STL namespace.
virtual void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray)
Definition: StdSegExp.cpp:519
StdSegExp()
Default constructor.
Definition: StdSegExp.cpp:49
Fourier Expansion .
Definition: BasisType.h:52
Chebyshev Polynomials .
Definition: BasisType.h:56
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLess > m_stdStaticCondMatrixManager
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:700
virtual bool v_IsBoundaryInteriorExpansion()
Definition: StdSegExp.cpp:89
virtual int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false)
Definition: StdSegExp.cpp:749
virtual int v_CalcNumberOfCoefficients(const std::vector< unsigned int > &nummodes, int &modes_offset)
Definition: StdSegExp.cpp:623
Class representing a segment element in reference space.
Definition: StdSegExp.h:54
The base class for all shapes.
Definition: StdExpansion.h:69
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
Definition: StdExpansion.h:826
virtual int v_GetNverts() const
Definition: StdSegExp.cpp:608
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
Definition: StdSegExp.cpp:637
Principle Modified Functions .
Definition: BasisType.h:50
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey)
Definition: StdSegExp.cpp:677
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
Principle Orthogonal Functions .
Definition: BasisType.h:46
Defines a specification for a set of points.
Definition: Points.h:58
double NekDouble
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Inner product of inarray over region with respect to the expansion basis (this)->m_base[0] and return...
Definition: StdSegExp.cpp:472
virtual NekDouble v_PhysEvaluate(const Array< OneD, const NekDouble > &Lcoords, const Array< OneD, const NekDouble > &physvals)
Definition: StdSegExp.cpp:531
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:329
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Evaluate the derivative at the physical quadrature points given by inarray and return in outarray...
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:240
virtual void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdSegExp.cpp:556
virtual void v_FwdTrans_BndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdSegExp.cpp:325
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:131
virtual void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdSegExp.cpp:579
virtual void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdSegExp.cpp:248
virtual LibUtilities::ShapeType v_DetShapeType() const
Return Shape of region, using ShapeType enum list. i.e. Segment.
Definition: StdSegExp.cpp:84
Lagrange for SEM basis .
Definition: BasisType.h:53
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:723
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:359
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
virtual void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdSegExp.cpp:538
Describes the specification for a Basis.
Definition: Basis.h:50
virtual NekDouble v_PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:50
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Evaluate the derivative at the physical quadrature points given by inarray and return in outarray...
Definition: StdSegExp.cpp:155
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
Definition: StdSegExp.cpp:487
virtual void v_GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdSegExp.cpp:716