Nektar++
StdExpansion.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: StdExpansion.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: Definition of methods in class StdExpansion which is
32// the base class to all expansion shapes
33//
34///////////////////////////////////////////////////////////////////////////////
35
36#include <LibUtilities/Foundations/ManagerAccess.h> // for BasisManager, etc
38
40{
41/** \brief Default constructor */
43{
44}
45
46StdExpansion::StdExpansion(const int numcoeffs, const int numbases,
47 const LibUtilities::BasisKey &Ba,
48 const LibUtilities::BasisKey &Bb,
49 const LibUtilities::BasisKey &Bc)
50 : m_base(numbases), m_elmt_id(0), m_ncoeffs(numcoeffs),
51 m_stdMatrixManager(std::bind(&StdExpansion::CreateStdMatrix, this,
52 std::placeholders::_1),
53 std::string("StdExpansionStdMatrix")),
54 m_stdStaticCondMatrixManager(
55 std::bind(&StdExpansion::CreateStdStaticCondMatrix, this,
56 std::placeholders::_1),
57 std::string("StdExpansionStdStaticCondMatrix"))
58{
59 switch (m_base.size())
60 {
61 case 3:
63 "NULL Basis attempting to be used.");
65 /* Falls through. */
66 case 2:
68 "NULL Basis attempting to be used.");
70 /* Falls through. */
71 case 1:
73 "NULL Basis attempting to be used.");
75 break;
76 default:
77 break;
78 // ASSERTL0(false, "numbases incorrectly specified");
79 };
80
81} // end constructor
82
84 : std::enable_shared_from_this<StdExpansion>(T), m_base(T.m_base),
85 m_elmt_id(T.m_elmt_id), m_ncoeffs(T.m_ncoeffs),
86 m_stdMatrixManager(T.m_stdMatrixManager),
87 m_stdStaticCondMatrixManager(T.m_stdStaticCondMatrixManager)
88{
89}
90
91// Destructor
93{
94}
95
98{
99 NekDouble val;
100 int ntot = GetTotPoints();
101 Array<OneD, NekDouble> wsp(ntot);
102
103 if (sol == NullNekDouble1DArray)
104 {
105 Vmath::Vabs(ntot, phys, 1, wsp, 1);
106 }
107 else
108 {
109 Vmath::Vsub(ntot, sol, 1, phys, 1, wsp, 1);
110 Vmath::Vabs(ntot, wsp, 1, wsp, 1);
111 }
112
113 val = Vmath::Vamax(ntot, wsp, 1);
114 return val;
115}
116
119{
120 NekDouble val;
121 int ntot = GetTotPoints();
122 Array<OneD, NekDouble> wsp(ntot);
123
124 if (sol.size() == 0)
125 {
126 Vmath::Vmul(ntot, phys, 1, phys, 1, wsp, 1);
127 }
128 else
129 {
130 Vmath::Vsub(ntot, sol, 1, phys, 1, wsp, 1);
131 Vmath::Vmul(ntot, wsp, 1, wsp, 1, wsp, 1);
132 }
133
134 val = v_Integral(wsp);
135
136 return (val < 0.0) ? 0.0 : sqrt(val);
137}
138
141{
142 int i;
143 NekDouble val;
144 int ntot = GetTotPoints();
145 int coordim = v_GetCoordim();
146 Array<OneD, NekDouble> wsp(3 * ntot);
147 Array<OneD, NekDouble> wsp_deriv = wsp + ntot;
148 Array<OneD, NekDouble> sum = wsp_deriv + ntot;
149
150 if (sol == NullNekDouble1DArray)
151 {
152 Vmath::Vcopy(ntot, phys, 1, wsp, 1);
153 Vmath::Vmul(ntot, phys, 1, phys, 1, sum, 1);
154 }
155 else
156 {
157 Vmath::Vsub(ntot, sol, 1, phys, 1, wsp, 1);
158 Vmath::Vmul(ntot, wsp, 1, wsp, 1, sum, 1);
159 }
160
161 for (i = 0; i < coordim; ++i)
162 {
163 v_PhysDeriv(i, wsp, wsp_deriv);
164 Vmath::Vvtvp(ntot, wsp_deriv, 1, wsp_deriv, 1, sum, 1, sum, 1);
165 }
166
167 val = sqrt(v_Integral(sum));
168
169 return val;
170}
171
173 const StdMatrixKey &mkey)
174{
175 DNekBlkMatSharedPtr returnval;
176
177 DNekMatSharedPtr mat = GetStdMatrix(mkey);
178 int nbdry = NumBndryCoeffs(); // also checks to see if this is a boundary
179 // interior decomposed expansion
180 int nint = m_ncoeffs - nbdry;
186
187 int i, j;
188
189 Array<OneD, unsigned int> bmap(nbdry);
190 Array<OneD, unsigned int> imap(nint);
191 GetBoundaryMap(bmap);
192 GetInteriorMap(imap);
193
194 for (i = 0; i < nbdry; ++i)
195 {
196 for (j = 0; j < nbdry; ++j)
197 {
198 (*A)(i, j) = (*mat)(bmap[i], bmap[j]);
199 }
200
201 for (j = 0; j < nint; ++j)
202 {
203 (*B)(i, j) = (*mat)(bmap[i], imap[j]);
204 }
205 }
206
207 for (i = 0; i < nint; ++i)
208 {
209 for (j = 0; j < nbdry; ++j)
210 {
211 (*C)(i, j) = (*mat)(imap[i], bmap[j]);
212 }
213
214 for (j = 0; j < nint; ++j)
215 {
216 (*D)(i, j) = (*mat)(imap[i], imap[j]);
217 }
218 }
219
220 // Calculate static condensed system
221 if (nint)
222 {
223 D->Invert();
224 (*B) = (*B) * (*D);
225 (*A) = (*A) - (*B) * (*C);
226 }
227
228 // set up block matrix system
229 Array<OneD, unsigned int> exp_size(2);
230 exp_size[0] = nbdry;
231 exp_size[1] = nint;
232 returnval =
234
235 returnval->SetBlock(0, 0, A);
236 returnval->SetBlock(0, 1, B);
237 returnval->SetBlock(1, 0, C);
238 returnval->SetBlock(1, 1, D);
239
240 return returnval;
241}
242
244{
245 int i;
246 DNekMatSharedPtr returnval;
247
248 switch (mkey.GetMatrixType())
249 {
250 case eInvMass:
251 {
252 StdMatrixKey masskey(eMass, mkey.GetShapeType(), *this,
254 mkey.GetNodalPointsType());
255 DNekMatSharedPtr mmat = GetStdMatrix(masskey);
256
258 *mmat); // Populate standard mass matrix.
259 returnval->Invert();
260 }
261 break;
262 case eInvNBasisTrans:
263 {
264 StdMatrixKey tmpkey(eNBasisTrans, mkey.GetShapeType(), *this,
266 mkey.GetNodalPointsType());
267 DNekMatSharedPtr tmpmat = GetStdMatrix(tmpkey);
269 *tmpmat); // Populate matrix.
270 returnval->Invert();
271 }
272 break;
273 case eBwdMat:
274 {
275 int nq = GetTotPoints();
277 Array<OneD, NekDouble> tmpout(nq);
278
279 returnval =
281 Array<OneD, NekDouble> Bwd_data = returnval->GetPtr();
282
284 this->DetShapeType(), *this);
285 DNekMatSharedPtr MatBwdTrans = GetStdMatrix(matkey);
286 Array<OneD, NekDouble> BwdTrans_data = MatBwdTrans->GetPtr();
287
288 for (i = 0; i < m_ncoeffs; ++i)
289 {
290 Array<OneD, NekDouble> tmpinn = BwdTrans_data + nq * i;
291 tmpout = Bwd_data + i;
292
293 Vmath::Vcopy(nq, tmpinn, 1, tmpout, m_ncoeffs);
294 }
295 }
296 break;
297 case eBwdTrans:
298 {
299 int nq = GetTotPoints();
301 Array<OneD, NekDouble> tmpout(nq);
302
303 returnval =
305
306 for (i = 0; i < m_ncoeffs; ++i)
307 {
308 Vmath::Zero(m_ncoeffs, tmpin, 1);
309 tmpin[i] = 1.0;
310
311 BwdTrans_SumFac(tmpin, tmpout);
312
313 Vmath::Vcopy(nq, tmpout.get(), 1,
314 returnval->GetRawPtr() + i * nq, 1);
315 }
316 }
317 break;
318 case eIProductWRTBase:
319 {
320 int nq = GetTotPoints();
321 Array<OneD, NekDouble> tmpin(nq);
323
324 returnval =
326
327 for (i = 0; i < nq; ++i)
328 {
329 Vmath::Zero(nq, tmpin, 1);
330 tmpin[i] = 1.0;
331
332 IProductWRTBase_SumFac(tmpin, tmpout);
333
334 Vmath::Vcopy(m_ncoeffs, tmpout.get(), 1,
335 returnval->GetRawPtr() + i * m_ncoeffs, 1);
336 }
337 }
338 break;
340 {
341 int nq = GetTotPoints();
342 Array<OneD, NekDouble> tmpin(nq);
344
345 returnval =
347
348 for (i = 0; i < nq; ++i)
349 {
350 Vmath::Zero(nq, tmpin, 1);
351 tmpin[i] = 1.0;
352
353 IProductWRTDerivBase_SumFac(0, tmpin, tmpout);
354
355 Vmath::Vcopy(m_ncoeffs, tmpout.get(), 1,
356 returnval->GetRawPtr() + i * m_ncoeffs, 1);
357 }
358 }
359 break;
361 {
362 int nq = GetTotPoints();
363 Array<OneD, NekDouble> tmpin(nq);
365
366 returnval =
368
369 for (i = 0; i < nq; ++i)
370 {
371 Vmath::Zero(nq, tmpin, 1);
372 tmpin[i] = 1.0;
373
374 IProductWRTDerivBase_SumFac(1, tmpin, tmpout);
375
376 Vmath::Vcopy(m_ncoeffs, tmpout.get(), 1,
377 returnval->GetRawPtr() + i * m_ncoeffs, 1);
378 }
379 }
380 break;
382 {
383 int nq = GetTotPoints();
384 Array<OneD, NekDouble> tmpin(nq);
386
387 returnval =
389
390 for (i = 0; i < nq; ++i)
391 {
392 Vmath::Zero(nq, tmpin, 1);
393 tmpin[i] = 1.0;
394
395 IProductWRTDerivBase_SumFac(2, tmpin, tmpout);
396
397 Vmath::Vcopy(m_ncoeffs, tmpout.get(), 1,
398 returnval->GetRawPtr() + i * m_ncoeffs, 1);
399 }
400 }
401 break;
402 case eDerivBase0:
403 {
404 int nq = GetTotPoints();
405 returnval =
407 GenStdMatBwdDeriv(0, returnval);
408 }
409 break;
410 case eDerivBase1:
411 {
412 int nq = GetTotPoints();
413 returnval =
415 GenStdMatBwdDeriv(1, returnval);
416 }
417 break;
418 case eDerivBase2:
419 {
420 int nq = GetTotPoints();
421 returnval =
423 GenStdMatBwdDeriv(2, returnval);
424 }
425 break;
427 {
428 // check to see if equispaced basis
429 int nummodes = m_base[0]->GetNumModes();
430 bool equispaced = true;
431 for (i = 1; i < m_base.size(); ++i)
432 {
433 if (m_base[i]->GetNumModes() != nummodes)
434 {
435 equispaced = false;
436 }
437 }
438
439 ASSERTL0(equispaced,
440 "Currently need to have same num modes in all "
441 "directionmodes to use EquiSpacedToCoeff method");
442
443 int ntot = GetTotPoints();
444 Array<OneD, NekDouble> qmode(ntot);
446
447 returnval =
449 for (i = 0; i < m_ncoeffs; ++i)
450 {
451 // Get mode at quadrature points
452 FillMode(i, qmode);
453
454 // interpolate to equi spaced
455 PhysInterpToSimplexEquiSpaced(qmode, emode, nummodes);
456
457 // fill matrix
458 Vmath::Vcopy(m_ncoeffs, &emode[0], 1,
459 returnval->GetRawPtr() + i * m_ncoeffs, 1);
460 }
461 // invert matrix
462 returnval->Invert();
463 }
464 break;
465 case eMass:
466 case eHelmholtz:
467 case eLaplacian:
468 case eLaplacian00:
469 case eLaplacian01:
470 case eLaplacian02:
471 case eLaplacian10:
472 case eLaplacian11:
473 case eLaplacian12:
474 case eLaplacian20:
475 case eLaplacian21:
476 case eLaplacian22:
477 case eWeakDeriv0:
478 case eWeakDeriv1:
479 case eWeakDeriv2:
482 case eLinearAdvection:
485 {
487 returnval =
489 DNekMat &Mat = *returnval;
490
491 for (i = 0; i < m_ncoeffs; ++i)
492 {
493 Vmath::Zero(m_ncoeffs, tmp, 1);
494 tmp[i] = 1.0;
495
496 GeneralMatrixOp_MatFree(tmp, tmp, mkey);
497
498 Vmath::Vcopy(m_ncoeffs, &tmp[0], 1,
499 &(Mat.GetPtr())[0] + i * m_ncoeffs, 1);
500 }
501 }
502 break;
503 default:
504 {
506 "This type of matrix, " +
507 static_cast<std::string>(
509 ", can not be created using a general approach");
510 }
511 break;
512 }
513
514 return returnval;
515}
516
518 Array<OneD, NekDouble> &outarray,
519 const StdMatrixKey &mkey)
520{
521 switch (mkey.GetMatrixType())
522 {
523 case eMass:
524 MassMatrixOp(inarray, outarray, mkey);
525 break;
526 case eWeakDeriv0:
527 WeakDerivMatrixOp(0, inarray, outarray, mkey);
528 break;
529 case eWeakDeriv1:
530 WeakDerivMatrixOp(1, inarray, outarray, mkey);
531 break;
532 case eWeakDeriv2:
533 WeakDerivMatrixOp(2, inarray, outarray, mkey);
534 break;
536 WeakDirectionalDerivMatrixOp(inarray, outarray, mkey);
537 break;
539 MassLevelCurvatureMatrixOp(inarray, outarray, mkey);
540 break;
541 case eLinearAdvection:
542 LinearAdvectionMatrixOp(inarray, outarray, mkey);
543 break;
545 LinearAdvectionDiffusionReactionMatrixOp(inarray, outarray, mkey,
546 false);
547 break;
549 LinearAdvectionDiffusionReactionMatrixOp(inarray, outarray, mkey);
550 break;
551 case eLaplacian:
552 LaplacianMatrixOp(inarray, outarray, mkey);
553 break;
554 case eLaplacian00:
555 LaplacianMatrixOp(0, 0, inarray, outarray, mkey);
556 break;
557 case eLaplacian01:
558 LaplacianMatrixOp(0, 1, inarray, outarray, mkey);
559 break;
560 case eLaplacian02:
561 LaplacianMatrixOp(0, 2, inarray, outarray, mkey);
562 break;
563 case eLaplacian10:
564 LaplacianMatrixOp(1, 0, inarray, outarray, mkey);
565 break;
566 case eLaplacian11:
567 LaplacianMatrixOp(1, 1, inarray, outarray, mkey);
568 break;
569 case eLaplacian12:
570 LaplacianMatrixOp(1, 2, inarray, outarray, mkey);
571 break;
572 case eLaplacian20:
573 LaplacianMatrixOp(2, 0, inarray, outarray, mkey);
574 break;
575 case eLaplacian21:
576 LaplacianMatrixOp(2, 1, inarray, outarray, mkey);
577 break;
578 case eLaplacian22:
579 LaplacianMatrixOp(2, 2, inarray, outarray, mkey);
580 break;
581 case eHelmholtz:
582 HelmholtzMatrixOp(inarray, outarray, mkey);
583 break;
584 default:
586 "This matrix does not have an operator");
587 break;
588 }
589}
590
592 const Array<OneD, const NekDouble> &inarray,
593 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
594{
595 switch (mkey.GetMatrixType())
596 {
597 case eMass:
598 MassMatrixOp_MatFree(inarray, outarray, mkey);
599 break;
600 case eWeakDeriv0:
601 WeakDerivMatrixOp_MatFree(0, inarray, outarray, mkey);
602 break;
603 case eWeakDeriv1:
604 WeakDerivMatrixOp_MatFree(1, inarray, outarray, mkey);
605 break;
606 case eWeakDeriv2:
607 WeakDerivMatrixOp_MatFree(2, inarray, outarray, mkey);
608 break;
610 WeakDirectionalDerivMatrixOp_MatFree(inarray, outarray, mkey);
611 break;
613 MassLevelCurvatureMatrixOp_MatFree(inarray, outarray, mkey);
614 break;
615 case eLinearAdvection:
616 LinearAdvectionMatrixOp_MatFree(inarray, outarray, mkey);
617 break;
620 mkey, false);
621 break;
624 mkey);
625 break;
626 case eLaplacian:
627 LaplacianMatrixOp_MatFree(inarray, outarray, mkey);
628 break;
629 case eLaplacian00:
630 LaplacianMatrixOp_MatFree(0, 0, inarray, outarray, mkey);
631 break;
632 case eLaplacian01:
633 LaplacianMatrixOp_MatFree(0, 1, inarray, outarray, mkey);
634 break;
635 case eLaplacian02:
636 LaplacianMatrixOp_MatFree(0, 2, inarray, outarray, mkey);
637 break;
638 case eLaplacian10:
639 LaplacianMatrixOp_MatFree(1, 0, inarray, outarray, mkey);
640 break;
641 case eLaplacian11:
642 LaplacianMatrixOp_MatFree(1, 1, inarray, outarray, mkey);
643 break;
644 case eLaplacian12:
645 LaplacianMatrixOp_MatFree(1, 2, inarray, outarray, mkey);
646 break;
647 case eLaplacian20:
648 LaplacianMatrixOp_MatFree(2, 0, inarray, outarray, mkey);
649 break;
650 case eLaplacian21:
651 LaplacianMatrixOp_MatFree(2, 1, inarray, outarray, mkey);
652 break;
653 case eLaplacian22:
654 LaplacianMatrixOp_MatFree(2, 2, inarray, outarray, mkey);
655 break;
656 case eHelmholtz:
657 HelmholtzMatrixOp_MatFree(inarray, outarray, mkey);
658 break;
659 default:
661 "This matrix does not have an operator");
662 break;
663 }
664}
665
667 const Array<OneD, const NekDouble> &inarray,
668 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
669{
670 int nq = GetTotPoints();
672
673 v_BwdTrans(inarray, tmp);
674
675 if (mkey.HasVarCoeff(eVarCoeffMass))
676 {
677 Vmath::Vmul(nq, mkey.GetVarCoeff(eVarCoeffMass), 1, tmp, 1, tmp, 1);
678 }
679
680 v_IProductWRTBase(tmp, outarray);
681}
682
684 const int k1, const int k2, const Array<OneD, const NekDouble> &inarray,
685 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
686{
687 ASSERTL1(k1 >= 0 && k1 < GetCoordim(), "invalid first argument");
688 ASSERTL1(k2 >= 0 && k2 < GetCoordim(), "invalid second argument");
689
690 int nq = GetTotPoints();
692 Array<OneD, NekDouble> dtmp(nq);
693 VarCoeffType varcoefftypes[3][3] = {
697
698 ConstFactorType constcoefftypes[3][3] = {
702
703 v_BwdTrans(inarray, tmp);
704 v_PhysDeriv(k2, tmp, dtmp);
706 {
707 if (k1 == k2)
708 {
709 // By default, k1 == k2 has \sigma = 1 (diagonal entries)
710 if (mkey.HasVarCoeff(varcoefftypes[k1][k1]))
711 {
712 Vmath::Vmul(nq, mkey.GetVarCoeff(varcoefftypes[k1][k1]), 1,
713 dtmp, 1, dtmp, 1);
714 }
715 v_IProductWRTDerivBase_SumFac(k1, dtmp, outarray);
716 }
717 else
718 {
719 // By default, k1 != k2 has \sigma = 0 (off-diagonal entries)
720 if (mkey.HasVarCoeff(varcoefftypes[k1][k2]))
721 {
722 Vmath::Vmul(nq, mkey.GetVarCoeff(varcoefftypes[k1][k2]), 1,
723 dtmp, 1, dtmp, 1);
724 v_IProductWRTDerivBase_SumFac(k1, dtmp, outarray);
725 }
726 else if (mkey.HasVarCoeff(
727 varcoefftypes[k2][k1])) // Check symmetric varcoeff
728 {
729 Vmath::Vmul(nq, mkey.GetVarCoeff(varcoefftypes[k2][k1]), 1,
730 dtmp, 1, dtmp, 1);
731 v_IProductWRTDerivBase_SumFac(k1, dtmp, outarray);
732 }
733 else
734 {
735 Vmath::Zero(GetNcoeffs(), outarray, 1);
736 }
737 }
738 }
739 else if (mkey.ConstFactorExists(eFactorCoeffD00) &&
741 {
742 if (k1 == k2)
743 {
744 // By default, k1 == k2 has \sigma = 1 (diagonal entries)
745 if (mkey.ConstFactorExists(constcoefftypes[k1][k1]))
746 {
747 Vmath::Smul(nq, mkey.GetConstFactor(constcoefftypes[k1][k1]),
748 dtmp, 1, dtmp, 1);
749 }
750 v_IProductWRTDerivBase(k1, dtmp, outarray);
751 }
752 else
753 {
754 // By default, k1 != k2 has \sigma = 0 (off-diagonal entries)
755 if (mkey.ConstFactorExists(constcoefftypes[k1][k2]))
756 {
757 Vmath::Smul(nq, mkey.GetConstFactor(constcoefftypes[k1][k2]),
758 dtmp, 1, dtmp, 1);
759 v_IProductWRTDerivBase(k1, dtmp, outarray);
760 }
761 else
762 {
763 Vmath::Zero(GetNcoeffs(), outarray, 1);
764 }
765 }
766 }
767 else
768 {
769 // Multiply by svv tensor
771 {
772 Vmath::Vcopy(nq, dtmp, 1, tmp, 1);
773 SVVLaplacianFilter(dtmp, mkey);
774 Vmath::Vadd(nq, tmp, 1, dtmp, 1, dtmp, 1);
775 }
776 v_IProductWRTDerivBase(k1, dtmp, outarray);
777 }
778}
779
781 const Array<OneD, const NekDouble> &inarray,
782 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
783{
784 const int dim = GetCoordim();
785
786 int i, j;
787
790
791 if ((mkey.GetNVarCoeff() == 0 &&
794 {
795 // just call diagonal matrix form of laplcian operator
796 for (i = 0; i < dim; ++i)
797 {
798 LaplacianMatrixOp(i, i, inarray, store, mkey);
799 Vmath::Vadd(m_ncoeffs, store, 1, store2, 1, store2, 1);
800 }
801 }
802 else
803 {
804 const MatrixType mtype[3][3] = {
809
810 for (i = 0; i < dim; i++)
811 {
812 for (j = 0; j < dim; j++)
813 {
815 mkey, mtype[i][j]);
816 LaplacianMatrixOp(i, j, inarray, store, *mkeyij);
817 Vmath::Vadd(m_ncoeffs, store, 1, store2, 1, store2, 1);
818 }
819 }
820 }
821
822 Vmath::Vcopy(m_ncoeffs, store2.get(), 1, outarray.get(), 1);
823}
824
826 const int k1, const Array<OneD, const NekDouble> &inarray,
827 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
828{
830 int nq = GetTotPoints();
831
832 v_BwdTrans(inarray, tmp);
833 v_PhysDeriv(k1, tmp, tmp);
834
836 if (mkey.HasVarCoeff(keys[k1]))
837 {
838 Vmath::Vmul(nq, &(mkey.GetVarCoeff(keys[k1]))[0], 1, &tmp[0], 1,
839 &tmp[0], 1);
840 }
841
842 v_IProductWRTBase(tmp, outarray);
843}
844
846 const Array<OneD, const NekDouble> &inarray,
847 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
848{
849 int nq = GetTotPoints();
850
851 Array<OneD, NekDouble> tmp(nq), Dtmp(nq);
852 Array<OneD, NekDouble> Mtmp(nq), Mout(m_ncoeffs);
853
854 v_BwdTrans(inarray, tmp);
856
857 v_IProductWRTBase(Dtmp, outarray);
858
859 // Compte M_{div tv}
860 Vmath::Vmul(nq, &(mkey.GetVarCoeff(eVarCoeffMFDiv))[0], 1, &tmp[0], 1,
861 &Mtmp[0], 1);
862
863 v_IProductWRTBase(Mtmp, Mout);
864
865 // Add D_tv + M_{div tv}
866 Vmath::Vadd(m_ncoeffs, &Mout[0], 1, &outarray[0], 1, &outarray[0], 1);
867}
868
870 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
871 [[maybe_unused]] Array<OneD, NekDouble> &outarray,
872 [[maybe_unused]] const StdMatrixKey &mkey)
873{
874}
875
877 const Array<OneD, const NekDouble> &inarray,
878 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
879{
880 int i, ndir = 0;
881
882 VarCoeffType varcoefftypes[] = {eVarCoeffVelX, eVarCoeffVelY,
884 // Count advection velocities
885 for (auto &x : varcoefftypes)
886 {
887 if (mkey.HasVarCoeff(x))
888 {
889 ndir++;
890 }
891 }
892
893 ASSERTL0(ndir, "Must define at least one advection velocity");
894 ASSERTL1(ndir <= GetCoordim(),
895 "Number of constants is larger than coordinate dimensions");
896
897 int totpts = GetTotPoints();
898 Array<OneD, NekDouble> tmp(3 * totpts);
899 Array<OneD, NekDouble> tmp_deriv = tmp + totpts;
900 Array<OneD, NekDouble> tmp_adv = tmp_deriv + totpts;
901
902 v_BwdTrans(inarray, tmp); // transform to PhysSpace
903
904 // Evaluate advection (u dx + v dy + w dz)
905 Vmath::Zero(totpts, tmp_adv, 1);
906 for (i = 0; i < ndir; ++i)
907 {
908 v_PhysDeriv(i, tmp, tmp_deriv);
909 Vmath::Vvtvp(totpts, mkey.GetVarCoeff(varcoefftypes[i]), 1, tmp_deriv,
910 1, tmp_adv, 1, tmp_adv, 1);
911 }
912
913 v_IProductWRTBase(tmp_adv, outarray);
914}
915
917 const Array<OneD, const NekDouble> &inarray,
918 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey,
919 bool addDiffusionTerm)
920{
921 int i, ndir = 0;
922
923 VarCoeffType varcoefftypes[] = {eVarCoeffVelX, eVarCoeffVelY,
925 // Count advection velocities
926 for (auto &x : varcoefftypes)
927 {
928 if (mkey.HasVarCoeff(x))
929 {
930 ndir++;
931 }
932 }
933
934 ASSERTL0(ndir, "Must define at least one advection velocity");
935 ASSERTL1(ndir <= GetCoordim(),
936 "Number of constants is larger than coordinate dimensions");
937
939 int totpts = GetTotPoints();
940 Array<OneD, NekDouble> tmp(3 * totpts);
941 Array<OneD, NekDouble> tmp_deriv = tmp + totpts;
942 Array<OneD, NekDouble> tmp_adv = tmp_deriv + totpts;
943
944 v_BwdTrans(inarray, tmp); // transform this mode \phi_i into PhysSpace
945
946 // calculate advection u dx + v dy + ..
947 Vmath::Zero(totpts, tmp_adv, 1);
948 for (i = 0; i < ndir; ++i)
949 {
950 v_PhysDeriv(i, tmp, tmp_deriv);
951 Vmath::Vvtvp(totpts, mkey.GetVarCoeff(varcoefftypes[i]), 1, tmp_deriv,
952 1, tmp_adv, 1, tmp_adv, 1);
953 }
954
955 // Add reaction term if lambda != 0.0
956 if (lambda)
957 {
958 // Add mass varcoeff
959 if (mkey.HasVarCoeff(eVarCoeffMass))
960 {
961 Vmath::Vmul(totpts, mkey.GetVarCoeff(eVarCoeffMass), 1, tmp, 1, tmp,
962 1);
963 }
964
965 Vmath::Svtvp(totpts, -lambda, tmp, 1, tmp_adv, 1, tmp_adv, 1);
966 }
967
968 // Create mass matrix = Advection - Reaction
969 v_IProductWRTBase(tmp_adv,
970 outarray); // Create mass matrix of Advection - Reaction
971
972 // Add Laplacian matrix
973 if (addDiffusionTerm)
974 {
976 StdMatrixKey mkeylap(eLaplacian, DetShapeType(), *this,
977 mkey.GetConstFactors(), mkey.GetVarCoeffs(),
978 mkey.GetNodalPointsType());
979 LaplacianMatrixOp(inarray, lap, mkeylap);
980
981 Vmath::Vadd(m_ncoeffs, lap, 1, outarray, 1, outarray,
982 1); // += Laplacian
983 }
984}
985
987 const Array<OneD, const NekDouble> &inarray,
988 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
989{
992 StdMatrixKey mkeymass(eMass, DetShapeType(), *this);
993 StdMatrixKey mkeylap(eLaplacian, DetShapeType(), *this,
994 mkey.GetConstFactors(), mkey.GetVarCoeffs(),
995 mkey.GetNodalPointsType());
996
997 MassMatrixOp(inarray, tmp, mkeymass);
998 LaplacianMatrixOp(inarray, outarray, mkeylap);
999
1000 Blas::Daxpy(m_ncoeffs, lambda, tmp, 1, outarray, 1);
1001}
1002
1003// VIRTUAL INLINE FUNCTIONS FROM HEADER FILE
1005 const Array<OneD, const NekDouble> &Lcoord,
1006 const Array<OneD, const NekDouble> &physvals)
1007{
1008 return v_StdPhysEvaluate(Lcoord, physvals);
1009}
1010
1012 [[maybe_unused]] const std::vector<unsigned int> &nummodes,
1013 [[maybe_unused]] int &modes_offset)
1014{
1015 NEKERROR(ErrorUtil::efatal, "This function is not defined for this class");
1016 return 0;
1017}
1018
1020 [[maybe_unused]] const Array<OneD, const NekDouble> &Fx,
1021 [[maybe_unused]] Array<OneD, NekDouble> &outarray)
1022{
1023 NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
1024}
1025
1027 [[maybe_unused]] const Array<OneD, const NekDouble> &Fx,
1028 [[maybe_unused]] const Array<OneD, const NekDouble> &Fy,
1029 [[maybe_unused]] Array<OneD, NekDouble> &outarray)
1030{
1031 NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
1032}
1033
1035 [[maybe_unused]] const Array<OneD, const NekDouble> &Fx,
1036 [[maybe_unused]] const Array<OneD, const NekDouble> &Fy,
1037 [[maybe_unused]] const Array<OneD, const NekDouble> &Fz,
1038 [[maybe_unused]] Array<OneD, NekDouble> &outarray)
1039{
1040 NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
1041}
1042
1044 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &Fvec,
1045 [[maybe_unused]] Array<OneD, NekDouble> &outarray)
1046{
1047 NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
1048}
1049
1051 [[maybe_unused]] const LocalRegions::MatrixKey &mkey)
1052{
1053 NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
1055}
1056
1058 [[maybe_unused]] const LocalRegions::MatrixKey &mkey)
1059{
1060 NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
1061}
1062
1064 [[maybe_unused]] StdRegions::Orientation dir,
1065 [[maybe_unused]] Array<OneD, const NekDouble> &inarray,
1066 [[maybe_unused]] Array<OneD, NekDouble> &outarray)
1067{
1068 NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1069}
1070
1072 [[maybe_unused]] Array<OneD, NekDouble> &coeffs,
1073 [[maybe_unused]] StdRegions::Orientation dir)
1074{
1075 NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1076}
1077
1079 [[maybe_unused]] const Array<OneD, const NekDouble> &Lcoord,
1080 [[maybe_unused]] const Array<OneD, const NekDouble> &physvals)
1081
1082{
1083 NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1084 return 0;
1085}
1086
1088 [[maybe_unused]] const Array<OneD, const NekDouble> &xi,
1089 [[maybe_unused]] Array<OneD, NekDouble> &eta)
1090{
1091 NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1092}
1093
1095 [[maybe_unused]] const Array<OneD, const NekDouble> &eta,
1096 [[maybe_unused]] Array<OneD, NekDouble> &xi)
1097{
1098 NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1099}
1100
1102 [[maybe_unused]] const int i, [[maybe_unused]] const int k) const
1103{
1104 ASSERTL0(false, "This function is not valid or not defined");
1106}
1107
1109 [[maybe_unused]] const int i, [[maybe_unused]] const int j) const
1110{
1111 ASSERTL0(false, "This function is not valid or not defined");
1113}
1114
1116{
1117 ASSERTL0(false, "This function is not valid or not defined");
1118
1120}
1121
1122std::shared_ptr<StdExpansion> StdExpansion::v_GetStdExp(void) const
1123{
1124 ASSERTL0(false, "This method is not defined for this expansion");
1125 StdExpansionSharedPtr returnval;
1126 return returnval;
1127}
1128
1129std::shared_ptr<StdExpansion> StdExpansion::v_GetLinStdExp(void) const
1130{
1131 ASSERTL0(false, "This method is not defined for this expansion");
1132 StdExpansionSharedPtr returnval;
1133 return returnval;
1134}
1135
1137{
1138 ASSERTL0(false, "This function has not been defined for this expansion");
1139 return false;
1140}
1141
1143{
1144 return false;
1145}
1146
1148 [[maybe_unused]] const int dir,
1149 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
1150 [[maybe_unused]] Array<OneD, NekDouble> &outarray)
1151{
1152 NEKERROR(ErrorUtil::efatal, "This method has not been defined");
1153}
1154
1155/**
1156 *
1157 */
1159 [[maybe_unused]] const Array<OneD, const NekDouble> &direction,
1160 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
1161 [[maybe_unused]] Array<OneD, NekDouble> &outarray)
1162{
1163 NEKERROR(ErrorUtil::efatal, "This method has not been defined");
1164}
1165
1166/**
1167 *
1168 */
1170 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
1171 [[maybe_unused]] Array<OneD, NekDouble> &outarray)
1172{
1173 NEKERROR(ErrorUtil::efatal, "This method has not been defined");
1174}
1175
1176/**
1177 * @brief Integrates the specified function over the domain.
1178 * @see StdRegions#StdExpansion#Integral.
1179 */
1181 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray)
1182{
1183 NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1184 "local expansions");
1185 return 0;
1186}
1187
1188/**
1189 * @brief Calculate the derivative of the physical points
1190 * @see StdRegions#StdExpansion#PhysDeriv
1191 */
1193 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
1194 [[maybe_unused]] Array<OneD, NekDouble> &out_d1,
1195 [[maybe_unused]] Array<OneD, NekDouble> &out_d2,
1196 [[maybe_unused]] Array<OneD, NekDouble> &out_d3)
1197{
1198 NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1199 "local expansions");
1200}
1201
1203 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
1204 [[maybe_unused]] Array<OneD, NekDouble> &out_ds)
1205{
1206 NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1207 "local expansions");
1208}
1210 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
1211 [[maybe_unused]] Array<OneD, NekDouble> &out_dn)
1212{
1213 NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1214 "local expansions");
1215}
1216
1217/**
1218 * @brief Calculate the derivative of the physical points in a
1219 * given direction
1220 * @see StdRegions#StdExpansion#PhysDeriv
1221 */
1223 [[maybe_unused]] const int dir,
1224 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
1225 [[maybe_unused]] Array<OneD, NekDouble> &out_d0)
1226
1227{
1228 NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1229 "specific element types");
1230}
1231
1232/**
1233 * @brief Physical derivative along a direction vector.
1234 * @see StdRegions#StdExpansion#PhysDirectionalDeriv
1235 */
1237 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
1238 [[maybe_unused]] const Array<OneD, const NekDouble> &direction,
1239 [[maybe_unused]] Array<OneD, NekDouble> &outarray)
1240{
1241 NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1242 "specific element types");
1243}
1244
1246 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
1247 [[maybe_unused]] Array<OneD, NekDouble> &out_d1,
1248 [[maybe_unused]] Array<OneD, NekDouble> &out_d2,
1249 [[maybe_unused]] Array<OneD, NekDouble> &out_d3)
1250{
1251 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1252}
1253
1255 [[maybe_unused]] const int dir,
1256 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
1257 [[maybe_unused]] Array<OneD, NekDouble> &outarray)
1258{
1259 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1260}
1261
1263 [[maybe_unused]] const Array<OneD, const NekDouble> &coords,
1264 [[maybe_unused]] const Array<OneD, const NekDouble> &physvals)
1265{
1266 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1267 return 0;
1268}
1269
1271 [[maybe_unused]] const Array<OneD, DNekMatSharedPtr> &I,
1272 [[maybe_unused]] const Array<OneD, const NekDouble> &physvals)
1273{
1274 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1275 return 0;
1276}
1277
1279 [[maybe_unused]] const Array<OneD, const NekDouble> &coords,
1280 [[maybe_unused]] int mode)
1281{
1282 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1283 return 0;
1284}
1285
1287 [[maybe_unused]] const Array<OneD, NekDouble> &coord,
1288 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
1289 [[maybe_unused]] std::array<NekDouble, 3> &firstOrderDerivs)
1290{
1292 "PhysEvaluate first order derivative method does not exist"
1293 " for this shape type: " +
1294 static_cast<std::string>(
1296 return 0;
1297}
1298
1300 [[maybe_unused]] const Array<OneD, NekDouble> &coord,
1301 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
1302 [[maybe_unused]] std::array<NekDouble, 3> &firstOrderDerivs,
1303 [[maybe_unused]] std::array<NekDouble, 6> &secondOrderDerivs)
1304{
1306 "PhysEvaluate second order derivative method does not exist"
1307 " for this shape type: " +
1308 static_cast<std::string>(
1310 return 0;
1311}
1312
1313void StdExpansion::v_FillMode([[maybe_unused]] const int mode,
1314 [[maybe_unused]] Array<OneD, NekDouble> &outarray)
1315{
1316 NEKERROR(ErrorUtil::efatal, "This function has not "
1317 "been defined for this shape");
1318}
1319
1321 [[maybe_unused]] const StdMatrixKey &mkey)
1322{
1323 NEKERROR(ErrorUtil::efatal, "This function has not "
1324 "been defined for this element");
1325 DNekMatSharedPtr returnval;
1326 return returnval;
1327}
1328
1330 [[maybe_unused]] const StdMatrixKey &mkey)
1331{
1332 NEKERROR(ErrorUtil::efatal, "This function has not "
1333 "been defined for this element");
1334 DNekMatSharedPtr returnval;
1335 return returnval;
1336}
1337
1339 [[maybe_unused]] Array<OneD, NekDouble> &coords_0,
1340 [[maybe_unused]] Array<OneD, NekDouble> &coords_1,
1341 [[maybe_unused]] Array<OneD, NekDouble> &coords_2)
1342{
1343 NEKERROR(ErrorUtil::efatal, "Write coordinate definition method");
1344}
1345
1347 [[maybe_unused]] const Array<OneD, const NekDouble> &Lcoord,
1348 [[maybe_unused]] Array<OneD, NekDouble> &coord)
1349{
1350 NEKERROR(ErrorUtil::efatal, "Write coordinate definition method");
1351}
1352
1354{
1355 return GetShapeDimension();
1356}
1357
1359 [[maybe_unused]] Array<OneD, unsigned int> &outarray)
1360{
1361 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1362}
1363
1365 [[maybe_unused]] Array<OneD, unsigned int> &outarray)
1366{
1367 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1368}
1369
1370int StdExpansion::v_GetVertexMap([[maybe_unused]] const int localVertexId,
1371 [[maybe_unused]] bool useCoeffPacking)
1372{
1373 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1374 return 0;
1375}
1376
1378 [[maybe_unused]] const int tid,
1379 [[maybe_unused]] Array<OneD, unsigned int> &maparray,
1380 [[maybe_unused]] Array<OneD, int> &signarray,
1381 [[maybe_unused]] Orientation traceOrient, [[maybe_unused]] int P,
1382 [[maybe_unused]] int Q)
1383{
1384 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1385}
1386
1388 [[maybe_unused]] const unsigned int traceid,
1389 [[maybe_unused]] Array<OneD, unsigned int> &maparray)
1390{
1391 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1392}
1393
1395 [[maybe_unused]] const unsigned int tid,
1396 [[maybe_unused]] Array<OneD, unsigned int> &maparray,
1397 [[maybe_unused]] Array<OneD, int> &signarray,
1398 [[maybe_unused]] Orientation traceOrient, [[maybe_unused]] int P,
1399 [[maybe_unused]] int Q)
1400{
1401 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1402}
1403
1405 [[maybe_unused]] const int tid,
1406 [[maybe_unused]] Array<OneD, unsigned int> &maparray,
1407 [[maybe_unused]] Array<OneD, int> &signarray,
1408 [[maybe_unused]] const Orientation traceOrient)
1409{
1410 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1411}
1412
1413void StdExpansion::v_GetTraceNumModes([[maybe_unused]] const int tid,
1414 [[maybe_unused]] int &numModes0,
1415 [[maybe_unused]] int &numModes1,
1416 [[maybe_unused]] Orientation traceOrient)
1417{
1418 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1419}
1420
1422 [[maybe_unused]] const int vertex,
1423 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
1424 [[maybe_unused]] NekDouble &outarray)
1425{
1426 NEKERROR(ErrorUtil::efatal, "Method does not exist for "
1427 "this shape or library");
1428}
1429
1431 const Array<OneD, const NekDouble> &inarray,
1432 Array<OneD, NekDouble> &outarray)
1433{
1434 v_MultiplyByStdQuadratureMetric(inarray, outarray);
1435}
1436
1438 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
1439 [[maybe_unused]] Array<OneD, NekDouble> &outarray)
1440{
1442 "Method does not exist for this shape or library");
1443}
1444
1446 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
1447 [[maybe_unused]] Array<OneD, NekDouble> &outarray)
1448{
1449 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1450}
1451
1453 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
1454 [[maybe_unused]] Array<OneD, NekDouble> &outarray,
1455 [[maybe_unused]] bool multiplybyweights)
1456{
1457 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1458}
1459
1460/**
1461 *
1462 */
1464 [[maybe_unused]] const Array<OneD, const NekDouble> &direction,
1465 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
1466 [[maybe_unused]] Array<OneD, NekDouble> &outarray)
1467{
1468 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1469}
1470
1472 [[maybe_unused]] const int dir,
1473 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
1474 [[maybe_unused]] Array<OneD, NekDouble> &outarray)
1475{
1476 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1477}
1478
1480 Array<OneD, NekDouble> &outarray,
1481 const StdMatrixKey &mkey)
1482{
1483 // If this function is not reimplemented on shape level, the function
1484 // below will be called
1485 MassMatrixOp_MatFree(inarray, outarray, mkey);
1486}
1487
1489 const Array<OneD, const NekDouble> &inarray,
1490 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1491{
1492 // If this function is not reimplemented on shape level, the function
1493 // below will be called
1494 LaplacianMatrixOp_MatFree(inarray, outarray, mkey);
1495}
1496
1498 [[maybe_unused]] Array<OneD, NekDouble> &array,
1499 [[maybe_unused]] const StdMatrixKey &mkey)
1500{
1501 ASSERTL0(false, "This function is not defined in StdExpansion.");
1502}
1503
1505 [[maybe_unused]] Array<OneD, NekDouble> &array,
1506 [[maybe_unused]] const NekDouble alpha,
1507 [[maybe_unused]] const NekDouble exponent,
1508 [[maybe_unused]] const NekDouble cutoff)
1509{
1510 ASSERTL0(false, "This function is not defined in StdExpansion.");
1511}
1512
1514 [[maybe_unused]] int numMin,
1515 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
1516 [[maybe_unused]] Array<OneD, NekDouble> &outarray)
1517{
1518 ASSERTL0(false, "This function is not defined in StdExpansion.");
1519}
1520
1522 const int k1, const int k2, const Array<OneD, const NekDouble> &inarray,
1523 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1524{
1525 // If this function is not reimplemented on shape level, the function
1526 // below will be called
1527 LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
1528}
1529
1531 const int i, const Array<OneD, const NekDouble> &inarray,
1532 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1533{
1534 // If this function is not reimplemented on shape level, the function
1535 // below will be called
1536 WeakDerivMatrixOp_MatFree(i, inarray, outarray, mkey);
1537}
1538
1540 const Array<OneD, const NekDouble> &inarray,
1541 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1542{
1543 // If this function is not reimplemented on shape level, the function
1544 // below will be called
1545 WeakDirectionalDerivMatrixOp_MatFree(inarray, outarray, mkey);
1546}
1547
1549 const Array<OneD, const NekDouble> &inarray,
1550 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1551{
1552 // If this function is not reimplemented on shape level, the function
1553 // below will be called
1554 MassLevelCurvatureMatrixOp_MatFree(inarray, outarray, mkey);
1555}
1556
1558 const Array<OneD, const NekDouble> &inarray,
1559 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1560{
1561 // If this function is not reimplemented on shape level, the function
1562 // below will be called
1563 LinearAdvectionMatrixOp_MatFree(inarray, outarray, mkey);
1564}
1565
1567 const Array<OneD, const NekDouble> &inarray,
1568 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey,
1569 bool addDiffusionTerm)
1570{
1571 // If this function is not reimplemented on shape level, the function
1572 // below will be called
1574 addDiffusionTerm);
1575}
1576
1578 const Array<OneD, const NekDouble> &inarray,
1579 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1580{
1581 // If this function is not reimplemented on shape level, the function
1582 // below will be called
1583 HelmholtzMatrixOp_MatFree(inarray, outarray, mkey);
1584}
1585
1587 const Array<OneD, const NekDouble> &inarray,
1588 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1589{
1590 // If this function is not reimplemented on shape level, the function
1591 // below will be called
1592 LaplacianMatrixOp_MatFree_GenericImpl(inarray, outarray, mkey);
1593}
1594
1596 [[maybe_unused]] const Array<OneD, const NekDouble> &inarray,
1597 [[maybe_unused]] Array<OneD, NekDouble> &outarray,
1598 [[maybe_unused]] Array<OneD, NekDouble> &wsp)
1599{
1600 ASSERTL0(false, "Not implemented.");
1601}
1602
1604 const Array<OneD, const NekDouble> &inarray,
1605 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1606{
1607 // If this function is not reimplemented on shape level, the function
1608 // below will be called
1609 HelmholtzMatrixOp_MatFree_GenericImpl(inarray, outarray, mkey);
1610}
1611
1613 [[maybe_unused]] const DNekScalMatSharedPtr &m_transformationmatrix)
1614{
1615 NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
1616 return NullDNekMatSharedPtr;
1617}
1618
1620 const Array<OneD, const NekDouble> &inarray,
1621 Array<OneD, NekDouble> &outarray, int npset)
1622{
1624 DNekMatSharedPtr intmat;
1625
1626 int nqtot = GetTotPoints();
1627 int np = 0;
1628 if (npset == -1) // use values from basis num points()
1629 {
1630 int nqbase;
1631 for (int i = 0; i < m_base.size(); ++i)
1632 {
1633 nqbase = m_base[i]->GetNumPoints();
1634 np = std::max(np, nqbase);
1635 }
1636
1637 StdMatrixKey Ikey(ePhysInterpToEquiSpaced, shape, *this);
1638 intmat = GetStdMatrix(Ikey);
1639 }
1640 else
1641 {
1642 np = npset;
1643
1644 ConstFactorMap cmap;
1645 cmap[eFactorConst] = np;
1646 StdMatrixKey Ikey(ePhysInterpToEquiSpaced, shape, *this, cmap);
1647 intmat = GetStdMatrix(Ikey);
1648 }
1649
1650 NekVector<NekDouble> in(nqtot, inarray, eWrapper);
1652 LibUtilities::GetNumberOfCoefficients(shape, np, np, np), outarray,
1653 eWrapper);
1654 out = (*intmat) * in;
1655}
1656
1658 [[maybe_unused]] Array<OneD, int> &conn, [[maybe_unused]] bool standard)
1659{
1661 "GetSimplexEquiSpacedConnectivity not"
1662 " implemented for " +
1663 static_cast<std::string>(
1665}
1666
1668 const Array<OneD, const NekDouble> &inarray,
1669 Array<OneD, NekDouble> &outarray)
1670{
1672
1673 // inarray has to be consistent with NumModes definition
1674 // There is also a check in GetStdMatrix to see if all
1675 // modes are of the same size
1676 ConstFactorMap cmap;
1677
1678 cmap[eFactorConst] = m_base[0]->GetNumModes();
1679 StdMatrixKey Ikey(eEquiSpacedToCoeffs, shape, *this, cmap);
1680 DNekMatSharedPtr intmat = GetStdMatrix(Ikey);
1681
1683 NekVector<NekDouble> out(m_ncoeffs, outarray, eWrapper);
1684 out = (*intmat) * in;
1685}
1686
1687} // namespace Nektar::StdRegions
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#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
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:265
Describes the specification for a Basis.
Definition: Basis.h:45
Defines a specification for a set of points.
Definition: Points.h:50
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
The base class for all shapes.
Definition: StdExpansion.h:65
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:669
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false)
virtual void v_LocCollapsedToLocCoord(const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi)
virtual ~StdExpansion()
Destructor.
virtual void v_PhysDeriv_n(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
virtual void v_PhysDirectionalDeriv(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
Physical derivative along a direction vector.
StdExpansion()
Default Constructor.
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:124
virtual void v_GetTraceToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:134
virtual void v_GetCoord(const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
void GeneralMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void FillMode(const int mode, Array< OneD, NekDouble > &outarray)
This function fills the array outarray with the mode-th mode of the expansion.
Definition: StdExpansion.h:491
virtual NekDouble v_StdPhysEvaluate(const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
void PhysInterpToSimplexEquiSpaced(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int npset=-1)
This function performs an interpolation from the physical space points provided at input into an arra...
void WeakDirectionalDerivMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:801
void WeakDerivMatrixOp_MatFree(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
DNekBlkMatSharedPtr CreateStdStaticCondMatrix(const StdMatrixKey &mkey)
Create the static condensation of a matrix when using a boundary interior decomposition.
void LaplacianMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_IProductWRTDirectionalDerivBase_SumFac(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void MassLevelCurvatureMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:808
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2)
void LinearAdvectionDiffusionReactionMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:603
virtual void v_GetElmtTraceToTraceMap(const unsigned int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
virtual void v_NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
virtual void v_GetSimplexEquiSpacedConnectivity(Array< OneD, int > &conn, bool standard=true)
virtual void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
virtual std::shared_ptr< StdExpansion > v_GetLinStdExp(void) const
void WeakDerivMatrixOp(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:793
virtual void v_GetInteriorMap(Array< OneD, unsigned int > &outarray)
virtual void v_SetCoeffsToOrientation(StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
NekDouble Linf(const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
Function to evaluate the discrete error where is given by the array sol.
void EquiSpacedToCoeffs(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs a projection/interpolation from the equispaced points sometimes used in post-p...
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:752
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_WeakDirectionalDerivMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual NekDouble v_PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
virtual DNekMatSharedPtr v_BuildInverseTransformationMatrix(const DNekScalMatSharedPtr &m_transformationmatrix)
NekDouble L2(const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
Function to evaluate the discrete error, where is given by the array sol.
virtual void v_GetVertexPhysVals(const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray)
void IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void LinearAdvectionMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:815
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
virtual void v_MassLevelCurvatureMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_GetTraceCoeffMap(const unsigned int traceid, Array< OneD, unsigned int > &maparray)
virtual bool v_IsBoundaryInteriorExpansion() const
virtual void v_MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:759
virtual void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray)
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_ExponentialFilter(Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
void MassLevelCurvatureMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:367
virtual void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:674
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey)
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_GetTraceNumModes(const int fid, int &numModes0, int &numModes1, Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2)
virtual void v_PhysDeriv_s(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
virtual void v_FwdTransBndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray)
void GenStdMatBwdDeriv(const int dir, DNekMatSharedPtr &mat)
virtual void v_MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual const LibUtilities::PointsKey v_GetNodalPointsKey() const
virtual void v_WeakDerivMatrixOp(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LinearAdvectionDiffusionReactionMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
Definition: StdExpansion.h:822
void IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
virtual void v_LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
virtual NekDouble v_PhysEvaluateBasis(const Array< OneD, const NekDouble > &coords, int mode)
virtual void v_GetTraceInteriorToElementMap(const int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eForwards)
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)=0
Calculates the inner product of a given function f with the different modes of the expansion.
virtual void v_LinearAdvectionDiffusionReactionMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
void HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:837
virtual void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_DropLocStaticCondMatrix(const LocalRegions::MatrixKey &mkey)
void WeakDirectionalDerivMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_LaplacianMatrixOp_MatFree_Kernel(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
NekDouble H1(const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
Function to evaluate the discrete error, where is given by the array sol.
virtual int v_CalcNumberOfCoefficients(const std::vector< unsigned int > &nummodes, int &modes_offset)
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
Integrates the specified function over the domain.
void SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
Definition: StdExpansion.h:773
virtual std::shared_ptr< StdExpansion > v_GetStdExp() const
void HelmholtzMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)=0
virtual LibUtilities::PointsKey v_GetTracePointsKey(const int i, const int j) const
Array< OneD, LibUtilities::BasisSharedPtr > m_base
virtual void v_LinearAdvectionMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
virtual void v_IProductWRTDirectionalDerivBase(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void LinearAdvectionMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2, Array< OneD, NekDouble > &out_d3)
Calculate the derivative of the physical points.
virtual const LibUtilities::BasisKey v_GetTraceBasisKey(const int i, const int k) const
virtual void v_StdPhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2, Array< OneD, NekDouble > &out_d3)
NekDouble StdPhysEvaluate(const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
void MassMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void GeneralMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix(const LocalRegions::MatrixKey &mkey)
LibUtilities::ShapeType GetShapeType() const
Definition: StdMatrixKey.h:88
const VarCoeffMap & GetVarCoeffs() const
Definition: StdMatrixKey.h:169
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:83
bool HasVarCoeff(const StdRegions::VarCoeffType &coeff) const
Definition: StdMatrixKey.h:174
LibUtilities::PointsType GetNodalPointsType() const
Definition: StdMatrixKey.h:93
const ConstFactorMap & GetConstFactors() const
Definition: StdMatrixKey.h:138
const Array< OneD, const NekDouble > & GetVarCoeff(const StdRegions::VarCoeffType &coeff) const
Definition: StdMatrixKey.h:148
NekDouble GetConstFactor(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:124
bool ConstFactorExists(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:133
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition: Blas.hpp:135
const char *const ShapeTypeMap[SIZE_ShapeType]
Definition: ShapeType.hpp:75
BasisManagerT & BasisManager(void)
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
int GetNumberOfCoefficients(ShapeType shape, std::vector< unsigned int > &modes, int offset=0)
Definition: ShapeType.hpp:305
@ P
Monomial polynomials .
Definition: BasisType.h:62
static const PointsKey NullPointsKey(0, eNoPointsType)
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
const char *const MatrixTypeMap[]
Definition: StdRegions.hpp:139
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:430
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:431
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:376
std::shared_ptr< StdMatrixKey > StdMatrixKeySharedPtr
Definition: StdMatrixKey.h:199
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:77
static DNekScalBlkMatSharedPtr NullDNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:85
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:79
static DNekMatSharedPtr NullDNekMatSharedPtr
Definition: NekTypeDefs.hpp:83
static Array< OneD, NekDouble > NullNekDouble1DArray
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble
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.hpp:72
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
void Vabs(int n, const T *x, const int incx, T *y, const int incy)
vabs: y = |x|
Definition: Vmath.hpp:352
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.hpp:366
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
T Vamax(int n, const T *x, const int incx)
Return the maximum absolute element in x called vamax to avoid conflict with max.
Definition: Vmath.hpp:685
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
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.hpp:220
STL namespace.
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294