Nektar++
Expansion2D.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: Expansion2D.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: File for Expansion2D routines
32//
33///////////////////////////////////////////////////////////////////////////////
34
42#include <LocalRegions/SegExp.h>
45
46using namespace std;
47
49{
51 : StdExpansion(), Expansion(pGeom), StdExpansion2D()
52{
53}
54
56{
57 DNekScalMatSharedPtr returnval;
59
61 "Geometric information is not set up");
62
63 switch (mkey.GetMatrixType())
64 {
66 {
67 if ((m_metricinfo->GetGtype() == SpatialDomains::eDeformed) ||
69 {
70 NekDouble one = 1.0;
71 DNekMatSharedPtr mat = GenMatrix(mkey);
72
73 returnval =
75 }
76 else
77 {
78 NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
80
81 returnval =
83 }
84 }
85 break;
87 {
88 MatrixKey masskey(mkey, StdRegions::eMass);
89 DNekScalMat &MassMat = *GetLocMatrix(masskey);
90
91 // Generate a local copy of traceMat
94
96 "Need to specify eFactorGJP to construct "
97 "a HelmholtzGJP matrix");
98
100
101 factor /= MassMat.Scale();
102
103 int ntot = MassMat.GetRows() * MassMat.GetColumns();
104
105 Vmath::Svtvp(ntot, factor, &NDTraceMat->GetPtr()[0], 1,
106 MassMat.GetRawPtr(), 1, &NDTraceMat->GetPtr()[0], 1);
107
109 MassMat.Scale(), NDTraceMat);
110 }
111 break;
113 {
114 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
115 {
116 NekDouble one = 1.0;
118 DetShapeType(), *this);
119 DNekMatSharedPtr mat = GenMatrix(masskey);
120 mat->Invert();
121
122 returnval =
124 }
125 else
126 {
127 NekDouble fac = 1.0 / (m_metricinfo->GetJac(ptsKeys))[0];
128 DNekMatSharedPtr mat = GetStdMatrix(mkey);
129
130 returnval =
132 }
133 }
134 break;
138 {
139 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed ||
147 {
148 NekDouble one = 1.0;
149 DNekMatSharedPtr mat = GenMatrix(mkey);
150
151 returnval =
153 }
154 else
155 {
156 NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
158 m_metricinfo->GetDerivFactors(ptsKeys);
159 int dir = 0;
161 {
162 dir = 0;
163 }
165 {
166 dir = 1;
167 }
169 {
170 dir = 2;
171 }
172
174 mkey.GetShapeType(), *this);
176 mkey.GetShapeType(), *this);
177
178 DNekMat &deriv0 = *GetStdMatrix(deriv0key);
179 DNekMat &deriv1 = *GetStdMatrix(deriv1key);
180
181 int rows = deriv0.GetRows();
182 int cols = deriv1.GetColumns();
183
184 DNekMatSharedPtr WeakDeriv =
186 (*WeakDeriv) =
187 df[2 * dir][0] * deriv0 + df[2 * dir + 1][0] * deriv1;
188
190 jac, WeakDeriv);
191 }
192 }
193 break;
195 {
196 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed ||
197 mkey.GetNVarCoeff())
198 {
199 NekDouble one = 1.0;
200 DNekMatSharedPtr mat = GenMatrix(mkey);
201
202 returnval =
204 }
205 else
206 {
207 int shapedim = 2;
208
209 // dfdirxi = tan_{xi_x} * d \xi/dx
210 // + tan_{xi_y} * d \xi/dy
211 // + tan_{xi_z} * d \xi/dz
212 // dfdireta = tan_{eta_x} * d \eta/dx
213 // + tan_{xi_y} * d \xi/dy
214 // + tan_{xi_z} * d \xi/dz
215 NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
217 m_metricinfo->GetDerivFactors(ptsKeys);
218
219 Array<OneD, NekDouble> direction =
221
222 // d / dx = df[0]*deriv0 + df[1]*deriv1
223 // d / dy = df[2]*deriv0 + df[3]*deriv1
224 // d / dz = df[4]*deriv0 + df[5]*deriv1
225
226 // dfdir[dir] = e \cdot (d/dx, d/dy, d/dz)
227 // = (e^0 * df[0] + e^1 * df[2]
228 // + e^2 * df[4]) * deriv0
229 // + (e^0 * df[1] + e^1 * df[3]
230 // + e^2 * df[5]) * deriv1
231 // dfdir[dir] = e^0 * df[2 * 0 + dir]
232 // + e^1 * df[2 * 1 + dir]
233 // + e^2 * df [ 2 * 2 + dir]
234 Array<OneD, Array<OneD, NekDouble>> dfdir(shapedim);
235 Expansion::ComputeGmatcdotMF(df, direction, dfdir);
236
239
240 dfdirxi[StdRegions::eVarCoeffWeakDeriv] = dfdir[0];
241 dfdireta[StdRegions::eVarCoeffWeakDeriv] = dfdir[1];
242
244 mkey.GetShapeType(), *this,
247 mkey.GetShapeType(), *this,
249
250 DNekMat &derivxi = *GetStdMatrix(derivxikey);
251 DNekMat &deriveta = *GetStdMatrix(derivetakey);
252
253 int rows = derivxi.GetRows();
254 int cols = deriveta.GetColumns();
255
256 DNekMatSharedPtr WeakDirDeriv =
258
259 (*WeakDirDeriv) = derivxi + deriveta;
260
261 // Add (\nabla \cdot e^k ) Mass
263 DiveMass[StdRegions::eVarCoeffMass] =
265 StdRegions::StdMatrixKey stdmasskey(
266 StdRegions::eMass, mkey.GetShapeType(), *this,
268
269 DNekMatSharedPtr DiveMassmat = GetStdMatrix(stdmasskey);
270
271 (*WeakDirDeriv) = (*WeakDirDeriv) + (*DiveMassmat);
272
274 jac, WeakDirDeriv);
275 }
276 break;
277 }
279 {
280 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed ||
292 {
293 NekDouble one = 1.0;
294 DNekMatSharedPtr mat = GenMatrix(mkey);
295
296 returnval =
298 }
299 else
300 {
302 mkey.GetShapeType(), *this);
304 mkey.GetShapeType(), *this);
306 mkey.GetShapeType(), *this);
307
308 DNekMat &lap00 = *GetStdMatrix(lap00key);
309 DNekMat &lap01 = *GetStdMatrix(lap01key);
310 DNekMat &lap11 = *GetStdMatrix(lap11key);
311
312 NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
314 m_metricinfo->GetGmat(ptsKeys);
315
316 int rows = lap00.GetRows();
317 int cols = lap00.GetColumns();
318
319 DNekMatSharedPtr lap =
321
322 (*lap) = gmat[0][0] * lap00 +
323 gmat[1][0] * (lap01 + Transpose(lap01)) +
324 gmat[3][0] * lap11;
325
326 returnval =
328 }
329 }
330 break;
332 {
333 DNekMatSharedPtr mat = GenMatrix(mkey);
334
336 }
337 break;
339 {
341
342 MatrixKey masskey(mkey, StdRegions::eMass);
343 DNekScalMat &MassMat = *GetLocMatrix(masskey);
344
345 MatrixKey lapkey(mkey, StdRegions::eLaplacian);
346 DNekScalMat &LapMat = *GetLocMatrix(lapkey);
347
348 int rows = LapMat.GetRows();
349 int cols = LapMat.GetColumns();
350
351 DNekMatSharedPtr helm =
353
354 NekDouble one = 1.0;
355 (*helm) = LapMat + factor * MassMat;
356
357 returnval =
359 }
360 break;
362 {
363 MatrixKey helmkey(mkey, StdRegions::eHelmholtz);
364 DNekScalMat &HelmMat = *GetLocMatrix(helmkey);
365
366 // Generate a local copy of traceMat
369
371 "Need to specify eFactorGJP to construct "
372 "a HelmholtzGJP matrix");
373
375
376 factor /= HelmMat.Scale();
377
378 int ntot = HelmMat.GetRows() * HelmMat.GetColumns();
379
380 Vmath::Svtvp(ntot, factor, &NDTraceMat->GetPtr()[0], 1,
381 HelmMat.GetRawPtr(), 1, &NDTraceMat->GetPtr()[0], 1);
382
384 HelmMat.Scale(), NDTraceMat);
385 }
386 break;
388 {
390
391 // Construct mass matrix
392 // Check for mass-specific varcoeffs to avoid unncessary
393 // re-computation of the elemental matrix every time step
396 {
397 massVarcoeffs[StdRegions::eVarCoeffMass] =
399 }
400 MatrixKey masskey(StdRegions::eMass, mkey.GetShapeType(), *this,
401 mkey.GetConstFactors(), massVarcoeffs);
402 DNekScalMat &MassMat = *GetLocMatrix(masskey);
403
404 // Construct laplacian matrix (Check for varcoeffs)
405 // Take all varcoeffs if one or more are detected
406 // TODO We might want to have a map
407 // from MatrixType to Vector of Varcoeffs and vice-versa
419 {
420 lapVarcoeffs = mkey.GetVarCoeffs();
421 }
422 MatrixKey lapkey(StdRegions::eLaplacian, mkey.GetShapeType(), *this,
423 mkey.GetConstFactors(), lapVarcoeffs);
424 DNekScalMat &LapMat = *GetLocMatrix(lapkey);
425
426 // Construct advection matrix
427 // Check for varcoeffs not required;
428 // assume advection velocity is always time-dependent
430 DNekScalMat &AdvMat = *GetLocMatrix(advkey);
431
432 int rows = LapMat.GetRows();
433 int cols = LapMat.GetColumns();
434
435 DNekMatSharedPtr adr =
437
438 NekDouble one = 1.0;
439 (*adr) = LapMat - lambda * MassMat + AdvMat;
440
442
443 // Clear memory (Repeat varcoeff checks)
444 DropLocMatrix(advkey);
445 DropLocMatrix(masskey);
446 DropLocMatrix(lapkey);
447 }
448 break;
450 {
451 // Copied mostly from ADR solve to have fine-grain control
452 // over updating only advection matrix, relevant for performance!
454
455 // Construct mass matrix (Check for varcoeffs)
458 {
459 massVarcoeffs[StdRegions::eVarCoeffMass] =
461 }
462 MatrixKey masskey(StdRegions::eMass, mkey.GetShapeType(), *this,
463 mkey.GetConstFactors(), massVarcoeffs);
464 DNekScalMat &MassMat = *GetLocMatrix(masskey);
465
466 // Construct laplacian matrix (Check for varcoeffs)
478 {
479 lapVarcoeffs = mkey.GetVarCoeffs();
480 }
481 MatrixKey lapkey(StdRegions::eLaplacian, mkey.GetShapeType(), *this,
482 mkey.GetConstFactors(), lapVarcoeffs);
483 DNekScalMat &LapMat = *GetLocMatrix(lapkey);
484
485 // Construct advection matrix
486 // (assume advection velocity defined and non-zero)
487 // Could check L2(AdvectionVelocity) or HasVarCoeff
489 DNekScalMat &AdvMat = *GetLocMatrix(advkey);
490
491 // Generate a local copy of traceMat
493 *this, mkey.GetConstFactors());
494 DNekScalMat &NDTraceMat = *GetLocMatrix(gjpkey);
495
498 "Need to specify eFactorGJP to construct "
499 "a LinearAdvectionDiffusionReactionGJP matrix");
500
501 int rows = LapMat.GetRows();
502 int cols = LapMat.GetColumns();
503
504 DNekMatSharedPtr adr =
506
507 NekDouble one = 1.0;
508 (*adr) =
509 LapMat - lambda * MassMat + AdvMat + gjpfactor * NDTraceMat;
510
512
513 // Clear memory
514 DropLocMatrix(advkey);
515 DropLocMatrix(masskey);
516 DropLocMatrix(lapkey);
517 }
518 break;
520 {
521 NekDouble one = 1.0;
523
525 }
526 break;
528 {
529 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
530 {
531 NekDouble one = 1.0;
532 DNekMatSharedPtr mat = GenMatrix(mkey);
533
534 returnval =
536 }
537 else
538 {
539 NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
540 DNekMatSharedPtr mat = GetStdMatrix(mkey);
541
542 returnval =
544 }
545 }
546 break;
550 {
551 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
552 {
553 NekDouble one = 1.0;
554 DNekMatSharedPtr mat = GenMatrix(mkey);
555
556 returnval =
558 }
559 else
560 {
561 NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
562
564 m_metricinfo->GetDerivFactors(ptsKeys);
565 int dir = 0;
567 {
568 dir = 0;
569 }
571 {
572 dir = 1;
573 }
575 {
576 dir = 2;
577 }
578
580 mkey.GetShapeType(), *this);
582 mkey.GetShapeType(), *this);
583
584 DNekMat &stdiprod0 = *GetStdMatrix(iProdDeriv0Key);
585 DNekMat &stdiprod1 = *GetStdMatrix(iProdDeriv0Key);
586
587 int rows = stdiprod0.GetRows();
588 int cols = stdiprod1.GetColumns();
589
590 DNekMatSharedPtr mat =
592 (*mat) =
593 df[2 * dir][0] * stdiprod0 + df[2 * dir + 1][0] * stdiprod1;
594
595 returnval =
597 }
598 }
599 break;
600
602 {
603 NekDouble one = 1.0;
604
606 *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
607
608 DNekMatSharedPtr mat = GenMatrix(hkey);
609
610 mat->Invert();
611
613 }
614 break;
616 {
618 "Matrix only setup for quad elements currently");
619 DNekMatSharedPtr m_Ix;
620 Array<OneD, NekDouble> coords(1, 0.0);
622 int edge = static_cast<int>(factors[StdRegions::eFactorGaussEdge]);
623
624 coords[0] = (edge == 0 || edge == 3) ? -1.0 : 1.0;
625
626 m_Ix = m_base[(edge + 1) % 2]->GetI(coords);
627
628 returnval =
630 }
631 break;
633 {
634 NekDouble one = 1.0;
636 *this, mkey.GetConstFactors(),
637 mkey.GetVarCoeffs());
638 DNekScalBlkMatSharedPtr helmStatCond =
639 GetLocStaticCondMatrix(helmkey);
640 DNekScalMatSharedPtr A = helmStatCond->GetBlock(0, 0);
642
644 }
645 break;
646 default:
647 {
648 NekDouble one = 1.0;
649 DNekMatSharedPtr mat = GenMatrix(mkey);
650
652 }
653 break;
654 }
655
656 return returnval;
657}
658
660 const int edge, const ExpansionSharedPtr &EdgeExp,
663{
664 ASSERTL1(GetCoordim() == 2, "Routine only set up for two-dimensions");
665
667 GetTraceNormal(edge);
668
669 if (m_requireNeg.size() == 0)
670 {
671 int nedges = GetNtraces();
672 m_requireNeg.resize(nedges);
673
674 for (int i = 0; i < nedges; ++i)
675 {
676 m_requireNeg[i] = false;
677
678 ExpansionSharedPtr edgeExp = m_traceExp[i].lock();
679
680 if (edgeExp->GetRightAdjacentElementExp())
681 {
682 if (edgeExp->GetRightAdjacentElementExp()
683 ->GetGeom()
684 ->GetGlobalID() == GetGeom()->GetGlobalID())
685 {
686 m_requireNeg[i] = true;
687 }
688 }
689 }
690 }
691
692 // We allow the case of mixed polynomial order by supporting only
693 // those modes on the edge common to both adjoining elements. This
694 // is enforced here by taking the minimum size and padding with
695 // zeros.
696 int nquad_e = min(EdgeExp->GetNumPoints(0), int(normals[0].size()));
697
698 int nEdgePts = EdgeExp->GetTotPoints();
699 Array<OneD, NekDouble> edgePhys(nEdgePts);
700 Vmath::Vmul(nquad_e, normals[0], 1, Fx, 1, edgePhys, 1);
701 Vmath::Vvtvp(nquad_e, normals[1], 1, Fy, 1, edgePhys, 1, edgePhys, 1);
702
703 Expansion1DSharedPtr locExp = EdgeExp->as<Expansion1D>();
704
705 if (m_requireNeg[edge])
706 {
707 if (locExp->GetRightAdjacentElementExp()->GetGeom()->GetGlobalID() ==
708 m_geom->GetGlobalID())
709 {
710 Vmath::Neg(nquad_e, edgePhys, 1);
711 }
712 }
713
714 AddEdgeNormBoundaryInt(edge, EdgeExp, edgePhys, outarray);
715}
716
718 const int edge, const ExpansionSharedPtr &EdgeExp,
720{
721 int i;
722
723 if (m_requireNeg.size() == 0)
724 {
725 int nedges = GetNtraces();
726 m_requireNeg.resize(nedges);
727
728 for (i = 0; i < nedges; ++i)
729 {
730 m_requireNeg[i] = false;
731
732 ExpansionSharedPtr edgeExp = m_traceExp[i].lock();
733
734 if (edgeExp->GetRightAdjacentElementExp())
735 {
736 if (edgeExp->GetRightAdjacentElementExp()
737 ->GetGeom()
738 ->GetGlobalID() == GetGeom()->GetGlobalID())
739 {
740 m_requireNeg[i] = true;
741 }
742 }
743 }
744 }
745
747 GetBasisNumModes(1), 0, edge, GetTraceOrient(edge));
748
750
751 // Order of the element
752 int order_e = map->size();
753 // Order of the trace
754 int n_coeffs = EdgeExp->GetNcoeffs();
755
756 Array<OneD, NekDouble> edgeCoeffs(n_coeffs);
757 if (n_coeffs != order_e) // Going to orthogonal space
758 {
759 EdgeExp->FwdTrans(Fn, edgeCoeffs);
760 Expansion1DSharedPtr locExp = EdgeExp->as<Expansion1D>();
761
762 if (m_requireNeg[edge])
763 {
764 Vmath::Neg(n_coeffs, edgeCoeffs, 1);
765 }
766
767 Array<OneD, NekDouble> coeff(n_coeffs, 0.0);
769 ((LibUtilities::BasisType)1); // 1-->Ortho_A
770 LibUtilities::BasisKey bkey_ortho(btype,
771 EdgeExp->GetBasis(0)->GetNumModes(),
772 EdgeExp->GetBasis(0)->GetPointsKey());
773 LibUtilities::BasisKey bkey(EdgeExp->GetBasis(0)->GetBasisType(),
774 EdgeExp->GetBasis(0)->GetNumModes(),
775 EdgeExp->GetBasis(0)->GetPointsKey());
776 LibUtilities::InterpCoeff1D(bkey, edgeCoeffs, bkey_ortho, coeff);
777
778 // Cutting high frequencies
779 for (i = order_e; i < n_coeffs; i++)
780 {
781 coeff[i] = 0.0;
782 }
783
784 LibUtilities::InterpCoeff1D(bkey_ortho, coeff, bkey, edgeCoeffs);
785
787 LibUtilities::eSegment, *EdgeExp);
788 EdgeExp->MassMatrixOp(edgeCoeffs, edgeCoeffs, masskey);
789 }
790 else
791 {
792 EdgeExp->IProductWRTBase(Fn, edgeCoeffs);
793
794 Expansion1DSharedPtr locExp = EdgeExp->as<Expansion1D>();
795
796 if (m_requireNeg[edge])
797 {
798 Vmath::Neg(n_coeffs, edgeCoeffs, 1);
799 }
800 }
801
802 // Implementation for all the basis except Gauss points
803 if (EdgeExp->GetBasis(0)->GetBasisType() != LibUtilities::eGauss_Lagrange)
804 {
805 // add data to outarray if forward edge normal is outwards
806 for (i = 0; i < order_e; ++i)
807 {
808 outarray[(*map)[i].index] += (*map)[i].sign * edgeCoeffs[i];
809 }
810 }
811 else
812 {
813 int nCoeffs0, nCoeffs1;
814 int j;
815
819 *this, factors);
820
821 DNekMatSharedPtr mat_gauss = m_stdMatrixManager[key];
822
823 switch (edge)
824 {
825 case 0:
826 {
827 nCoeffs1 = m_base[1]->GetNumModes();
828
829 for (i = 0; i < order_e; ++i)
830 {
831 for (j = 0; j < nCoeffs1; j++)
832 {
833 outarray[(*map)[i].index + j * order_e] +=
834 mat_gauss->GetPtr()[j] * (*map)[i].sign *
835 edgeCoeffs[i];
836 }
837 }
838 break;
839 }
840 case 1:
841 {
842 nCoeffs0 = m_base[0]->GetNumModes();
843
844 for (i = 0; i < order_e; ++i)
845 {
846 for (j = 0; j < nCoeffs0; j++)
847 {
848 outarray[(*map)[i].index - j] +=
849 mat_gauss->GetPtr()[order_e - 1 - j] *
850 (*map)[i].sign * edgeCoeffs[i];
851 }
852 }
853 break;
854 }
855 case 2:
856 {
857 nCoeffs1 = m_base[1]->GetNumModes();
858
859 for (i = 0; i < order_e; ++i)
860 {
861 for (j = 0; j < nCoeffs1; j++)
862 {
863 outarray[(*map)[i].index - j * order_e] +=
864 mat_gauss->GetPtr()[order_e - 1 - j] *
865 (*map)[i].sign * edgeCoeffs[i];
866 }
867 }
868 break;
869 }
870 case 3:
871 {
872 nCoeffs0 = m_base[0]->GetNumModes();
873
874 for (i = 0; i < order_e; ++i)
875 {
876 for (j = 0; j < nCoeffs0; j++)
877 {
878 outarray[(*map)[i].index + j] +=
879 mat_gauss->GetPtr()[j] * (*map)[i].sign *
880 edgeCoeffs[i];
881 }
882 }
883 break;
884 }
885 default:
886 ASSERTL0(false, "edge value (< 3) is out of range");
887 break;
888 }
889 }
890}
891
894{
895 int i, cnt = 0;
896 int nedges = GetNtraces();
898
899 for (i = 0; i < nedges; ++i)
900 {
901 EdgeExp[i]->SetCoeffsToOrientation(
902 GetTraceOrient(i), e_tmp = inout + cnt, e_tmp = inout + cnt);
903 cnt += GetTraceNcoeffs(i);
904 }
905}
906
907/**
908 * Computes the C matrix entries due to the presence of the identity
909 * matrix in Eqn. 32.
910 */
914 Array<OneD, NekDouble> &outarray,
915 const StdRegions::VarCoeffMap &varcoeffs)
916{
917 int i, e, cnt;
918 int order_e, nquad_e;
919 int nedges = GetNtraces();
920
921 cnt = 0;
922 for (e = 0; e < nedges; ++e)
923 {
924 order_e = EdgeExp[e]->GetNcoeffs();
925 nquad_e = EdgeExp[e]->GetNumPoints(0);
926
929 Array<OneD, NekDouble> edgeCoeffs(order_e);
930 Array<OneD, NekDouble> edgePhys(nquad_e);
931
932 for (i = 0; i < order_e; ++i)
933 {
934 edgeCoeffs[i] = inarray[i + cnt];
935 }
936 cnt += order_e;
937
938 EdgeExp[e]->BwdTrans(edgeCoeffs, edgePhys);
939
940 // Multiply by variable coefficient
941 /// @TODO: Document this
942 // StdRegions::VarCoeffType VarCoeff[3] = {StdRegions::eVarCoeffD00,
943 // StdRegions::eVarCoeffD11,
944 // StdRegions::eVarCoeffD22};
945 // StdRegions::VarCoeffMap::const_iterator x;
946 // Array<OneD, NekDouble> varcoeff_work(nquad_e);
947
948 // if ((x = varcoeffs.find(VarCoeff[dir])) != varcoeffs.end())
949 // {
950 // GetPhysEdgeVarCoeffsFromElement(e,EdgeExp[e],x->second,varcoeff_work);
951 // Vmath::Vmul(nquad_e,varcoeff_work,1,EdgeExp[e]->GetPhys(),1,EdgeExp[e]->UpdatePhys(),1);
952 // }
953
954 if (varcoeffs.find(StdRegions::eVarCoeffMF1x) != varcoeffs.end())
955 {
956 // MMF case
957 Array<OneD, NekDouble> ncdotMF_e =
958 GetnEdgecdotMF(dir, e, EdgeExp[e], normals, varcoeffs);
959
960 Vmath::Vmul(nquad_e, ncdotMF_e, 1, edgePhys, 1, edgePhys, 1);
961 }
962 else
963 {
964 Vmath::Vmul(nquad_e, normals[dir], 1, edgePhys, 1, edgePhys, 1);
965 }
966
967 AddEdgeBoundaryInt(e, EdgeExp[e], edgePhys, outarray, varcoeffs);
968 }
969}
970
972 const int dir, Array<OneD, ExpansionSharedPtr> &EdgeExp,
973 Array<OneD, Array<OneD, NekDouble>> &edgeCoeffs,
974 Array<OneD, NekDouble> &outarray)
975{
976 int e;
977 int nquad_e;
978 int nedges = GetNtraces();
979
980 for (e = 0; e < nedges; ++e)
981 {
982 nquad_e = EdgeExp[e]->GetNumPoints(0);
983
984 Array<OneD, NekDouble> edgePhys(nquad_e);
987
988 EdgeExp[e]->BwdTrans(edgeCoeffs[e], edgePhys);
989
990 Vmath::Vmul(nquad_e, normals[dir], 1, edgePhys, 1, edgePhys, 1);
991
992 AddEdgeBoundaryInt(e, EdgeExp[e], edgePhys, outarray);
993 }
994}
995
996/**
997 * For a given edge add the \tilde{F}_1j contributions
998 */
1000 ExpansionSharedPtr &EdgeExp,
1001 Array<OneD, NekDouble> &edgePhys,
1002 Array<OneD, NekDouble> &outarray,
1003 const StdRegions::VarCoeffMap &varcoeffs)
1004{
1005 int i;
1006 int order_e = EdgeExp->GetNcoeffs();
1007 int nquad_e = EdgeExp->GetNumPoints(0);
1010 Array<OneD, NekDouble> coeff(order_e);
1011
1012 GetTraceToElementMap(edge, map, sign, v_GetTraceOrient(edge));
1013
1017 StdRegions::VarCoeffMap::const_iterator x;
1018
1019 /// @TODO Variable coeffs
1020 if ((x = varcoeffs.find(VarCoeff[0])) != varcoeffs.end())
1021 {
1022 Array<OneD, NekDouble> work(nquad_e);
1023 GetPhysEdgeVarCoeffsFromElement(edge, EdgeExp, x->second.GetValue(),
1024 work);
1025 Vmath::Vmul(nquad_e, work, 1, edgePhys, 1, edgePhys, 1);
1026 }
1027
1028 EdgeExp->IProductWRTBase(edgePhys, coeff);
1029
1030 // add data to out array
1031 for (i = 0; i < order_e; ++i)
1032 {
1033 outarray[map[i]] += sign[i] * coeff[i];
1034 }
1035}
1036
1037// This method assumes that data in EdgeExp is orientated according to
1038// elemental counter clockwise format AddHDGHelmholtzTraceTerms with
1039// directions
1041 const NekDouble tau, const Array<OneD, const NekDouble> &inarray,
1043 const StdRegions::VarCoeffMap &dirForcing, Array<OneD, NekDouble> &outarray)
1044{
1045 ASSERTL0(&inarray[0] != &outarray[0],
1046 "Input and output arrays use the same memory");
1047
1048 int e, cnt, order_e, nedges = GetNtraces();
1050
1051 cnt = 0;
1052
1053 for (e = 0; e < nedges; ++e)
1054 {
1055 order_e = EdgeExp[e]->GetNcoeffs();
1056 Array<OneD, NekDouble> edgeCoeffs(order_e);
1057 Array<OneD, NekDouble> edgePhys(EdgeExp[e]->GetTotPoints());
1058
1059 Vmath::Vcopy(order_e, tmp = inarray + cnt, 1, edgeCoeffs, 1);
1060 EdgeExp[e]->BwdTrans(edgeCoeffs, edgePhys);
1061 AddHDGHelmholtzEdgeTerms(tau, e, EdgeExp, edgePhys, dirForcing,
1062 outarray);
1063
1064 cnt += order_e;
1065 }
1066}
1067
1068// evaluate additional terms in HDG edges. Not that this assumes that
1069// edges are unpacked into local cartesian order.
1071 const NekDouble tau, const int edge,
1073 const StdRegions::VarCoeffMap &varcoeffs, Array<OneD, NekDouble> &outarray)
1074{
1075 bool mmf = (varcoeffs.find(StdRegions::eVarCoeffMF1x) != varcoeffs.end());
1076 int i, j, n;
1077 int nquad_e = EdgeExp[edge]->GetNumPoints(0);
1078 int order_e = EdgeExp[edge]->GetNcoeffs();
1079 int coordim = mmf ? 2 : GetCoordim();
1080 int ncoeffs = GetNcoeffs();
1081
1082 Array<OneD, NekDouble> inval(nquad_e);
1083 Array<OneD, NekDouble> outcoeff(order_e);
1084 Array<OneD, NekDouble> tmpcoeff(ncoeffs);
1085
1087 GetTraceNormal(edge);
1088
1091
1093
1094 DNekVec Coeffs(ncoeffs, outarray, eWrapper);
1095 DNekVec Tmpcoeff(ncoeffs, tmpcoeff, eWrapper);
1096
1097 GetTraceToElementMap(edge, emap, sign, edgedir);
1098
1102
1106
1107 StdRegions::VarCoeffMap::const_iterator x;
1108 /// @TODO: What direction to use here??
1109 if ((x = varcoeffs.find(VarCoeff[0])) != varcoeffs.end())
1110 {
1111 Array<OneD, NekDouble> work(nquad_e);
1112 GetPhysEdgeVarCoeffsFromElement(edge, EdgeExp[edge],
1113 x->second.GetValue(), work);
1114 Vmath::Vmul(nquad_e, work, 1, edgePhys, 1, edgePhys, 1);
1115 }
1116
1117 //================================================================
1118 // Add F = \tau <phi_i,in_phys>
1119 // Fill edge and take inner product
1120 EdgeExp[edge]->IProductWRTBase(edgePhys, outcoeff);
1121 // add data to out array
1122 for (i = 0; i < order_e; ++i)
1123 {
1124 outarray[emap[i]] += sign[i] * tau * outcoeff[i];
1125 }
1126 //================================================================
1127
1128 //===============================================================
1129 // Add -\sum_i D_i^T M^{-1} G_i + E_i M^{-1} G_i =
1130 // \sum_i D_i M^{-1} G_i term
1131
1132 // Two independent direction
1133 DNekScalMatSharedPtr invMass;
1134 for (n = 0; n < coordim; ++n)
1135 {
1136 if (mmf)
1137 {
1139 Weight[StdRegions::eVarCoeffMass] = GetMFMag(n, varcoeffs);
1140
1141 MatrixKey invMasskey(StdRegions::eInvMass, DetShapeType(), *this,
1143
1144 invMass = GetLocMatrix(invMasskey);
1145
1146 Array<OneD, NekDouble> ncdotMF_e =
1147 GetnEdgecdotMF(n, edge, EdgeExp[edge], normals, varcoeffs);
1148
1149 Vmath::Vmul(nquad_e, ncdotMF_e, 1, edgePhys, 1, inval, 1);
1150 }
1151 else
1152 {
1153 Vmath::Vmul(nquad_e, normals[n], 1, edgePhys, 1, inval, 1);
1155 }
1156
1157 // Multiply by variable coefficient
1158 /// @TODO: Document this (probably not needed)
1159 // StdRegions::VarCoeffMap::const_iterator x;
1160 // if ((x = varcoeffs.find(VarCoeff[n])) !=
1161 // varcoeffs.end())
1162 // {
1163 // GetPhysEdgeVarCoeffsFromElement(edge,EdgeExp[edge],x->second,varcoeff_work);
1164 // Vmath::Vmul(nquad_e,varcoeff_work,1,EdgeExp[edge]->GetPhys(),1,EdgeExp[edge]->UpdatePhys(),1);
1165 // }
1166
1167 EdgeExp[edge]->IProductWRTBase(inval, outcoeff);
1168
1169 // M^{-1} G
1170 for (i = 0; i < ncoeffs; ++i)
1171 {
1172 tmpcoeff[i] = 0;
1173 for (j = 0; j < order_e; ++j)
1174 {
1175 tmpcoeff[i] += (*invMass)(i, emap[j]) * sign[j] * outcoeff[j];
1176 }
1177 }
1178
1179 if (mmf)
1180 {
1181 StdRegions::VarCoeffMap VarCoeffDirDeriv;
1182 VarCoeffDirDeriv[StdRegions::eVarCoeffMF] =
1183 GetMF(n, coordim, varcoeffs);
1184 VarCoeffDirDeriv[StdRegions::eVarCoeffMFDiv] =
1185 GetMFDiv(n, varcoeffs);
1186
1189 VarCoeffDirDeriv);
1190
1191 DNekScalMat &Dmat = *GetLocMatrix(Dmatkey);
1192
1193 Coeffs = Coeffs + Dmat * Tmpcoeff;
1194 }
1195 else
1196 {
1197 if (varcoeffs.find(VarCoeff[n]) != varcoeffs.end())
1198 {
1199 MatrixKey mkey(DerivType[n], DetShapeType(), *this,
1201
1202 DNekScalMat &Dmat = *GetLocMatrix(mkey);
1203 Coeffs = Coeffs + Dmat * Tmpcoeff;
1204 }
1205 else
1206 {
1207 DNekScalMat &Dmat = *GetLocMatrix(DerivType[n]);
1208 Coeffs = Coeffs + Dmat * Tmpcoeff;
1209 }
1210 }
1211 }
1212}
1213
1214/**
1215 * Extracts the variable coefficients along an edge
1216 */
1218 const int edge, ExpansionSharedPtr &EdgeExp,
1219 const Array<OneD, const NekDouble> &varcoeff,
1220 Array<OneD, NekDouble> &outarray)
1221{
1223 Array<OneD, NekDouble> edgetmp(EdgeExp->GetNcoeffs());
1224
1225 // FwdTrans varcoeffs
1226 FwdTrans(varcoeff, tmp);
1227
1228 // Map to edge
1232 GetTraceToElementMap(edge, emap, sign, edgedir);
1233
1234 for (unsigned int i = 0; i < EdgeExp->GetNcoeffs(); ++i)
1235 {
1236 edgetmp[i] = tmp[emap[i]];
1237 }
1238
1239 // BwdTrans
1240 EdgeExp->BwdTrans(edgetmp, outarray);
1241}
1242
1243/**
1244 * Computes matrices needed for the HDG formulation. References to
1245 * equations relate to the following paper:
1246 * R. M. Kirby, S. J. Sherwin, B. Cockburn, To CG or to HDG: A
1247 * Comparative Study, J. Sci. Comp P1-30
1248 * DOI 10.1007/s10915-011-9501-7
1249 */
1251{
1252 DNekMatSharedPtr returnval;
1253
1254 switch (mkey.GetMatrixType())
1255 {
1256 // (Z^e)^{-1} (Eqn. 33, P22)
1258 {
1260 "HybridDGHelmholtz matrix not set up "
1261 "for non boundary-interior expansions");
1262
1263 int i, j, k;
1264 NekDouble lambdaval =
1267 int ncoeffs = GetNcoeffs();
1268 int nedges = GetNtraces();
1269 int shapedim = 2;
1270 const StdRegions::VarCoeffMap &varcoeffs = mkey.GetVarCoeffs();
1271 bool mmf =
1272 (varcoeffs.find(StdRegions::eVarCoeffMF1x) != varcoeffs.end());
1273
1277 ExpansionSharedPtr EdgeExp;
1278
1279 int order_e, coordim = GetCoordim();
1284 DNekMat LocMat(ncoeffs, ncoeffs);
1285
1286 returnval =
1288 DNekMat &Mat = *returnval;
1289 Vmath::Zero(ncoeffs * ncoeffs, Mat.GetPtr(), 1);
1290
1294
1295 StdRegions::VarCoeffMap::const_iterator x;
1296
1297 for (i = 0; i < coordim; ++i)
1298 {
1299 if (mmf)
1300 {
1301 if (i < shapedim)
1302 {
1303 StdRegions::VarCoeffMap VarCoeffDirDeriv;
1304 VarCoeffDirDeriv[StdRegions::eVarCoeffMF] =
1305 GetMF(i, shapedim, varcoeffs);
1306 VarCoeffDirDeriv[StdRegions::eVarCoeffMFDiv] =
1307 GetMFDiv(i, varcoeffs);
1308
1310 DetShapeType(), *this,
1312 VarCoeffDirDeriv);
1313
1314 DNekScalMat &Dmat = *GetLocMatrix(Dmatkey);
1315
1318 GetMFMag(i, mkey.GetVarCoeffs());
1319
1320 MatrixKey invMasskey(
1323
1324 DNekScalMat &invMass = *GetLocMatrix(invMasskey);
1325
1326 Mat = Mat + Dmat * invMass * Transpose(Dmat);
1327 }
1328 }
1329 else if (mkey.HasVarCoeff(Coeffs[i]))
1330 {
1331 MatrixKey DmatkeyL(DerivType[i], DetShapeType(), *this,
1333 mkey.GetVarCoeffAsMap(Coeffs[i]));
1334
1335 MatrixKey DmatkeyR(DerivType[i], DetShapeType(), *this);
1336
1337 DNekScalMat &DmatL = *GetLocMatrix(DmatkeyL);
1338 DNekScalMat &DmatR = *GetLocMatrix(DmatkeyR);
1339 Mat = Mat + DmatL * invMass * Transpose(DmatR);
1340 }
1341 else
1342 {
1343 DNekScalMat &Dmat = *GetLocMatrix(DerivType[i]);
1344 Mat = Mat + Dmat * invMass * Transpose(Dmat);
1345 }
1346 }
1347
1348 // Add Mass Matrix Contribution for Helmholtz problem
1350 Mat = Mat + lambdaval * Mass;
1351
1352 // Add tau*E_l using elemental mass matrices on each edge
1353 for (i = 0; i < nedges; ++i)
1354 {
1355 EdgeExp = GetTraceExp(i);
1356 order_e = EdgeExp->GetNcoeffs();
1357
1358 int nq = EdgeExp->GetNumPoints(0);
1359 GetTraceToElementMap(i, emap, sign, edgedir);
1360
1361 // @TODO: Document
1362 StdRegions::VarCoeffMap edgeVarCoeffs;
1364 {
1367 i, EdgeExp, mkey.GetVarCoeff(StdRegions::eVarCoeffD00),
1368 mu);
1369 edgeVarCoeffs[StdRegions::eVarCoeffMass] = mu;
1370 }
1371 DNekScalMat &eMass = *EdgeExp->GetLocMatrix(
1373 edgeVarCoeffs);
1374 // DNekScalMat &eMass =
1375 // *EdgeExp->GetLocMatrix(StdRegions::eMass);
1376
1377 for (j = 0; j < order_e; ++j)
1378 {
1379 for (k = 0; k < order_e; ++k)
1380 {
1381 Mat(emap[j], emap[k]) =
1382 Mat(emap[j], emap[k]) +
1383 tau * sign[j] * sign[k] * eMass(j, k);
1384 }
1385 }
1386 }
1387 }
1388 break;
1389 // U^e (P22)
1391 {
1392 int i, j, k;
1393 int nbndry = NumDGBndryCoeffs();
1394 int ncoeffs = GetNcoeffs();
1395 int nedges = GetNtraces();
1397
1398 Array<OneD, NekDouble> lambda(nbndry);
1399 DNekVec Lambda(nbndry, lambda, eWrapper);
1400 Array<OneD, NekDouble> ulam(ncoeffs);
1401 DNekVec Ulam(ncoeffs, ulam, eWrapper);
1402 Array<OneD, NekDouble> f(ncoeffs);
1403 DNekVec F(ncoeffs, f, eWrapper);
1404
1405 Array<OneD, ExpansionSharedPtr> EdgeExp(nedges);
1406 // declare matrix space
1407 returnval =
1409 DNekMat &Umat = *returnval;
1410
1411 // Z^e matrix
1413 *this, mkey.GetConstFactors(),
1414 mkey.GetVarCoeffs());
1415 DNekScalMat &invHmat = *GetLocMatrix(newkey);
1416
1419
1420 for (i = 0; i < nedges; ++i)
1421 {
1422 EdgeExp[i] = GetTraceExp(i);
1423 }
1424
1425 // for each degree of freedom of the lambda space
1426 // calculate Umat entry
1427 // Generate Lambda to U_lambda matrix
1428 for (j = 0; j < nbndry; ++j)
1429 {
1430 // standard basis vectors e_j
1431 Vmath::Zero(nbndry, &lambda[0], 1);
1432 Vmath::Zero(ncoeffs, &f[0], 1);
1433 lambda[j] = 1.0;
1434
1435 SetTraceToGeomOrientation(EdgeExp, lambda);
1436
1437 // Compute F = [I D_1 M^{-1} D_2 M^{-1}] C e_j
1438 AddHDGHelmholtzTraceTerms(tau, lambda, EdgeExp,
1439 mkey.GetVarCoeffs(), f);
1440
1441 // Compute U^e_j
1442 Ulam = invHmat * F; // generate Ulam from lambda
1443
1444 // fill column of matrix
1445 for (k = 0; k < ncoeffs; ++k)
1446 {
1447 Umat(k, j) = Ulam[k];
1448 }
1449 }
1450 }
1451 break;
1452 // Q_0, Q_1, Q_2 matrices (P23)
1453 // Each are a product of a row of Eqn 32 with the C matrix.
1454 // Rather than explicitly computing all of Eqn 32, we note each
1455 // row is almost a multiple of U^e, so use that as our starting
1456 // point.
1460 {
1461 int i = 0;
1462 int j = 0;
1463 int k = 0;
1464 int dir = 0;
1465 int nbndry = NumDGBndryCoeffs();
1466 int ncoeffs = GetNcoeffs();
1467 int nedges = GetNtraces();
1468 int shapedim = 2;
1469
1470 Array<OneD, NekDouble> lambda(nbndry);
1471 DNekVec Lambda(nbndry, lambda, eWrapper);
1472 Array<OneD, ExpansionSharedPtr> EdgeExp(nedges);
1473
1474 Array<OneD, NekDouble> ulam(ncoeffs);
1475 DNekVec Ulam(ncoeffs, ulam, eWrapper);
1476 Array<OneD, NekDouble> f(ncoeffs);
1477 DNekVec F(ncoeffs, f, eWrapper);
1478
1479 // declare matrix space
1480 returnval =
1482 DNekMat &Qmat = *returnval;
1483
1484 // Lambda to U matrix
1486 *this, mkey.GetConstFactors(),
1487 mkey.GetVarCoeffs());
1488 DNekScalMat &lamToU = *GetLocMatrix(lamToUkey);
1489
1490 // Inverse mass matrix
1492
1493 for (i = 0; i < nedges; ++i)
1494 {
1495 EdgeExp[i] = GetTraceExp(i);
1496 }
1497
1498 // Weak Derivative matrix
1500 switch (mkey.GetMatrixType())
1501 {
1503 dir = 0;
1505 break;
1507 dir = 1;
1509 break;
1511 dir = 2;
1513 break;
1514 default:
1515 ASSERTL0(false, "Direction not known");
1516 break;
1517 }
1518
1519 const StdRegions::VarCoeffMap &varcoeffs = mkey.GetVarCoeffs();
1520 if (varcoeffs.find(StdRegions::eVarCoeffMF1x) != varcoeffs.end())
1521 {
1522 StdRegions::VarCoeffMap VarCoeffDirDeriv;
1523 VarCoeffDirDeriv[StdRegions::eVarCoeffMF] =
1524 GetMF(dir, shapedim, varcoeffs);
1525 VarCoeffDirDeriv[StdRegions::eVarCoeffMFDiv] =
1526 GetMFDiv(dir, varcoeffs);
1527
1528 MatrixKey Dmatkey(
1530 StdRegions::NullConstFactorMap, VarCoeffDirDeriv);
1531
1532 Dmat = GetLocMatrix(Dmatkey);
1533
1536 GetMFMag(dir, mkey.GetVarCoeffs());
1537
1540 Weight);
1541
1542 invMass = *GetLocMatrix(invMasskey);
1543 }
1544 else
1545 {
1549
1550 Dmat = GetLocMatrix(DerivType[dir]);
1551
1553 *this);
1554 invMass = *GetLocMatrix(invMasskey);
1555 }
1556
1557 // for each degree of freedom of the lambda space
1558 // calculate Qmat entry
1559 // Generate Lambda to Q_lambda matrix
1560 for (j = 0; j < nbndry; ++j)
1561 {
1562 Vmath::Zero(nbndry, &lambda[0], 1);
1563 lambda[j] = 1.0;
1564
1565 // for lambda[j] = 1 this is the solution to ulam
1566 for (k = 0; k < ncoeffs; ++k)
1567 {
1568 Ulam[k] = lamToU(k, j);
1569 }
1570
1571 // -D^T ulam
1572 Vmath::Neg(ncoeffs, &ulam[0], 1);
1573 F = Transpose(*Dmat) * Ulam;
1574
1575 SetTraceToGeomOrientation(EdgeExp, lambda);
1576
1577 // Add the C terms resulting from the I's on the
1578 // diagonals of Eqn 32
1579 AddNormTraceInt(dir, lambda, EdgeExp, f, mkey.GetVarCoeffs());
1580
1581 // finally multiply by inverse mass matrix
1582 Ulam = invMass * F;
1583
1584 // fill column of matrix (Qmat is in column major format)
1585 Vmath::Vcopy(ncoeffs, &ulam[0], 1,
1586 &(Qmat.GetPtr())[0] + j * ncoeffs, 1);
1587 }
1588 }
1589 break;
1590 // Matrix K (P23)
1592 {
1593 int i, j, e, cnt;
1594 int order_e, nquad_e;
1595 int nbndry = NumDGBndryCoeffs();
1596 int coordim = GetCoordim();
1597 int nedges = GetNtraces();
1599 StdRegions::VarCoeffMap::const_iterator x;
1600 const StdRegions::VarCoeffMap &varcoeffs = mkey.GetVarCoeffs();
1601 bool mmf =
1602 (varcoeffs.find(StdRegions::eVarCoeffMF1x) != varcoeffs.end());
1603
1604 Array<OneD, NekDouble> work, varcoeff_work;
1606 Array<OneD, ExpansionSharedPtr> EdgeExp(nedges);
1607 Array<OneD, NekDouble> lam(nbndry);
1608
1612
1613 // declare matrix space
1614 returnval =
1616 DNekMat &BndMat = *returnval;
1617
1618 DNekScalMatSharedPtr LamToQ[3];
1619
1620 // Matrix to map Lambda to U
1622 *this, mkey.GetConstFactors(),
1623 mkey.GetVarCoeffs());
1624 DNekScalMat &LamToU = *GetLocMatrix(LamToUkey);
1625
1626 // Matrix to map Lambda to Q0
1628 *this, mkey.GetConstFactors(),
1629 mkey.GetVarCoeffs());
1630 LamToQ[0] = GetLocMatrix(LamToQ0key);
1631
1632 // Matrix to map Lambda to Q1
1634 *this, mkey.GetConstFactors(),
1635 mkey.GetVarCoeffs());
1636 LamToQ[1] = GetLocMatrix(LamToQ1key);
1637
1638 // Matrix to map Lambda to Q2 for 3D coordinates
1639 if (coordim == 3)
1640 {
1641 MatrixKey LamToQ2key(
1643 mkey.GetConstFactors(), mkey.GetVarCoeffs());
1644 LamToQ[2] = GetLocMatrix(LamToQ2key);
1645 }
1646
1647 // Set up edge segment expansions from local geom info
1648 for (i = 0; i < nedges; ++i)
1649 {
1650 EdgeExp[i] = GetTraceExp(i);
1651 }
1652
1653 // Set up matrix derived from <mu, Q_lam.n - \tau (U_lam - Lam) >
1654 for (i = 0; i < nbndry; ++i)
1655 {
1656 cnt = 0;
1657
1658 Vmath::Zero(nbndry, lam, 1);
1659 lam[i] = 1.0;
1660 SetTraceToGeomOrientation(EdgeExp, lam);
1661
1662 for (e = 0; e < nedges; ++e)
1663 {
1664 order_e = EdgeExp[e]->GetNcoeffs();
1665 nquad_e = EdgeExp[e]->GetNumPoints(0);
1666
1667 normals = GetTraceNormal(e);
1668 edgedir = GetTraceOrient(e);
1669
1670 work = Array<OneD, NekDouble>(nquad_e);
1671 varcoeff_work = Array<OneD, NekDouble>(nquad_e);
1672
1673 GetTraceToElementMap(e, emap, sign, edgedir);
1674
1675 StdRegions::VarCoeffType VarCoeff[3] = {
1678
1679 // Q0 * n0 (BQ_0 terms)
1680 Array<OneD, NekDouble> edgeCoeffs(order_e);
1681 Array<OneD, NekDouble> edgePhys(nquad_e);
1682 for (j = 0; j < order_e; ++j)
1683 {
1684 edgeCoeffs[j] = sign[j] * (*LamToQ[0])(emap[j], i);
1685 }
1686
1687 EdgeExp[e]->BwdTrans(edgeCoeffs, edgePhys);
1688 // @TODO Var coeffs
1689 // Multiply by variable coefficient
1690 // if ((x =
1691 // varcoeffs.find(VarCoeff[0]))
1692 // != varcoeffs.end())
1693 // {
1694 // GetPhysEdgeVarCoeffsFromElement(e,EdgeExp[e],x->second,varcoeff_work);
1695 // Vmath::Vmul(nquad_e,varcoeff_work,1,EdgeExp[e]->GetPhys(),1,EdgeExp[e]->UpdatePhys(),1);
1696 // }
1697 if (mmf)
1698 {
1700 0, e, EdgeExp[e], normals, varcoeffs);
1701 Vmath::Vmul(nquad_e, ncdotMF, 1, edgePhys, 1, work, 1);
1702 }
1703 else
1704 {
1705 Vmath::Vmul(nquad_e, normals[0], 1, edgePhys, 1, work,
1706 1);
1707 }
1708
1709 // Q1 * n1 (BQ_1 terms)
1710 for (j = 0; j < order_e; ++j)
1711 {
1712 edgeCoeffs[j] = sign[j] * (*LamToQ[1])(emap[j], i);
1713 }
1714
1715 EdgeExp[e]->BwdTrans(edgeCoeffs, edgePhys);
1716
1717 // @TODO var coeffs
1718 // Multiply by variable coefficients
1719 // if ((x =
1720 // varcoeffs.find(VarCoeff[1]))
1721 // != varcoeffs.end())
1722 // {
1723 // GetPhysEdgeVarCoeffsFromElement(e,EdgeExp[e],x->second,varcoeff_work);
1724 // Vmath::Vmul(nquad_e,varcoeff_work,1,EdgeExp[e]->GetPhys(),1,EdgeExp[e]->UpdatePhys(),1);
1725 // }
1726
1727 if (mmf)
1728 {
1730 1, e, EdgeExp[e], normals, varcoeffs);
1731 Vmath::Vvtvp(nquad_e, ncdotMF, 1, edgePhys, 1, work, 1,
1732 work, 1);
1733 }
1734 else
1735 {
1736 Vmath::Vvtvp(nquad_e, normals[1], 1, edgePhys, 1, work,
1737 1, work, 1);
1738 }
1739
1740 // Q2 * n2 (BQ_2 terms)
1741 if (coordim == 3)
1742 {
1743 for (j = 0; j < order_e; ++j)
1744 {
1745 edgeCoeffs[j] = sign[j] * (*LamToQ[2])(emap[j], i);
1746 }
1747
1748 EdgeExp[e]->BwdTrans(edgeCoeffs, edgePhys);
1749 // @TODO var coeffs
1750 // Multiply by variable coefficients
1751 // if ((x =
1752 // varcoeffs.find(VarCoeff[2]))
1753 // != varcoeffs.end())
1754 // {
1755 // GetPhysEdgeVarCoeffsFromElement(e,EdgeExp[e],x->second,varcoeff_work);
1756 // Vmath::Vmul(nquad_e,varcoeff_work,1,EdgeExp[e]->GetPhys(),1,EdgeExp[e]->UpdatePhys(),1);
1757 // }
1758
1759 Vmath::Vvtvp(nquad_e, normals[2], 1, edgePhys, 1, work,
1760 1, work, 1);
1761 }
1762
1763 // - tau (ulam - lam)
1764 // Corresponds to the G and BU terms.
1765 for (j = 0; j < order_e; ++j)
1766 {
1767 edgeCoeffs[j] =
1768 sign[j] * LamToU(emap[j], i) - lam[cnt + j];
1769 }
1770
1771 EdgeExp[e]->BwdTrans(edgeCoeffs, edgePhys);
1772
1773 // Multiply by variable coefficients
1774 if ((x = varcoeffs.find(VarCoeff[0])) != varcoeffs.end())
1775 {
1777 e, EdgeExp[e], x->second.GetValue(), varcoeff_work);
1778 Vmath::Vmul(nquad_e, varcoeff_work, 1, edgePhys, 1,
1779 edgePhys, 1);
1780 }
1781
1782 Vmath::Svtvp(nquad_e, -tau, edgePhys, 1, work, 1, work, 1);
1783 /// TODO: Add variable coeffs
1784 EdgeExp[e]->IProductWRTBase(work, edgeCoeffs);
1785
1786 EdgeExp[e]->SetCoeffsToOrientation(edgedir, edgeCoeffs,
1787 edgeCoeffs);
1788
1789 for (j = 0; j < order_e; ++j)
1790 {
1791 BndMat(cnt + j, i) = edgeCoeffs[j];
1792 }
1793
1794 cnt += order_e;
1795 }
1796 }
1797 }
1798 break;
1799 // HDG postprocessing
1801 {
1803 mkey.GetConstFactors(), mkey.GetVarCoeffs());
1804 DNekScalMat &LapMat = *GetLocMatrix(lapkey);
1805
1807 LapMat.GetRows(), LapMat.GetColumns());
1808 DNekMatSharedPtr lmat = returnval;
1809
1810 (*lmat) = LapMat;
1811
1812 // replace first column with inner product wrt 1
1813 int nq = GetTotPoints();
1814 Array<OneD, NekDouble> tmp(nq);
1816 Vmath::Fill(nq, 1.0, tmp, 1);
1817 IProductWRTBase(tmp, outarray);
1818
1819 Vmath::Vcopy(m_ncoeffs, &outarray[0], 1, &(lmat->GetPtr())[0], 1);
1820 lmat->Invert();
1821 }
1822 break;
1824 {
1825 int ntraces = GetNtraces();
1826 int ncoords = GetCoordim();
1827 int nphys = GetTotPoints();
1829 Array<OneD, NekDouble> phys(nphys);
1830 returnval =
1832 DNekMat &Mat = *returnval;
1833 Vmath::Zero(m_ncoeffs * m_ncoeffs, Mat.GetPtr(), 1);
1834
1836
1837 for (int d = 0; d < ncoords; ++d)
1838 {
1839 Deriv[d] = Array<OneD, NekDouble>(nphys);
1840 }
1841
1842 Array<OneD, int> tracepts(ntraces);
1843 Array<OneD, ExpansionSharedPtr> traceExp(ntraces);
1844 int maxtpts = 0;
1845 for (int t = 0; t < ntraces; ++t)
1846 {
1847 traceExp[t] = GetTraceExp(t);
1848 tracepts[t] = traceExp[t]->GetTotPoints();
1849 maxtpts = (maxtpts > tracepts[t]) ? maxtpts : tracepts[t];
1850 }
1851
1852 Array<OneD, NekDouble> val(maxtpts), tmp, tmp1;
1853
1854 Array<OneD, Array<OneD, NekDouble>> dphidn(ntraces);
1855 for (int t = 0; t < ntraces; ++t)
1856 {
1857 dphidn[t] =
1858 Array<OneD, NekDouble>(m_ncoeffs * tracepts[t], 0.0);
1859 }
1860
1861 for (int i = 0; i < m_ncoeffs; ++i)
1862 {
1863 FillMode(i, phys);
1864 PhysDeriv(phys, Deriv[0], Deriv[1], Deriv[2]);
1865
1866 for (int t = 0; t < ntraces; ++t)
1867 {
1868 const NormalVector norm = GetTraceNormal(t);
1869
1872 traceExp[t]->GetBasis(0)->GetBasisKey();
1873 bool DoInterp = (fromkey != tokey);
1874
1875 Array<OneD, NekDouble> n(tracepts[t]);
1876 ;
1877 for (int d = 0; d < ncoords; ++d)
1878 {
1879 // if variable p may need to interpolate
1880 if (DoInterp)
1881 {
1882 LibUtilities::Interp1D(fromkey, norm[d], tokey, n);
1883 }
1884 else
1885 {
1886 n = norm[d];
1887 }
1888
1889 GetTracePhysVals(t, traceExp[t], Deriv[d], val,
1890 v_GetTraceOrient(t));
1891
1892 Vmath::Vvtvp(tracepts[t], n, 1, val, 1,
1893 tmp = dphidn[t] + i * tracepts[t], 1,
1894 tmp1 = dphidn[t] + i * tracepts[t], 1);
1895 }
1896 }
1897 }
1898
1899 for (int t = 0; t < ntraces; ++t)
1900 {
1901 int nt = tracepts[t];
1902 NekDouble h, p;
1903 TraceNormLen(t, h, p);
1904
1905 // scaling from GJP paper
1906 NekDouble scale =
1907 (p == 1) ? 0.02 * h * h : 0.8 * pow(p + 1, -4.0) * h * h;
1908
1909 for (int i = 0; i < m_ncoeffs; ++i)
1910 {
1911 for (int j = i; j < m_ncoeffs; ++j)
1912 {
1913 Vmath::Vmul(nt, dphidn[t] + i * nt, 1,
1914 dphidn[t] + j * nt, 1, val, 1);
1915 Mat(i, j) =
1916 Mat(i, j) + scale * traceExp[t]->Integral(val);
1917 }
1918 }
1919 }
1920
1921 // fill in symmetric components.
1922 for (int i = 0; i < m_ncoeffs; ++i)
1923 {
1924 for (int j = 0; j < i; ++j)
1925 {
1926 Mat(i, j) = Mat(j, i);
1927 }
1928 }
1929 }
1930 break;
1931 default:
1932 ASSERTL0(false,
1933 "This matrix type cannot be generated from this class");
1934 break;
1935 }
1936
1937 return returnval;
1938}
1939
1940// Evaluate Coefficients of weak deriviative in the direction dir
1941// given the input coefficicents incoeffs and the imposed
1942// boundary values in EdgeExp (which will have its phys space updated);
1944 const Array<OneD, const NekDouble> &incoeffs,
1946 Array<OneD, Array<OneD, NekDouble>> &edgeCoeffs,
1948{
1952
1953 int ncoeffs = GetNcoeffs();
1954
1956 DNekScalMat &Dmat = *GetLocMatrix(DerivType[dir]);
1957
1958 Array<OneD, NekDouble> coeffs = incoeffs;
1959 DNekVec Coeffs(ncoeffs, coeffs, eWrapper);
1960
1961 Coeffs = Transpose(Dmat) * Coeffs;
1962 Vmath::Neg(ncoeffs, coeffs, 1);
1963
1964 // Add the boundary integral including the relevant part of
1965 // the normal
1966 AddNormTraceInt(dir, EdgeExp, edgeCoeffs, coeffs);
1967
1968 DNekVec Out_d(ncoeffs, out_d, eWrapper);
1969
1970 Out_d = InvMass * Coeffs;
1971}
1972
1974{
1979
1981 const int edge, const Array<OneD, const NekDouble> &primCoeffs,
1982 DNekMatSharedPtr &inoutmat)
1983{
1985 "Not set up for non boundary-interior expansions");
1986 ASSERTL1(inoutmat->GetRows() == inoutmat->GetColumns(),
1987 "Assuming that input matrix was square");
1988 int i, j;
1989 int id1, id2;
1990 ExpansionSharedPtr edgeExp = m_traceExp[edge].lock();
1991 int order_e = edgeExp->GetNcoeffs();
1992
1995
1996 StdRegions::VarCoeffMap varcoeffs;
1997 varcoeffs[StdRegions::eVarCoeffMass] = primCoeffs;
1998
2001 varcoeffs);
2002 DNekScalMat &edgemat = *edgeExp->GetLocMatrix(mkey);
2003
2004 // Now need to identify a map which takes the local edge
2005 // mass matrix to the matrix stored in inoutmat;
2006 // This can currently be deduced from the size of the matrix
2007
2008 // - if inoutmat.m_rows() == v_NCoeffs() it is a full
2009 // matrix system
2010
2011 // - if inoutmat.m_rows() == v_NumBndCoeffs() it is a
2012 // boundary CG system
2013
2014 // - if inoutmat.m_rows() == v_NumDGBndCoeffs() it is a
2015 // trace DG system
2016 int rows = inoutmat->GetRows();
2017
2018 if (rows == GetNcoeffs())
2019 {
2020 GetTraceToElementMap(edge, map, sign, v_GetTraceOrient(edge));
2021 }
2022 else if (rows == NumBndryCoeffs())
2023 {
2024 int nbndry = NumBndryCoeffs();
2025 Array<OneD, unsigned int> bmap(nbndry);
2026
2027 GetTraceToElementMap(edge, map, sign, v_GetTraceOrient(edge));
2028
2029 GetBoundaryMap(bmap);
2030
2031 for (i = 0; i < order_e; ++i)
2032 {
2033 for (j = 0; j < nbndry; ++j)
2034 {
2035 if (map[i] == bmap[j])
2036 {
2037 map[i] = j;
2038 break;
2039 }
2040 }
2041 ASSERTL1(j != nbndry, "Did not find number in map");
2042 }
2043 }
2044 else if (rows == NumDGBndryCoeffs())
2045 {
2046 // possibly this should be a separate method
2047 int cnt = 0;
2048 map = Array<OneD, unsigned int>(order_e);
2049 sign = Array<OneD, int>(order_e, 1);
2050
2051 for (i = 0; i < edge; ++i)
2052 {
2053 cnt += GetTraceNcoeffs(i);
2054 }
2055
2056 for (i = 0; i < order_e; ++i)
2057 {
2058 map[i] = cnt++;
2059 }
2060 // check for mapping reversal
2062 {
2063 switch (edgeExp->GetBasis(0)->GetBasisType())
2064 {
2066 reverse(map.get(), map.get() + order_e);
2067 break;
2069 reverse(map.get(), map.get() + order_e);
2070 break;
2072 {
2073 swap(map[0], map[1]);
2074 for (i = 3; i < order_e; i += 2)
2075 {
2076 sign[i] = -1;
2077 }
2078 }
2079 break;
2080 default:
2081 ASSERTL0(false,
2082 "Edge boundary type not valid for this method");
2083 }
2084 }
2085 }
2086 else
2087 {
2088 ASSERTL0(false, "Could not identify matrix type from dimension");
2089 }
2090
2091 for (i = 0; i < order_e; ++i)
2092 {
2093 id1 = map[i];
2094 for (j = 0; j < order_e; ++j)
2095 {
2096 id2 = map[j];
2097 (*inoutmat)(id1, id2) += edgemat(i, j) * sign[i] * sign[j];
2098 }
2099 }
2100}
2101
2102/**
2103 * Given an edge and vector of element coefficients:
2104 * - maps those elemental coefficients corresponding to the edge into
2105 * an edge-vector.
2106 * - resets the element coefficients
2107 * - multiplies the edge vector by the edge mass matrix
2108 * - maps the edge coefficients back onto the elemental coefficients
2109 */
2111 const int edgeid, const Array<OneD, const NekDouble> &primCoeffs,
2112 const Array<OneD, NekDouble> &incoeffs, Array<OneD, NekDouble> &coeffs)
2113{
2115 "Not set up for non boundary-interior expansions");
2116 int i;
2117 ExpansionSharedPtr edgeExp = m_traceExp[edgeid].lock();
2118 int order_e = edgeExp->GetNcoeffs();
2119
2122
2123 StdRegions::VarCoeffMap varcoeffs;
2124 varcoeffs[StdRegions::eVarCoeffMass] = primCoeffs;
2125
2128 varcoeffs);
2129 DNekScalMat &edgemat = *edgeExp->GetLocMatrix(mkey);
2130
2131 NekVector<NekDouble> vEdgeCoeffs(order_e);
2132
2133 GetTraceToElementMap(edgeid, map, sign, v_GetTraceOrient(edgeid));
2134
2135 for (i = 0; i < order_e; ++i)
2136 {
2137 vEdgeCoeffs[i] = incoeffs[map[i]] * sign[i];
2138 }
2139
2140 vEdgeCoeffs = edgemat * vEdgeCoeffs;
2141
2142 for (i = 0; i < order_e; ++i)
2143 {
2144 coeffs[map[i]] += vEdgeCoeffs[i] * sign[i];
2145 }
2146}
2147
2149 const DNekScalMatSharedPtr &r_bnd)
2150{
2151 MatrixStorage storage = eFULL;
2152 DNekMatSharedPtr m_vertexmatrix;
2153
2154 int nVerts, vid1, vid2, vMap1, vMap2;
2155 NekDouble VertexValue;
2156
2157 nVerts = GetNverts();
2158
2159 m_vertexmatrix =
2160 MemoryManager<DNekMat>::AllocateSharedPtr(nVerts, nVerts, 0.0, storage);
2161 DNekMat &VertexMat = (*m_vertexmatrix);
2162
2163 for (vid1 = 0; vid1 < nVerts; ++vid1)
2164 {
2165 vMap1 = GetVertexMap(vid1);
2166
2167 for (vid2 = 0; vid2 < nVerts; ++vid2)
2168 {
2169 vMap2 = GetVertexMap(vid2);
2170 VertexValue = (*r_bnd)(vMap1, vMap2);
2171 VertexMat.SetValue(vid1, vid2, VertexValue);
2172 }
2173 }
2174
2175 return m_vertexmatrix;
2176}
2177
2179{
2181 GetTraceBasisKey(traceid), m_geom->GetEdge(traceid));
2182}
2183
2185{
2186 int n, j;
2187 int nEdgeCoeffs;
2188 int nBndCoeffs = NumBndryCoeffs();
2189
2190 Array<OneD, unsigned int> bmap(nBndCoeffs);
2191 GetBoundaryMap(bmap);
2192
2193 // Map from full system to statically condensed system (i.e reverse
2194 // GetBoundaryMap)
2195 map<int, int> invmap;
2196 for (j = 0; j < nBndCoeffs; ++j)
2197 {
2198 invmap[bmap[j]] = j;
2199 }
2200
2201 // Number of interior edge coefficients
2202 nEdgeCoeffs = GetTraceNcoeffs(eid) - 2;
2203
2205
2206 Array<OneD, unsigned int> edgemaparray(nEdgeCoeffs);
2207 Array<OneD, unsigned int> maparray(nEdgeCoeffs);
2208 Array<OneD, int> signarray(nEdgeCoeffs, 1);
2209 StdRegions::Orientation eOrient = geom->GetEorient(eid);
2210
2211 // maparray is the location of the edge within the matrix
2212 GetTraceInteriorToElementMap(eid, maparray, signarray, eOrient);
2213
2214 for (n = 0; n < nEdgeCoeffs; ++n)
2215 {
2216 edgemaparray[n] = invmap[maparray[n]];
2217 }
2218
2219 return edgemaparray;
2220}
2221
2223{
2225}
2226
2228 Array<OneD, int> &idmap, const int nq0,
2229 [[maybe_unused]] const int nq1)
2230{
2231 if (idmap.size() != nq0)
2232 {
2233 idmap = Array<OneD, int>(nq0);
2234 }
2235 switch (orient)
2236 {
2238 // Fwd
2239 for (int i = 0; i < nq0; ++i)
2240 {
2241 idmap[i] = i;
2242 }
2243 break;
2245 {
2246 // Bwd
2247 for (int i = 0; i < nq0; ++i)
2248 {
2249 idmap[i] = nq0 - 1 - i;
2250 }
2251 }
2252 break;
2253 default:
2254 ASSERTL0(false, "Unknown orientation");
2255 break;
2256 }
2257}
2258
2259// Compute edgenormal \cdot vector
2261 const int dir, const int edge, ExpansionSharedPtr &EdgeExp_e,
2262 const Array<OneD, const Array<OneD, NekDouble>> &normals,
2263 const StdRegions::VarCoeffMap &varcoeffs)
2264{
2265 int nquad_e = EdgeExp_e->GetNumPoints(0);
2266 int coordim = GetCoordim();
2267 int nquad0 = m_base[0]->GetNumPoints();
2268 int nquad1 = m_base[1]->GetNumPoints();
2269 int nqtot = nquad0 * nquad1;
2270
2271 StdRegions::VarCoeffType MMFCoeffs[15] = {
2280
2281 StdRegions::VarCoeffMap::const_iterator MFdir;
2282
2283 Array<OneD, NekDouble> ncdotMF(nqtot, 0.0);
2284 Array<OneD, NekDouble> tmp(nqtot);
2285 Array<OneD, NekDouble> tmp_e(nquad_e);
2286 for (int k = 0; k < coordim; k++)
2287 {
2288 MFdir = varcoeffs.find(MMFCoeffs[dir * 5 + k]);
2289 tmp = MFdir->second.GetValue();
2290
2291 GetPhysEdgeVarCoeffsFromElement(edge, EdgeExp_e, tmp, tmp_e);
2292
2293 Vmath::Vvtvp(nquad_e, &tmp_e[0], 1, &normals[k][0], 1, &ncdotMF[0], 1,
2294 &ncdotMF[0], 1);
2295 }
2296 return ncdotMF;
2297}
2298
2300 const Array<OneD, Array<OneD, NekDouble>> &vec)
2301{
2303 GetLeftAdjacentElementExp()->GetTraceNormal(
2305
2306 int nq = GetTotPoints();
2308 Vmath::Vmul(nq, &vec[0][0], 1, &normals[0][0], 1, &Fn[0], 1);
2309 Vmath::Vvtvp(nq, &vec[1][0], 1, &normals[1][0], 1, &Fn[0], 1, &Fn[0], 1);
2310 Vmath::Vvtvp(nq, &vec[2][0], 1, &normals[2][0], 1, &Fn[0], 1, &Fn[0], 1);
2311
2312 return StdExpansion::Integral(Fn);
2313}
2314
2316{
2318
2319 int nverts = geom->GetNumVerts();
2320
2321 // vertices on edges
2322 SpatialDomains::PointGeom ev0 = *geom->GetVertex(traceid);
2323 SpatialDomains::PointGeom ev1 = *geom->GetVertex((traceid + 1) % nverts);
2324
2325 // vertex on adjacent edge to ev0
2327 *geom->GetVertex((traceid + (nverts - 1)) % nverts);
2328
2329 // calculate perpendicular distance of normal length
2330 // from first vertex
2331 NekDouble h1 = ev0.dist(vadj);
2333
2334 Dx.Sub(ev1, ev0);
2335 Dx1.Sub(vadj, ev0);
2336
2337 NekDouble d1 = Dx.dot(Dx1);
2338 NekDouble lenDx = Dx.dot(Dx);
2339 h = sqrt(h1 * h1 - d1 * d1 / lenDx);
2340
2341 // perpendicular distanace from second vertex
2342 SpatialDomains::PointGeom vadj1 = *geom->GetVertex((traceid + 2) % nverts);
2343
2344 h1 = ev1.dist(vadj1);
2345 Dx1.Sub(vadj1, ev1);
2346 d1 = Dx.dot(Dx1);
2347
2348 h = (h + sqrt(h1 * h1 - d1 * d1 / lenDx)) * 0.5;
2349
2350 int dirn = (geom->GetDir(traceid) == 0) ? 1 : 0;
2351
2352 p = (NekDouble)(GetBasisNumModes(dirn) - 1);
2353}
2354} // namespace Nektar::LocalRegions
#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
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:47
Describes the specification for a Basis.
Definition: Basis.h:45
void v_AddEdgeNormBoundaryInt(const int edge, const ExpansionSharedPtr &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray) override
DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey) override
std::vector< bool > m_requireNeg
Definition: Expansion2D.h:112
void v_AddRobinTraceContribution(const int traceid, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs) override
void v_SetUpPhysNormals(const int edge) override
void v_AddRobinMassMatrix(const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat) override
Array< OneD, NekDouble > GetnEdgecdotMF(const int dir, const int edge, ExpansionSharedPtr &EdgeExp_e, const Array< OneD, const Array< OneD, NekDouble > > &normals, const StdRegions::VarCoeffMap &varcoeffs)
void SetTraceToGeomOrientation(Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, NekDouble > &inout)
DNekMatSharedPtr v_BuildVertexMatrix(const DNekScalMatSharedPtr &r_bnd) override
NekDouble v_VectorFlux(const Array< OneD, Array< OneD, NekDouble > > &vec) override
DNekScalMatSharedPtr CreateMatrix(const MatrixKey &mkey)
Definition: Expansion2D.cpp:55
Array< OneD, unsigned int > GetTraceInverseBoundaryMap(int eid)
void v_ReOrientTracePhysMap(const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1) override
Expansion2D(SpatialDomains::Geometry2DSharedPtr pGeom)
Definition: Expansion2D.cpp:50
void AddNormTraceInt(const int dir, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &edgeCoeffs, Array< OneD, NekDouble > &outarray)
SpatialDomains::Geometry2DSharedPtr GetGeom2D() const
Definition: Expansion2D.h:164
void AddHDGHelmholtzEdgeTerms(const NekDouble tau, const int edge, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, NekDouble > &edgePhys, const StdRegions::VarCoeffMap &dirForcing, Array< OneD, NekDouble > &outarray)
void v_DGDeriv(const int dir, const Array< OneD, const NekDouble > &incoeffs, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &edgeCoeffs, Array< OneD, NekDouble > &out_d) override
void v_TraceNormLen(const int traceid, NekDouble &h, NekDouble &p) override
void GetPhysEdgeVarCoeffsFromElement(const int edge, ExpansionSharedPtr &EdgeExp, const Array< OneD, const NekDouble > &varcoeff, Array< OneD, NekDouble > &outarray)
void AddEdgeBoundaryInt(const int edge, ExpansionSharedPtr &EdgeExp, Array< OneD, NekDouble > &edgePhys, Array< OneD, NekDouble > &outarray, const StdRegions::VarCoeffMap &varcoeffs=StdRegions::NullVarCoeffMap)
void AddHDGHelmholtzTraceTerms(const NekDouble tau, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, const StdRegions::VarCoeffMap &dirForcing, Array< OneD, NekDouble > &outarray)
void v_GenTraceExp(const int traceid, ExpansionSharedPtr &exp) override
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:167
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:273
void DropLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:90
ExpansionSharedPtr GetLeftAdjacentElementExp() const
Definition: Expansion.h:441
DNekMatSharedPtr BuildVertexMatrix(const DNekScalMatSharedPtr &r_bnd)
Definition: Expansion.cpp:101
void ComputeGmatcdotMF(const Array< TwoD, const NekDouble > &df, const Array< OneD, const NekDouble > &direction, Array< OneD, Array< OneD, NekDouble > > &dfdir)
Definition: Expansion.cpp:603
virtual void v_ComputeTraceNormal(const int id)
Definition: Expansion.cpp:868
Array< OneD, NekDouble > GetMFMag(const int dir, const StdRegions::VarCoeffMap &varcoeffs)
Definition: Expansion.cpp:709
int GetLeftAdjacentElementTrace() const
Definition: Expansion.h:454
void GetTracePhysVals(const int trace, const StdRegions::StdExpansionSharedPtr &TraceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=StdRegions::eNoOrientation)
Definition: Expansion.h:200
std::map< int, ExpansionWeakPtr > m_traceExp
Definition: Expansion.h:272
void AddEdgeNormBoundaryInt(const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:114
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:274
virtual StdRegions::Orientation v_GetTraceOrient(int trace)
Definition: Expansion.cpp:813
ExpansionSharedPtr GetTraceExp(const int traceid)
Definition: Expansion.h:414
StdRegions::Orientation GetTraceOrient(int trace)
Definition: Expansion.h:168
IndexMapValuesSharedPtr GetIndexMap(const IndexMapKey &ikey)
Definition: Expansion.h:146
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:84
Array< OneD, NekDouble > GetMFDiv(const int dir, const StdRegions::VarCoeffMap &varcoeffs)
Definition: Expansion.cpp:686
void TraceNormLen(const int traceid, NekDouble &h, NekDouble &p)
Definition: Expansion.h:250
const NormalVector & GetTraceNormal(const int id)
Definition: Expansion.cpp:251
Array< OneD, NekDouble > GetMF(const int dir, const int shapedim, const StdRegions::VarCoeffMap &varcoeffs)
Definition: Expansion.cpp:633
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
PointGeomSharedPtr GetVertex(int i) const
Returns vertex i of this object.
Definition: Geometry.h:355
void Sub(PointGeom &a, PointGeom &b)
Definition: PointGeom.cpp:122
NekDouble dot(PointGeom &a)
retun the dot product between this and input a
Definition: PointGeom.cpp:185
NekDouble dist(PointGeom &a)
return distance between this and input a
Definition: PointGeom.cpp:178
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:669
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:124
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:134
void FillMode(const int mode, Array< OneD, NekDouble > &outarray)
This function fills the array outarray with the mode-th mode of the expansion.
Definition: StdExpansion.h:491
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:603
const LibUtilities::PointsKeyVector GetPointsKeys() const
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix(const LocalRegions::MatrixKey &mkey)
Definition: StdExpansion.h:641
const LibUtilities::BasisKey GetTraceBasisKey(const int i, int k=-1) const
This function returns the basis key belonging to the i-th trace.
Definition: StdExpansion.h:299
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
Definition: StdExpansion.h:679
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
this function calculates the inner product of a given function f with the different modes of the expa...
Definition: StdExpansion.h:528
int GetNtraces() const
Returns the number of trace elements connected to this element.
Definition: StdExpansion.h:351
int GetNverts() const
This function returns the number of vertices of the expansion domain.
Definition: StdExpansion.h:246
void GetTraceToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
Definition: StdExpansion.h:684
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:367
void GetTraceInteriorToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eForwards)
Definition: StdExpansion.h:708
int GetTraceNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th trace.
Definition: StdExpansion.h:261
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:844
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, DNekMat, StdMatrixKey::opLess > m_stdMatrixManager
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:169
void PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Definition: StdExpansion.h:849
Array< OneD, LibUtilities::BasisSharedPtr > m_base
LibUtilities::ShapeType GetShapeType() const
Definition: StdMatrixKey.h:88
const VarCoeffMap & GetVarCoeffs() const
Definition: StdMatrixKey.h:169
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:83
bool HasVarCoeff(const StdRegions::VarCoeffType &coeff) const
Definition: StdMatrixKey.h:174
const ConstFactorMap & GetConstFactors() const
Definition: StdMatrixKey.h:138
const Array< OneD, const NekDouble > & GetVarCoeff(const StdRegions::VarCoeffType &coeff) const
Definition: StdMatrixKey.h:148
NekDouble GetConstFactor(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:124
bool ConstFactorExists(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:133
const VarCoeffMap GetVarCoeffAsMap(const VarCoeffType &coeff) const
Definition: StdMatrixKey.h:158
void Interp1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
this function interpolates a 1D function evaluated at the quadrature points of the basis fbasis0 to ...
Definition: Interp.cpp:47
void InterpCoeff1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
Definition: InterpCoeff.cpp:42
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:231
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points.
Definition: BasisType.h:57
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:56
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:66
std::shared_ptr< IndexMapValues > IndexMapValuesSharedPtr
Definition: IndexMapKey.h:126
std::shared_ptr< Expansion1D > Expansion1DSharedPtr
Definition: Expansion1D.h:50
@ eNoGeomType
No type defined.
@ eDeformed
Geometry is curved or has non-constant factors.
std::shared_ptr< Geometry2D > Geometry2DSharedPtr
Definition: Geometry.h:62
std::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:51
@ eLinearAdvectionDiffusionReactionGJP
Definition: StdRegions.hpp:105
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:402
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:403
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:347
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:346
std::vector< double > d(NPUPPER *NPUPPER)
StdRegions::ConstFactorMap factors
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:79
static Array< OneD, NekDouble > NullNekDouble1DArray
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvp (scalar times vector plus vector): z = alpha*x + y.
Definition: Vmath.hpp:396
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.hpp:366
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.hpp:54
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294