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 <boost/core/ignore_unused.hpp>
37
38#include <LibUtilities/Foundations/ManagerAccess.h> // for BasisManager, etc
40
41namespace Nektar
42{
43namespace StdRegions
44{
45/** \brief Default constructor */
47{
48}
49
50StdExpansion::StdExpansion(const int numcoeffs, const int numbases,
51 const LibUtilities::BasisKey &Ba,
52 const LibUtilities::BasisKey &Bb,
53 const LibUtilities::BasisKey &Bc)
54 : m_base(numbases), m_elmt_id(0), m_ncoeffs(numcoeffs),
55 m_stdMatrixManager(std::bind(&StdExpansion::CreateStdMatrix, this,
56 std::placeholders::_1),
57 std::string("StdExpansionStdMatrix")),
58 m_stdStaticCondMatrixManager(
59 std::bind(&StdExpansion::CreateStdStaticCondMatrix, this,
60 std::placeholders::_1),
61 std::string("StdExpansionStdStaticCondMatrix"))
62{
63 switch (m_base.size())
64 {
65 case 3:
67 "NULL Basis attempting to be used.");
69 /* Falls through. */
70 case 2:
72 "NULL Basis attempting to be used.");
74 /* Falls through. */
75 case 1:
77 "NULL Basis attempting to be used.");
79 break;
80 default:
81 break;
82 // ASSERTL0(false, "numbases incorrectly specified");
83 };
84
85} // end constructor
86
88 : std::enable_shared_from_this<StdExpansion>(T), m_base(T.m_base),
89 m_elmt_id(T.m_elmt_id), m_ncoeffs(T.m_ncoeffs),
90 m_stdMatrixManager(T.m_stdMatrixManager),
91 m_stdStaticCondMatrixManager(T.m_stdStaticCondMatrixManager)
92{
93}
94
95// Destructor
97{
98}
99
102{
103 NekDouble val;
104 int ntot = GetTotPoints();
105 Array<OneD, NekDouble> wsp(ntot);
106
107 if (sol == NullNekDouble1DArray)
108 {
109 Vmath::Vabs(ntot, phys, 1, wsp, 1);
110 }
111 else
112 {
113 Vmath::Vsub(ntot, sol, 1, phys, 1, wsp, 1);
114 Vmath::Vabs(ntot, wsp, 1, wsp, 1);
115 }
116
117 val = Vmath::Vamax(ntot, wsp, 1);
118 return val;
119}
120
123{
124 NekDouble val;
125 int ntot = GetTotPoints();
126 Array<OneD, NekDouble> wsp(ntot);
127
128 if (sol.size() == 0)
129 {
130 Vmath::Vmul(ntot, phys, 1, phys, 1, wsp, 1);
131 }
132 else
133 {
134 Vmath::Vsub(ntot, sol, 1, phys, 1, wsp, 1);
135 Vmath::Vmul(ntot, wsp, 1, wsp, 1, wsp, 1);
136 }
137
138 val = v_Integral(wsp);
139
140 return (val < 0.0) ? 0.0 : sqrt(val);
141}
142
145{
146 int i;
147 NekDouble val;
148 int ntot = GetTotPoints();
149 int coordim = v_GetCoordim();
150 Array<OneD, NekDouble> wsp(3 * ntot);
151 Array<OneD, NekDouble> wsp_deriv = wsp + ntot;
152 Array<OneD, NekDouble> sum = wsp_deriv + ntot;
153
154 if (sol == NullNekDouble1DArray)
155 {
156 Vmath::Vcopy(ntot, phys, 1, wsp, 1);
157 Vmath::Vmul(ntot, phys, 1, phys, 1, sum, 1);
158 }
159 else
160 {
161 Vmath::Vsub(ntot, sol, 1, phys, 1, wsp, 1);
162 Vmath::Vmul(ntot, wsp, 1, wsp, 1, sum, 1);
163 }
164
165 for (i = 0; i < coordim; ++i)
166 {
167 v_PhysDeriv(i, wsp, wsp_deriv);
168 Vmath::Vvtvp(ntot, wsp_deriv, 1, wsp_deriv, 1, sum, 1, sum, 1);
169 }
170
171 val = sqrt(v_Integral(sum));
172
173 return val;
174}
175
177 const StdMatrixKey &mkey)
178{
179 DNekBlkMatSharedPtr returnval;
180
181 DNekMatSharedPtr mat = GetStdMatrix(mkey);
182 int nbdry = NumBndryCoeffs(); // also checks to see if this is a boundary
183 // interior decomposed expansion
184 int nint = m_ncoeffs - nbdry;
190
191 int i, j;
192
193 Array<OneD, unsigned int> bmap(nbdry);
194 Array<OneD, unsigned int> imap(nint);
195 GetBoundaryMap(bmap);
196 GetInteriorMap(imap);
197
198 for (i = 0; i < nbdry; ++i)
199 {
200 for (j = 0; j < nbdry; ++j)
201 {
202 (*A)(i, j) = (*mat)(bmap[i], bmap[j]);
203 }
204
205 for (j = 0; j < nint; ++j)
206 {
207 (*B)(i, j) = (*mat)(bmap[i], imap[j]);
208 }
209 }
210
211 for (i = 0; i < nint; ++i)
212 {
213 for (j = 0; j < nbdry; ++j)
214 {
215 (*C)(i, j) = (*mat)(imap[i], bmap[j]);
216 }
217
218 for (j = 0; j < nint; ++j)
219 {
220 (*D)(i, j) = (*mat)(imap[i], imap[j]);
221 }
222 }
223
224 // Calculate static condensed system
225 if (nint)
226 {
227 D->Invert();
228 (*B) = (*B) * (*D);
229 (*A) = (*A) - (*B) * (*C);
230 }
231
232 // set up block matrix system
233 Array<OneD, unsigned int> exp_size(2);
234 exp_size[0] = nbdry;
235 exp_size[1] = nint;
236 returnval =
238
239 returnval->SetBlock(0, 0, A);
240 returnval->SetBlock(0, 1, B);
241 returnval->SetBlock(1, 0, C);
242 returnval->SetBlock(1, 1, D);
243
244 return returnval;
245}
246
248{
249 int i;
250 DNekMatSharedPtr returnval;
251
252 switch (mkey.GetMatrixType())
253 {
254 case eInvMass:
255 {
256 StdMatrixKey masskey(eMass, mkey.GetShapeType(), *this,
258 mkey.GetNodalPointsType());
259 DNekMatSharedPtr mmat = GetStdMatrix(masskey);
260
262 *mmat); // Populate standard mass matrix.
263 returnval->Invert();
264 }
265 break;
266 case eInvNBasisTrans:
267 {
268 StdMatrixKey tmpkey(eNBasisTrans, mkey.GetShapeType(), *this,
270 mkey.GetNodalPointsType());
271 DNekMatSharedPtr tmpmat = GetStdMatrix(tmpkey);
273 *tmpmat); // Populate matrix.
274 returnval->Invert();
275 }
276 break;
277 case eBwdMat:
278 {
279 int nq = GetTotPoints();
281 Array<OneD, NekDouble> tmpout(nq);
282
283 returnval =
285 Array<OneD, NekDouble> Bwd_data = returnval->GetPtr();
286
288 this->DetShapeType(), *this);
289 DNekMatSharedPtr MatBwdTrans = GetStdMatrix(matkey);
290 Array<OneD, NekDouble> BwdTrans_data = MatBwdTrans->GetPtr();
291
292 for (int i = 0; i < m_ncoeffs; ++i)
293 {
294 Array<OneD, NekDouble> tmpinn = BwdTrans_data + nq * i;
295 Array<OneD, NekDouble> tmpout = Bwd_data + i;
296
297 Vmath::Vcopy(nq, tmpinn, 1, tmpout, m_ncoeffs);
298 }
299 }
300 break;
301 case eBwdTrans:
302 {
303 int nq = GetTotPoints();
305 Array<OneD, NekDouble> tmpout(nq);
306
307 returnval =
309
310 for (int i = 0; i < m_ncoeffs; ++i)
311 {
312 Vmath::Zero(m_ncoeffs, tmpin, 1);
313 tmpin[i] = 1.0;
314
315 BwdTrans_SumFac(tmpin, tmpout);
316
317 Vmath::Vcopy(nq, tmpout.get(), 1,
318 returnval->GetRawPtr() + i * nq, 1);
319 }
320 }
321 break;
322 case eIProductWRTBase:
323 {
324 int nq = GetTotPoints();
325 Array<OneD, NekDouble> tmpin(nq);
327
328 returnval =
330
331 for (i = 0; i < nq; ++i)
332 {
333 Vmath::Zero(nq, tmpin, 1);
334 tmpin[i] = 1.0;
335
336 IProductWRTBase_SumFac(tmpin, tmpout);
337
338 Vmath::Vcopy(m_ncoeffs, tmpout.get(), 1,
339 returnval->GetRawPtr() + i * m_ncoeffs, 1);
340 }
341 }
342 break;
344 {
345 int nq = GetTotPoints();
346 Array<OneD, NekDouble> tmpin(nq);
348
349 returnval =
351
352 for (i = 0; i < nq; ++i)
353 {
354 Vmath::Zero(nq, tmpin, 1);
355 tmpin[i] = 1.0;
356
357 IProductWRTDerivBase_SumFac(0, tmpin, tmpout);
358
359 Vmath::Vcopy(m_ncoeffs, tmpout.get(), 1,
360 returnval->GetRawPtr() + i * m_ncoeffs, 1);
361 }
362 }
363 break;
365 {
366 int nq = GetTotPoints();
367 Array<OneD, NekDouble> tmpin(nq);
369
370 returnval =
372
373 for (i = 0; i < nq; ++i)
374 {
375 Vmath::Zero(nq, tmpin, 1);
376 tmpin[i] = 1.0;
377
378 IProductWRTDerivBase_SumFac(1, tmpin, tmpout);
379
380 Vmath::Vcopy(m_ncoeffs, tmpout.get(), 1,
381 returnval->GetRawPtr() + i * m_ncoeffs, 1);
382 }
383 }
384 break;
386 {
387 int nq = GetTotPoints();
388 Array<OneD, NekDouble> tmpin(nq);
390
391 returnval =
393
394 for (i = 0; i < nq; ++i)
395 {
396 Vmath::Zero(nq, tmpin, 1);
397 tmpin[i] = 1.0;
398
399 IProductWRTDerivBase_SumFac(2, tmpin, tmpout);
400
401 Vmath::Vcopy(m_ncoeffs, tmpout.get(), 1,
402 returnval->GetRawPtr() + i * m_ncoeffs, 1);
403 }
404 }
405 break;
406 case eDerivBase0:
407 {
408 int nq = GetTotPoints();
409 returnval =
411 GenStdMatBwdDeriv(0, returnval);
412 }
413 break;
414 case eDerivBase1:
415 {
416 int nq = GetTotPoints();
417 returnval =
419 GenStdMatBwdDeriv(1, returnval);
420 }
421 break;
422 case eDerivBase2:
423 {
424 int nq = GetTotPoints();
425 returnval =
427 GenStdMatBwdDeriv(2, returnval);
428 }
429 break;
431 {
432 // check to see if equispaced basis
433 int nummodes = m_base[0]->GetNumModes();
434 bool equispaced = true;
435 for (int i = 1; i < m_base.size(); ++i)
436 {
437 if (m_base[i]->GetNumModes() != nummodes)
438 {
439 equispaced = false;
440 }
441 }
442
443 ASSERTL0(equispaced,
444 "Currently need to have same num modes in all "
445 "directionmodes to use EquiSpacedToCoeff method");
446
447 int ntot = GetTotPoints();
448 Array<OneD, NekDouble> qmode(ntot);
450
451 returnval =
453 for (int i = 0; i < m_ncoeffs; ++i)
454 {
455 // Get mode at quadrature points
456 FillMode(i, qmode);
457
458 // interpolate to equi spaced
459 PhysInterpToSimplexEquiSpaced(qmode, emode, nummodes);
460
461 // fill matrix
462 Vmath::Vcopy(m_ncoeffs, &emode[0], 1,
463 returnval->GetRawPtr() + i * m_ncoeffs, 1);
464 }
465 // invert matrix
466 returnval->Invert();
467 }
468 break;
469 case eMass:
470 case eHelmholtz:
471 case eLaplacian:
472 case eLaplacian00:
473 case eLaplacian01:
474 case eLaplacian02:
475 case eLaplacian10:
476 case eLaplacian11:
477 case eLaplacian12:
478 case eLaplacian20:
479 case eLaplacian21:
480 case eLaplacian22:
481 case eWeakDeriv0:
482 case eWeakDeriv1:
483 case eWeakDeriv2:
486 case eLinearAdvection:
489 {
491 returnval =
493 DNekMat &Mat = *returnval;
494
495 for (i = 0; i < m_ncoeffs; ++i)
496 {
497 Vmath::Zero(m_ncoeffs, tmp, 1);
498 tmp[i] = 1.0;
499
500 GeneralMatrixOp_MatFree(tmp, tmp, mkey);
501
502 Vmath::Vcopy(m_ncoeffs, &tmp[0], 1,
503 &(Mat.GetPtr())[0] + i * m_ncoeffs, 1);
504 }
505 }
506 break;
507 default:
508 {
510 "This type of matrix, " +
511 static_cast<std::string>(
513 ", can not be created using a general approach");
514 }
515 break;
516 }
517
518 return returnval;
519}
520
522 Array<OneD, NekDouble> &outarray,
523 const StdMatrixKey &mkey)
524{
525 switch (mkey.GetMatrixType())
526 {
527 case eMass:
528 MassMatrixOp(inarray, outarray, mkey);
529 break;
530 case eWeakDeriv0:
531 WeakDerivMatrixOp(0, inarray, outarray, mkey);
532 break;
533 case eWeakDeriv1:
534 WeakDerivMatrixOp(1, inarray, outarray, mkey);
535 break;
536 case eWeakDeriv2:
537 WeakDerivMatrixOp(2, inarray, outarray, mkey);
538 break;
540 WeakDirectionalDerivMatrixOp(inarray, outarray, mkey);
541 break;
543 MassLevelCurvatureMatrixOp(inarray, outarray, mkey);
544 break;
545 case eLinearAdvection:
546 LinearAdvectionMatrixOp(inarray, outarray, mkey);
547 break;
549 LinearAdvectionDiffusionReactionMatrixOp(inarray, outarray, mkey,
550 false);
551 break;
553 LinearAdvectionDiffusionReactionMatrixOp(inarray, outarray, mkey);
554 break;
555 case eLaplacian:
556 LaplacianMatrixOp(inarray, outarray, mkey);
557 break;
558 case eLaplacian00:
559 LaplacianMatrixOp(0, 0, inarray, outarray, mkey);
560 break;
561 case eLaplacian01:
562 LaplacianMatrixOp(0, 1, inarray, outarray, mkey);
563 break;
564 case eLaplacian02:
565 LaplacianMatrixOp(0, 2, inarray, outarray, mkey);
566 break;
567 case eLaplacian10:
568 LaplacianMatrixOp(1, 0, inarray, outarray, mkey);
569 break;
570 case eLaplacian11:
571 LaplacianMatrixOp(1, 1, inarray, outarray, mkey);
572 break;
573 case eLaplacian12:
574 LaplacianMatrixOp(1, 2, inarray, outarray, mkey);
575 break;
576 case eLaplacian20:
577 LaplacianMatrixOp(2, 0, inarray, outarray, mkey);
578 break;
579 case eLaplacian21:
580 LaplacianMatrixOp(2, 1, inarray, outarray, mkey);
581 break;
582 case eLaplacian22:
583 LaplacianMatrixOp(2, 2, inarray, outarray, mkey);
584 break;
585 case eHelmholtz:
586 HelmholtzMatrixOp(inarray, outarray, mkey);
587 break;
588 default:
590 "This matrix does not have an operator");
591 break;
592 }
593}
594
596 const Array<OneD, const NekDouble> &inarray,
597 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
598{
599 switch (mkey.GetMatrixType())
600 {
601 case eMass:
602 MassMatrixOp_MatFree(inarray, outarray, mkey);
603 break;
604 case eWeakDeriv0:
605 WeakDerivMatrixOp_MatFree(0, inarray, outarray, mkey);
606 break;
607 case eWeakDeriv1:
608 WeakDerivMatrixOp_MatFree(1, inarray, outarray, mkey);
609 break;
610 case eWeakDeriv2:
611 WeakDerivMatrixOp_MatFree(2, inarray, outarray, mkey);
612 break;
614 WeakDirectionalDerivMatrixOp_MatFree(inarray, outarray, mkey);
615 break;
617 MassLevelCurvatureMatrixOp_MatFree(inarray, outarray, mkey);
618 break;
619 case eLinearAdvection:
620 LinearAdvectionMatrixOp_MatFree(inarray, outarray, mkey);
621 break;
624 mkey, false);
625 break;
628 mkey);
629 break;
630 case eLaplacian:
631 LaplacianMatrixOp_MatFree(inarray, outarray, mkey);
632 break;
633 case eLaplacian00:
634 LaplacianMatrixOp_MatFree(0, 0, inarray, outarray, mkey);
635 break;
636 case eLaplacian01:
637 LaplacianMatrixOp_MatFree(0, 1, inarray, outarray, mkey);
638 break;
639 case eLaplacian02:
640 LaplacianMatrixOp_MatFree(0, 2, inarray, outarray, mkey);
641 break;
642 case eLaplacian10:
643 LaplacianMatrixOp_MatFree(1, 0, inarray, outarray, mkey);
644 break;
645 case eLaplacian11:
646 LaplacianMatrixOp_MatFree(1, 1, inarray, outarray, mkey);
647 break;
648 case eLaplacian12:
649 LaplacianMatrixOp_MatFree(1, 2, inarray, outarray, mkey);
650 break;
651 case eLaplacian20:
652 LaplacianMatrixOp_MatFree(2, 0, inarray, outarray, mkey);
653 break;
654 case eLaplacian21:
655 LaplacianMatrixOp_MatFree(2, 1, inarray, outarray, mkey);
656 break;
657 case eLaplacian22:
658 LaplacianMatrixOp_MatFree(2, 2, inarray, outarray, mkey);
659 break;
660 case eHelmholtz:
661 HelmholtzMatrixOp_MatFree(inarray, outarray, mkey);
662 break;
663 default:
665 "This matrix does not have an operator");
666 break;
667 }
668}
669
671 const Array<OneD, const NekDouble> &inarray,
672 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
673{
674 int nq = GetTotPoints();
676
677 v_BwdTrans(inarray, tmp);
678
679 if (mkey.HasVarCoeff(eVarCoeffMass))
680 {
681 Vmath::Vmul(nq, mkey.GetVarCoeff(eVarCoeffMass), 1, tmp, 1, tmp, 1);
682 }
683
684 v_IProductWRTBase(tmp, outarray);
685}
686
688 const int k1, const int k2, const Array<OneD, const NekDouble> &inarray,
689 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
690{
691 ASSERTL1(k1 >= 0 && k1 < GetCoordim(), "invalid first argument");
692 ASSERTL1(k2 >= 0 && k2 < GetCoordim(), "invalid second argument");
693
694 int nq = GetTotPoints();
696 Array<OneD, NekDouble> dtmp(nq);
697 VarCoeffType varcoefftypes[3][3] = {
701
702 ConstFactorType constcoefftypes[3][3] = {
706
707 v_BwdTrans(inarray, tmp);
708 v_PhysDeriv(k2, tmp, dtmp);
710 {
711 if (k1 == k2)
712 {
713 // By default, k1 == k2 has \sigma = 1 (diagonal entries)
714 if (mkey.HasVarCoeff(varcoefftypes[k1][k1]))
715 {
716 Vmath::Vmul(nq, mkey.GetVarCoeff(varcoefftypes[k1][k1]), 1,
717 dtmp, 1, dtmp, 1);
718 }
719 v_IProductWRTDerivBase_SumFac(k1, dtmp, outarray);
720 }
721 else
722 {
723 // By default, k1 != k2 has \sigma = 0 (off-diagonal entries)
724 if (mkey.HasVarCoeff(varcoefftypes[k1][k2]))
725 {
726 Vmath::Vmul(nq, mkey.GetVarCoeff(varcoefftypes[k1][k2]), 1,
727 dtmp, 1, dtmp, 1);
728 v_IProductWRTDerivBase_SumFac(k1, dtmp, outarray);
729 }
730 else if (mkey.HasVarCoeff(
731 varcoefftypes[k2][k1])) // Check symmetric varcoeff
732 {
733 Vmath::Vmul(nq, mkey.GetVarCoeff(varcoefftypes[k2][k1]), 1,
734 dtmp, 1, dtmp, 1);
735 v_IProductWRTDerivBase_SumFac(k1, dtmp, outarray);
736 }
737 else
738 {
739 Vmath::Zero(GetNcoeffs(), outarray, 1);
740 }
741 }
742 }
743 else if (mkey.ConstFactorExists(eFactorCoeffD00) &&
745 {
746 if (k1 == k2)
747 {
748 // By default, k1 == k2 has \sigma = 1 (diagonal entries)
749 if (mkey.ConstFactorExists(constcoefftypes[k1][k1]))
750 {
751 Vmath::Smul(nq, mkey.GetConstFactor(constcoefftypes[k1][k1]),
752 dtmp, 1, dtmp, 1);
753 }
754 v_IProductWRTDerivBase(k1, dtmp, outarray);
755 }
756 else
757 {
758 // By default, k1 != k2 has \sigma = 0 (off-diagonal entries)
759 if (mkey.ConstFactorExists(constcoefftypes[k1][k2]))
760 {
761 Vmath::Smul(nq, mkey.GetConstFactor(constcoefftypes[k1][k2]),
762 dtmp, 1, dtmp, 1);
763 v_IProductWRTDerivBase(k1, dtmp, outarray);
764 }
765 else
766 {
767 Vmath::Zero(GetNcoeffs(), outarray, 1);
768 }
769 }
770 }
771 else
772 {
773 // Multiply by svv tensor
775 {
776 Vmath::Vcopy(nq, dtmp, 1, tmp, 1);
777 SVVLaplacianFilter(dtmp, mkey);
778 Vmath::Vadd(nq, tmp, 1, dtmp, 1, dtmp, 1);
779 }
780 v_IProductWRTDerivBase(k1, dtmp, outarray);
781 }
782}
783
785 const Array<OneD, const NekDouble> &inarray,
786 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
787{
788 const int dim = GetCoordim();
789
790 int i, j;
791
794
795 if ((mkey.GetNVarCoeff() == 0 &&
798 {
799 // just call diagonal matrix form of laplcian operator
800 for (i = 0; i < dim; ++i)
801 {
802 LaplacianMatrixOp(i, i, inarray, store, mkey);
803 Vmath::Vadd(m_ncoeffs, store, 1, store2, 1, store2, 1);
804 }
805 }
806 else
807 {
808 const MatrixType mtype[3][3] = {
813
814 for (i = 0; i < dim; i++)
815 {
816 for (j = 0; j < dim; j++)
817 {
819 mkey, mtype[i][j]);
820 LaplacianMatrixOp(i, j, inarray, store, *mkeyij);
821 Vmath::Vadd(m_ncoeffs, store, 1, store2, 1, store2, 1);
822 }
823 }
824 }
825
826 Vmath::Vcopy(m_ncoeffs, store2.get(), 1, outarray.get(), 1);
827}
828
830 const int k1, const Array<OneD, const NekDouble> &inarray,
831 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
832{
834 int nq = GetTotPoints();
835
836 v_BwdTrans(inarray, tmp);
837 v_PhysDeriv(k1, tmp, tmp);
838
840 if (mkey.HasVarCoeff(keys[k1]))
841 {
842 Vmath::Vmul(nq, &(mkey.GetVarCoeff(keys[k1]))[0], 1, &tmp[0], 1,
843 &tmp[0], 1);
844 }
845
846 v_IProductWRTBase(tmp, outarray);
847}
848
850 const Array<OneD, const NekDouble> &inarray,
851 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
852{
853 int nq = GetTotPoints();
854
855 Array<OneD, NekDouble> tmp(nq), Dtmp(nq);
856 Array<OneD, NekDouble> Mtmp(nq), Mout(m_ncoeffs);
857
858 v_BwdTrans(inarray, tmp);
860
861 v_IProductWRTBase(Dtmp, outarray);
862
863 // Compte M_{div tv}
864 Vmath::Vmul(nq, &(mkey.GetVarCoeff(eVarCoeffMFDiv))[0], 1, &tmp[0], 1,
865 &Mtmp[0], 1);
866
867 v_IProductWRTBase(Mtmp, Mout);
868
869 // Add D_tv + M_{div tv}
870 Vmath::Vadd(m_ncoeffs, &Mout[0], 1, &outarray[0], 1, &outarray[0], 1);
871}
872
874 const Array<OneD, const NekDouble> &inarray,
875 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
876{
877 boost::ignore_unused(inarray, outarray, mkey);
878}
879
881 const Array<OneD, const NekDouble> &inarray,
882 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
883{
884 int i, ndir = 0;
885
886 VarCoeffType varcoefftypes[] = {eVarCoeffVelX, eVarCoeffVelY,
888 // Count advection velocities
889 for (auto &x : varcoefftypes)
890 {
891 if (mkey.HasVarCoeff(x))
892 {
893 ndir++;
894 }
895 }
896
897 ASSERTL0(ndir, "Must define at least one advection velocity");
898 ASSERTL1(ndir <= GetCoordim(),
899 "Number of constants is larger than coordinate dimensions");
900
901 int totpts = GetTotPoints();
902 Array<OneD, NekDouble> tmp(3 * totpts);
903 Array<OneD, NekDouble> tmp_deriv = tmp + totpts;
904 Array<OneD, NekDouble> tmp_adv = tmp_deriv + totpts;
905
906 v_BwdTrans(inarray, tmp); // transform to PhysSpace
907
908 // Evaluate advection (u dx + v dy + w dz)
909 Vmath::Zero(totpts, tmp_adv, 1);
910 for (i = 0; i < ndir; ++i)
911 {
912 v_PhysDeriv(i, tmp, tmp_deriv);
913 Vmath::Vvtvp(totpts, mkey.GetVarCoeff(varcoefftypes[i]), 1, tmp_deriv,
914 1, tmp_adv, 1, tmp_adv, 1);
915 }
916
917 v_IProductWRTBase(tmp_adv, outarray);
918}
919
921 const Array<OneD, const NekDouble> &inarray,
922 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey,
923 bool addDiffusionTerm)
924{
925 int i, ndir = 0;
926
927 VarCoeffType varcoefftypes[] = {eVarCoeffVelX, eVarCoeffVelY,
929 // Count advection velocities
930 for (auto &x : varcoefftypes)
931 {
932 if (mkey.HasVarCoeff(x))
933 {
934 ndir++;
935 }
936 }
937
938 ASSERTL0(ndir, "Must define at least one advection velocity");
939 ASSERTL1(ndir <= GetCoordim(),
940 "Number of constants is larger than coordinate dimensions");
941
943 int totpts = GetTotPoints();
944 Array<OneD, NekDouble> tmp(3 * totpts);
945 Array<OneD, NekDouble> tmp_deriv = tmp + totpts;
946 Array<OneD, NekDouble> tmp_adv = tmp_deriv + totpts;
947
948 v_BwdTrans(inarray, tmp); // transform this mode \phi_i into PhysSpace
949
950 // calculate advection u dx + v dy + ..
951 Vmath::Zero(totpts, tmp_adv, 1);
952 for (i = 0; i < ndir; ++i)
953 {
954 v_PhysDeriv(i, tmp, tmp_deriv);
955 Vmath::Vvtvp(totpts, mkey.GetVarCoeff(varcoefftypes[i]), 1, tmp_deriv,
956 1, tmp_adv, 1, tmp_adv, 1);
957 }
958
959 // Add reaction term if lambda != 0.0
960 if (lambda)
961 {
962 Vmath::Svtvp(totpts, -lambda, tmp, 1, tmp_adv, 1, tmp_adv, 1);
963 }
964
965 // Create mass matrix = Advection - Reaction
966 v_IProductWRTBase(tmp_adv,
967 outarray); // Create mass matrix of Advection - Reaction
968
969 // Add Laplacian matrix
970 if (addDiffusionTerm)
971 {
973 StdMatrixKey mkeylap(eLaplacian, DetShapeType(), *this,
974 mkey.GetConstFactors(), mkey.GetVarCoeffs(),
975 mkey.GetNodalPointsType());
976 LaplacianMatrixOp(inarray, lap, mkeylap);
977
978 Vmath::Vadd(m_ncoeffs, lap, 1, outarray, 1, outarray,
979 1); // += Laplacian
980 }
981}
982
984 const Array<OneD, const NekDouble> &inarray,
985 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
986{
989 StdMatrixKey mkeymass(eMass, DetShapeType(), *this);
990 StdMatrixKey mkeylap(eLaplacian, DetShapeType(), *this,
991 mkey.GetConstFactors(), mkey.GetVarCoeffs(),
992 mkey.GetNodalPointsType());
993
994 MassMatrixOp(inarray, tmp, mkeymass);
995 LaplacianMatrixOp(inarray, outarray, mkeylap);
996
997 Blas::Daxpy(m_ncoeffs, lambda, tmp, 1, outarray, 1);
998}
999
1000// VIRTUAL INLINE FUNCTIONS FROM HEADER FILE
1002 const Array<OneD, const NekDouble> &Lcoord,
1003 const Array<OneD, const NekDouble> &physvals)
1004{
1005 return v_StdPhysEvaluate(Lcoord, physvals);
1006}
1007
1009 const std::vector<unsigned int> &nummodes, int &modes_offset)
1010{
1011 boost::ignore_unused(nummodes, modes_offset);
1012 NEKERROR(ErrorUtil::efatal, "This function is not defined for this class");
1013 return 0;
1014}
1015
1018{
1019 boost::ignore_unused(Fx, outarray);
1020 NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
1021}
1022
1026{
1027 boost::ignore_unused(Fx, Fy, outarray);
1028 NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
1029}
1030
1035{
1036 boost::ignore_unused(Fx, Fy, Fz, outarray);
1037 NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
1038}
1039
1041 const Array<OneD, const Array<OneD, NekDouble>> &Fvec,
1042 Array<OneD, NekDouble> &outarray)
1043{
1044 boost::ignore_unused(Fvec, outarray);
1045 NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
1046}
1047
1049 const LocalRegions::MatrixKey &mkey)
1050{
1051 boost::ignore_unused(mkey);
1052 NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
1054}
1055
1057 const LocalRegions::MatrixKey &mkey)
1058{
1059 boost::ignore_unused(mkey);
1060 NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
1061}
1062
1065 Array<OneD, NekDouble> &outarray)
1066{
1067 boost::ignore_unused(dir, inarray, outarray);
1068 NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1069}
1070
1073{
1074 boost::ignore_unused(coeffs, dir);
1075 NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1076}
1077
1079 const Array<OneD, const NekDouble> &Lcoord,
1080 const Array<OneD, const NekDouble> &physvals)
1081
1082{
1083 boost::ignore_unused(Lcoord, physvals);
1084 NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1085 return 0;
1086}
1087
1090{
1091 boost::ignore_unused(xi, eta);
1092 NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1093}
1094
1097{
1098 boost::ignore_unused(eta, xi);
1099 NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1100}
1101
1103 const int k) const
1104{
1105 boost::ignore_unused(i, k);
1106 ASSERTL0(false, "This function is not valid or not defined");
1108}
1109
1111 const int j) const
1112{
1113 boost::ignore_unused(i, j);
1114 ASSERTL0(false, "This function is not valid or not defined");
1116}
1117
1119{
1120 ASSERTL0(false, "This function is not valid or not defined");
1121
1123}
1124
1125std::shared_ptr<StdExpansion> StdExpansion::v_GetStdExp(void) const
1126{
1127 ASSERTL0(false, "This method is not defined for this expansion");
1128 StdExpansionSharedPtr returnval;
1129 return returnval;
1130}
1131
1132std::shared_ptr<StdExpansion> StdExpansion::v_GetLinStdExp(void) const
1133{
1134 ASSERTL0(false, "This method is not defined for this expansion");
1135 StdExpansionSharedPtr returnval;
1136 return returnval;
1137}
1138
1140{
1141 ASSERTL0(false, "This function has not been defined for this expansion");
1142 return false;
1143}
1144
1146{
1147 return false;
1148}
1149
1151 const int dir, const Array<OneD, const NekDouble> &inarray,
1152 Array<OneD, NekDouble> &outarray)
1153{
1154 boost::ignore_unused(dir, inarray, outarray);
1155 NEKERROR(ErrorUtil::efatal, "This method has not been defined");
1156}
1157
1158/**
1159 *
1160 */
1162 const Array<OneD, const NekDouble> &direction,
1163 const Array<OneD, const NekDouble> &inarray,
1164 Array<OneD, NekDouble> &outarray)
1165{
1166 boost::ignore_unused(direction, inarray, outarray);
1167 NEKERROR(ErrorUtil::efatal, "This method has not been defined");
1168}
1169
1170/**
1171 *
1172 */
1174 const Array<OneD, const NekDouble> &inarray,
1175 Array<OneD, NekDouble> &outarray)
1176{
1177 boost::ignore_unused(inarray, outarray);
1178 NEKERROR(ErrorUtil::efatal, "This method has not been defined");
1179}
1180
1181/**
1182 * @brief Integrates the specified function over the domain.
1183 * @see StdRegions#StdExpansion#Integral.
1184 */
1186{
1187 boost::ignore_unused(inarray);
1188 NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1189 "local expansions");
1190 return 0;
1191}
1192
1193/**
1194 * @brief Calculate the derivative of the physical points
1195 * @see StdRegions#StdExpansion#PhysDeriv
1196 */
1198 Array<OneD, NekDouble> &out_d1,
1199 Array<OneD, NekDouble> &out_d2,
1200 Array<OneD, NekDouble> &out_d3)
1201{
1202 boost::ignore_unused(inarray, out_d1, out_d2, out_d3);
1203 NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1204 "local expansions");
1205}
1206
1208 Array<OneD, NekDouble> &out_ds)
1209{
1210 boost::ignore_unused(inarray, out_ds);
1211 NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1212 "local expansions");
1213}
1215 Array<OneD, NekDouble> &out_dn)
1216{
1217 boost::ignore_unused(inarray, out_dn);
1218 NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1219 "local expansions");
1220}
1221
1222/**
1223 * @brief Calculate the derivative of the physical points in a
1224 * given direction
1225 * @see StdRegions#StdExpansion#PhysDeriv
1226 */
1227void StdExpansion::v_PhysDeriv(const int dir,
1228 const Array<OneD, const NekDouble> &inarray,
1229 Array<OneD, NekDouble> &out_d0)
1230
1231{
1232 boost::ignore_unused(dir, inarray, out_d0);
1233 NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1234 "specific element types");
1235}
1236
1237/**
1238 * @brief Physical derivative along a direction vector.
1239 * @see StdRegions#StdExpansion#PhysDirectionalDeriv
1240 */
1242 const Array<OneD, const NekDouble> &inarray,
1243 const Array<OneD, const NekDouble> &direction,
1244 Array<OneD, NekDouble> &outarray)
1245{
1246 boost::ignore_unused(inarray, direction, outarray);
1247 NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1248 "specific element types");
1249}
1250
1252 Array<OneD, NekDouble> &out_d1,
1253 Array<OneD, NekDouble> &out_d2,
1254 Array<OneD, NekDouble> &out_d3)
1255{
1256 boost::ignore_unused(inarray, out_d1, out_d2, out_d3);
1257 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1258}
1259
1261 const Array<OneD, const NekDouble> &inarray,
1262 Array<OneD, NekDouble> &outarray)
1263{
1264 boost::ignore_unused(dir, inarray, outarray);
1265 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1266}
1267
1269 const Array<OneD, const NekDouble> &coords,
1270 const Array<OneD, const NekDouble> &physvals)
1271{
1272 boost::ignore_unused(coords, physvals);
1273 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1274 return 0;
1275}
1276
1279 const Array<OneD, const NekDouble> &physvals)
1280{
1281 boost::ignore_unused(I, physvals);
1282 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1283 return 0;
1284}
1285
1287 const Array<OneD, const NekDouble> &coords, int mode)
1288{
1289 boost::ignore_unused(coords, mode);
1290 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1291 return 0;
1292}
1293
1295 const Array<OneD, NekDouble> &coord,
1296 const Array<OneD, const NekDouble> &inarray,
1297 std::array<NekDouble, 3> &firstOrderDerivs)
1298{
1299 boost::ignore_unused(coord, inarray, firstOrderDerivs);
1301 "PhysEvaluate first order derivative method does not exist"
1302 " for this shape type: " +
1303 static_cast<std::string>(
1305 return 0;
1306}
1307
1309 const Array<OneD, NekDouble> &coord,
1310 const Array<OneD, const NekDouble> &inarray,
1311 std::array<NekDouble, 3> &firstOrderDerivs,
1312 std::array<NekDouble, 6> &secondOrderDerivs)
1313{
1314 boost::ignore_unused(coord, inarray, firstOrderDerivs, secondOrderDerivs);
1316 "PhysEvaluate second order derivative method does not exist"
1317 " for this shape type: " +
1318 static_cast<std::string>(
1320 return 0;
1321}
1322
1324{
1325 boost::ignore_unused(mode, outarray);
1326 NEKERROR(ErrorUtil::efatal, "This function has not "
1327 "been defined for this shape");
1328}
1329
1331{
1332 boost::ignore_unused(mkey);
1333 NEKERROR(ErrorUtil::efatal, "This function has not "
1334 "been defined for this element");
1335 DNekMatSharedPtr returnval;
1336 return returnval;
1337}
1338
1340{
1341 boost::ignore_unused(mkey);
1342 NEKERROR(ErrorUtil::efatal, "This function has not "
1343 "been defined for this element");
1344 DNekMatSharedPtr returnval;
1345 return returnval;
1346}
1347
1349 Array<OneD, NekDouble> &coords_1,
1350 Array<OneD, NekDouble> &coords_2)
1351{
1352 boost::ignore_unused(coords_0, coords_1, coords_2);
1353 NEKERROR(ErrorUtil::efatal, "Write coordinate definition method");
1354}
1355
1358{
1359 boost::ignore_unused(Lcoord, coord);
1360 NEKERROR(ErrorUtil::efatal, "Write coordinate definition method");
1361}
1362
1364{
1365 return GetShapeDimension();
1366}
1367
1369{
1370 boost::ignore_unused(outarray);
1371 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1372}
1373
1375{
1376 boost::ignore_unused(outarray);
1377 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1378}
1379
1380int StdExpansion::v_GetVertexMap(const int localVertexId, bool useCoeffPacking)
1381{
1382 boost::ignore_unused(localVertexId, useCoeffPacking);
1383 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1384 return 0;
1385}
1386
1388 Array<OneD, unsigned int> &maparray,
1389 Array<OneD, int> &signarray,
1390 Orientation traceOrient, int P, int Q)
1391{
1392 boost::ignore_unused(tid, maparray, signarray, traceOrient, P, Q);
1393 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1394}
1395
1396void StdExpansion::v_GetTraceCoeffMap(const unsigned int traceid,
1397 Array<OneD, unsigned int> &maparray)
1398{
1399 boost::ignore_unused(traceid, maparray);
1400 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1401}
1402
1403void StdExpansion::v_GetElmtTraceToTraceMap(const unsigned int tid,
1404 Array<OneD, unsigned int> &maparray,
1405 Array<OneD, int> &signarray,
1406 Orientation traceOrient, int P,
1407 int Q)
1408{
1409 boost::ignore_unused(tid, maparray, signarray, traceOrient, P, Q);
1410 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1411}
1412
1414 const int tid, Array<OneD, unsigned int> &maparray,
1415 Array<OneD, int> &signarray, const Orientation traceOrient)
1416{
1417 boost::ignore_unused(tid, maparray, signarray, traceOrient);
1418 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1419}
1420
1421void StdExpansion::v_GetTraceNumModes(const int tid, int &numModes0,
1422 int &numModes1, Orientation traceOrient)
1423{
1424 boost::ignore_unused(tid, traceOrient, numModes0, numModes1);
1425 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1426}
1427
1429 const int vertex, const Array<OneD, const NekDouble> &inarray,
1430 NekDouble &outarray)
1431{
1432 boost::ignore_unused(vertex, inarray, outarray);
1433 NEKERROR(ErrorUtil::efatal, "Method does not exist for "
1434 "this shape or library");
1435}
1436
1438 const Array<OneD, const NekDouble> &inarray,
1439 Array<OneD, NekDouble> &outarray)
1440{
1441 v_MultiplyByStdQuadratureMetric(inarray, outarray);
1442}
1443
1445 const Array<OneD, const NekDouble> &inarray,
1446 Array<OneD, NekDouble> &outarray)
1447{
1448 boost::ignore_unused(inarray, outarray);
1450 "Method does not exist for this shape or library");
1451}
1452
1454 const Array<OneD, const NekDouble> &inarray,
1455 Array<OneD, NekDouble> &outarray)
1456{
1457 boost::ignore_unused(inarray, outarray);
1458 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1459}
1460
1462 const Array<OneD, const NekDouble> &inarray,
1463 Array<OneD, NekDouble> &outarray, bool multiplybyweights)
1464{
1465 boost::ignore_unused(inarray, outarray, multiplybyweights);
1466 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1467}
1468
1469/**
1470 *
1471 */
1473 const Array<OneD, const NekDouble> &direction,
1474 const Array<OneD, const NekDouble> &inarray,
1475 Array<OneD, NekDouble> &outarray)
1476{
1477 boost::ignore_unused(direction, inarray, outarray);
1478 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1479}
1480
1482 const int dir, const Array<OneD, const NekDouble> &inarray,
1483 Array<OneD, NekDouble> &outarray)
1484{
1485 boost::ignore_unused(dir, inarray, outarray);
1486 NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1487}
1488
1490 Array<OneD, NekDouble> &outarray,
1491 const StdMatrixKey &mkey)
1492{
1493 // If this function is not reimplemented on shape level, the function
1494 // below will be called
1495 MassMatrixOp_MatFree(inarray, outarray, mkey);
1496}
1497
1499 const Array<OneD, const NekDouble> &inarray,
1500 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1501{
1502 // If this function is not reimplemented on shape level, the function
1503 // below will be called
1504 LaplacianMatrixOp_MatFree(inarray, outarray, mkey);
1505}
1506
1508 const StdMatrixKey &mkey)
1509{
1510 boost::ignore_unused(array, mkey);
1511 ASSERTL0(false, "This function is not defined in StdExpansion.");
1512}
1513
1515 const NekDouble alpha,
1516 const NekDouble exponent,
1517 const NekDouble cutoff)
1518{
1519 boost::ignore_unused(array, alpha, exponent, cutoff);
1520 ASSERTL0(false, "This function is not defined in StdExpansion.");
1521}
1522
1524 int numMin, const Array<OneD, const NekDouble> &inarray,
1525 Array<OneD, NekDouble> &outarray)
1526{
1527 boost::ignore_unused(numMin, inarray, outarray);
1528 ASSERTL0(false, "This function is not defined in StdExpansion.");
1529}
1530
1532 const int k1, const int k2, const Array<OneD, const NekDouble> &inarray,
1533 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1534{
1535 // If this function is not reimplemented on shape level, the function
1536 // below will be called
1537 LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
1538}
1539
1541 const int i, const Array<OneD, const NekDouble> &inarray,
1542 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1543{
1544 // If this function is not reimplemented on shape level, the function
1545 // below will be called
1546 WeakDerivMatrixOp_MatFree(i, inarray, outarray, mkey);
1547}
1548
1550 const Array<OneD, const NekDouble> &inarray,
1551 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1552{
1553 // If this function is not reimplemented on shape level, the function
1554 // below will be called
1555 WeakDirectionalDerivMatrixOp_MatFree(inarray, outarray, mkey);
1556}
1557
1559 const Array<OneD, const NekDouble> &inarray,
1560 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1561{
1562 // If this function is not reimplemented on shape level, the function
1563 // below will be called
1564 MassLevelCurvatureMatrixOp_MatFree(inarray, outarray, mkey);
1565}
1566
1568 const Array<OneD, const NekDouble> &inarray,
1569 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1570{
1571 // If this function is not reimplemented on shape level, the function
1572 // below will be called
1573 LinearAdvectionMatrixOp_MatFree(inarray, outarray, mkey);
1574}
1575
1577 const Array<OneD, const NekDouble> &inarray,
1578 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey,
1579 bool addDiffusionTerm)
1580{
1581 // If this function is not reimplemented on shape level, the function
1582 // below will be called
1584 addDiffusionTerm);
1585}
1586
1588 const Array<OneD, const NekDouble> &inarray,
1589 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1590{
1591 // If this function is not reimplemented on shape level, the function
1592 // below will be called
1593 HelmholtzMatrixOp_MatFree(inarray, outarray, mkey);
1594}
1595
1597 const Array<OneD, const NekDouble> &inarray,
1598 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1599{
1600 // If this function is not reimplemented on shape level, the function
1601 // below will be called
1602 LaplacianMatrixOp_MatFree_GenericImpl(inarray, outarray, mkey);
1603}
1604
1606 const Array<OneD, const NekDouble> &inarray,
1608{
1609 boost::ignore_unused(inarray, outarray, wsp);
1610 ASSERTL0(false, "Not implemented.");
1611}
1612
1614 const Array<OneD, const NekDouble> &inarray,
1615 Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1616{
1617 // If this function is not reimplemented on shape level, the function
1618 // below will be called
1619 HelmholtzMatrixOp_MatFree_GenericImpl(inarray, outarray, mkey);
1620}
1621
1623 const DNekScalMatSharedPtr &m_transformationmatrix)
1624{
1625 boost::ignore_unused(m_transformationmatrix);
1626 NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
1627 return NullDNekMatSharedPtr;
1628}
1629
1631 const Array<OneD, const NekDouble> &inarray,
1632 Array<OneD, NekDouble> &outarray, int npset)
1633{
1635 DNekMatSharedPtr intmat;
1636
1637 int nqtot = GetTotPoints();
1638 int np = 0;
1639 if (npset == -1) // use values from basis num points()
1640 {
1641 int nqbase;
1642 for (int i = 0; i < m_base.size(); ++i)
1643 {
1644 nqbase = m_base[i]->GetNumPoints();
1645 np = std::max(np, nqbase);
1646 }
1647
1648 StdMatrixKey Ikey(ePhysInterpToEquiSpaced, shape, *this);
1649 intmat = GetStdMatrix(Ikey);
1650 }
1651 else
1652 {
1653 np = npset;
1654
1655 ConstFactorMap cmap;
1656 cmap[eFactorConst] = np;
1657 StdMatrixKey Ikey(ePhysInterpToEquiSpaced, shape, *this, cmap);
1658 intmat = GetStdMatrix(Ikey);
1659 }
1660
1661 NekVector<NekDouble> in(nqtot, inarray, eWrapper);
1663 LibUtilities::GetNumberOfCoefficients(shape, np, np, np), outarray,
1664 eWrapper);
1665 out = (*intmat) * in;
1666}
1667
1669 bool standard)
1670{
1671 boost::ignore_unused(conn, standard);
1673 "GetSimplexEquiSpacedConnectivity not"
1674 " implemented for " +
1675 static_cast<std::string>(
1677}
1678
1680 const Array<OneD, const NekDouble> &inarray,
1681 Array<OneD, NekDouble> &outarray)
1682{
1684
1685 // inarray has to be consistent with NumModes definition
1686 // There is also a check in GetStdMatrix to see if all
1687 // modes are of the same size
1688 ConstFactorMap cmap;
1689
1690 cmap[eFactorConst] = m_base[0]->GetNumModes();
1691 StdMatrixKey Ikey(eEquiSpacedToCoeffs, shape, *this, cmap);
1692 DNekMatSharedPtr intmat = GetStdMatrix(Ikey);
1693
1695 NekVector<NekDouble> out(m_ncoeffs, outarray, eWrapper);
1696 out = (*intmat) * in;
1697}
1698
1699} // namespace StdRegions
1700} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:272
Describes the specification for a Basis.
Definition: Basis.h:47
Defines a specification for a set of points.
Definition: Points.h:55
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:71
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:675
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:130
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:140
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:497
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:807
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:814
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:609
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:799
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:758
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:821
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:765
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:373
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:680
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:828
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:843
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:779
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:92
const VarCoeffMap & GetVarCoeffs() const
Definition: StdMatrixKey.h:173
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:87
bool HasVarCoeff(const StdRegions::VarCoeffType &coeff) const
Definition: StdMatrixKey.h:178
LibUtilities::PointsType GetNodalPointsType() const
Definition: StdMatrixKey.h:97
const ConstFactorMap & GetConstFactors() const
Definition: StdMatrixKey.h:142
const Array< OneD, const NekDouble > & GetVarCoeff(const StdRegions::VarCoeffType &coeff) const
Definition: StdMatrixKey.h:152
NekDouble GetConstFactor(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:128
bool ConstFactorExists(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:137
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:137
const char *const ShapeTypeMap[SIZE_ShapeType]
Definition: ShapeType.hpp:79
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:310
@ P
Monomial polynomials .
Definition: BasisType.h:64
static const PointsKey NullPointsKey(0, eNoPointsType)
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
const char *const MatrixTypeMap[]
Definition: StdRegions.hpp:145
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:408
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:409
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:353
std::shared_ptr< StdMatrixKey > StdMatrixKeySharedPtr
Definition: StdMatrixKey.h:203
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
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.cpp:207
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:617
void Vabs(int n, const T *x, const int incx, T *y, const int incy)
vabs: y = |x|
Definition: Vmath.cpp:548
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.cpp:569
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.cpp:354
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.cpp:245
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487
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.cpp:994
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191
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.cpp:414
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294