Nektar++
StdSegExp.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: StdSegExp.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: Routines within Standard Segment Expansions
32//
33///////////////////////////////////////////////////////////////////////////////
34
37
38using namespace std;
39
40namespace Nektar::StdRegions
41{
42/** \brief Constructor using BasisKey class for quadrature points and
43 * order definition
44 *
45 * \param Ba BasisKey class definition containing order and quadrature
46 * points.
47 */
49 : StdExpansion(Ba.GetNumModes(), 1, Ba),
50 StdExpansion1D(Ba.GetNumModes(), Ba)
51{
52}
53
54/** \brief Return Shape of region, using ShapeType enum list.
55 * i.e. Segment
56 */
58{
60}
61
63{
64
65 bool returnval = false;
66
68 {
69 returnval = true;
70 }
71
73 {
74 returnval = true;
75 }
76
77 return returnval;
78}
79
80//---------------------------------------------------------------------
81// Integration Methods
82//---------------------------------------------------------------------
83
84/** \brief Integrate the physical point list \a inarray over region
85 * and return the value
86 *
87 * \param inarray definition of function to be integrated evauluated at
88 * quadrature point of expansion.
89 * \return returns \f$\int^1_{-1} u(\xi_1)d \xi_1 \f$ where \f$inarray[i]
90 * = u(\xi_{1i}) \f$
91 */
93{
94 NekDouble Int = 0.0;
95 int nquad0 = m_base[0]->GetNumPoints();
96 Array<OneD, NekDouble> tmp(nquad0);
99
100 // multiply by integration constants
101 Vmath::Vmul(nquad0, inarray, 1, w0, 1, tmp, 1);
102
103 Int = Vmath::Vsum(nquad0, tmp, 1);
104
105 return Int;
106}
107
108//---------------------------------------------------------------------
109// Differentiation Methods
110//---------------------------------------------------------------------
111
112/** \brief Evaluate the derivative \f$ d/d{\xi_1} \f$ at the physical
113 * quadrature points given by \a inarray and return in \a outarray.
114 *
115 * This is a wrapper around StdExpansion1D::Tensor_Deriv
116 * \param inarray array of a function evaluated at the quadrature points
117 * \param outarray the resulting array of the derivative \f$
118 * du/d_{\xi_1}|_{\xi_{1i}} \f$ will be stored in the array \a outarra
119 */
120
123 [[maybe_unused]] Array<OneD, NekDouble> &out_d1,
124 [[maybe_unused]] Array<OneD, NekDouble> &out_d2)
125{
126 PhysTensorDeriv(inarray, out_d0);
127}
128
129void StdSegExp::v_PhysDeriv([[maybe_unused]] const int dir,
130 const Array<OneD, const NekDouble> &inarray,
131 Array<OneD, NekDouble> &outarray)
132{
133 ASSERTL1(dir == 0, "input dir is out of range");
134 PhysTensorDeriv(inarray, outarray);
135 // PhysDeriv(inarray, outarray);
136}
137
140 [[maybe_unused]] Array<OneD, NekDouble> &out_d1,
141 [[maybe_unused]] Array<OneD, NekDouble> &out_d2)
142{
143 PhysTensorDeriv(inarray, out_d0);
144}
145
146void StdSegExp::v_StdPhysDeriv([[maybe_unused]] const int dir,
147 const Array<OneD, const NekDouble> &inarray,
148 Array<OneD, NekDouble> &outarray)
149{
150 ASSERTL1(dir == 0, "input dir is out of range");
151 PhysTensorDeriv(inarray, outarray);
152}
153
154//---------------------------------------------------------------------
155// Transforms
156//---------------------------------------------------------------------
157
158/** \brief Backward transform from coefficient space given
159 * in \a inarray and evaluate at the physical quadrature
160 * points \a outarray
161 *
162 * Operation can be evaluated as \f$ u(\xi_{1i}) =
163 * \sum_{p=0}^{order-1} \hat{u}_p \phi_p(\xi_{1i}) \f$ or equivalently
164 * \f$ {\bf u} = {\bf B}^T {\bf \hat{u}} \f$ where
165 * \f${\bf B}[i][j] = \phi_i(\xi_{1j}), \mbox{\_coeffs}[p] = {\bf
166 * \hat{u}}[p] \f$
167 *
168 * The function takes the coefficient array \a inarray as
169 * input for the transformation
170 *
171 * \param inarray: the coeffficients of the expansion
172 *
173 * \param outarray: the resulting array of the values of the function at
174 * the physical quadrature points will be stored in the array \a outarray
175 */
176
178 Array<OneD, NekDouble> &outarray)
179{
180 int nquad = m_base[0]->GetNumPoints();
181
182 if (m_base[0]->Collocation())
183 {
184 Vmath::Vcopy(nquad, inarray, 1, outarray, 1);
185 }
186 else
187 {
188 Blas::Dgemv('N', nquad, m_base[0]->GetNumModes(), 1.0,
189 (m_base[0]->GetBdata()).get(), nquad, &inarray[0], 1, 0.0,
190 &outarray[0], 1);
191 }
192}
193
195 const Array<OneD, const NekDouble> &inarray,
196 Array<OneD, NekDouble> &outarray)
197{
198 int n_coeffs = inarray.size();
199
200 Array<OneD, NekDouble> coeff(n_coeffs);
201 Array<OneD, NekDouble> coeff_tmp(n_coeffs, 0.0);
204
205 int nmodes0 = m_base[0]->GetNumModes();
206
207 Vmath::Vcopy(n_coeffs, inarray, 1, coeff_tmp, 1);
208
209 const LibUtilities::PointsKey Pkey0(nmodes0,
211
212 LibUtilities::BasisKey b0(m_base[0]->GetBasisType(), nmodes0, Pkey0);
213
214 LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A, nmodes0, Pkey0);
215
216 LibUtilities::InterpCoeff1D(b0, coeff_tmp, bortho0, coeff);
217
218 Vmath::Zero(n_coeffs, coeff_tmp, 1);
219
220 Vmath::Vcopy(numMin, tmp = coeff, 1, tmp2 = coeff_tmp, 1);
221
222 LibUtilities::InterpCoeff1D(bortho0, coeff_tmp, b0, outarray);
223}
224
225/**
226 * \brief Forward transform from physical quadrature space stored in \a
227 * inarray and evaluate the expansion coefficients and store in \a
228 * outarray
229 *
230 * Perform a forward transform using a Galerkin projection by taking the
231 * inner product of the physical points and multiplying by the inverse
232 * of the mass matrix using the Solve method of the standard matrix
233 * container holding the local mass matrix, i.e. \f$ {\bf \hat{u}} =
234 * {\bf M}^{-1} {\bf I} \f$ where \f$ {\bf I}[p] = \int^1_{-1}
235 * \phi_p(\xi_1) u(\xi_1) d\xi_1 \f$
236 *
237 * This function stores the expansion coefficients calculated by the
238 * transformation in the coefficient space array \a outarray
239 *
240 * \param inarray: array of physical quadrature points to be transformed
241 * \param outarray: the coeffficients of the expansion
242 */
244 Array<OneD, NekDouble> &outarray)
245{
246 if (m_base[0]->Collocation())
247 {
248 Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
249 }
250 else
251 {
252 v_IProductWRTBase(inarray, outarray);
253
254 // get Mass matrix inverse
255 StdMatrixKey masskey(eInvMass, v_DetShapeType(), *this);
256 DNekMatSharedPtr matsys = GetStdMatrix(masskey);
257
258 NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
260
261 out = (*matsys) * in;
262 }
263}
264
266 const Array<OneD, const NekDouble> &inarray,
267 Array<OneD, NekDouble> &outarray)
268{
269 if (m_base[0]->Collocation())
270 {
271 Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
272 }
273 else
274 {
275 int nInteriorDofs = m_ncoeffs - 2;
276 int offset = 0;
277
278 switch (m_base[0]->GetBasisType())
279 {
281 {
282 offset = 1;
283 }
284 break;
286 {
287 nInteriorDofs = m_ncoeffs;
288 offset = 0;
289 }
290 break;
293 {
294 offset = 2;
295 }
296 break;
297 default:
298 ASSERTL0(false, "This type of FwdTrans is not defined for this "
299 "expansion type");
300 }
301
302 fill(outarray.get(), outarray.get() + m_ncoeffs, 0.0);
303
305 {
306 outarray[GetVertexMap(0)] = inarray[0];
307 outarray[GetVertexMap(1)] = inarray[m_base[0]->GetNumPoints() - 1];
308
309 if (m_ncoeffs > 2)
310 {
311 // ideally, we would like to have tmp0 to be replaced by
312 // outarray (currently MassMatrixOp does not allow aliasing)
315
316 StdMatrixKey masskey(eMass, v_DetShapeType(), *this);
317 MassMatrixOp(outarray, tmp0, masskey);
318 v_IProductWRTBase(inarray, tmp1);
319
320 Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
321
322 // get Mass matrix inverse (only of interior DOF)
323 DNekMatSharedPtr matsys =
324 (m_stdStaticCondMatrixManager[masskey])->GetBlock(1, 1);
325
326 Blas::Dgemv('N', nInteriorDofs, nInteriorDofs, 1.0,
327 &(matsys->GetPtr())[0], nInteriorDofs,
328 tmp1.get() + offset, 1, 0.0,
329 outarray.get() + offset, 1);
330 }
331 }
332 else
333 {
334 StdSegExp::v_FwdTrans(inarray, outarray);
335 }
336 }
337}
338
340 Array<OneD, NekDouble> &outarray)
341{
342 v_BwdTrans(inarray, outarray);
343}
344
345//---------------------------------------------------------------------
346// Inner product functions
347//---------------------------------------------------------------------
348
349/** \brief Inner product of \a inarray over region with respect to
350 * expansion basis \a base and return in \a outarray
351 *
352 * Calculate \f$ I[p] = \int^{1}_{-1} \phi_p(\xi_1) u(\xi_1) d\xi_1
353 * = \sum_{i=0}^{nq-1} \phi_p(\xi_{1i}) u(\xi_{1i}) w_i \f$ where
354 * \f$ outarray[p] = I[p], inarray[i] = u(\xi_{1i}), base[p*nq+i] =
355 * \phi_p(\xi_{1i}) \f$.
356 *
357 * \param base an array defining the local basis for the inner product
358 * usually passed from Basis->GetBdata() or Basis->GetDbdata()
359 * \param inarray: physical point array of function to be integrated
360 * \f$ u(\xi_1) \f$
361 * \param coll_check flag to identify when a Basis->Collocation() call
362 * should be performed to see if this is a GLL_Lagrange basis with a
363 * collocation property. (should be set to 0 if taking the inner
364 * product with respect to the derivative of basis)
365 * \param outarray the values of the inner product with respect to
366 * each basis over region will be stored in the array \a outarray as
367 * output of the function
368 */
369
371 const Array<OneD, const NekDouble> &inarray,
372 Array<OneD, NekDouble> &outarray,
373 [[maybe_unused]] int coll_check)
374{
375 int nquad = m_base[0]->GetNumPoints();
376 Array<OneD, NekDouble> tmp(nquad);
378
379 Vmath::Vmul(nquad, inarray, 1, w, 1, tmp, 1);
380
381 /* Comment below was a bug for collocated basis
382 if(coll_check&&m_base[0]->Collocation())
383 {
384 Vmath::Vcopy(nquad, tmp, 1, outarray, 1);
385 }
386 else
387 {
388 Blas::Dgemv('T',nquad,m_ncoeffs,1.0,base.get(),nquad,
389 &tmp[0],1,0.0,outarray.get(),1);
390 }*/
391
392 // Correct implementation
393 Blas::Dgemv('T', nquad, m_ncoeffs, 1.0, base.get(), nquad, &tmp[0], 1, 0.0,
394 outarray.get(), 1);
395}
396
397/** \brief Inner product of \a inarray over region with respect to the
398 * expansion basis (this)->m_base[0] and return in \a outarray
399 *
400 * Wrapper call to StdSegExp::IProductWRTBase
401 * \param inarray array of function values evaluated at the physical
402 * collocation points
403 * \param outarray the values of the inner product with respect to
404 * each basis over region will be stored in the array \a outarray as
405 * output of the function
406 */
408 Array<OneD, NekDouble> &outarray)
409{
410 v_IProductWRTBase(m_base[0]->GetBdata(), inarray, outarray, 1);
411}
412
414 const int dir, const Array<OneD, const NekDouble> &inarray,
415 Array<OneD, NekDouble> &outarray)
416{
417 StdSegExp::IProductWRTDerivBase_SumFac(dir, inarray, outarray);
418}
419
421 [[maybe_unused]] const int dir, const Array<OneD, const NekDouble> &inarray,
422 Array<OneD, NekDouble> &outarray)
423{
424 ASSERTL1(dir == 0, "input dir is out of range");
425 StdSegExp::v_IProductWRTBase(m_base[0]->GetDbdata(), inarray, outarray, 1);
426}
427
429 const Array<OneD, const NekDouble> &inarray,
430 Array<OneD, NekDouble> &outarray, bool multiplybyweights)
431{
432 int nquad = m_base[0]->GetNumPoints();
433 Array<OneD, NekDouble> tmp(nquad);
435 Array<OneD, const NekDouble> base = m_base[0]->GetBdata();
436
437 if (multiplybyweights)
438 {
439 Vmath::Vmul(nquad, inarray, 1, w, 1, tmp, 1);
440
441 Blas::Dgemv('T', nquad, m_ncoeffs, 1.0, base.get(), nquad, &tmp[0], 1,
442 0.0, outarray.get(), 1);
443 }
444 else
445 {
446 Blas::Dgemv('T', nquad, m_ncoeffs, 1.0, base.get(), nquad, &inarray[0],
447 1, 0.0, outarray.get(), 1);
448 }
449}
450
451//----------------------------
452// Evaluation
453//----------------------------
456{
457 eta[0] = xi[0];
458}
459
462{
463 xi[0] = eta[0];
464}
465
466void StdSegExp::v_FillMode(const int mode, Array<OneD, NekDouble> &outarray)
467{
468 int nquad = m_base[0]->GetNumPoints();
469 const NekDouble *base = m_base[0]->GetBdata().get();
470
471 ASSERTL2(mode <= m_ncoeffs,
472 "calling argument mode is larger than total expansion order");
473
474 Vmath::Vcopy(nquad, (NekDouble *)base + mode * nquad, 1, &outarray[0], 1);
475}
476
478 const Array<OneD, const NekDouble> &coords, int mode)
479{
480 return StdExpansion::BaryEvaluateBasis<0>(coords[0], mode);
481}
482
484 Array<OneD, NekDouble> &outarray,
485 [[maybe_unused]] const StdMatrixKey &mkey)
486{
487 int nquad = m_base[0]->GetNumPoints();
488
489 Array<OneD, NekDouble> physValues(nquad);
490 Array<OneD, NekDouble> dPhysValuesdx(nquad);
491
492 v_BwdTrans(inarray, physValues);
493
494 // Laplacian matrix operation
495 v_PhysDeriv(physValues, dPhysValuesdx);
496 v_IProductWRTBase(m_base[0]->GetDbdata(), dPhysValuesdx, outarray, 1);
497}
498
500 Array<OneD, NekDouble> &outarray,
501 const StdMatrixKey &mkey)
502{
503 int nquad = m_base[0]->GetNumPoints();
504
505 Array<OneD, NekDouble> physValues(nquad);
506 Array<OneD, NekDouble> dPhysValuesdx(nquad);
508
509 v_BwdTrans(inarray, physValues);
510
511 // mass matrix operation
512 v_IProductWRTBase((m_base[0]->GetBdata()), physValues, wsp, 1);
513
514 // Laplacian matrix operation
515 v_PhysDeriv(physValues, dPhysValuesdx);
516 v_IProductWRTBase(m_base[0]->GetDbdata(), dPhysValuesdx, outarray, 1);
518 outarray.get(), 1);
519}
520
522 const StdMatrixKey &mkey)
523{
524 // Generate an orthogonal expansion
525 int nq = m_base[0]->GetNumPoints();
526 int nmodes = m_base[0]->GetNumModes();
527 // Declare orthogonal basis.
529
531 StdSegExp OrthoExp(B);
532
533 // SVV parameters loaded from the .xml case file
534 NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDiffCoeff);
535 int cutoff = (int)(mkey.GetConstFactor(eFactorSVVCutoffRatio)) * nmodes;
536
537 Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
538
539 // project onto modal space.
540 OrthoExp.FwdTrans(array, orthocoeffs);
541
542 //
543 for (int j = 0; j < nmodes; ++j)
544 {
545 if (j >= cutoff) // to filter out only the "high-modes"
546 {
547 orthocoeffs[j] *=
548 (SvvDiffCoeff *
549 exp(-(j - nmodes) * (j - nmodes) /
550 ((NekDouble)((j - cutoff + 1) * (j - cutoff + 1)))));
551 }
552 else
553 {
554 orthocoeffs[j] *= 0.0;
555 }
556 }
557
558 // backward transform to physical space
559 OrthoExp.BwdTrans(orthocoeffs, array);
560}
561
563 const NekDouble alpha,
564 const NekDouble exponent,
565 const NekDouble cutoff)
566{
567 // Generate an orthogonal expansion
568 int nq = m_base[0]->GetNumPoints();
569 int nmodes = m_base[0]->GetNumModes();
570 int P = nmodes - 1;
571 // Declare orthogonal basis.
573
575 StdSegExp OrthoExp(B);
576
577 // Cutoff
578 int Pcut = cutoff * P;
579
580 // Project onto orthogonal space.
581 Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
582 OrthoExp.FwdTrans(array, orthocoeffs);
583
584 //
585 NekDouble fac;
586 for (int j = 0; j < nmodes; ++j)
587 {
588 // to filter out only the "high-modes"
589 if (j > Pcut)
590 {
591 fac = (NekDouble)(j - Pcut) / ((NekDouble)(P - Pcut));
592 fac = pow(fac, exponent);
593 orthocoeffs[j] *= exp(-alpha * fac);
594 }
595 }
596
597 // backward transform to physical space
598 OrthoExp.BwdTrans(orthocoeffs, array);
599}
600
601// up to here
603 const Array<OneD, const NekDouble> &inarray,
604 Array<OneD, NekDouble> &outarray)
605{
606 int nquad0 = m_base[0]->GetNumPoints();
607
608 const Array<OneD, const NekDouble> &w0 = m_base[0]->GetW();
609
610 Vmath::Vmul(nquad0, inarray.get(), 1, w0.get(), 1, outarray.get(), 1);
611}
612
614 [[maybe_unused]] Array<OneD, NekDouble> &coords_1,
615 [[maybe_unused]] Array<OneD, NekDouble> &coords_2)
616{
617 Blas::Dcopy(GetNumPoints(0), (m_base[0]->GetZ()).get(), 1, &coords_0[0], 1);
618}
619
620//---------------------------------------------------------------------
621// Helper functions
622//---------------------------------------------------------------------
623
625{
626 return 2;
627}
628
630{
631 return 2;
632}
633
634int StdSegExp::v_GetTraceNcoeffs([[maybe_unused]] const int i) const
635{
636 return 1;
637}
638
639int StdSegExp::v_GetTraceIntNcoeffs([[maybe_unused]] const int i) const
640{
641 return 0;
642}
643
644int StdSegExp::v_GetTraceNumPoints([[maybe_unused]] const int i) const
645{
646 return 1;
647}
648
650{
651 return 2;
652}
653
655{
656 return 2;
657}
658
660 const std::vector<unsigned int> &nummodes, int &modes_offset)
661{
662 int nmodes = nummodes[modes_offset];
663 modes_offset += 1;
664
665 return nmodes;
666}
667
668//---------------------------------------------------------------------
669// Wrapper functions
670//---------------------------------------------------------------------
671
673{
675 MatrixType mattype;
676
677 switch (mattype = mkey.GetMatrixType())
678 {
680 {
681 int nq = m_base[0]->GetNumPoints();
682
683 // take definition from key
685 {
686 nq = (int)mkey.GetConstFactor(eFactorConst);
687 }
688
690 Array<OneD, NekDouble> coords(1);
693
694 for (int i = 0; i < neq; ++i)
695 {
696 coords[0] = -1.0 + 2 * i / (NekDouble)(neq - 1);
697 I = m_base[0]->GetI(coords);
698 Vmath::Vcopy(nq, I->GetRawPtr(), 1, Mat->GetRawPtr() + i, neq);
699 }
700 }
701 break;
702 case eFwdTrans:
703 {
704 Mat =
706 StdMatrixKey iprodkey(eIProductWRTBase, v_DetShapeType(), *this);
707 DNekMat &Iprod = *GetStdMatrix(iprodkey);
708 StdMatrixKey imasskey(eInvMass, v_DetShapeType(), *this);
709 DNekMat &Imass = *GetStdMatrix(imasskey);
710
711 (*Mat) = Imass * Iprod;
712 }
713 break;
714 default:
715 {
717
718 if (mattype == eMass)
719 {
720 // For Fourier basis set the imaginary component
721 // of mean mode to have a unit diagonal component
722 // in mass matrix
724 {
725 (*Mat)(1, 1) = 1.0;
726 }
727 }
728 }
729 break;
730 }
731
732 return Mat;
733}
734
736{
737 return v_GenMatrix(mkey);
738}
739
740//---------------------------------------------------------------------
741// Mappings
742//---------------------------------------------------------------------
743
745{
746 if (outarray.size() != NumBndryCoeffs())
747 {
749 }
750 const LibUtilities::BasisType Btype = GetBasisType(0);
751 int nummodes = m_base[0]->GetNumModes();
752
753 outarray[0] = 0;
754
755 switch (Btype)
756 {
761 outarray[1] = nummodes - 1;
762 break;
765 outarray[1] = 1;
766 break;
767 default:
768 ASSERTL0(0, "Mapping array is not defined for this expansion");
769 break;
770 }
771}
772
774{
775 int i;
776 if (outarray.size() != GetNcoeffs() - NumBndryCoeffs())
777 {
779 }
780 const LibUtilities::BasisType Btype = GetBasisType(0);
781
782 switch (Btype)
783 {
788 for (i = 0; i < GetNcoeffs() - 2; i++)
789 {
790 outarray[i] = i + 1;
791 }
792 break;
795 for (i = 0; i < GetNcoeffs() - 2; i++)
796 {
797 outarray[i] = i + 2;
798 }
799 break;
800 default:
801 ASSERTL0(0, "Mapping array is not defined for this expansion");
802 break;
803 }
804}
805
806int StdSegExp::v_GetVertexMap(int localVertexId,
807 [[maybe_unused]] bool useCoeffPacking)
808{
809 ASSERTL0((localVertexId == 0) || (localVertexId == 1),
810 "local vertex id"
811 "must be between 0 or 1");
812
813 int localDOF = localVertexId;
814
816 (localVertexId == 1))
817 {
818 localDOF = m_base[0]->GetNumModes() - 1;
819 }
820 return localDOF;
821}
822
824 Array<OneD, int> &conn, [[maybe_unused]] bool standard)
825{
826 int np = m_base[0]->GetNumPoints();
827
828 conn = Array<OneD, int>(2 * (np - 1));
829 int cnt = 0;
830 for (int i = 0; i < np - 1; ++i)
831 {
832 conn[cnt++] = i;
833 conn[cnt++] = i + 1;
834 }
835}
836
837/** \brief Get the map of the coefficient location to teh
838 * local trace coefficients
839 */
840
841void StdSegExp::v_GetTraceCoeffMap(const unsigned int traceid,
843{
844 int order0 = m_base[0]->GetNumModes();
845
846 ASSERTL0(traceid < 2, "eid must be between 0 and 1");
847
848 if (maparray.size() != 1)
849 {
850 maparray = Array<OneD, unsigned int>(1);
851 }
852
853 const LibUtilities::BasisType bType = GetBasisType(0);
854
855 if (bType == LibUtilities::eModified_A)
856 {
857 maparray[0] = (traceid == 0) ? 0 : 1;
858 }
859 else if (bType == LibUtilities::eGLL_Lagrange ||
861 {
862 maparray[0] = (traceid == 0) ? 0 : order0 - 1;
863 }
864 else
865 {
866 ASSERTL0(false, "Unknown Basis");
867 }
868}
869
872 Array<OneD, int> &signarray,
873 [[maybe_unused]] Orientation orient,
874 [[maybe_unused]] int P,
875 [[maybe_unused]] int Q)
876{
877 v_GetTraceCoeffMap(tid, maparray);
878
879 if (signarray.size() != 1)
880 {
881 signarray = Array<OneD, int>(1, 1);
882 }
883 else
884 {
885 signarray[0] = 1;
886 }
887}
888
890 [[maybe_unused]] const unsigned int eid,
891 Array<OneD, unsigned int> &maparray, Array<OneD, int> &signarray,
892 [[maybe_unused]] Orientation orient, [[maybe_unused]] int P,
893 [[maybe_unused]] int Q)
894{
895 // parameters for higher dimnesion traces
896 if (maparray.size() != 1)
897 {
898 maparray = Array<OneD, unsigned int>(1);
899 }
900
901 maparray[0] = 0;
902
903 if (signarray.size() != 1)
904 {
905 signarray = Array<OneD, int>(1, 1);
906 }
907 else
908 {
909 signarray[0] = 1;
910 }
911}
912
913} // namespace Nektar::StdRegions
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#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.
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Evaluate the derivative at the physical quadrature points given by inarray and return in outarray.
The base class for all shapes.
Definition: StdExpansion.h:65
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:124
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:156
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:603
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:752
void IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
Definition: StdExpansion.h:679
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space.
Definition: StdExpansion.h:424
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:205
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space.
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLess > m_stdStaticCondMatrixManager
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:218
Array< OneD, LibUtilities::BasisSharedPtr > m_base
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:83
NekDouble GetConstFactor(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:124
bool ConstFactorExists(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:133
Class representing a segment element in reference space All interface of this class sits in StdExpans...
Definition: StdSegExp.h:45
void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Inner product of inarray over region with respect to the expansion basis (this)->m_base[0] and return...
Definition: StdSegExp.cpp:407
int v_NumDGBndryCoeffs() const override
Definition: StdSegExp.cpp:654
void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
Definition: StdSegExp.cpp:483
void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdSegExp.cpp:194
void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey) override
Definition: StdSegExp.cpp:521
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
Definition: StdSegExp.cpp:428
void v_GetInteriorMap(Array< OneD, unsigned int > &outarray) override
Definition: StdSegExp.cpp:773
NekDouble v_PhysEvaluateBasis(const Array< OneD, const NekDouble > &coords, int mode) final
Definition: StdSegExp.cpp:477
LibUtilities::ShapeType v_DetShapeType() const override
Return Shape of region, using ShapeType enum list. i.e. Segment.
Definition: StdSegExp.cpp:57
StdSegExp(const LibUtilities::BasisKey &Ba)
Constructor using BasisKey class for quadrature points and order definition.
Definition: StdSegExp.cpp:48
int v_NumBndryCoeffs() const override
Definition: StdSegExp.cpp:649
void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
Definition: StdSegExp.cpp:499
void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdSegExp.cpp:339
void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdSegExp.cpp:413
void v_LocCollapsedToLocCoord(const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi) override
Definition: StdSegExp.cpp:460
void v_GetTraceCoeffMap(const unsigned int traceid, Array< OneD, unsigned int > &maparray) override
Get the map of the coefficient location to teh local trace coefficients.
Definition: StdSegExp.cpp:841
int v_GetNverts() const final
Definition: StdSegExp.cpp:624
void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Forward transform from physical quadrature space stored in inarray and evaluate the expansion coeffic...
Definition: StdSegExp.cpp:243
int v_GetTraceNumPoints(const int i) const final
Definition: StdSegExp.cpp:644
void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray) override
Definition: StdSegExp.cpp:466
int v_GetTraceNcoeffs(const int i) const final
Definition: StdSegExp.cpp:634
void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray) override
Evaluate the derivative at the physical quadrature points given by inarray and return in outarray.
Definition: StdSegExp.cpp:121
void v_GetSimplexEquiSpacedConnectivity(Array< OneD, int > &conn, bool standard=true) override
Definition: StdSegExp.cpp:823
void v_GetCoords(Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2) override
Definition: StdSegExp.cpp:613
int v_GetTraceIntNcoeffs(const int i) const final
Definition: StdSegExp.cpp:639
NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray) override
Integrate the physical point list inarray over region and return the value.
Definition: StdSegExp.cpp:92
int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false) override
Definition: StdSegExp.cpp:806
void v_ExponentialFilter(Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff) override
Definition: StdSegExp.cpp:562
DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey) override
Definition: StdSegExp.cpp:735
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdSegExp.cpp:420
void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdSegExp.cpp:602
void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Backward transform from coefficient space given in inarray and evaluate at the physical quadrature po...
Definition: StdSegExp.cpp:177
int v_GetNtraces() const final
Definition: StdSegExp.cpp:629
void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray) override
Definition: StdSegExp.cpp:744
void v_GetTraceToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient, int P, int Q) override
Definition: StdSegExp.cpp:870
void v_GetElmtTraceToTraceMap(const unsigned int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient, int P, int Q) override
Definition: StdSegExp.cpp:889
void v_StdPhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray) override
Definition: StdSegExp.cpp:138
void v_LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta) override
Definition: StdSegExp.cpp:454
bool v_IsBoundaryInteriorExpansion() const override
Definition: StdSegExp.cpp:62
int v_CalcNumberOfCoefficients(const std::vector< unsigned int > &nummodes, int &modes_offset) override
Definition: StdSegExp.cpp:659
DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
Definition: StdSegExp.cpp:672
void v_FwdTransBndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdSegExp.cpp:265
static void Dgemv(const char &trans, const int &m, const int &n, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = alpha A x plus beta y where A[m x n].
Definition: Blas.hpp:211
static void Dcopy(const int &n, const double *x, const int &incx, double *y, const int &incy)
BLAS level 1: Copy x to y.
Definition: Blas.hpp:128
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
void InterpCoeff1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
Definition: InterpCoeff.cpp:42
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:49
@ P
Monomial polynomials .
Definition: BasisType.h:62
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points.
Definition: BasisType.h:57
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:42
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:56
@ eChebyshev
Chebyshev Polynomials.
Definition: BasisType.h:61
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48
@ eFourier
Fourier Expansion .
Definition: BasisType.h:55
std::vector< double > w(NPUPPER)
std::vector< double > z(NPUPPER)
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
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.hpp:608
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
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.