Nektar++
Loading...
Searching...
No Matches
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 */
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
146//---------------------------------------------------------------------
147// Transforms
148//---------------------------------------------------------------------
149
150/** \brief Backward transform from coefficient space given
151 * in \a inarray and evaluate at the physical quadrature
152 * points \a outarray
153 *
154 * Operation can be evaluated as \f$ u(\xi_{1i}) =
155 * \sum_{p=0}^{order-1} \hat{u}_p \phi_p(\xi_{1i}) \f$ or equivalently
156 * \f$ {\bf u} = {\bf B}^T {\bf \hat{u}} \f$ where
157 * \f${\bf B}[i][j] = \phi_i(\xi_{1j}), \mbox{\_coeffs}[p] = {\bf
158 * \hat{u}}[p] \f$
159 *
160 * The function takes the coefficient array \a inarray as
161 * input for the transformation
162 *
163 * \param inarray: the coeffficients of the expansion
164 *
165 * \param outarray: the resulting array of the values of the function at
166 * the physical quadrature points will be stored in the array \a outarray
167 */
168
170 Array<OneD, NekDouble> &outarray)
171{
172 int nquad = m_base[0]->GetNumPoints();
173
174 if (m_base[0]->Collocation())
175 {
176 Vmath::Vcopy(nquad, inarray, 1, outarray, 1);
177 }
178 else
179 {
180 Blas::Dgemv('N', nquad, m_base[0]->GetNumModes(), 1.0,
181 (m_base[0]->GetBdata()).data(), nquad, &inarray[0], 1, 0.0,
182 &outarray[0], 1);
183 }
184}
185
187 const Array<OneD, const NekDouble> &inarray,
188 Array<OneD, NekDouble> &outarray)
189{
190 int n_coeffs = inarray.size();
191
192 Array<OneD, NekDouble> coeff(n_coeffs);
193 Array<OneD, NekDouble> coeff_tmp(n_coeffs, 0.0);
196
197 int nmodes0 = m_base[0]->GetNumModes();
198
199 Vmath::Vcopy(n_coeffs, inarray, 1, coeff_tmp, 1);
200
201 const LibUtilities::PointsKey Pkey0(nmodes0,
203
204 LibUtilities::BasisKey b0(m_base[0]->GetBasisType(), nmodes0, Pkey0);
205
206 LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A, nmodes0, Pkey0);
207
208 LibUtilities::InterpCoeff1D(b0, coeff_tmp, bortho0, coeff);
209
210 Vmath::Zero(n_coeffs, coeff_tmp, 1);
211
212 Vmath::Vcopy(numMin, tmp = coeff, 1, tmp2 = coeff_tmp, 1);
213
214 LibUtilities::InterpCoeff1D(bortho0, coeff_tmp, b0, outarray);
215}
216
217/**
218 * \brief Forward transform from physical quadrature space stored in \a
219 * inarray and evaluate the expansion coefficients and store in \a
220 * outarray
221 *
222 * Perform a forward transform using a Galerkin projection by taking the
223 * inner product of the physical points and multiplying by the inverse
224 * of the mass matrix using the Solve method of the standard matrix
225 * container holding the local mass matrix, i.e. \f$ {\bf \hat{u}} =
226 * {\bf M}^{-1} {\bf I} \f$ where \f$ {\bf I}[p] = \int^1_{-1}
227 * \phi_p(\xi_1) u(\xi_1) d\xi_1 \f$
228 *
229 * This function stores the expansion coefficients calculated by the
230 * transformation in the coefficient space array \a outarray
231 *
232 * \param inarray: array of physical quadrature points to be transformed
233 * \param outarray: the coeffficients of the expansion
234 */
236 Array<OneD, NekDouble> &outarray)
237{
238 if (m_base[0]->Collocation())
239 {
240 Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
241 }
242 else
243 {
244 v_IProductWRTBase(inarray, outarray);
245
246 // get Mass matrix inverse
247 StdMatrixKey masskey(eInvMass, v_DetShapeType(), *this);
248 DNekMatSharedPtr matsys = GetStdMatrix(masskey);
249
250 NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
252
253 out = (*matsys) * in;
254 }
255}
256
258 const Array<OneD, const NekDouble> &inarray,
259 Array<OneD, NekDouble> &outarray)
260{
261 if (m_base[0]->Collocation())
262 {
263 Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
264 }
265 else
266 {
267 int nInteriorDofs = m_ncoeffs - 2;
268 int offset = 0;
269
270 switch (m_base[0]->GetBasisType())
271 {
273 {
274 offset = 1;
275 }
276 break;
278 {
279 nInteriorDofs = m_ncoeffs;
280 offset = 0;
281 }
282 break;
285 {
286 offset = 2;
287 }
288 break;
289 default:
290 ASSERTL0(false, "This type of FwdTrans is not defined for this "
291 "expansion type");
292 }
293
294 fill(outarray.data(), outarray.data() + m_ncoeffs, 0.0);
295
297 {
298 outarray[GetVertexMap(0)] = inarray[0];
299 outarray[GetVertexMap(1)] = inarray[m_base[0]->GetNumPoints() - 1];
300
301 if (m_ncoeffs > 2)
302 {
303 // ideally, we would like to have tmp0 to be replaced by
304 // outarray (currently MassMatrixOp does not allow aliasing)
307
308 StdMatrixKey masskey(eMass, v_DetShapeType(), *this);
309 MassMatrixOp(outarray, tmp0, masskey);
310 v_IProductWRTBase(inarray, tmp1);
311
312 Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
313
314 // get Mass matrix inverse (only of interior DOF)
315 DNekMatSharedPtr matsys =
316 (m_stdStaticCondMatrixManager[masskey])->GetBlock(1, 1);
317
318 Blas::Dgemv('N', nInteriorDofs, nInteriorDofs, 1.0,
319 &(matsys->GetPtr())[0], nInteriorDofs,
320 tmp1.data() + offset, 1, 0.0,
321 outarray.data() + offset, 1);
322 }
323 }
324 else
325 {
326 StdSegExp::v_FwdTrans(inarray, outarray);
327 }
328 }
329}
330
332 Array<OneD, NekDouble> &outarray)
333{
334 v_BwdTrans(inarray, outarray);
335}
336
337//---------------------------------------------------------------------
338// Inner product functions
339//---------------------------------------------------------------------
340
341/** \brief Inner product of \a inarray over region with respect to
342 * expansion basis \a base and return in \a outarray
343 *
344 * Calculate \f$ I[p] = \int^{1}_{-1} \phi_p(\xi_1) u(\xi_1) d\xi_1
345 * = \sum_{i=0}^{nq-1} \phi_p(\xi_{1i}) u(\xi_{1i}) w_i \f$ where
346 * \f$ outarray[p] = I[p], inarray[i] = u(\xi_{1i}), base[p*nq+i] =
347 * \phi_p(\xi_{1i}) \f$.
348 *
349 * \param base an array defining the local basis for the inner product
350 * usually passed from Basis->GetBdata() or Basis->GetDbdata()
351 * \param inarray: physical point array of function to be integrated
352 * \f$ u(\xi_1) \f$
353 * \param coll_check flag to identify when a Basis->Collocation() call
354 * should be performed to see if this is a GLL_Lagrange basis with a
355 * collocation property. (should be set to 0 if taking the inner
356 * product with respect to the derivative of basis)
357 * \param outarray the values of the inner product with respect to
358 * each basis over region will be stored in the array \a outarray as
359 * output of the function
360 */
361
363 const Array<OneD, const NekDouble> &inarray,
364 Array<OneD, NekDouble> &outarray,
365 [[maybe_unused]] int coll_check)
366{
367 int nquad = m_base[0]->GetNumPoints();
368 Array<OneD, NekDouble> tmp(nquad);
370
371 Vmath::Vmul(nquad, inarray, 1, w, 1, tmp, 1);
372
373 /* Comment below was a bug for collocated basis
374 if(coll_check&&m_base[0]->Collocation())
375 {
376 Vmath::Vcopy(nquad, tmp, 1, outarray, 1);
377 }
378 else
379 {
380 Blas::Dgemv('T',nquad,m_ncoeffs,1.0,base.data(),nquad,
381 &tmp[0],1,0.0,outarray.data(),1);
382 }*/
383
384 // Correct implementation
385 Blas::Dgemv('T', nquad, m_ncoeffs, 1.0, base.data(), nquad, &tmp[0], 1, 0.0,
386 outarray.data(), 1);
387}
388
389/** \brief Inner product of \a inarray over region with respect to the
390 * expansion basis (this)->m_base[0] and return in \a outarray
391 *
392 * Wrapper call to StdSegExp::IProductWRTBase
393 * \param inarray array of function values evaluated at the physical
394 * collocation points
395 * \param outarray the values of the inner product with respect to
396 * each basis over region will be stored in the array \a outarray as
397 * output of the function
398 */
400 Array<OneD, NekDouble> &outarray)
401{
402 v_IProductWRTBase(m_base[0]->GetBdata(), inarray, outarray, 1);
403}
404
406 const int dir, const Array<OneD, const NekDouble> &inarray,
407 Array<OneD, NekDouble> &outarray)
408{
409 StdSegExp::IProductWRTDerivBase_SumFac(dir, inarray, outarray);
410}
411
413 [[maybe_unused]] const int dir, const Array<OneD, const NekDouble> &inarray,
414 Array<OneD, NekDouble> &outarray)
415{
416 ASSERTL1(dir == 0, "input dir is out of range");
417 StdSegExp::v_IProductWRTBase(m_base[0]->GetDbdata(), inarray, outarray, 1);
418}
419
421 const Array<OneD, const NekDouble> &inarray,
422 Array<OneD, NekDouble> &outarray, bool multiplybyweights)
423{
424 int nquad = m_base[0]->GetNumPoints();
425 Array<OneD, NekDouble> tmp(nquad);
427 Array<OneD, const NekDouble> base = m_base[0]->GetBdata();
428
429 if (multiplybyweights)
430 {
431 Vmath::Vmul(nquad, inarray, 1, w, 1, tmp, 1);
432
433 Blas::Dgemv('T', nquad, m_ncoeffs, 1.0, base.data(), nquad, &tmp[0], 1,
434 0.0, outarray.data(), 1);
435 }
436 else
437 {
438 Blas::Dgemv('T', nquad, m_ncoeffs, 1.0, base.data(), nquad, &inarray[0],
439 1, 0.0, outarray.data(), 1);
440 }
441}
442
443//----------------------------
444// Evaluation
445//----------------------------
451
454{
455 xi[0] = eta[0];
456}
457
458void StdSegExp::v_FillMode(const int mode, Array<OneD, NekDouble> &outarray)
459{
460 int nquad = m_base[0]->GetNumPoints();
461 const NekDouble *base = m_base[0]->GetBdata().data();
462
463 ASSERTL2(mode <= m_ncoeffs,
464 "calling argument mode is larger than total expansion order");
465
466 Vmath::Vcopy(nquad, (NekDouble *)base + mode * nquad, 1, &outarray[0], 1);
467}
468
470 const Array<OneD, const NekDouble> &coords, int mode)
471{
472 return StdExpansion::BaryEvaluateBasis<0>(coords[0], mode);
473}
474
476 const Array<OneD, NekDouble> &coord,
477 const Array<OneD, const NekDouble> &inarray,
478 std::array<NekDouble, 3> &firstOrderDerivs)
479{
480 return StdExpansion1D::BaryTensorDeriv(coord, inarray, firstOrderDerivs);
481}
482
484 const Array<OneD, NekDouble> &coord,
485 const Array<OneD, const NekDouble> &inarray,
486 std::array<NekDouble, 3> &firstOrderDerivs,
487 std::array<NekDouble, 6> &secondOrderDerivs)
488{
489 return StdExpansion1D::BaryTensorDeriv(coord, inarray, firstOrderDerivs,
490 secondOrderDerivs);
491}
493 Array<OneD, NekDouble> &outarray,
494 [[maybe_unused]] const StdMatrixKey &mkey)
495{
496 int nquad = m_base[0]->GetNumPoints();
497
498 Array<OneD, NekDouble> physValues(nquad);
499 Array<OneD, NekDouble> dPhysValuesdx(nquad);
500
501 v_BwdTrans(inarray, physValues);
502
503 // Laplacian matrix operation
504 v_PhysDeriv(physValues, dPhysValuesdx);
505 v_IProductWRTBase(m_base[0]->GetDbdata(), dPhysValuesdx, outarray, 1);
506}
507
508void StdSegExp::v_LaplacianMatrixOp(const int k1, const int k2,
509 const Array<OneD, const NekDouble> &inarray,
510 Array<OneD, NekDouble> &outarray,
511 const StdMatrixKey &mkey)
512{
513 StdExpansion::LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
514}
515
517 Array<OneD, NekDouble> &outarray,
518 const StdMatrixKey &mkey)
519{
520 int nquad = m_base[0]->GetNumPoints();
521
522 Array<OneD, NekDouble> physValues(nquad);
523 Array<OneD, NekDouble> dPhysValuesdx(nquad);
525
526 v_BwdTrans(inarray, physValues);
527
528 // mass matrix operation
529 v_IProductWRTBase((m_base[0]->GetBdata()), physValues, wsp, 1);
530
531 // Laplacian matrix operation
532 v_PhysDeriv(physValues, dPhysValuesdx);
533 v_IProductWRTBase(m_base[0]->GetDbdata(), dPhysValuesdx, outarray, 1);
535 outarray.data(), 1);
536}
537
539 const StdMatrixKey &mkey)
540{
541 // Generate an orthogonal expansion
542 int nq = m_base[0]->GetNumPoints();
543 int nmodes = m_base[0]->GetNumModes();
544 // Declare orthogonal basis.
546
548 StdSegExp OrthoExp(B);
549
550 // SVV parameters loaded from the .xml case file
551 NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDiffCoeff);
552 int cutoff = (int)(mkey.GetConstFactor(eFactorSVVCutoffRatio)) * nmodes;
553
554 Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
555
556 // project onto modal space.
557 OrthoExp.FwdTrans(array, orthocoeffs);
558
559 //
560 for (int j = 0; j < nmodes; ++j)
561 {
562 if (j >= cutoff) // to filter out only the "high-modes"
563 {
564 orthocoeffs[j] *=
565 (SvvDiffCoeff *
566 exp(-(j - nmodes) * (j - nmodes) /
567 ((NekDouble)((j - cutoff + 1) * (j - cutoff + 1)))));
568 }
569 else
570 {
571 orthocoeffs[j] *= 0.0;
572 }
573 }
574
575 // backward transform to physical space
576 OrthoExp.BwdTrans(orthocoeffs, array);
577}
578
580 const NekDouble alpha,
581 const NekDouble exponent,
582 const NekDouble cutoff)
583{
584 // Generate an orthogonal expansion
585 int nq = m_base[0]->GetNumPoints();
586 int nmodes = m_base[0]->GetNumModes();
587 int P = nmodes - 1;
588 // Declare orthogonal basis.
590
592 StdSegExp OrthoExp(B);
593
594 // Cutoff
595 int Pcut = cutoff * P;
596
597 // Project onto orthogonal space.
598 Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
599 OrthoExp.FwdTrans(array, orthocoeffs);
600
601 //
602 NekDouble fac;
603 for (int j = 0; j < nmodes; ++j)
604 {
605 // to filter out only the "high-modes"
606 if (j > Pcut)
607 {
608 fac = (NekDouble)(j - Pcut) / ((NekDouble)(P - Pcut));
609 fac = pow(fac, exponent);
610 orthocoeffs[j] *= exp(-alpha * fac);
611 }
612 }
613
614 // backward transform to physical space
615 OrthoExp.BwdTrans(orthocoeffs, array);
616}
617
618// up to here
620 const Array<OneD, const NekDouble> &inarray,
621 Array<OneD, NekDouble> &outarray)
622{
623 int nquad0 = m_base[0]->GetNumPoints();
624
625 const Array<OneD, const NekDouble> &w0 = m_base[0]->GetW();
626
627 Vmath::Vmul(nquad0, inarray.data(), 1, w0.data(), 1, outarray.data(), 1);
628}
629
631 [[maybe_unused]] Array<OneD, NekDouble> &coords_1,
632 [[maybe_unused]] Array<OneD, NekDouble> &coords_2)
633{
634 Blas::Dcopy(GetNumPoints(0), (m_base[0]->GetZ()).data(), 1, &coords_0[0],
635 1);
636}
637
638//---------------------------------------------------------------------
639// Helper functions
640//---------------------------------------------------------------------
641
643{
644 return 2;
645}
646
648{
649 return 2;
650}
651
652int StdSegExp::v_GetTraceNcoeffs([[maybe_unused]] const int i) const
653{
654 return 1;
655}
656
657int StdSegExp::v_GetTraceIntNcoeffs([[maybe_unused]] const int i) const
658{
659 return 0;
660}
661
662int StdSegExp::v_GetTraceNumPoints([[maybe_unused]] const int i) const
663{
664 return 1;
665}
666
668{
669 return 2;
670}
671
673{
674 return 2;
675}
676
678 const std::vector<unsigned int> &nummodes, int &modes_offset)
679{
680 int nmodes = nummodes[modes_offset];
681 modes_offset += 1;
682
683 return nmodes;
684}
685
686//---------------------------------------------------------------------
687// Wrapper functions
688//---------------------------------------------------------------------
689
691{
693 MatrixType mattype;
694
695 switch (mattype = mkey.GetMatrixType())
696 {
698 {
699 int nq = m_base[0]->GetNumPoints();
700
701 // take definition from key
703 {
704 nq = (int)mkey.GetConstFactor(eFactorConst);
705 }
706
708 Array<OneD, NekDouble> coords(1);
711
712 for (int i = 0; i < neq; ++i)
713 {
714 coords[0] = -1.0 + 2 * i / (NekDouble)(neq - 1);
715 I = m_base[0]->GetI(coords);
716 Vmath::Vcopy(nq, I->GetRawPtr(), 1, Mat->GetRawPtr() + i, neq);
717 }
718 }
719 break;
720 case eFwdTrans:
721 {
722 Mat =
724 StdMatrixKey iprodkey(eIProductWRTBase, v_DetShapeType(), *this);
725 DNekMat &Iprod = *GetStdMatrix(iprodkey);
726 StdMatrixKey imasskey(eInvMass, v_DetShapeType(), *this);
727 DNekMat &Imass = *GetStdMatrix(imasskey);
728
729 (*Mat) = Imass * Iprod;
730 }
731 break;
732 default:
733 {
735
736 if (mattype == eMass)
737 {
738 // For Fourier basis set the imaginary component
739 // of mean mode to have a unit diagonal component
740 // in mass matrix
742 {
743 (*Mat)(1, 1) = 1.0;
744 }
745 }
746 }
747 break;
748 }
749
750 return Mat;
751}
752
757
758//---------------------------------------------------------------------
759// Mappings
760//---------------------------------------------------------------------
761
763{
764 if (outarray.size() != NumBndryCoeffs())
765 {
767 }
768 const LibUtilities::BasisType Btype = GetBasisType(0);
769 int nummodes = m_base[0]->GetNumModes();
770
771 outarray[0] = 0;
772
773 switch (Btype)
774 {
779 outarray[1] = nummodes - 1;
780 break;
783 outarray[1] = 1;
784 break;
785 default:
786 ASSERTL0(0, "Mapping array is not defined for this expansion");
787 break;
788 }
789}
790
792{
793 int i;
794 if (outarray.size() != GetNcoeffs() - NumBndryCoeffs())
795 {
797 }
798 const LibUtilities::BasisType Btype = GetBasisType(0);
799
800 switch (Btype)
801 {
806 for (i = 0; i < GetNcoeffs() - 2; i++)
807 {
808 outarray[i] = i + 1;
809 }
810 break;
813 for (i = 0; i < GetNcoeffs() - 2; i++)
814 {
815 outarray[i] = i + 2;
816 }
817 break;
818 default:
819 ASSERTL0(0, "Mapping array is not defined for this expansion");
820 break;
821 }
822}
823
824int StdSegExp::v_GetVertexMap(int localVertexId,
825 [[maybe_unused]] bool useCoeffPacking)
826{
827 ASSERTL0((localVertexId == 0) || (localVertexId == 1),
828 "local vertex id"
829 "must be between 0 or 1");
830
831 int localDOF = localVertexId;
832
834 (localVertexId == 1))
835 {
836 localDOF = m_base[0]->GetNumModes() - 1;
837 }
838 return localDOF;
839}
840
842 Array<OneD, int> &conn, [[maybe_unused]] bool standard)
843{
844 int np = m_base[0]->GetNumPoints();
845
846 conn = Array<OneD, int>(2 * (np - 1));
847 int cnt = 0;
848 for (int i = 0; i < np - 1; ++i)
849 {
850 conn[cnt++] = i;
851 conn[cnt++] = i + 1;
852 }
853}
854
855/** \brief Get the map of the coefficient location to teh
856 * local trace coefficients
857 */
858
859void StdSegExp::v_GetTraceCoeffMap(const unsigned int traceid,
861{
862 int order0 = m_base[0]->GetNumModes();
863
864 ASSERTL0(traceid < 2, "eid must be between 0 and 1");
865
866 if (maparray.size() != 1)
867 {
868 maparray = Array<OneD, unsigned int>(1);
869 }
870
871 const LibUtilities::BasisType bType = GetBasisType(0);
872
873 if (bType == LibUtilities::eModified_A)
874 {
875 maparray[0] = (traceid == 0) ? 0 : 1;
876 }
877 else if (bType == LibUtilities::eGLL_Lagrange ||
879 {
880 maparray[0] = (traceid == 0) ? 0 : order0 - 1;
881 }
882 else
883 {
884 ASSERTL0(false, "Unknown Basis");
885 }
886}
887
890 Array<OneD, int> &signarray,
891 [[maybe_unused]] Orientation orient,
892 [[maybe_unused]] int P,
893 [[maybe_unused]] int Q)
894{
895 v_GetTraceCoeffMap(tid, maparray);
896
897 if (signarray.size() != 1)
898 {
899 signarray = Array<OneD, int>(1, 1);
900 }
901 else
902 {
903 signarray[0] = 1;
904 }
905}
906
908 [[maybe_unused]] const unsigned int eid,
909 Array<OneD, unsigned int> &maparray, Array<OneD, int> &signarray,
910 [[maybe_unused]] Orientation orient, [[maybe_unused]] int P,
911 [[maybe_unused]] int Q)
912{
913 // parameters for higher dimnesion traces
914 if (maparray.size() != 1)
915 {
916 maparray = Array<OneD, unsigned int>(1);
917 }
918
919 maparray[0] = 0;
920
921 if (signarray.size() != 1)
922 {
923 signarray = Array<OneD, int>(1, 1);
924 }
925 else
926 {
927 signarray[0] = 1;
928 }
929}
930
931} // namespace Nektar::StdRegions
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
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.
NekDouble BaryTensorDeriv(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
The base class for all shapes.
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space.
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
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.
Array< OneD, LibUtilities::BasisSharedPtr > m_base
MatrixType GetMatrixType() const
NekDouble GetConstFactor(const ConstFactorType &factor) const
bool ConstFactorExists(const ConstFactorType &factor) const
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...
int v_NumDGBndryCoeffs() const override
NekDouble v_PhysEvalFirstDeriv(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey) override
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
void v_GetInteriorMap(Array< OneD, unsigned int > &outarray) override
NekDouble v_PhysEvaluateBasis(const Array< OneD, const NekDouble > &coords, int mode) final
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
NekDouble v_PhysEvalFirstSecondDeriv(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs, std::array< NekDouble, 6 > &secondOrderDerivs) override
int v_NumBndryCoeffs() const override
void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_LocCollapsedToLocCoord(const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi) override
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.
int v_GetNverts() const final
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...
int v_GetTraceNumPoints(const int i) const final
void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray) override
int v_GetTraceNcoeffs(const int i) const final
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.
void v_GetSimplexEquiSpacedConnectivity(Array< OneD, int > &conn, bool standard=true) override
void v_GetCoords(Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2) override
int v_GetTraceIntNcoeffs(const int i) const final
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
void v_ExponentialFilter(Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff) override
DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey) override
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
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...
int v_GetNtraces() const final
void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray) override
void v_GetTraceToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient, int P, int Q) override
void v_GetElmtTraceToTraceMap(const unsigned int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient, int P, int Q) override
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
void v_LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta) override
bool v_IsBoundaryInteriorExpansion() const override
Definition StdSegExp.cpp:62
int v_CalcNumberOfCoefficients(const std::vector< unsigned int > &nummodes, int &modes_offset) override
DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
void v_FwdTransBndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
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)
@ 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::shared_ptr< DNekMat > DNekMatSharedPtr
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.