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 [[maybe_unused]] int numcoeffs,
48 [[maybe_unused]] const LibUtilities::BasisKey &Ba,
49 [[maybe_unused]] const LibUtilities::BasisKey &Bb,
50 [[maybe_unused]] const LibUtilities::BasisKey &Bc)
51{
52}
53
57{
58 const int nquad0 = m_base[0]->GetNumPoints();
59 const int nquad1 = m_base[1]->GetNumPoints();
60 const int nquad2 = m_base[2]->GetNumPoints();
61
62 Array<OneD, NekDouble> wsp(nquad0 * nquad1 * nquad2);
63
64 // copy inarray to wsp in case inarray is used as outarray
65 Vmath::Vcopy(nquad0 * nquad1 * nquad2, &inarray[0], 1, &wsp[0], 1);
66
67 if (out_dx.size() > 0)
68 {
69 NekDouble *D0 = &((m_base[0]->GetD())->GetPtr())[0];
70
71 Blas::Dgemm('N', 'N', nquad0, nquad1 * nquad2, nquad0, 1.0, D0, nquad0,
72 &wsp[0], nquad0, 0.0, &out_dx[0], nquad0);
73 }
74
75 if (out_dy.size() > 0)
76 {
77 NekDouble *D1 = &((m_base[1]->GetD())->GetPtr())[0];
78 for (int j = 0; j < nquad2; ++j)
79 {
80 Blas::Dgemm('N', 'T', nquad0, nquad1, nquad1, 1.0,
81 &wsp[j * nquad0 * nquad1], nquad0, D1, nquad1, 0.0,
82 &out_dy[j * nquad0 * nquad1], nquad0);
83 }
84 }
85
86 if (out_dz.size() > 0)
87 {
88 NekDouble *D2 = &((m_base[2]->GetD())->GetPtr())[0];
89
90 Blas::Dgemm('N', 'T', nquad0 * nquad1, nquad2, nquad2, 1.0, &wsp[0],
91 nquad0 * nquad1, D2, nquad2, 0.0, &out_dz[0],
92 nquad0 * nquad1);
93 }
94}
95
100 const Array<OneD, const NekDouble> &inarray,
102 bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
103{
104 v_BwdTrans_SumFacKernel(base0, base1, base2, inarray, outarray, wsp,
105 doCheckCollDir0, doCheckCollDir1, doCheckCollDir2);
106}
107
109 const Array<OneD, const NekDouble> &base0,
110 const Array<OneD, const NekDouble> &base1,
111 const Array<OneD, const NekDouble> &base2,
112 const Array<OneD, const NekDouble> &inarray,
114 bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
115{
116 v_IProductWRTBase_SumFacKernel(base0, base1, base2, inarray, outarray, wsp,
117 doCheckCollDir0, doCheckCollDir1,
118 doCheckCollDir2);
119}
120
122{
123 ASSERTL1((dir == 0) || (dir == 1) || (dir == 2), "Invalid direction.");
124
125 const int nq0 = m_base[0]->GetNumPoints();
126 const int nq1 = m_base[1]->GetNumPoints();
127 const int nq2 = m_base[2]->GetNumPoints();
128 const int nq = nq0 * nq1 * nq2;
129 const int nm0 = m_base[0]->GetNumModes();
130 const int nm1 = m_base[1]->GetNumModes();
131
132 Array<OneD, NekDouble> alloc(4 * nq + m_ncoeffs + nm0 * nq2 * (nq1 + nm1),
133 0.0);
134 Array<OneD, NekDouble> tmp1(alloc); // Quad metric
135 Array<OneD, NekDouble> tmp2(alloc + nq); // Dir1 metric
136 Array<OneD, NekDouble> tmp3(alloc + 2 * nq); // Dir2 metric
137 Array<OneD, NekDouble> tmp4(alloc + 3 * nq); // Dir3 metric
138 Array<OneD, NekDouble> tmp5(alloc + 4 * nq); // iprod tmp
139 Array<OneD, NekDouble> wsp(tmp5 + m_ncoeffs); // Wsp
140 switch (dir)
141 {
142 case 0:
143 for (int i = 0; i < nq; i++)
144 {
145 tmp2[i] = 1.0;
147 m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
148 m_base[2]->GetBdata(), tmp2, tmp5, wsp, false, true, true);
149
150 tmp2[i] = 0.0;
151
152 for (int j = 0; j < m_ncoeffs; j++)
153 {
154 (*mat)(j, i) = tmp5[j];
155 }
156 }
157 break;
158 case 1:
159 for (int i = 0; i < nq; i++)
160 {
161 tmp2[i] = 1.0;
163 m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
164 m_base[2]->GetBdata(), tmp2, tmp5, wsp, true, false, true);
165
166 tmp2[i] = 0.0;
167
168 for (int j = 0; j < m_ncoeffs; j++)
169 {
170 (*mat)(j, i) = tmp5[j];
171 }
172 }
173 break;
174 case 2:
175 for (int i = 0; i < nq; i++)
176 {
177 tmp2[i] = 1.0;
179 m_base[0]->GetBdata(), m_base[1]->GetBdata(),
180 m_base[2]->GetDbdata(), tmp2, tmp5, wsp, true, true, false);
181 tmp2[i] = 0.0;
182
183 for (int j = 0; j < m_ncoeffs; j++)
184 {
185 (*mat)(j, i) = tmp5[j];
186 }
187 }
188 break;
189 default:
190 NEKERROR(ErrorUtil::efatal, "Not a 2D expansion.");
191 break;
192 }
193}
194
196 const Array<OneD, const NekDouble> &coords,
197 const Array<OneD, const NekDouble> &physvals)
198{
200
201 WARNINGL2(coords[0] >= -1 - NekConstants::kNekZeroTol, "coord[0] < -1");
202 WARNINGL2(coords[0] <= 1 + NekConstants::kNekZeroTol, "coord[0] > 1");
203 WARNINGL2(coords[1] >= -1 - NekConstants::kNekZeroTol, "coord[1] < -1");
204 WARNINGL2(coords[1] <= 1 + NekConstants::kNekZeroTol, "coord[1] > 1");
205 WARNINGL2(coords[2] >= -1 - NekConstants::kNekZeroTol, "coord[2] < -1");
206 WARNINGL2(coords[2] <= 1 + NekConstants::kNekZeroTol, "coord[2] > 1");
207
208 // Obtain local collapsed coordinate from Cartesian coordinate.
209 LocCoordToLocCollapsed(coords, eta);
210
211 const int nq0 = m_base[0]->GetNumPoints();
212 const int nq1 = m_base[1]->GetNumPoints();
213 const int nq2 = m_base[2]->GetNumPoints();
214
215 Array<OneD, NekDouble> wsp1(nq1 * nq2), wsp2(nq2);
216
217 // Construct the 2D square...
218 const NekDouble *ptr = &physvals[0];
219 for (int i = 0; i < nq1 * nq2; ++i, ptr += nq0)
220 {
221 wsp1[i] = StdExpansion::BaryEvaluate<0>(eta[0], ptr);
222 }
223
224 for (int i = 0; i < nq2; ++i)
225 {
226 wsp2[i] = StdExpansion::BaryEvaluate<1>(eta[1], &wsp1[i * nq1]);
227 }
228
229 return StdExpansion::BaryEvaluate<2>(eta[2], &wsp2[0]);
230}
231
234 const Array<OneD, const NekDouble> &physvals)
235{
236 NekDouble value;
237
238 int Qx = m_base[0]->GetNumPoints();
239 int Qy = m_base[1]->GetNumPoints();
240 int Qz = m_base[2]->GetNumPoints();
241
242 Array<OneD, NekDouble> sumFactorization_qr =
243 Array<OneD, NekDouble>(Qy * Qz);
244 Array<OneD, NekDouble> sumFactorization_r = Array<OneD, NekDouble>(Qz);
245
246 // Lagrangian interpolation matrix
247 NekDouble *interpolatingNodes = nullptr;
248
249 // Interpolate first coordinate direction
250 interpolatingNodes = &I[0]->GetPtr()[0];
251
252 Blas::Dgemv('T', Qx, Qy * Qz, 1.0, &physvals[0], Qx, &interpolatingNodes[0],
253 1, 0.0, &sumFactorization_qr[0], 1);
254
255 // Interpolate in second coordinate direction
256 interpolatingNodes = &I[1]->GetPtr()[0];
257
258 Blas::Dgemv('T', Qy, Qz, 1.0, &sumFactorization_qr[0], Qy,
259 &interpolatingNodes[0], 1, 0.0, &sumFactorization_r[0], 1);
260
261 // Interpolate in third coordinate direction
262 interpolatingNodes = &I[2]->GetPtr()[0];
263 value = Blas::Ddot(Qz, interpolatingNodes, 1, &sumFactorization_r[0], 1);
264
265 return value;
266}
267
269 [[maybe_unused]] const Array<OneD, NekDouble> &coord,
270 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
271 [[maybe_unused]] std::array<NekDouble, 3> &firstOrderDerivs)
272{
273 return 0;
274}
275
276/**
277 * @param inarray Input coefficients.
278 * @param output Output coefficients.
279 * @param mkey Matrix key
280 */
282 const Array<OneD, const NekDouble> &inarray,
284{
285 if (mkey.GetNVarCoeff() == 0 &&
288 {
289 // This implementation is only valid when there are no
290 // coefficients associated to the Laplacian operator
291 int nqtot = GetTotPoints();
292
293 const Array<OneD, const NekDouble> &base0 = m_base[0]->GetBdata();
294 const Array<OneD, const NekDouble> &base1 = m_base[1]->GetBdata();
295 const Array<OneD, const NekDouble> &base2 = m_base[2]->GetBdata();
296
297 // Allocate temporary storage
298 Array<OneD, NekDouble> wsp0(7 * nqtot);
299 Array<OneD, NekDouble> wsp1(wsp0 + nqtot);
300
301 if (!(m_base[0]->Collocation() && m_base[1]->Collocation() &&
302 m_base[2]->Collocation()))
303 {
304 // LAPLACIAN MATRIX OPERATION
305 // wsp0 = u = B * u_hat
306 // wsp1 = du_dxi1 = D_xi1 * wsp0 = D_xi1 * u
307 // wsp2 = du_dxi2 = D_xi2 * wsp0 = D_xi2 * u
308 BwdTrans_SumFacKernel(base0, base1, base2, inarray, wsp0, wsp1,
309 true, true, true);
310 LaplacianMatrixOp_MatFree_Kernel(wsp0, outarray, wsp1);
311 }
312 else
313 {
314 LaplacianMatrixOp_MatFree_Kernel(inarray, outarray, wsp1);
315 }
316 }
317 else
318 {
320 mkey);
321 }
322}
323
325 const Array<OneD, const NekDouble> &inarray,
327{
328 if (mkey.GetNVarCoeff() == 0 &&
330 {
331 using std::max;
332
333 int nquad0 = m_base[0]->GetNumPoints();
334 int nquad1 = m_base[1]->GetNumPoints();
335 int nquad2 = m_base[2]->GetNumPoints();
336 int nmodes0 = m_base[0]->GetNumModes();
337 int nmodes1 = m_base[1]->GetNumModes();
338 int nmodes2 = m_base[2]->GetNumModes();
339 int wspsize = max(nquad0 * nmodes2 * (nmodes1 + nquad1),
340 nquad0 * nquad1 * (nquad2 + nmodes0) +
341 nmodes0 * nmodes1 * nquad2);
342
344
345 const Array<OneD, const NekDouble> &base0 = m_base[0]->GetBdata();
346 const Array<OneD, const NekDouble> &base1 = m_base[1]->GetBdata();
347 const Array<OneD, const NekDouble> &base2 = m_base[2]->GetBdata();
348 Array<OneD, NekDouble> wsp0(8 * wspsize);
349 Array<OneD, NekDouble> wsp1(wsp0 + 1 * wspsize);
350 Array<OneD, NekDouble> wsp2(wsp0 + 2 * wspsize);
351
352 if (!(m_base[0]->Collocation() && m_base[1]->Collocation() &&
353 m_base[2]->Collocation()))
354 {
355 // MASS MATRIX OPERATION
356 // The following is being calculated:
357 // wsp0 = B * u_hat = u
358 // wsp1 = W * wsp0
359 // outarray = B^T * wsp1 = B^T * W * B * u_hat = M * u_hat
360 BwdTrans_SumFacKernel(base0, base1, base2, inarray, wsp0, wsp2,
361 true, true, true);
362 MultiplyByQuadratureMetric(wsp0, wsp1);
363 IProductWRTBase_SumFacKernel(base0, base1, base2, wsp1, outarray,
364 wsp2, true, true, true);
365 LaplacianMatrixOp_MatFree_Kernel(wsp0, wsp1, wsp2);
366 }
367 else
368 {
369 // specialised implementation for the classical spectral
370 // element method
371 MultiplyByQuadratureMetric(inarray, outarray);
372 LaplacianMatrixOp_MatFree_Kernel(inarray, wsp1, wsp2);
373 }
374
375 // outarray = lambda * outarray + wsp1
376 // = (lambda * M + L ) * u_hat
377 Vmath::Svtvp(m_ncoeffs, lambda, &outarray[0], 1, &wsp1[0], 1,
378 &outarray[0], 1);
379 }
380 else
381 {
383 mkey);
384 }
385}
386
388 const Array<OneD, const NekDouble> &inarray)
389{
390 const int nqtot = GetTotPoints();
393 return Vmath::Vsum(nqtot, tmp, 1);
394}
395
397{
398 NEKERROR(ErrorUtil::efatal, "This function is not valid or not defined");
399 return 0;
400}
401
402int StdExpansion3D::v_GetEdgeNcoeffs([[maybe_unused]] const int i) const
403{
404 NEKERROR(ErrorUtil::efatal, "This function is not valid or not defined");
405 return 0;
406}
407
409 [[maybe_unused]] const int tid,
410 [[maybe_unused]] Array<OneD, unsigned int> &maparray,
411 [[maybe_unused]] Array<OneD, int> &signarray,
412 [[maybe_unused]] Orientation traceOrient)
413{
414 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
415}
416
419 Array<OneD, int> &signarray,
420 Orientation traceOrient, int P,
421 int Q)
422{
423 Array<OneD, unsigned int> map1, map2;
424 GetTraceCoeffMap(tid, map1);
425 GetElmtTraceToTraceMap(tid, map2, signarray, traceOrient, P, Q);
426
427 if (maparray.size() != map2.size())
428 {
429 maparray = Array<OneD, unsigned int>(map2.size());
430 }
431
432 for (int i = 0; i < map2.size(); ++i)
433 {
434 maparray[i] = map1[map2[i]];
435 }
436}
437
439 [[maybe_unused]] const int facedir,
440 const LibUtilities::BasisType faceDirBasisType, const int numpoints,
441 const int nummodes)
442{
443
444 switch (faceDirBasisType)
445 {
447 {
448 const LibUtilities::PointsKey pkey(
451 pkey);
452 }
455 {
456 const LibUtilities::PointsKey pkey(
459 pkey);
460 }
462 {
463 const LibUtilities::PointsKey pkey(
466 pkey);
467 }
469 {
470 const LibUtilities::PointsKey pkey(
473 pkey);
474 }
477 {
478 const LibUtilities::PointsKey pkey(
481 pkey);
482 }
483 default:
484 {
485 NEKERROR(ErrorUtil::efatal, "expansion type unknown");
486 break;
487 }
488 }
489
490 // Keep things happy by returning a value.
492}
493
495 const int facedir, const LibUtilities::BasisType faceDirBasisType,
496 const int numpoints, const int nummodes)
497{
498 switch (faceDirBasisType)
499 {
501 {
502 const LibUtilities::PointsKey pkey(
505 pkey);
506 }
510 {
511 switch (facedir)
512 {
513 case 0:
514 {
515 const LibUtilities::PointsKey pkey(
518 nummodes, pkey);
519 }
520 case 1:
521 {
522 const LibUtilities::PointsKey pkey(
523 numpoints, LibUtilities::eGaussRadauMAlpha1Beta0);
525 nummodes, pkey);
526 }
527 default:
528 {
529
530 NEKERROR(ErrorUtil::efatal, "invalid value to flag");
531 break;
532 }
533 }
534 break;
535 }
536
538 {
539 switch (facedir)
540 {
541 case 0:
542 {
543 const LibUtilities::PointsKey pkey(
546 nummodes, pkey);
547 }
548 case 1:
549 {
550 const LibUtilities::PointsKey pkey(
551 numpoints, LibUtilities::eGaussRadauMAlpha1Beta0);
553 nummodes, pkey);
554 }
555 default:
556 {
557 NEKERROR(ErrorUtil::efatal, "invalid value to flag");
558 break;
559 }
560 }
561 break;
562 }
563
568 {
569 switch (facedir)
570 {
571 case 0:
572 {
573 const LibUtilities::PointsKey pkey(
576 nummodes, pkey);
577 }
578 case 1:
579 {
580 const LibUtilities::PointsKey pkey(
581 numpoints, LibUtilities::eGaussRadauMAlpha1Beta0);
583 nummodes, pkey);
584 }
585 default:
586 {
587 NEKERROR(ErrorUtil::efatal, "invalid value to flag");
588 break;
589 }
590 }
591 break;
592 }
593 default:
594 {
595 NEKERROR(ErrorUtil::efatal, "expansion type unknown");
596 break;
597 }
598 }
599
600 // Keep things happy by returning a value.
602}
603} // 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.
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