Nektar++
StdExpansion3D.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File StdExpansion3D.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 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Daughter of StdExpansion. This class contains routine
32 // which are common to 3D expansion. Typically this inolves physiocal
33 // space operations.
34 //
35 ///////////////////////////////////////////////////////////////////////////////
36 
37 #include <boost/core/ignore_unused.hpp>
38 
40 
41 #ifdef max
42 #undef max
43 #endif
44 
45 namespace Nektar
46 {
47  namespace StdRegions
48  {
50  {
51  }
52 
54  const LibUtilities::BasisKey &Ba,
55  const LibUtilities::BasisKey &Bb,
56  const LibUtilities::BasisKey &Bc) :
57  StdExpansion(numcoeffs,3,Ba,Bb,Bc)
58  {
59  }
60 
62  StdExpansion(T)
63  {
64  }
65 
67  {
68  }
70  const Array<OneD, const NekDouble> &inarray,
71  Array<OneD, NekDouble> &out_dx,
72  Array<OneD, NekDouble> &out_dy,
73  Array<OneD, NekDouble> &out_dz)
74  {
75  const int nquad0 = m_base[0]->GetNumPoints();
76  const int nquad1 = m_base[1]->GetNumPoints();
77  const int nquad2 = m_base[2]->GetNumPoints();
78 
79  Array<OneD, NekDouble> wsp(nquad0*nquad1*nquad2);
80 
81  // copy inarray to wsp in case inarray is used as outarray
82  Vmath::Vcopy(nquad0*nquad1*nquad2, &inarray[0], 1, &wsp[0], 1);
83 
84  if (out_dx.size() > 0)
85  {
86  NekDouble *D0 = &((m_base[0]->GetD())->GetPtr())[0];
87 
88  Blas::Dgemm('N','N', nquad0,nquad1*nquad2,nquad0,1.0,
89  D0,nquad0,&wsp[0],nquad0,0.0,&out_dx[0],nquad0);
90  }
91 
92  if (out_dy.size() > 0)
93  {
94  NekDouble *D1 = &((m_base[1]->GetD())->GetPtr())[0];
95  for (int j = 0; j < nquad2; ++j)
96  {
97  Blas::Dgemm('N', 'T', nquad0, nquad1, nquad1,
98  1.0, &wsp[j*nquad0*nquad1], nquad0,
99  D1, nquad1,
100  0.0, &out_dy[j*nquad0*nquad1], nquad0);
101  }
102  }
103 
104  if (out_dz.size() > 0)
105  {
106  NekDouble *D2 = &((m_base[2]->GetD())->GetPtr())[0];
107 
108  Blas::Dgemm('N','T',nquad0*nquad1,nquad2,nquad2,1.0,
109  &wsp[0],nquad0*nquad1,D2,nquad2,0.0,&out_dz[0],
110  nquad0*nquad1);
111  }
112  }
113 
115  const Array<OneD, const NekDouble>& base0,
116  const Array<OneD, const NekDouble>& base1,
117  const Array<OneD, const NekDouble>& base2,
118  const Array<OneD, const NekDouble>& inarray,
119  Array<OneD, NekDouble>& outarray,
121  bool doCheckCollDir0,
122  bool doCheckCollDir1,
123  bool doCheckCollDir2)
124  {
125  v_BwdTrans_SumFacKernel(base0, base1, base2, inarray, outarray, wsp, doCheckCollDir0, doCheckCollDir1, doCheckCollDir2);
126  }
127 
129  const Array<OneD, const NekDouble>& base0,
130  const Array<OneD, const NekDouble>& base1,
131  const Array<OneD, const NekDouble>& base2,
132  const Array<OneD, const NekDouble>& inarray,
133  Array<OneD, NekDouble> &outarray,
135  bool doCheckCollDir0,
136  bool doCheckCollDir1,
137  bool doCheckCollDir2)
138  {
139  v_IProductWRTBase_SumFacKernel(base0, base1, base2, inarray, outarray, wsp, doCheckCollDir0, doCheckCollDir1, doCheckCollDir2);
140  }
141 
143  const int dir,
144  DNekMatSharedPtr &mat)
145  {
146  ASSERTL1((dir==0)||(dir==1)||(dir==2),"Invalid direction.");
147 
148  const int nq0 = m_base[0]->GetNumPoints();
149  const int nq1 = m_base[1]->GetNumPoints();
150  const int nq2 = m_base[2]->GetNumPoints();
151  const int nq = nq0*nq1*nq2;
152  const int nm0 = m_base[0]->GetNumModes();
153  const int nm1 = m_base[1]->GetNumModes();
154 
155  Array<OneD, NekDouble> alloc(4*nq + m_ncoeffs + nm0*nq2*(nq1+nm1),0.0);
156  Array<OneD, NekDouble> tmp1 (alloc); // Quad metric
157  Array<OneD, NekDouble> tmp2 (alloc + nq); // Dir1 metric
158  Array<OneD, NekDouble> tmp3 (alloc + 2*nq); // Dir2 metric
159  Array<OneD, NekDouble> tmp4 (alloc + 3*nq); // Dir3 metric
160  Array<OneD, NekDouble> tmp5 (alloc + 4*nq); // iprod tmp
161  Array<OneD, NekDouble> wsp (tmp5 + m_ncoeffs); // Wsp
162  switch(dir)
163  {
164  case 0:
165  for(int i=0; i<nq;i++)
166  {
167  tmp2[i] = 1.0;
168  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
169  m_base[1]->GetBdata(),
170  m_base[2]->GetBdata(),
171  tmp2,tmp5,wsp,
172  false,true,true);
173 
174  tmp2[i] = 0.0;
175 
176  for(int j=0; j<m_ncoeffs;j++)
177  {
178  (*mat)(j,i) = tmp5[j];
179  }
180  }
181  break;
182  case 1:
183  for(int i=0; i<nq;i++)
184  {
185  tmp2[i] = 1.0;
186  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
187  m_base[1]->GetDbdata(),
188  m_base[2]->GetBdata(),
189  tmp2,tmp5,wsp,
190  true,false,true);
191 
192  tmp2[i] = 0.0;
193 
194  for(int j=0; j<m_ncoeffs;j++)
195  {
196  (*mat)(j,i) = tmp5[j];
197  }
198  }
199  break;
200  case 2:
201  for(int i=0; i<nq;i++)
202  {
203  tmp2[i] = 1.0;
204  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
205  m_base[1]->GetBdata(),
206  m_base[2]->GetDbdata(),
207  tmp2,tmp5,wsp,
208  true,true,false);
209  tmp2[i] = 0.0;
210 
211  for(int j=0; j<m_ncoeffs;j++)
212  {
213  (*mat)(j,i) = tmp5[j];
214  }
215  }
216  break;
217  default:
218  NEKERROR(ErrorUtil::efatal, "Not a 2D expansion.");
219  break;
220  }
221  }
222 
224  const Array<OneD, const NekDouble> &coords,
225  const Array<OneD, const NekDouble> &physvals)
226  {
227  Array<OneD, NekDouble> eta(3);
228 
229  WARNINGL2(coords[0] >= -1 - NekConstants::kNekZeroTol,
230  "coord[0] < -1");
231  WARNINGL2(coords[0] <= 1 + NekConstants::kNekZeroTol,
232  "coord[0] > 1");
233  WARNINGL2(coords[1] >= -1 - NekConstants::kNekZeroTol,
234  "coord[1] < -1");
235  WARNINGL2(coords[1] <= 1 + NekConstants::kNekZeroTol,
236  "coord[1] > 1");
237  WARNINGL2(coords[2] >= -1 - NekConstants::kNekZeroTol,
238  "coord[2] < -1");
239  WARNINGL2(coords[2] <= 1 + NekConstants::kNekZeroTol,
240  "coord[2] > 1");
241 
242  // Obtain local collapsed corodinate from Cartesian coordinate.
243  LocCoordToLocCollapsed(coords, eta);
244 
245  const int nq0 = m_base[0]->GetNumPoints();
246  const int nq1 = m_base[1]->GetNumPoints();
247  const int nq2 = m_base[2]->GetNumPoints();
248 
249  Array<OneD, NekDouble> wsp1(nq1 * nq2), wsp2(nq2);
250 
251  // Construct the 2D square...
252  const NekDouble *ptr = &physvals[0];
253  for (int i = 0; i < nq1 * nq2; ++i, ptr += nq0)
254  {
255  wsp1[i] = StdExpansion::BaryEvaluate<0>(eta[0], ptr);
256  }
257 
258  for (int i = 0; i < nq2; ++i)
259  {
260  wsp2[i] = StdExpansion::BaryEvaluate<1>(eta[1], &wsp1[i * nq1]);
261  }
262 
263  return StdExpansion::BaryEvaluate<2>(eta[2], &wsp2[0]);
264  }
265 
268  const Array<OneD, const NekDouble> &physvals)
269  {
270  NekDouble value;
271 
272  int Qx = m_base[0]->GetNumPoints();
273  int Qy = m_base[1]->GetNumPoints();
274  int Qz = m_base[2]->GetNumPoints();
275 
276  Array<OneD, NekDouble> sumFactorization_qr = Array<OneD, NekDouble>(Qy*Qz);
277  Array<OneD, NekDouble> sumFactorization_r = Array<OneD, NekDouble>(Qz);
278 
279  // Lagrangian interpolation matrix
280  NekDouble *interpolatingNodes = 0;
281 
282  // Interpolate first coordinate direction
283  interpolatingNodes = &I[0]->GetPtr()[0];
284 
285  Blas::Dgemv('T',Qx,Qy*Qz,1.0,&physvals[0],Qx,&interpolatingNodes[0], 1, 0.0, &sumFactorization_qr[0], 1);
286 
287  // Interpolate in second coordinate direction
288  interpolatingNodes = &I[1]->GetPtr()[0];
289 
290  Blas::Dgemv('T',Qy,Qz,1.0,&sumFactorization_qr[0],Qy,&interpolatingNodes[0],1,0.0,&sumFactorization_r[0], 1);
291 
292  // Interpolate in third coordinate direction
293  interpolatingNodes = &I[2]->GetPtr()[0];
294  value = Blas::Ddot(Qz, interpolatingNodes, 1, &sumFactorization_r[0], 1);
295 
296  return value;
297  }
298 
299 
300  /**
301  * @param inarray Input coefficients.
302  * @param output Output coefficients.
303  * @param mkey Matrix key
304  */
306  const Array<OneD, const NekDouble> &inarray,
307  Array<OneD,NekDouble> &outarray,
308  const StdRegions::StdMatrixKey &mkey)
309  {
310  if ( mkey.GetNVarCoeff() == 0 && !mkey.ConstFactorExists(StdRegions::eFactorCoeffD00) &&
312  {
313  // This implementation is only valid when there are no
314  // coefficients associated to the Laplacian operator
315  int nqtot = GetTotPoints();
316 
317  const Array<OneD, const NekDouble>& base0 = m_base[0]->GetBdata();
318  const Array<OneD, const NekDouble>& base1 = m_base[1]->GetBdata();
319  const Array<OneD, const NekDouble>& base2 = m_base[2]->GetBdata();
320 
321  // Allocate temporary storage
322  Array<OneD,NekDouble> wsp0(7*nqtot);
323  Array<OneD,NekDouble> wsp1(wsp0+nqtot);
324 
325  if(!(m_base[0]->Collocation() && m_base[1]->Collocation() &&
326  m_base[2]->Collocation()))
327  {
328  // LAPLACIAN MATRIX OPERATION
329  // wsp0 = u = B * u_hat
330  // wsp1 = du_dxi1 = D_xi1 * wsp0 = D_xi1 * u
331  // wsp2 = du_dxi2 = D_xi2 * wsp0 = D_xi2 * u
332  BwdTrans_SumFacKernel(base0,base1,base2,inarray,wsp0,wsp1,true,true,true);
333  LaplacianMatrixOp_MatFree_Kernel(wsp0,outarray,wsp1);
334  }
335  else
336  {
337  LaplacianMatrixOp_MatFree_Kernel(inarray,outarray,wsp1);
338  }
339  }
340  else
341  {
343  }
344  }
345 
346 
348  const Array<OneD, const NekDouble> &inarray,
349  Array<OneD,NekDouble> &outarray,
350  const StdRegions::StdMatrixKey &mkey)
351  {
353  {
354  using std::max;
355 
356  int nquad0 = m_base[0]->GetNumPoints();
357  int nquad1 = m_base[1]->GetNumPoints();
358  int nquad2 = m_base[2]->GetNumPoints();
359  int nmodes0 = m_base[0]->GetNumModes();
360  int nmodes1 = m_base[1]->GetNumModes();
361  int nmodes2 = m_base[2]->GetNumModes();
362  int wspsize = max(nquad0*nmodes2*(nmodes1+nquad1),
363  nquad0*nquad1*(nquad2+nmodes0)+
364  nmodes0*nmodes1*nquad2);
365 
367 
368  const Array<OneD, const NekDouble>& base0 = m_base[0]->GetBdata ();
369  const Array<OneD, const NekDouble>& base1 = m_base[1]->GetBdata ();
370  const Array<OneD, const NekDouble>& base2 = m_base[2]->GetBdata ();
371  Array<OneD,NekDouble> wsp0(8*wspsize);
372  Array<OneD,NekDouble> wsp1(wsp0+1*wspsize);
373  Array<OneD,NekDouble> wsp2(wsp0+2*wspsize);
374 
375  if(!(m_base[0]->Collocation() && m_base[1]->Collocation() &&
376  m_base[2]->Collocation()))
377  {
378  // MASS MATRIX OPERATION
379  // The following is being calculated:
380  // wsp0 = B * u_hat = u
381  // wsp1 = W * wsp0
382  // outarray = B^T * wsp1 = B^T * W * B * u_hat = M * u_hat
383  BwdTrans_SumFacKernel (base0,base1,base2,inarray,
384  wsp0,wsp2,true,true,true);
385  MultiplyByQuadratureMetric (wsp0,wsp1);
386  IProductWRTBase_SumFacKernel (base0,base1,base2,wsp1,
387  outarray,wsp2,true,true,true);
388  LaplacianMatrixOp_MatFree_Kernel(wsp0,wsp1,wsp2);
389  }
390  else
391  {
392  // specialised implementation for the classical spectral
393  // element method
394  MultiplyByQuadratureMetric (inarray,outarray);
395  LaplacianMatrixOp_MatFree_Kernel(inarray,wsp1,wsp2);
396  }
397 
398  // outarray = lambda * outarray + wsp1
399  // = (lambda * M + L ) * u_hat
400  Vmath::Svtvp(m_ncoeffs,lambda,&outarray[0],1,&wsp1[0],1,
401  &outarray[0],1);
402  }
403  else
404  {
406  }
407  }
408 
410  const Array<OneD, const NekDouble>& inarray)
411  {
412  const int nqtot = GetTotPoints();
414  v_MultiplyByStdQuadratureMetric(inarray, tmp);
415  return Vmath::Vsum(nqtot, tmp, 1);
416  }
417 
419  {
420  NEKERROR(ErrorUtil::efatal, "This function is not valid or not defined");
421  return 0;
422  }
423 
424  int StdExpansion3D::v_GetEdgeNcoeffs(const int i) const
425  {
426  boost::ignore_unused(i);
427  NEKERROR(ErrorUtil::efatal, "This function is not valid or not defined");
428  return 0;
429  }
430 
432  const int tid,
433  Array<OneD, unsigned int> &maparray,
434  Array<OneD, int> &signarray,
435  Orientation traceOrient)
436  {
437  boost::ignore_unused(tid,maparray,signarray,traceOrient);
438  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
439  }
440 
442  const int facedir,
443  const LibUtilities::BasisType faceDirBasisType,
444  const int numpoints,
445  const int nummodes)
446  {
447  boost::ignore_unused(facedir);
448 
449  switch(faceDirBasisType)
450  {
452  {
453  const LibUtilities::PointsKey pkey(
455  return LibUtilities::BasisKey(
456  LibUtilities::eModified_A, nummodes, pkey);
457  }
460  {
461  const LibUtilities::PointsKey pkey(
463  return LibUtilities::BasisKey(
464  LibUtilities::eModified_A, nummodes, pkey);
465  }
467  {
468  const LibUtilities::PointsKey pkey(
470  return LibUtilities::BasisKey(
471  LibUtilities::eGLL_Lagrange, nummodes, pkey);
472  }
474  {
475  const LibUtilities::PointsKey pkey(
477  return LibUtilities::BasisKey(
478  LibUtilities::eOrtho_A, nummodes, pkey);
479  }
482  {
483  const LibUtilities::PointsKey pkey(
485  return LibUtilities::BasisKey(
486  LibUtilities::eOrtho_A, nummodes, pkey);
487  }
488  default:
489  {
490  NEKERROR(ErrorUtil::efatal, "expansion type unknown");
491  break;
492  }
493  }
494 
495  // Keep things happy by returning a value.
497  }
498 
500  const int facedir,
501  const LibUtilities::BasisType faceDirBasisType,
502  const int numpoints,
503  const int nummodes)
504  {
505  switch(faceDirBasisType)
506  {
508  {
509  const LibUtilities::PointsKey pkey(
511  return LibUtilities::BasisKey(
512  LibUtilities::eModified_A, nummodes, pkey);
513  }
517  {
518  switch (facedir)
519  {
520  case 0:
521  {
522  const LibUtilities::PointsKey pkey(
523  numpoints+1,
525  return LibUtilities::BasisKey(
526  LibUtilities::eModified_A, nummodes, pkey);
527  }
528  case 1:
529  {
530 // const LibUtilities::PointsKey pkey(
531 // numpoints+1,
532 // LibUtilities::eGaussLobattoLegendre);
533  const LibUtilities::PointsKey pkey(
534  numpoints,
536  return LibUtilities::BasisKey(
537  LibUtilities::eModified_B, nummodes, pkey);
538  }
539  default:
540  {
541 
542  NEKERROR(ErrorUtil::efatal,"invalid value to flag");
543  break;
544  }
545  }
546  break;
547  }
548 
550  {
551  switch (facedir)
552  {
553  case 0:
554  {
555  const LibUtilities::PointsKey pkey(
556  numpoints,
558  return LibUtilities::BasisKey(
559  LibUtilities::eOrtho_A, nummodes, pkey);
560  }
561  case 1:
562  {
563  const LibUtilities::PointsKey pkey(
564  numpoints,
566  return LibUtilities::BasisKey(
567  LibUtilities::eOrtho_B, nummodes, pkey);
568  }
569  default:
570  {
571  NEKERROR(ErrorUtil::efatal,"invalid value to flag");
572  break;
573  }
574  }
575  break;
576  }
577 
582  {
583  switch (facedir)
584  {
585  case 0:
586  {
587  const LibUtilities::PointsKey pkey(
588  numpoints,
590  return LibUtilities::BasisKey(
591  LibUtilities::eOrtho_A, nummodes, pkey);
592  }
593  case 1:
594  {
595  const LibUtilities::PointsKey pkey(
596  numpoints,
598  return LibUtilities::BasisKey(
599  LibUtilities::eOrtho_B, nummodes, pkey);
600  }
601  default:
602  {
603  NEKERROR(ErrorUtil::efatal,"invalid value to flag");
604  break;
605  }
606  }
607  break;
608  }
609  default:
610  {
611  NEKERROR(ErrorUtil::efatal,"expansion type unknown");
612  break;
613  }
614  }
615 
616  // Keep things happy by returning a value.
618  }
619  }//end namespace
620 }//end namespace
#define WARNINGL2(condition, msg)
Definition: ErrorUtil.hpp:275
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
Describes the specification for a Basis.
Definition: Basis.h:50
Defines a specification for a set of points.
Definition: Points.h:60
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
Integrates the specified function over the domain.
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
virtual int v_GetNedges(void) const
void BwdTrans_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
virtual void v_GetEdgeInteriorToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards)
void IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
virtual void v_BwdTrans_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)=0
virtual int v_GetEdgeNcoeffs(const int i) const
virtual NekDouble v_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_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
virtual void v_IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)=0
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d1, Array< OneD, NekDouble > &outarray_d2, Array< OneD, NekDouble > &outarray_d3)
Calculate the 3D derivative in the local tensor/collapsed coordinate at the physical points.
virtual void v_GenStdMatBwdDeriv(const int dir, DNekMatSharedPtr &mat)
The base class for all shapes.
Definition: StdExpansion.h:63
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:134
void LaplacianMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
Definition: StdExpansion.h:982
virtual void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:733
void HelmholtzMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void LaplacianMatrixOp_MatFree_Kernel(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
NekDouble GetConstFactor(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:121
bool ConstFactorExists(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:130
static void Dgemv(const char &trans, const int &m, const int &n, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = A x where A[m x n].
Definition: Blas.hpp:265
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
Definition: Blas.hpp:197
static void Dgemm(const char &transa, const char &transb, const int &m, const int &n, const int &k, const double &alpha, const double *a, const int &lda, const double *b, const int &ldb, const double &beta, double *c, const int &ldc)
BLAS level 3: Matrix-matrix multiply C = A x B where op(A)[m x k], op(B)[k x n], C[m x n] DGEMM perfo...
Definition: Blas.hpp:394
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ eGaussRadauMAlpha1Beta0
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:58
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:49
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:45
@ eModified_C
Principle Modified Functions .
Definition: BasisType.h:50
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:54
@ eOrtho_C
Principle Orthogonal Functions .
Definition: BasisType.h:47
@ eModifiedPyr_C
Principle Modified Functions .
Definition: BasisType.h:52
@ eOrtho_B
Principle Orthogonal Functions .
Definition: BasisType.h:46
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48
@ eOrthoPyr_C
Principle Orthogonal Functions .
Definition: BasisType.h:51
static const NekDouble kNekZeroTol
LibUtilities::BasisKey EvaluateTriFaceBasisKey(const int facedir, const LibUtilities::BasisType faceDirBasisType, const int numpoints, const int nummodes)
LibUtilities::BasisKey EvaluateQuadFaceBasisKey(const int facedir, const LibUtilities::BasisType faceDirBasisType, const int numpoints, const int nummodes)
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
double NekDouble
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:565
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:846
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199