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 &Bb,
55  const LibUtilities::BasisKey &Bc)
56  : StdExpansion(numcoeffs, 3, Ba, Bb, Bc)
57 {
58 }
59 
61 {
62 }
63 
65 {
66 }
70 {
71  const int nquad0 = m_base[0]->GetNumPoints();
72  const int nquad1 = m_base[1]->GetNumPoints();
73  const int nquad2 = m_base[2]->GetNumPoints();
74 
75  Array<OneD, NekDouble> wsp(nquad0 * nquad1 * nquad2);
76 
77  // copy inarray to wsp in case inarray is used as outarray
78  Vmath::Vcopy(nquad0 * nquad1 * nquad2, &inarray[0], 1, &wsp[0], 1);
79 
80  if (out_dx.size() > 0)
81  {
82  NekDouble *D0 = &((m_base[0]->GetD())->GetPtr())[0];
83 
84  Blas::Dgemm('N', 'N', nquad0, nquad1 * nquad2, nquad0, 1.0, D0, nquad0,
85  &wsp[0], nquad0, 0.0, &out_dx[0], nquad0);
86  }
87 
88  if (out_dy.size() > 0)
89  {
90  NekDouble *D1 = &((m_base[1]->GetD())->GetPtr())[0];
91  for (int j = 0; j < nquad2; ++j)
92  {
93  Blas::Dgemm('N', 'T', nquad0, nquad1, nquad1, 1.0,
94  &wsp[j * nquad0 * nquad1], nquad0, D1, nquad1, 0.0,
95  &out_dy[j * nquad0 * nquad1], nquad0);
96  }
97  }
98 
99  if (out_dz.size() > 0)
100  {
101  NekDouble *D2 = &((m_base[2]->GetD())->GetPtr())[0];
102 
103  Blas::Dgemm('N', 'T', nquad0 * nquad1, nquad2, nquad2, 1.0, &wsp[0],
104  nquad0 * nquad1, D2, nquad2, 0.0, &out_dz[0],
105  nquad0 * nquad1);
106  }
107 }
108 
110  const Array<OneD, const NekDouble> &base0,
111  const Array<OneD, const NekDouble> &base1,
112  const Array<OneD, const NekDouble> &base2,
113  const Array<OneD, const NekDouble> &inarray,
115  bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
116 {
117  v_BwdTrans_SumFacKernel(base0, base1, base2, inarray, outarray, wsp,
118  doCheckCollDir0, doCheckCollDir1, doCheckCollDir2);
119 }
120 
122  const Array<OneD, const NekDouble> &base0,
123  const Array<OneD, const NekDouble> &base1,
124  const Array<OneD, const NekDouble> &base2,
125  const Array<OneD, const NekDouble> &inarray,
127  bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
128 {
129  v_IProductWRTBase_SumFacKernel(base0, base1, base2, inarray, outarray, wsp,
130  doCheckCollDir0, doCheckCollDir1,
131  doCheckCollDir2);
132 }
133 
135 {
136  ASSERTL1((dir == 0) || (dir == 1) || (dir == 2), "Invalid direction.");
137 
138  const int nq0 = m_base[0]->GetNumPoints();
139  const int nq1 = m_base[1]->GetNumPoints();
140  const int nq2 = m_base[2]->GetNumPoints();
141  const int nq = nq0 * nq1 * nq2;
142  const int nm0 = m_base[0]->GetNumModes();
143  const int nm1 = m_base[1]->GetNumModes();
144 
145  Array<OneD, NekDouble> alloc(4 * nq + m_ncoeffs + nm0 * nq2 * (nq1 + nm1),
146  0.0);
147  Array<OneD, NekDouble> tmp1(alloc); // Quad metric
148  Array<OneD, NekDouble> tmp2(alloc + nq); // Dir1 metric
149  Array<OneD, NekDouble> tmp3(alloc + 2 * nq); // Dir2 metric
150  Array<OneD, NekDouble> tmp4(alloc + 3 * nq); // Dir3 metric
151  Array<OneD, NekDouble> tmp5(alloc + 4 * nq); // iprod tmp
152  Array<OneD, NekDouble> wsp(tmp5 + m_ncoeffs); // Wsp
153  switch (dir)
154  {
155  case 0:
156  for (int i = 0; i < nq; i++)
157  {
158  tmp2[i] = 1.0;
160  m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
161  m_base[2]->GetBdata(), tmp2, tmp5, wsp, false, true, true);
162 
163  tmp2[i] = 0.0;
164 
165  for (int j = 0; j < m_ncoeffs; j++)
166  {
167  (*mat)(j, i) = tmp5[j];
168  }
169  }
170  break;
171  case 1:
172  for (int i = 0; i < nq; i++)
173  {
174  tmp2[i] = 1.0;
176  m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
177  m_base[2]->GetBdata(), tmp2, tmp5, wsp, true, false, true);
178 
179  tmp2[i] = 0.0;
180 
181  for (int j = 0; j < m_ncoeffs; j++)
182  {
183  (*mat)(j, i) = tmp5[j];
184  }
185  }
186  break;
187  case 2:
188  for (int i = 0; i < nq; i++)
189  {
190  tmp2[i] = 1.0;
192  m_base[0]->GetBdata(), m_base[1]->GetBdata(),
193  m_base[2]->GetDbdata(), tmp2, tmp5, wsp, true, true, false);
194  tmp2[i] = 0.0;
195 
196  for (int j = 0; j < m_ncoeffs; j++)
197  {
198  (*mat)(j, i) = tmp5[j];
199  }
200  }
201  break;
202  default:
203  NEKERROR(ErrorUtil::efatal, "Not a 2D expansion.");
204  break;
205  }
206 }
207 
209  const Array<OneD, const NekDouble> &coords,
210  const Array<OneD, const NekDouble> &physvals)
211 {
212  Array<OneD, NekDouble> eta(3);
213 
214  WARNINGL2(coords[0] >= -1 - NekConstants::kNekZeroTol, "coord[0] < -1");
215  WARNINGL2(coords[0] <= 1 + NekConstants::kNekZeroTol, "coord[0] > 1");
216  WARNINGL2(coords[1] >= -1 - NekConstants::kNekZeroTol, "coord[1] < -1");
217  WARNINGL2(coords[1] <= 1 + NekConstants::kNekZeroTol, "coord[1] > 1");
218  WARNINGL2(coords[2] >= -1 - NekConstants::kNekZeroTol, "coord[2] < -1");
219  WARNINGL2(coords[2] <= 1 + NekConstants::kNekZeroTol, "coord[2] > 1");
220 
221  // Obtain local collapsed corodinate from Cartesian coordinate.
222  LocCoordToLocCollapsed(coords, eta);
223 
224  const int nq0 = m_base[0]->GetNumPoints();
225  const int nq1 = m_base[1]->GetNumPoints();
226  const int nq2 = m_base[2]->GetNumPoints();
227 
228  Array<OneD, NekDouble> wsp1(nq1 * nq2), wsp2(nq2);
229 
230  // Construct the 2D square...
231  const NekDouble *ptr = &physvals[0];
232  for (int i = 0; i < nq1 * nq2; ++i, ptr += nq0)
233  {
234  wsp1[i] = StdExpansion::BaryEvaluate<0>(eta[0], ptr);
235  }
236 
237  for (int i = 0; i < nq2; ++i)
238  {
239  wsp2[i] = StdExpansion::BaryEvaluate<1>(eta[1], &wsp1[i * nq1]);
240  }
241 
242  return StdExpansion::BaryEvaluate<2>(eta[2], &wsp2[0]);
243 }
244 
247  const Array<OneD, const NekDouble> &physvals)
248 {
249  NekDouble value;
250 
251  int Qx = m_base[0]->GetNumPoints();
252  int Qy = m_base[1]->GetNumPoints();
253  int Qz = m_base[2]->GetNumPoints();
254 
255  Array<OneD, NekDouble> sumFactorization_qr =
256  Array<OneD, NekDouble>(Qy * Qz);
257  Array<OneD, NekDouble> sumFactorization_r = Array<OneD, NekDouble>(Qz);
258 
259  // Lagrangian interpolation matrix
260  NekDouble *interpolatingNodes = 0;
261 
262  // Interpolate first coordinate direction
263  interpolatingNodes = &I[0]->GetPtr()[0];
264 
265  Blas::Dgemv('T', Qx, Qy * Qz, 1.0, &physvals[0], Qx, &interpolatingNodes[0],
266  1, 0.0, &sumFactorization_qr[0], 1);
267 
268  // Interpolate in second coordinate direction
269  interpolatingNodes = &I[1]->GetPtr()[0];
270 
271  Blas::Dgemv('T', Qy, Qz, 1.0, &sumFactorization_qr[0], Qy,
272  &interpolatingNodes[0], 1, 0.0, &sumFactorization_r[0], 1);
273 
274  // Interpolate in third coordinate direction
275  interpolatingNodes = &I[2]->GetPtr()[0];
276  value = Blas::Ddot(Qz, interpolatingNodes, 1, &sumFactorization_r[0], 1);
277 
278  return value;
279 }
280 
281 /**
282  * @param inarray Input coefficients.
283  * @param output Output coefficients.
284  * @param mkey Matrix key
285  */
287  const Array<OneD, const NekDouble> &inarray,
288  Array<OneD, NekDouble> &outarray, const StdRegions::StdMatrixKey &mkey)
289 {
290  if (mkey.GetNVarCoeff() == 0 &&
293  {
294  // This implementation is only valid when there are no
295  // coefficients associated to the Laplacian operator
296  int nqtot = GetTotPoints();
297 
298  const Array<OneD, const NekDouble> &base0 = m_base[0]->GetBdata();
299  const Array<OneD, const NekDouble> &base1 = m_base[1]->GetBdata();
300  const Array<OneD, const NekDouble> &base2 = m_base[2]->GetBdata();
301 
302  // Allocate temporary storage
303  Array<OneD, NekDouble> wsp0(7 * nqtot);
304  Array<OneD, NekDouble> wsp1(wsp0 + nqtot);
305 
306  if (!(m_base[0]->Collocation() && m_base[1]->Collocation() &&
307  m_base[2]->Collocation()))
308  {
309  // LAPLACIAN MATRIX OPERATION
310  // wsp0 = u = B * u_hat
311  // wsp1 = du_dxi1 = D_xi1 * wsp0 = D_xi1 * u
312  // wsp2 = du_dxi2 = D_xi2 * wsp0 = D_xi2 * u
313  BwdTrans_SumFacKernel(base0, base1, base2, inarray, wsp0, wsp1,
314  true, true, true);
315  LaplacianMatrixOp_MatFree_Kernel(wsp0, outarray, wsp1);
316  }
317  else
318  {
319  LaplacianMatrixOp_MatFree_Kernel(inarray, outarray, wsp1);
320  }
321  }
322  else
323  {
325  mkey);
326  }
327 }
328 
330  const Array<OneD, const NekDouble> &inarray,
331  Array<OneD, NekDouble> &outarray, const StdRegions::StdMatrixKey &mkey)
332 {
333  if (mkey.GetNVarCoeff() == 0 &&
335  {
336  using std::max;
337 
338  int nquad0 = m_base[0]->GetNumPoints();
339  int nquad1 = m_base[1]->GetNumPoints();
340  int nquad2 = m_base[2]->GetNumPoints();
341  int nmodes0 = m_base[0]->GetNumModes();
342  int nmodes1 = m_base[1]->GetNumModes();
343  int nmodes2 = m_base[2]->GetNumModes();
344  int wspsize = max(nquad0 * nmodes2 * (nmodes1 + nquad1),
345  nquad0 * nquad1 * (nquad2 + nmodes0) +
346  nmodes0 * nmodes1 * nquad2);
347 
349 
350  const Array<OneD, const NekDouble> &base0 = m_base[0]->GetBdata();
351  const Array<OneD, const NekDouble> &base1 = m_base[1]->GetBdata();
352  const Array<OneD, const NekDouble> &base2 = m_base[2]->GetBdata();
353  Array<OneD, NekDouble> wsp0(8 * wspsize);
354  Array<OneD, NekDouble> wsp1(wsp0 + 1 * wspsize);
355  Array<OneD, NekDouble> wsp2(wsp0 + 2 * wspsize);
356 
357  if (!(m_base[0]->Collocation() && m_base[1]->Collocation() &&
358  m_base[2]->Collocation()))
359  {
360  // MASS MATRIX OPERATION
361  // The following is being calculated:
362  // wsp0 = B * u_hat = u
363  // wsp1 = W * wsp0
364  // outarray = B^T * wsp1 = B^T * W * B * u_hat = M * u_hat
365  BwdTrans_SumFacKernel(base0, base1, base2, inarray, wsp0, wsp2,
366  true, true, true);
367  MultiplyByQuadratureMetric(wsp0, wsp1);
368  IProductWRTBase_SumFacKernel(base0, base1, base2, wsp1, outarray,
369  wsp2, true, true, true);
370  LaplacianMatrixOp_MatFree_Kernel(wsp0, wsp1, wsp2);
371  }
372  else
373  {
374  // specialised implementation for the classical spectral
375  // element method
376  MultiplyByQuadratureMetric(inarray, outarray);
377  LaplacianMatrixOp_MatFree_Kernel(inarray, wsp1, wsp2);
378  }
379 
380  // outarray = lambda * outarray + wsp1
381  // = (lambda * M + L ) * u_hat
382  Vmath::Svtvp(m_ncoeffs, lambda, &outarray[0], 1, &wsp1[0], 1,
383  &outarray[0], 1);
384  }
385  else
386  {
388  mkey);
389  }
390 }
391 
393  const Array<OneD, const NekDouble> &inarray)
394 {
395  const int nqtot = GetTotPoints();
397  v_MultiplyByStdQuadratureMetric(inarray, tmp);
398  return Vmath::Vsum(nqtot, tmp, 1);
399 }
400 
402 {
403  NEKERROR(ErrorUtil::efatal, "This function is not valid or not defined");
404  return 0;
405 }
406 
407 int StdExpansion3D::v_GetEdgeNcoeffs(const int i) const
408 {
409  boost::ignore_unused(i);
410  NEKERROR(ErrorUtil::efatal, "This function is not valid or not defined");
411  return 0;
412 }
413 
415  const int tid, Array<OneD, unsigned int> &maparray,
416  Array<OneD, int> &signarray, Orientation traceOrient)
417 {
418  boost::ignore_unused(tid, maparray, signarray, traceOrient);
419  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
420 }
421 
423  Array<OneD, unsigned int> &maparray,
424  Array<OneD, int> &signarray,
425  Orientation traceOrient, int P,
426  int Q)
427 {
428  Array<OneD, unsigned int> map1, map2;
429  GetTraceCoeffMap(tid, map1);
430  GetElmtTraceToTraceMap(tid, map2, signarray, traceOrient, P, Q);
431 
432  if (maparray.size() != map2.size())
433  {
434  maparray = Array<OneD, unsigned int>(map2.size());
435  }
436 
437  for (int i = 0; i < map2.size(); ++i)
438  {
439  maparray[i] = map1[map2[i]];
440  }
441 }
442 
444  const int facedir, const LibUtilities::BasisType faceDirBasisType,
445  const int numpoints, const int nummodes)
446 {
447  boost::ignore_unused(facedir);
448 
449  switch (faceDirBasisType)
450  {
452  {
453  const LibUtilities::PointsKey pkey(
456  pkey);
457  }
460  {
461  const LibUtilities::PointsKey pkey(
462  numpoints + 1, LibUtilities::eGaussLobattoLegendre);
464  pkey);
465  }
467  {
468  const LibUtilities::PointsKey pkey(
471  pkey);
472  }
474  {
475  const LibUtilities::PointsKey pkey(
478  pkey);
479  }
482  {
483  const LibUtilities::PointsKey pkey(
484  numpoints + 1, LibUtilities::eGaussLobattoLegendre);
486  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, const LibUtilities::BasisType faceDirBasisType,
501  const int numpoints, const int nummodes)
502 {
503  switch (faceDirBasisType)
504  {
506  {
507  const LibUtilities::PointsKey pkey(
510  pkey);
511  }
515  {
516  switch (facedir)
517  {
518  case 0:
519  {
520  const LibUtilities::PointsKey pkey(
521  numpoints + 1, LibUtilities::eGaussLobattoLegendre);
523  nummodes, pkey);
524  }
525  case 1:
526  {
527  // const LibUtilities::PointsKey pkey(
528  // numpoints+1,
529  // LibUtilities::eGaussLobattoLegendre);
530  const LibUtilities::PointsKey pkey(
531  numpoints, LibUtilities::eGaussRadauMAlpha1Beta0);
533  nummodes, pkey);
534  }
535  default:
536  {
537 
538  NEKERROR(ErrorUtil::efatal, "invalid value to flag");
539  break;
540  }
541  }
542  break;
543  }
544 
546  {
547  switch (facedir)
548  {
549  case 0:
550  {
551  const LibUtilities::PointsKey pkey(
554  nummodes, pkey);
555  }
556  case 1:
557  {
558  const LibUtilities::PointsKey pkey(
559  numpoints, LibUtilities::eGaussRadauMAlpha1Beta0);
561  nummodes, pkey);
562  }
563  default:
564  {
565  NEKERROR(ErrorUtil::efatal, "invalid value to flag");
566  break;
567  }
568  }
569  break;
570  }
571 
576  {
577  switch (facedir)
578  {
579  case 0:
580  {
581  const LibUtilities::PointsKey pkey(
584  nummodes, pkey);
585  }
586  case 1:
587  {
588  const LibUtilities::PointsKey pkey(
589  numpoints, LibUtilities::eGaussRadauMAlpha1Beta0);
591  nummodes, pkey);
592  }
593  default:
594  {
595  NEKERROR(ErrorUtil::efatal, "invalid value to flag");
596  break;
597  }
598  }
599  break;
600  }
601  default:
602  {
603  NEKERROR(ErrorUtil::efatal, "expansion type unknown");
604  break;
605  }
606  }
607 
608  // Keep things happy by returning a value.
610 }
611 } // namespace StdRegions
612 } // namespace Nektar
#define WARNINGL2(condition, msg)
Definition: ErrorUtil.hpp:273
#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:249
Describes the specification for a Basis.
Definition: Basis.h:50
Defines a specification for a set of points.
Definition: Points.h:59
virtual void v_GetTraceToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient, int P, int Q)
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:71
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
void GetElmtTraceToTraceMap(const unsigned int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
Definition: StdExpansion.h:707
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:974
void GetTraceCoeffMap(const unsigned int traceid, Array< OneD, unsigned int > &maparray)
Definition: StdExpansion.h:701
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:731
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:126
bool ConstFactorExists(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:135
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:246
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:182
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:368
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:53
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:51
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:44
@ eModified_C
Principle Modified Functions .
Definition: BasisType.h:52
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:58
@ eOrtho_C
Principle Orthogonal Functions .
Definition: BasisType.h:48
@ eModifiedPyr_C
Principle Modified Functions.
Definition: BasisType.h:55
@ eOrtho_B
Principle Orthogonal Functions .
Definition: BasisType.h:46
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50
@ eOrthoPyr_C
Principle Orthogonal Functions .
Definition: BasisType.h:53
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:75
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:622
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:895
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255