35 #include <boost/core/ignore_unused.hpp> 50 StdQuadExp::StdQuadExp()
59 StdExpansion (Ba.GetNumModes()*Bb.GetNumModes(),2,Ba,Bb),
104 boost::ignore_unused(out_d2);
126 ASSERTL1(
false,
"input dir is out of range");
137 boost::ignore_unused(out_d2);
159 if(
m_base[0]->Collocation() &&
m_base[1]->Collocation())
162 inarray, 1, outarray, 1);
175 m_base[1]->GetNumModes());
179 inarray,outarray,wsp,
true,
true);
198 bool doCheckCollDir0,
199 bool doCheckCollDir1)
201 int nquad0 =
m_base[0]->GetNumPoints();
202 int nquad1 =
m_base[1]->GetNumPoints();
203 int nmodes0 =
m_base[0]->GetNumModes();
204 int nmodes1 =
m_base[1]->GetNumModes();
206 bool colldir0 = doCheckCollDir0?(
m_base[0]->Collocation()):
false;
207 bool colldir1 = doCheckCollDir1?(
m_base[1]->Collocation()):
false;
209 if(colldir0 && colldir1)
215 Blas::Dgemm(
'N',
'T', nquad0, nquad1,nmodes1, 1.0, &inarray[0], nquad0,
216 base1.get(), nquad1, 0.0, &outarray[0], nquad0);
220 Blas::Dgemm(
'N',
'N', nquad0,nmodes1,nmodes0,1.0, base0.get(),
221 nquad0, &inarray[0], nmodes0,0.0,&outarray[0], nquad0);
225 ASSERTL1(wsp.num_elements()>=nquad0*nmodes1,
"Workspace size is not sufficient");
229 Blas::Dgemm(
'N',
'N', nquad0,nmodes1,nmodes0,1.0, base0.get(),
230 nquad0, &inarray[0], nmodes0,0.0,&wsp[0], nquad0);
231 Blas::Dgemm(
'N',
'T', nquad0, nquad1,nmodes1, 1.0, &wsp[0], nquad0,
232 base1.get(), nquad1, 0.0, &outarray[0], nquad0);
240 if((
m_base[0]->Collocation())&&(
m_base[1]->Collocation()))
264 if((
m_base[0]->Collocation())&&(
m_base[1]->Collocation()))
271 int npoints[2] = {
m_base[0]->GetNumPoints(),
272 m_base[1]->GetNumPoints()};
273 int nmodes[2] = {
m_base[0]->GetNumModes(),
274 m_base[1]->GetNumModes()};
276 fill(outarray.get(), outarray.get()+
m_ncoeffs, 0.0 );
280 for(i = 0; i < 4; i++)
286 for(i = 0; i < npoints[0]; i++)
288 physEdge[0][i] = inarray[i];
289 physEdge[2][i] = inarray[npoints[0]*npoints[1]-1-i];
292 for(i = 0; i < npoints[1]; i++)
294 physEdge[1][i] = inarray[npoints[0]-1+i*npoints[0]];
295 physEdge[3][i] = inarray[(npoints[1]-1)*npoints[0]-i*npoints[0]];
305 for(i = 0; i < 4; i++)
307 segexp[i%2]->FwdTrans_BndConstrained(physEdge[i],coeffEdge[i]);
310 for(j=0; j < nmodes[i%2]; j++)
313 outarray[ mapArray[j] ] = sign * coeffEdge[i][j];
332 int nInteriorDofs =
m_ncoeffs - nBoundaryDofs;
340 for(i = 0; i < nInteriorDofs; i++)
342 rhs[i] = tmp1[ mapArray[i] ];
345 Blas::Dgemv(
'N',nInteriorDofs,nInteriorDofs,1.0, &(matsys->GetPtr())[0],
346 nInteriorDofs,rhs.get(),1,0.0,result.get(),1);
348 for(i = 0; i < nInteriorDofs; i++)
350 outarray[ mapArray[i] ] = result[i];
389 if(
m_base[0]->Collocation() &&
m_base[1]->Collocation())
402 bool multiplybyweights)
404 int nquad0 =
m_base[0]->GetNumPoints();
405 int nquad1 =
m_base[1]->GetNumPoints();
406 int order0 =
m_base[0]->GetNumModes();
408 if(multiplybyweights)
417 tmp,outarray,wsp,
true,
true);
424 inarray,outarray,wsp,
true,
true);
437 m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
451 ASSERTL0((dir==0)||(dir==1),
"input dir is out of range");
453 int nquad0 =
m_base[0]->GetNumPoints();
454 int nquad1 =
m_base[1]->GetNumPoints();
455 int nqtot = nquad0*nquad1;
456 int order0 =
m_base[0]->GetNumModes();
468 tmp,outarray,wsp,
true,
false);
474 tmp,outarray,wsp,
false,
true);
482 ASSERTL0((dir==0)||(dir==1),
"input dir is out of range");
500 m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
519 bool doCheckCollDir0,
520 bool doCheckCollDir1)
522 int nquad0 =
m_base[0]->GetNumPoints();
523 int nquad1 =
m_base[1]->GetNumPoints();
524 int nmodes0 =
m_base[0]->GetNumModes();
525 int nmodes1 =
m_base[1]->GetNumModes();
527 bool colldir0 = doCheckCollDir0?(
m_base[0]->Collocation()):
false;
528 bool colldir1 = doCheckCollDir1?(
m_base[1]->Collocation()):
false;
530 if(colldir0 && colldir1)
536 Blas::Dgemm(
'N',
'N',nmodes0,nmodes1, nquad1,1.0, inarray.get(),
537 nmodes0, base1.get(), nquad1, 0.0,outarray.get(),nmodes0);
541 Blas::Dgemm(
'T',
'N',nmodes0,nquad1,nquad0,1.0,base0.get(),
542 nquad0,inarray.get(),nquad0,0.0,outarray.get(),nmodes0);
546 ASSERTL1(wsp.num_elements()>=nquad1*nmodes0,
"Workspace size is not sufficient");
549 Blas::Dgemm(
'T',
'N',nmodes0,nquad1,nquad0,1.0,base0.get(),
550 nquad0,inarray.get(),nquad0,0.0,wsp.get(),nmodes0);
553 for(
int i = 0; i < nmodes0; ++i)
555 for(
int j = 0; j < nquad1; ++j)
558 base0.get()+i*nquad0,1,
559 inarray.get()+j*nquad0,1);
563 Blas::Dgemm(
'N',
'N',nmodes0,nmodes1, nquad1,1.0, wsp.get(),
564 nmodes0, base1.get(), nquad1, 0.0,outarray.get(),nmodes0);
590 int nquad0 =
m_base[0]->GetNumPoints();
591 int nquad1 =
m_base[1]->GetNumPoints();
594 int btmp0 =
m_base[0]->GetNumModes();
595 int mode0 = mode%btmp0;
596 int mode1 = mode/btmp0;
599 ASSERTL2(mode1 == (
int)floor((1.0*mode)/btmp0),
600 "Integer Truncation not Equiv to Floor");
603 "calling argument mode is larger than total expansion order");
605 for(i = 0; i < nquad1; ++i)
608 1, &outarray[0]+i*nquad0,1);
611 for(i = 0; i < nquad0; ++i)
614 &outarray[0]+i,nquad0,&outarray[0]+i,nquad0);
634 ASSERTL2((i >= 0)&&(i <= 3),
"edge id is out of range");
636 if((i == 0)||(i == 2))
648 ASSERTL2((i >= 0)&&(i <= 3),
"edge id is out of range");
650 if((i == 0)||(i == 2))
662 ASSERTL2((i >= 0)&&(i <= 3),
"edge id is out of range");
664 if((i == 0)||(i == 2))
676 ASSERTL2((edge >= 0)&&(edge <= 3),
"edge id is out of range");
678 if((edge == 0)||(edge == 2))
690 ASSERTL2((i >= 0)&&(i <= 3),
"edge id is out of range");
692 if((i == 0)||(i == 2))
714 "BasisType is not a boundary interior form");
718 "BasisType is not a boundary interior form");
728 "BasisType is not a boundary interior form");
732 "BasisType is not a boundary interior form");
738 const std::vector<unsigned int> &nummodes,
741 int nmodes = nummodes[modes_offset]*nummodes[modes_offset+1];
749 bool returnval =
false;
768 boost::ignore_unused(coords_2);
775 for(i = 0; i < nq1; ++i)
777 Blas::Dcopy(nq0,z0.get(), 1,&coords_0[0] + i*nq0,1);
790 int nummodes0, nummodes1;
791 int value1 = 0, value2 = 0;
797 nummodes0 =
m_base[0]->GetNumModes();
798 nummodes1 =
m_base[1]->GetNumModes();
810 value1 = 2*nummodes0;
813 ASSERTL0(0,
"Mapping array is not defined for this expansion");
817 for(i = 0; i < value1; i++)
827 value2 = value1+nummodes0-1;
833 ASSERTL0(0,
"Mapping array is not defined for this expansion");
837 for(i = 0; i < nummodes1-2; i++)
839 outarray[cnt++]=value1+i*nummodes0;
840 outarray[cnt++]=value2+i*nummodes0;
846 for(i = nummodes0*(nummodes1-1);i <
GetNcoeffs(); i++)
857 int nummodes0, nummodes1;
864 nummodes0 =
m_base[0]->GetNumModes();
865 nummodes1 =
m_base[1]->GetNumModes();
873 startvalue = nummodes0;
876 startvalue = 2*nummodes0;
879 ASSERTL0(0,
"Mapping array is not defined for this expansion");
892 ASSERTL0(0,
"Mapping array is not defined for this expansion");
896 for(i = 0; i < nummodes1-2; i++)
898 for(j = 0; j < nummodes0-2; j++)
900 outarray[cnt++]=startvalue+j;
902 startvalue+=nummodes0;
910 if(useCoeffPacking ==
true)
912 switch(localVertexId)
923 localDOF =
m_base[0]->GetNumModes()-1;
935 localDOF =
m_base[0]->GetNumModes() * (
m_base[1]->GetNumModes()-1);
939 localDOF =
m_base[0]->GetNumModes();
947 localDOF =
m_base[0]->GetNumModes()*
m_base[1]->GetNumModes()-1;
951 localDOF =
m_base[0]->GetNumModes()+1;
956 ASSERTL0(
false,
"eid must be between 0 and 3");
963 switch(localVertexId)
974 localDOF =
m_base[0]->GetNumModes()-1;
986 localDOF =
m_base[0]->GetNumModes()*
m_base[1]->GetNumModes()-1;
990 localDOF =
m_base[0]->GetNumModes()+1;
998 localDOF =
m_base[0]->GetNumModes() * (
m_base[1]->GetNumModes()-1);
1002 localDOF =
m_base[0]->GetNumModes();
1007 ASSERTL0(
false,
"eid must be between 0 and 3");
1020 const int nummodes0 =
m_base[0]->GetNumModes();
1021 const int nummodes1 =
m_base[1]->GetNumModes();
1025 if(maparray.num_elements() != nEdgeIntCoeffs)
1030 if(signarray.num_elements() != nEdgeIntCoeffs)
1036 fill( signarray.get() , signarray.get()+nEdgeIntCoeffs, 1 );
1045 for(i = 0; i < nEdgeIntCoeffs; i++)
1052 for(i = 1; i < nEdgeIntCoeffs; i+=2)
1061 for(i = 0; i < nEdgeIntCoeffs; i++)
1063 maparray[i] = (i+2)*nummodes0 + 1;
1068 for(i = 1; i < nEdgeIntCoeffs; i+=2)
1077 for(i = 0; i < nEdgeIntCoeffs; i++)
1079 maparray[i] = nummodes0+i+2;
1084 for(i = 1; i < nEdgeIntCoeffs; i+=2)
1093 for(i = 0; i < nEdgeIntCoeffs; i++)
1095 maparray[i] = (i+2)*nummodes0;
1100 for(i = 1; i < nEdgeIntCoeffs; i+=2)
1108 ASSERTL0(
false,
"eid must be between 0 and 3");
1118 for(i = 0; i < nEdgeIntCoeffs; i++)
1126 for(i = 0; i < nEdgeIntCoeffs; i++)
1128 maparray[i] = (i+2)*nummodes0 - 1;
1134 for(i = 0; i < nEdgeIntCoeffs; i++)
1136 maparray[i] = nummodes0*(nummodes1-1) + i + 1;
1142 for(i = 0; i < nEdgeIntCoeffs; i++)
1144 maparray[i] = nummodes0 * (i+1);
1149 ASSERTL0(
false,
"eid must be between 0 and 3");
1154 reverse( maparray.get() , maparray.get()+nEdgeIntCoeffs );
1159 ASSERTL0(
false,
"Mapping not defined for this type of basis");
1173 int order0 =
m_base[0]->GetNumModes();
1174 int order1 =
m_base[1]->GetNumModes();
1187 ASSERTL0(
false,
"eid must be between 0 and 3");
1190 bool checkForZeroedModes =
false;
1195 else if(P != numModes)
1197 checkForZeroedModes =
true;
1202 if (maparray.num_elements() !=
P)
1207 if(signarray.num_elements() !=
P)
1213 fill(signarray.get(), signarray.get()+
P, 1);
1222 for (i = 0; i <
P; i++)
1229 swap(maparray[0], maparray[1]);
1231 for(i = 3; i <
P; i+=2)
1240 for (i = 0; i <
P; i++)
1242 maparray[i] = i*order0 + 1;
1247 swap(maparray[0], maparray[1]);
1249 for(i = 3; i <
P; i+=2)
1258 for (i = 0; i <
P; i++)
1260 maparray[i] = order0+i;
1265 swap(maparray[0], maparray[1]);
1267 for (i = 3; i <
P; i+=2)
1276 for (i = 0; i <
P; i++)
1278 maparray[i] = i*order0;
1283 swap(maparray[0], maparray[1]);
1285 for (i = 3; i <
P; i+=2)
1293 ASSERTL0(
false,
"eid must be between 0 and 3");
1304 for (i = 0; i <
P; i++)
1312 for (i = 0; i <
P; i++)
1314 maparray[i] = (i+1)*order0 - 1;
1320 for (i = 0; i <
P; i++)
1322 maparray[i] = order0*(order1-1) + i;
1328 for (i = 0; i <
P; i++)
1330 maparray[i] = order0*i;
1335 ASSERTL0(
false,
"eid must be between 0 and 3");
1340 reverse(maparray.get(), maparray.get()+
P);
1345 ASSERTL0(
false,
"Mapping not defined for this type of basis");
1348 if (checkForZeroedModes)
1354 for (
int j = numModes; j <
P; j++)
1357 maparray[j] = maparray[0];
1362 ASSERTL0(
false,
"Different trace space edge dimension " 1363 "and element edge dimension not possible " 1364 "for GLL-Lagrange bases");
1388 int nq0 =
m_base[0]->GetNumPoints();
1389 int nq1 =
m_base[1]->GetNumPoints();
1413 for(
int i = 0; i < nq; ++i)
1415 for(
int j = 0; j < nq; ++j,++cnt)
1418 coords[cnt][0] = -1.0 + 2*j/(
NekDouble)(nq-1);
1419 coords[cnt][1] = -1.0 + 2*i/(
NekDouble)(nq-1);
1423 for(
int i = 0; i < neq; ++i)
1427 I[0] =
m_base[0]->GetI(coll);
1428 I[1] =
m_base[1]->GetI(coll+1);
1431 for (
int j = 0; j < nq1; ++j)
1437 Mat->GetRawPtr()+j*nq0*neq+i,neq);
1450 for(i = 0; i < order1; ++i)
1452 (*Mat)(order0*i+1,i*order0+1) = 1.0;
1458 for(i = 0; i < order0; ++i)
1460 (*Mat)(order0+i ,order0+i) = 1.0;
1473 (*Mat) = Imass*Iprod;
1481 int dir = (edge + 1) % 2;
1482 int nCoeffs =
m_base[dir]->GetNumModes();
1490 coords[0] = (edge == 0 || edge == 3) ? -1.0 : 1.0;
1498 Vmath::Vcopy(nCoeffs, m_Ix->GetPtr(), 1, Mat->GetPtr(), 1);
1528 if(inarray.get() == outarray.get())
1534 m_ncoeffs, tmp.get(), 1, 0.0, outarray.get(), 1);
1539 m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
1548 int qa =
m_base[0]->GetNumPoints();
1549 int qb =
m_base[1]->GetNumPoints();
1550 int nmodes_a =
m_base[0]->GetNumModes();
1551 int nmodes_b =
m_base[1]->GetNumModes();
1552 int nmodes = min(nmodes_a,nmodes_b);
1565 OrthoExp.
FwdTrans(array,orthocoeffs);
1574 for(
int j = 0; j < nmodes_a; ++j)
1576 for(
int k = 0; k < nmodes_b; ++k)
1580 pow((1.0*j)/(nmodes_a-1),cutoff*nmodes_a),
1581 pow((1.0*k)/(nmodes_b-1),cutoff*nmodes_b));
1583 orthocoeffs[j*nmodes_b+k] *= SvvDiffCoeff * fac;
1594 max_ab = max(max_ab,0);
1597 for(
int j = 0; j < nmodes_a; ++j)
1599 for(
int k = 0; k < nmodes_b; ++k)
1601 int maxjk = max(j,k);
1604 orthocoeffs[j*nmodes_b+k] *= SvvDiffCoeff *
1614 min(nmodes_a,nmodes_b));
1620 for(
int j = 0; j < nmodes_a; ++j)
1622 for(
int k = 0; k < nmodes_b; ++k)
1626 orthocoeffs[j*nmodes_b+k] *=
1627 (SvvDiffCoeff*exp(-(j+k-nmodes)*(j+k-nmodes)/
1633 orthocoeffs[j*nmodes_b+k] *= 0.0;
1641 OrthoExp.
BwdTrans(orthocoeffs,array);
1651 int qa =
m_base[0]->GetNumPoints();
1652 int qb =
m_base[1]->GetNumPoints();
1653 int nmodesA =
m_base[0]->GetNumModes();
1654 int nmodesB =
m_base[1]->GetNumModes();
1655 int P = nmodesA - 1;
1656 int Q = nmodesB - 1;
1667 int Pcut = cutoff*
P;
1668 int Qcut = cutoff*Q;
1672 OrthoExp.
FwdTrans(array,orthocoeffs);
1676 for(
int i = 0; i < nmodesA; ++i)
1678 for(
int j = 0; j < nmodesB; ++j)
1681 if(i > Pcut || j > Qcut)
1685 fac = max(fac1, fac2);
1686 fac = pow(fac, exponent);
1687 orthocoeffs[i*nmodesB+j] *= exp(-alpha*fac);
1693 OrthoExp.
BwdTrans(orthocoeffs,array);
1701 int n_coeffs = inarray.num_elements();
1709 int nmodes0 =
m_base[0]->GetNumModes();
1710 int nmodes1 =
m_base[1]->GetNumModes();
1711 int numMax = nmodes0;
1727 b0, b1, coeff_tmp, bortho0, bortho1, coeff);
1732 for (
int i = 0; i < numMin+1; ++i)
1736 tmp2 = coeff_tmp+cnt,1);
1742 bortho0, bortho1, coeff_tmp, b0, b1, outarray);
1791 int nquad0 =
m_base[0]->GetNumPoints();
1792 int nquad1 =
m_base[1]->GetNumPoints();
1798 for(i = 0; i < nquad1; ++i)
1801 w0.get(),1,outarray.get()+i*nquad0,1);
1804 for(i = 0; i < nquad0; ++i)
1806 Vmath::Vmul(nquad1,outarray.get()+i,nquad0,w1.get(),1,
1807 outarray.get()+i,nquad0);
1815 boost::ignore_unused(standard);
1817 int np1 =
m_base[0]->GetNumPoints();
1818 int np2 =
m_base[1]->GetNumPoints();
1819 int np = max(np1,np2);
1826 for(
int i = 0; i < np-1; ++i)
1829 for(
int j = 0; j < np-1; ++j)
1831 conn[cnt++] = row +j;
1832 conn[cnt++] = row +j+1;
1833 conn[cnt++] = rowp1 +j;
1835 conn[cnt++] = rowp1 +j+1;
1836 conn[cnt++] = rowp1 +j;
1837 conn[cnt++] = row +j+1;
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
virtual int v_GetNverts() const
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)
#define ASSERTL0(condition, msg)
virtual LibUtilities::BasisType v_GetEdgeBasisType(const int i) const
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
const int kSVVDGFiltermodesmax
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
static Array< OneD, NekDouble > NullNekDouble1DArray
void BwdTrans_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
#define sign(a, b)
return the sign(b)*a
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Calculate the derivative of the physical points.
void IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int NumBndryCoeffs(void) const
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey)
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
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...
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
virtual void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Principle Modified Functions .
virtual bool v_IsBoundaryInteriorExpansion()
Lagrange Polynomials using the Gauss points .
virtual int v_GetEdgeNumPoints(const int i) const
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
std::map< ConstFactorType, NekDouble > ConstFactorMap
NekDouble Integral(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &w0, const Array< OneD, const NekDouble > &w1)
virtual void v_GetInteriorMap(Array< OneD, unsigned int > &outarray)
virtual void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
const int kSVVDGFiltermodesmin
std::shared_ptr< DNekMat > DNekMatSharedPtr
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLess > m_stdStaticCondMatrixManager
MatrixType GetMatrixType() const
std::shared_ptr< Basis > BasisSharedPtr
static void Dgemm(const char &transa, const char &transb, const int &m, const int &n, const int &k, const double &alpha, const double *a, const int &lda, const double *b, const int &ldb, const double &beta, double *c, const int &ldc)
BLAS level 3: Matrix-matrix multiply C = A x B where A[m x n], B[n x k], C[m x k].
virtual void v_StdPhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d0, Array< OneD, NekDouble > &outarray_d1)
Calculate the 2D derivative in the local tensor/collapsed coordinate at the physical points...
1D Gauss-Gauss-Legendre quadrature points
virtual void v_FwdTrans_BndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Transform a given function from physical quadrature space to coefficient space.
void v_GeneralMatrixOp_MatOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
BasisManagerT & BasisManager(void)
virtual void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray)
virtual LibUtilities::ShapeType v_DetShapeType() const
void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1)
The base class for all shapes.
int GetEdgeNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th edge.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
virtual void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
const LibUtilities::BasisSharedPtr & GetBasis(int dir) const
This function gets the shared point to basis in the dir direction.
virtual void v_IProductWRTBase_MatOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual int v_DetCartesianDirOfEdge(const int edge)
virtual void v_ExponentialFilter(Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
void WeakDerivMatrixOp_MatFree(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
virtual void v_GetEdgeInteriorMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
Principle Orthogonal Functions .
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void v_GetEdgeToElementMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-1)
virtual void v_LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Defines a specification for a set of points.
virtual void v_WeakDerivMatrixOp(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
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 = A x where A[m x n].
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2)
NekDouble GetConstFactor(const ConstFactorType &factor) const
virtual void v_BwdTrans_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1)
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
virtual int v_NumBndryCoeffs() const
const NekDouble kSVVDGFilter[9][11]
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.
virtual void v_IProductWRTDerivBase_MatOp(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
const ConstFactorMap & GetConstFactors() const
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
virtual int v_CalcNumberOfCoefficients(const std::vector< unsigned int > &nummodes, int &modes_offset)
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
void GetEdgeToElementMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-1)
bool ConstFactorExists(const ConstFactorType &factor) const
int getNumberOfCoefficients(int Na, int Nb)
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Calculate the inner product of inarray with respect to the basis B=base0*base1 and put into outarray...
void MassMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false)
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space...
virtual int v_GetNedges() const
virtual const LibUtilities::BasisKey v_DetEdgeBasisKey(const int i) const
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
LibUtilities::NekManager< StdMatrixKey, DNekMat, StdMatrixKey::opLess > m_stdMatrixManager
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs
NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
Integrates the specified function over the domain.
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
virtual int v_GetEdgeNcoeffs(const int i) const
virtual void v_FillMode(const int mode, Array< OneD, NekDouble > &array)
Fill outarray with mode mode of expansion.
void IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
virtual void v_MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
std::shared_ptr< StdSegExp > StdSegExpSharedPtr
void Zero(int n, T *x, const int incx)
Zero vector.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void InterpCoeff2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
static void Dcopy(const int &n, const double *x, const int &incx, double *y, const int &incy)
BLAS level 1: Copy x to y.
LibUtilities::BasisType GetEdgeBasisType(const int i) const
This function returns the type of expansion basis on the i-th edge.
Describes the specification for a Basis.
1D Gauss-Lobatto-Legendre quadrature points
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.
virtual void v_GetSimplexEquiSpacedConnectivity(Array< OneD, int > &conn, bool standard=true)
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
virtual int v_NumDGBndryCoeffs() const
virtual void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space...