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
38
39#ifdef max
40#undef max
41#endif
42
43namespace Nektar::StdRegions
44{
45
47 const LibUtilities::BasisKey &Bb,
48 const LibUtilities::BasisKey &Bc)
49 : StdExpansion(numcoeffs, 3, Ba, Bb, Bc)
50{
51}
52
56{
57 const int nquad0 = m_base[0]->GetNumPoints();
58 const int nquad1 = m_base[1]->GetNumPoints();
59 const int nquad2 = m_base[2]->GetNumPoints();
60
61 Array<OneD, NekDouble> wsp(nquad0 * nquad1 * nquad2);
62
63 // copy inarray to wsp in case inarray is used as outarray
64 Vmath::Vcopy(nquad0 * nquad1 * nquad2, &inarray[0], 1, &wsp[0], 1);
65
66 if (out_dx.size() > 0)
67 {
68 NekDouble *D0 = &((m_base[0]->GetD())->GetPtr())[0];
69
70 Blas::Dgemm('N', 'N', nquad0, nquad1 * nquad2, nquad0, 1.0, D0, nquad0,
71 &wsp[0], nquad0, 0.0, &out_dx[0], nquad0);
72 }
73
74 if (out_dy.size() > 0)
75 {
76 NekDouble *D1 = &((m_base[1]->GetD())->GetPtr())[0];
77 for (int j = 0; j < nquad2; ++j)
78 {
79 Blas::Dgemm('N', 'T', nquad0, nquad1, nquad1, 1.0,
80 &wsp[j * nquad0 * nquad1], nquad0, D1, nquad1, 0.0,
81 &out_dy[j * nquad0 * nquad1], nquad0);
82 }
83 }
84
85 if (out_dz.size() > 0)
86 {
87 NekDouble *D2 = &((m_base[2]->GetD())->GetPtr())[0];
88
89 Blas::Dgemm('N', 'T', nquad0 * nquad1, nquad2, nquad2, 1.0, &wsp[0],
90 nquad0 * nquad1, D2, nquad2, 0.0, &out_dz[0],
91 nquad0 * nquad1);
92 }
93}
94
99 const Array<OneD, const NekDouble> &inarray,
101 bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
102{
103 v_BwdTrans_SumFacKernel(base0, base1, base2, inarray, outarray, wsp,
104 doCheckCollDir0, doCheckCollDir1, doCheckCollDir2);
105}
106
108 const Array<OneD, const NekDouble> &base0,
109 const Array<OneD, const NekDouble> &base1,
110 const Array<OneD, const NekDouble> &base2,
111 const Array<OneD, const NekDouble> &inarray,
113 bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
114{
115 v_IProductWRTBase_SumFacKernel(base0, base1, base2, inarray, outarray, wsp,
116 doCheckCollDir0, doCheckCollDir1,
117 doCheckCollDir2);
118}
119
121{
122 ASSERTL1((dir == 0) || (dir == 1) || (dir == 2), "Invalid direction.");
123
124 const int nq0 = m_base[0]->GetNumPoints();
125 const int nq1 = m_base[1]->GetNumPoints();
126 const int nq2 = m_base[2]->GetNumPoints();
127 const int nq = nq0 * nq1 * nq2;
128 const int nm0 = m_base[0]->GetNumModes();
129 const int nm1 = m_base[1]->GetNumModes();
130
131 Array<OneD, NekDouble> alloc(4 * nq + m_ncoeffs + nm0 * nq2 * (nq1 + nm1),
132 0.0);
133 Array<OneD, NekDouble> tmp1(alloc); // Quad metric
134 Array<OneD, NekDouble> tmp2(alloc + nq); // Dir1 metric
135 Array<OneD, NekDouble> tmp3(alloc + 2 * nq); // Dir2 metric
136 Array<OneD, NekDouble> tmp4(alloc + 3 * nq); // Dir3 metric
137 Array<OneD, NekDouble> tmp5(alloc + 4 * nq); // iprod tmp
138 Array<OneD, NekDouble> wsp(tmp5 + m_ncoeffs); // Wsp
139 switch (dir)
140 {
141 case 0:
142 for (int i = 0; i < nq; i++)
143 {
144 tmp2[i] = 1.0;
146 m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
147 m_base[2]->GetBdata(), tmp2, tmp5, wsp, false, true, true);
148
149 tmp2[i] = 0.0;
150
151 for (int j = 0; j < m_ncoeffs; j++)
152 {
153 (*mat)(j, i) = tmp5[j];
154 }
155 }
156 break;
157 case 1:
158 for (int i = 0; i < nq; i++)
159 {
160 tmp2[i] = 1.0;
162 m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
163 m_base[2]->GetBdata(), tmp2, tmp5, wsp, true, false, true);
164
165 tmp2[i] = 0.0;
166
167 for (int j = 0; j < m_ncoeffs; j++)
168 {
169 (*mat)(j, i) = tmp5[j];
170 }
171 }
172 break;
173 case 2:
174 for (int i = 0; i < nq; i++)
175 {
176 tmp2[i] = 1.0;
178 m_base[0]->GetBdata(), m_base[1]->GetBdata(),
179 m_base[2]->GetDbdata(), tmp2, tmp5, wsp, true, true, false);
180 tmp2[i] = 0.0;
181
182 for (int j = 0; j < m_ncoeffs; j++)
183 {
184 (*mat)(j, i) = tmp5[j];
185 }
186 }
187 break;
188 default:
189 NEKERROR(ErrorUtil::efatal, "Not a 2D expansion.");
190 break;
191 }
192}
193
195 const Array<OneD, const NekDouble> &coords,
196 const Array<OneD, const NekDouble> &physvals)
197{
199
200 WARNINGL2(coords[0] >= -1 - NekConstants::kNekZeroTol, "coord[0] < -1");
201 WARNINGL2(coords[0] <= 1 + NekConstants::kNekZeroTol, "coord[0] > 1");
202 WARNINGL2(coords[1] >= -1 - NekConstants::kNekZeroTol, "coord[1] < -1");
203 WARNINGL2(coords[1] <= 1 + NekConstants::kNekZeroTol, "coord[1] > 1");
204 WARNINGL2(coords[2] >= -1 - NekConstants::kNekZeroTol, "coord[2] < -1");
205 WARNINGL2(coords[2] <= 1 + NekConstants::kNekZeroTol, "coord[2] > 1");
206
207 // Obtain local collapsed coordinate from Cartesian coordinate.
208 LocCoordToLocCollapsed(coords, eta);
209
210 const int nq0 = m_base[0]->GetNumPoints();
211 const int nq1 = m_base[1]->GetNumPoints();
212 const int nq2 = m_base[2]->GetNumPoints();
213
214 Array<OneD, NekDouble> wsp1(nq1 * nq2), wsp2(nq2);
215
216 // Construct the 2D square...
217 const NekDouble *ptr = &physvals[0];
218 for (int i = 0; i < nq1 * nq2; ++i, ptr += nq0)
219 {
220 wsp1[i] = StdExpansion::BaryEvaluate<0>(eta[0], ptr);
221 }
222
223 for (int i = 0; i < nq2; ++i)
224 {
225 wsp2[i] = StdExpansion::BaryEvaluate<1>(eta[1], &wsp1[i * nq1]);
226 }
227
228 return StdExpansion::BaryEvaluate<2>(eta[2], &wsp2[0]);
229}
230
233 const Array<OneD, const NekDouble> &physvals)
234{
235 NekDouble value;
236
237 int Qx = m_base[0]->GetNumPoints();
238 int Qy = m_base[1]->GetNumPoints();
239 int Qz = m_base[2]->GetNumPoints();
240
241 Array<OneD, NekDouble> sumFactorization_qr =
242 Array<OneD, NekDouble>(Qy * Qz);
243 Array<OneD, NekDouble> sumFactorization_r = Array<OneD, NekDouble>(Qz);
244
245 // Lagrangian interpolation matrix
246 NekDouble *interpolatingNodes = nullptr;
247
248 // Interpolate first coordinate direction
249 interpolatingNodes = &I[0]->GetPtr()[0];
250
251 Blas::Dgemv('T', Qx, Qy * Qz, 1.0, &physvals[0], Qx, &interpolatingNodes[0],
252 1, 0.0, &sumFactorization_qr[0], 1);
253
254 // Interpolate in second coordinate direction
255 interpolatingNodes = &I[1]->GetPtr()[0];
256
257 Blas::Dgemv('T', Qy, Qz, 1.0, &sumFactorization_qr[0], Qy,
258 &interpolatingNodes[0], 1, 0.0, &sumFactorization_r[0], 1);
259
260 // Interpolate in third coordinate direction
261 interpolatingNodes = &I[2]->GetPtr()[0];
262 value = Blas::Ddot(Qz, interpolatingNodes, 1, &sumFactorization_r[0], 1);
263
264 return value;
265}
266
268 [[maybe_unused]] const Array<OneD, NekDouble> &coord,
269 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
270 [[maybe_unused]] std::array<NekDouble, 3> &firstOrderDerivs)
271{
272 return 0;
273}
274
275/**
276 * @param inarray Input coefficients.
277 * @param output Output coefficients.
278 * @param mkey Matrix key
279 */
281 const Array<OneD, const NekDouble> &inarray,
283{
284 if (mkey.GetNVarCoeff() == 0 &&
287 {
288 // This implementation is only valid when there are no
289 // coefficients associated to the Laplacian operator
290 int nqtot = GetTotPoints();
291
292 const Array<OneD, const NekDouble> &base0 = m_base[0]->GetBdata();
293 const Array<OneD, const NekDouble> &base1 = m_base[1]->GetBdata();
294 const Array<OneD, const NekDouble> &base2 = m_base[2]->GetBdata();
295
296 // Allocate temporary storage
297 Array<OneD, NekDouble> wsp0(7 * nqtot);
298 Array<OneD, NekDouble> wsp1(wsp0 + nqtot);
299
300 if (!(m_base[0]->Collocation() && m_base[1]->Collocation() &&
301 m_base[2]->Collocation()))
302 {
303 // LAPLACIAN MATRIX OPERATION
304 // wsp0 = u = B * u_hat
305 // wsp1 = du_dxi1 = D_xi1 * wsp0 = D_xi1 * u
306 // wsp2 = du_dxi2 = D_xi2 * wsp0 = D_xi2 * u
307 BwdTrans_SumFacKernel(base0, base1, base2, inarray, wsp0, wsp1,
308 true, true, true);
309 LaplacianMatrixOp_MatFree_Kernel(wsp0, outarray, wsp1);
310 }
311 else
312 {
313 LaplacianMatrixOp_MatFree_Kernel(inarray, outarray, wsp1);
314 }
315 }
316 else
317 {
319 mkey);
320 }
321}
322
324 const Array<OneD, const NekDouble> &inarray,
326{
327 if (mkey.GetNVarCoeff() == 0 &&
329 {
330 using std::max;
331
332 int nquad0 = m_base[0]->GetNumPoints();
333 int nquad1 = m_base[1]->GetNumPoints();
334 int nquad2 = m_base[2]->GetNumPoints();
335 int nmodes0 = m_base[0]->GetNumModes();
336 int nmodes1 = m_base[1]->GetNumModes();
337 int nmodes2 = m_base[2]->GetNumModes();
338 int wspsize = max(nquad0 * nmodes2 * (nmodes1 + nquad1),
339 nquad0 * nquad1 * (nquad2 + nmodes0) +
340 nmodes0 * nmodes1 * nquad2);
341
343
344 const Array<OneD, const NekDouble> &base0 = m_base[0]->GetBdata();
345 const Array<OneD, const NekDouble> &base1 = m_base[1]->GetBdata();
346 const Array<OneD, const NekDouble> &base2 = m_base[2]->GetBdata();
347 Array<OneD, NekDouble> wsp0(8 * wspsize);
348 Array<OneD, NekDouble> wsp1(wsp0 + 1 * wspsize);
349 Array<OneD, NekDouble> wsp2(wsp0 + 2 * wspsize);
350
351 if (!(m_base[0]->Collocation() && m_base[1]->Collocation() &&
352 m_base[2]->Collocation()))
353 {
354 // MASS MATRIX OPERATION
355 // The following is being calculated:
356 // wsp0 = B * u_hat = u
357 // wsp1 = W * wsp0
358 // outarray = B^T * wsp1 = B^T * W * B * u_hat = M * u_hat
359 BwdTrans_SumFacKernel(base0, base1, base2, inarray, wsp0, wsp2,
360 true, true, true);
361 MultiplyByQuadratureMetric(wsp0, wsp1);
362 IProductWRTBase_SumFacKernel(base0, base1, base2, wsp1, outarray,
363 wsp2, true, true, true);
364 LaplacianMatrixOp_MatFree_Kernel(wsp0, wsp1, wsp2);
365 }
366 else
367 {
368 // specialised implementation for the classical spectral
369 // element method
370 MultiplyByQuadratureMetric(inarray, outarray);
371 LaplacianMatrixOp_MatFree_Kernel(inarray, wsp1, wsp2);
372 }
373
374 // outarray = lambda * outarray + wsp1
375 // = (lambda * M + L ) * u_hat
376 Vmath::Svtvp(m_ncoeffs, lambda, &outarray[0], 1, &wsp1[0], 1,
377 &outarray[0], 1);
378 }
379 else
380 {
382 mkey);
383 }
384}
385
387 const Array<OneD, const NekDouble> &inarray)
388{
389 const int nqtot = GetTotPoints();
392 return Vmath::Vsum(nqtot, tmp, 1);
393}
394
396{
397 NEKERROR(ErrorUtil::efatal, "This function is not valid or not defined");
398 return 0;
399}
400
401int StdExpansion3D::v_GetEdgeNcoeffs([[maybe_unused]] const int i) const
402{
403 NEKERROR(ErrorUtil::efatal, "This function is not valid or not defined");
404 return 0;
405}
406
408 [[maybe_unused]] const int tid,
409 [[maybe_unused]] Array<OneD, unsigned int> &maparray,
410 [[maybe_unused]] Array<OneD, int> &signarray,
411 [[maybe_unused]] Orientation traceOrient)
412{
413 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
414}
415
418 Array<OneD, int> &signarray,
419 Orientation traceOrient, int P,
420 int Q)
421{
422 Array<OneD, unsigned int> map1, map2;
423 GetTraceCoeffMap(tid, map1);
424 GetElmtTraceToTraceMap(tid, map2, signarray, traceOrient, P, Q);
425
426 if (maparray.size() != map2.size())
427 {
428 maparray = Array<OneD, unsigned int>(map2.size());
429 }
430
431 for (int i = 0; i < map2.size(); ++i)
432 {
433 maparray[i] = map1[map2[i]];
434 }
435}
436
438 [[maybe_unused]] const int facedir,
439 const LibUtilities::BasisType faceDirBasisType, const int numpoints,
440 const int nummodes)
441{
442
443 switch (faceDirBasisType)
444 {
446 {
447 const LibUtilities::PointsKey pkey(
450 pkey);
451 }
454 {
455 const LibUtilities::PointsKey pkey(
458 pkey);
459 }
461 {
462 const LibUtilities::PointsKey pkey(
465 pkey);
466 }
468 {
469 const LibUtilities::PointsKey pkey(
472 pkey);
473 }
476 {
477 const LibUtilities::PointsKey pkey(
480 pkey);
481 }
482 default:
483 {
484 NEKERROR(ErrorUtil::efatal, "expansion type unknown");
485 break;
486 }
487 }
488
489 // Keep things happy by returning a value.
491}
492
494 const int facedir, const LibUtilities::BasisType faceDirBasisType,
495 const int numpoints, const int nummodes)
496{
497 switch (faceDirBasisType)
498 {
500 {
501 const LibUtilities::PointsKey pkey(
504 pkey);
505 }
509 {
510 switch (facedir)
511 {
512 case 0:
513 {
514 const LibUtilities::PointsKey pkey(
517 nummodes, pkey);
518 }
519 case 1:
520 {
521 const LibUtilities::PointsKey pkey(
522 numpoints, LibUtilities::eGaussRadauMAlpha1Beta0);
524 nummodes, pkey);
525 }
526 default:
527 {
528
529 NEKERROR(ErrorUtil::efatal, "invalid value to flag");
530 break;
531 }
532 }
533 break;
534 }
535
537 {
538 switch (facedir)
539 {
540 case 0:
541 {
542 const LibUtilities::PointsKey pkey(
545 nummodes, pkey);
546 }
547 case 1:
548 {
549 const LibUtilities::PointsKey pkey(
550 numpoints, LibUtilities::eGaussRadauMAlpha1Beta0);
552 nummodes, pkey);
553 }
554 default:
555 {
556 NEKERROR(ErrorUtil::efatal, "invalid value to flag");
557 break;
558 }
559 }
560 break;
561 }
562
567 {
568 switch (facedir)
569 {
570 case 0:
571 {
572 const LibUtilities::PointsKey pkey(
575 nummodes, pkey);
576 }
577 case 1:
578 {
579 const LibUtilities::PointsKey pkey(
580 numpoints, LibUtilities::eGaussRadauMAlpha1Beta0);
582 nummodes, pkey);
583 }
584 default:
585 {
586 NEKERROR(ErrorUtil::efatal, "invalid value to flag");
587 break;
588 }
589 }
590 break;
591 }
592 default:
593 {
594 NEKERROR(ErrorUtil::efatal, "expansion type unknown");
595 break;
596 }
597 }
598
599 // Keep things happy by returning a value.
601}
602} // namespace Nektar::StdRegions
#define WARNINGL2(condition, msg)
Definition: ErrorUtil.hpp:266
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
Describes the specification for a Basis.
Definition: Basis.h:45
Defines a specification for a set of points.
Definition: Points.h:50
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
void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
virtual int v_GetEdgeNcoeffs(const int i) const
NekDouble v_PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals) override
This function evaluates the expansion at a single (arbitrary) point of the domain.
NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray) override
Integrates the specified function over the domain.
void v_GenStdMatBwdDeriv(const int dir, DNekMatSharedPtr &mat) override
void v_GetTraceToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient, int P, int Q) override
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.
The base class for all shapes.
Definition: StdExpansion.h:65
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:134
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:699
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.
void GetTraceCoeffMap(const unsigned int traceid, Array< OneD, unsigned int > &maparray)
Definition: StdExpansion.h:693
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:723
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:124
bool ConstFactorExists(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:133
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 = alpha A x plus beta y where A[m x n].
Definition: Blas.hpp:211
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:163
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:383
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
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:49
@ P
Monomial polynomials .
Definition: BasisType.h:62
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:42
@ eModified_C
Principle Modified Functions .
Definition: BasisType.h:50
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:56
@ eOrtho_C
Principle Orthogonal Functions .
Definition: BasisType.h:46
@ eModifiedPyr_C
Principle Modified Functions.
Definition: BasisType.h:53
@ eOrtho_B
Principle Orthogonal Functions .
Definition: BasisType.h:44
@ 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)
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.hpp:396
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.hpp:608
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825