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.size()>=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.size()>=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);
597 int nquad0 =
m_base[0]->GetNumPoints();
598 int nquad1 =
m_base[1]->GetNumPoints();
601 int btmp0 =
m_base[0]->GetNumModes();
602 int mode0 = mode%btmp0;
603 int mode1 = mode/btmp0;
605 ASSERTL2(mode1 == (
int)floor((1.0*mode)/btmp0),
606 "Integer Truncation not Equiv to Floor");
609 "calling argument mode is larger than total expansion order");
611 for(i = 0; i < nquad1; ++i)
614 1, &outarray[0]+i*nquad0,1);
617 for(i = 0; i < nquad0; ++i)
620 &outarray[0]+i,nquad0,&outarray[0]+i,nquad0);
640 ASSERTL2((i >= 0)&&(i <= 3),
"edge id is out of range");
642 if((i == 0)||(i == 2))
654 ASSERTL2((i >= 0)&&(i <= 3),
"edge id is out of range");
656 if((i == 0)||(i == 2))
669 boost::ignore_unused(j);
670 ASSERTL2((i >= 0)&&(i <= 3),
"edge id is out of range");
672 if((i == 0)||(i == 2))
694 "BasisType is not a boundary interior form");
698 "BasisType is not a boundary interior form");
708 "BasisType is not a boundary interior form");
712 "BasisType is not a boundary interior form");
718 const std::vector<unsigned int> &nummodes,
721 int nmodes = nummodes[modes_offset]*nummodes[modes_offset+1];
729 bool returnval =
false;
748 boost::ignore_unused(coords_2);
755 for(i = 0; i < nq1; ++i)
757 Blas::Dcopy(nq0,z0.get(), 1,&coords_0[0] + i*nq0,1);
787 const int nm0 =
m_base[0]->GetNumModes();
788 const int nm1 =
m_base[1]->GetNumModes();
790 return StdExpansion::BaryEvaluateBasis<0>(coords[0], mode % nm1) *
791 StdExpansion::BaryEvaluateBasis<1>(coords[1], mode / nm0);
802 int nummodes0, nummodes1;
803 int value1 = 0, value2 = 0;
809 nummodes0 =
m_base[0]->GetNumModes();
810 nummodes1 =
m_base[1]->GetNumModes();
822 value1 = 2*nummodes0;
825 ASSERTL0(0,
"Mapping array is not defined for this expansion");
829 for(i = 0; i < value1; i++)
839 value2 = value1+nummodes0-1;
845 ASSERTL0(0,
"Mapping array is not defined for this expansion");
849 for(i = 0; i < nummodes1-2; i++)
851 outarray[cnt++]=value1+i*nummodes0;
852 outarray[cnt++]=value2+i*nummodes0;
858 for(i = nummodes0*(nummodes1-1);i <
GetNcoeffs(); i++)
869 int nummodes0, nummodes1;
876 nummodes0 =
m_base[0]->GetNumModes();
877 nummodes1 =
m_base[1]->GetNumModes();
885 startvalue = nummodes0;
888 startvalue = 2*nummodes0;
891 ASSERTL0(0,
"Mapping array is not defined for this expansion");
904 ASSERTL0(0,
"Mapping array is not defined for this expansion");
908 for(i = 0; i < nummodes1-2; i++)
910 for(j = 0; j < nummodes0-2; j++)
912 outarray[cnt++]=startvalue+j;
914 startvalue+=nummodes0;
922 if(useCoeffPacking ==
true)
924 switch(localVertexId)
935 localDOF =
m_base[0]->GetNumModes()-1;
947 localDOF =
m_base[0]->GetNumModes() * (
m_base[1]->GetNumModes()-1);
951 localDOF =
m_base[0]->GetNumModes();
959 localDOF =
m_base[0]->GetNumModes()*
m_base[1]->GetNumModes()-1;
963 localDOF =
m_base[0]->GetNumModes()+1;
968 ASSERTL0(
false,
"eid must be between 0 and 3");
975 switch(localVertexId)
986 localDOF =
m_base[0]->GetNumModes()-1;
998 localDOF =
m_base[0]->GetNumModes()*
m_base[1]->GetNumModes()-1;
1002 localDOF =
m_base[0]->GetNumModes()+1;
1010 localDOF =
m_base[0]->GetNumModes() * (
m_base[1]->GetNumModes()-1);
1014 localDOF =
m_base[0]->GetNumModes();
1019 ASSERTL0(
false,
"eid must be between 0 and 3");
1033 boost::ignore_unused(Q);
1036 int order0 =
m_base[0]->GetNumModes();
1037 int order1 =
m_base[1]->GetNumModes();
1050 ASSERTL0(
false,
"eid must be between 0 and 3");
1053 bool checkForZeroedModes =
false;
1058 else if(
P != numModes)
1060 checkForZeroedModes =
true;
1064 if (maparray.size() !=
P)
1069 if(signarray.size() !=
P)
1075 fill(signarray.get(), signarray.get()+
P, 1);
1086 for (i = 0; i <
P; i++)
1094 for (i = 0; i <
P; i++)
1096 maparray[i] = i*order0 + 1;
1102 for (i = 0; i <
P; i++)
1104 maparray[i] = order0+i;
1110 for (i = 0; i <
P; i++)
1112 maparray[i] = i*order0;
1117 ASSERTL0(
false,
"eid must be between 0 and 3");
1123 swap(maparray[0], maparray[1]);
1125 for(i = 3; i <
P; i+=2)
1138 for (i = 0; i <
P; i++)
1146 for (i = 0; i <
P; i++)
1148 maparray[i] = (i+1)*order0 - 1;
1154 for (i = 0; i <
P; i++)
1156 maparray[i] = order0*(order1-1) + i;
1162 for (i = 0; i <
P; i++)
1164 maparray[i] = order0*i;
1169 ASSERTL0(
false,
"eid must be between 0 and 3");
1175 reverse(maparray.get(), maparray.get()+
P);
1180 ASSERTL0(
false,
"Mapping not defined for this type of basis");
1183 if (checkForZeroedModes)
1189 for (
int j = numModes; j <
P; j++)
1192 maparray[j] = maparray[0];
1197 ASSERTL0(
false,
"Different trace space edge dimension "
1198 "and element edge dimension not possible "
1199 "for GLL-Lagrange bases");
1211 const int nummodes0 =
m_base[0]->GetNumModes();
1212 const int nummodes1 =
m_base[1]->GetNumModes();
1216 if(maparray.size() != nEdgeIntCoeffs)
1221 if(signarray.size() != nEdgeIntCoeffs)
1227 fill( signarray.get() , signarray.get()+nEdgeIntCoeffs, 1 );
1236 for(i = 0; i < nEdgeIntCoeffs; i++)
1244 for(i = 0; i < nEdgeIntCoeffs; i++)
1246 maparray[i] = (i+2)*nummodes0 + 1;
1252 for(i = 0; i < nEdgeIntCoeffs; i++)
1254 maparray[i] = nummodes0+i+2;
1260 for(i = 0; i < nEdgeIntCoeffs; i++)
1262 maparray[i] = (i+2)*nummodes0;
1267 ASSERTL0(
false,
"eid must be between 0 and 3");
1273 for(i = 1; i < nEdgeIntCoeffs; i+=2)
1285 for(i = 0; i < nEdgeIntCoeffs; i++)
1293 for(i = 0; i < nEdgeIntCoeffs; i++)
1295 maparray[i] = (i+2)*nummodes0 - 1;
1301 for(i = 0; i < nEdgeIntCoeffs; i++)
1303 maparray[i] = nummodes0*(nummodes1-1) + i + 1;
1309 for(i = 0; i < nEdgeIntCoeffs; i++)
1311 maparray[i] = nummodes0 * (i+1);
1316 ASSERTL0(
false,
"eid must be between 0 and 3");
1321 reverse( maparray.get() , maparray.get()+nEdgeIntCoeffs );
1326 ASSERTL0(
false,
"Mapping not defined for this type of basis");
1348 int nq0 =
m_base[0]->GetNumPoints();
1349 int nq1 =
m_base[1]->GetNumPoints();
1373 for(
int i = 0; i < nq; ++i)
1375 for(
int j = 0; j < nq; ++j,++cnt)
1378 coords[cnt][0] = -1.0 + 2*j/(
NekDouble)(nq-1);
1379 coords[cnt][1] = -1.0 + 2*i/(
NekDouble)(nq-1);
1383 for(
int i = 0; i < neq; ++i)
1387 I[0] =
m_base[0]->GetI(coll);
1388 I[1] =
m_base[1]->GetI(coll+1);
1391 for (
int j = 0; j < nq1; ++j)
1397 Mat->GetRawPtr()+j*nq0*neq+i,neq);
1410 for(i = 0; i < order1; ++i)
1412 (*Mat)(order0*i+1,i*order0+1) = 1.0;
1418 for(i = 0; i < order0; ++i)
1420 (*Mat)(order0+i ,order0+i) = 1.0;
1434 (*Mat) = Imass*Iprod;
1442 int dir = (edge + 1) % 2;
1443 int nCoeffs =
m_base[dir]->GetNumModes();
1451 coords[0] = (edge == 0 || edge == 3) ? -1.0 : 1.0;
1459 Vmath::Vcopy(nCoeffs, m_Ix->GetPtr(), 1, Mat->GetPtr(), 1);
1489 if(inarray.get() == outarray.get())
1495 m_ncoeffs, tmp.get(), 1, 0.0, outarray.get(), 1);
1500 m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
1509 int qa =
m_base[0]->GetNumPoints();
1510 int qb =
m_base[1]->GetNumPoints();
1511 int nmodes_a =
m_base[0]->GetNumModes();
1512 int nmodes_b =
m_base[1]->GetNumModes();
1513 int nmodes = min(nmodes_a,nmodes_b);
1526 OrthoExp.
FwdTrans(array,orthocoeffs);
1535 for(
int j = 0; j < nmodes_a; ++j)
1537 for(
int k = 0; k < nmodes_b; ++k)
1541 pow((1.0*j)/(nmodes_a-1),cutoff*nmodes_a),
1542 pow((1.0*k)/(nmodes_b-1),cutoff*nmodes_b));
1544 orthocoeffs[j*nmodes_b+k] *= SvvDiffCoeff * fac;
1555 max_ab = max(max_ab,0);
1558 for(
int j = 0; j < nmodes_a; ++j)
1560 for(
int k = 0; k < nmodes_b; ++k)
1562 int maxjk = max(j,k);
1565 orthocoeffs[j*nmodes_b+k] *= SvvDiffCoeff *
1575 min(nmodes_a,nmodes_b));
1581 for(
int j = 0; j < nmodes_a; ++j)
1583 for(
int k = 0; k < nmodes_b; ++k)
1587 orthocoeffs[j*nmodes_b+k] *=
1588 (SvvDiffCoeff*exp(-(j+k-nmodes)*(j+k-nmodes)/
1594 orthocoeffs[j*nmodes_b+k] *= 0.0;
1602 OrthoExp.
BwdTrans(orthocoeffs,array);
1612 int qa =
m_base[0]->GetNumPoints();
1613 int qb =
m_base[1]->GetNumPoints();
1614 int nmodesA =
m_base[0]->GetNumModes();
1615 int nmodesB =
m_base[1]->GetNumModes();
1616 int P = nmodesA - 1;
1617 int Q = nmodesB - 1;
1628 int Pcut = cutoff*
P;
1629 int Qcut = cutoff*Q;
1633 OrthoExp.
FwdTrans(array,orthocoeffs);
1637 for(
int i = 0; i < nmodesA; ++i)
1639 for(
int j = 0; j < nmodesB; ++j)
1642 if(i > Pcut || j > Qcut)
1646 fac = max(fac1, fac2);
1647 fac = pow(fac, exponent);
1648 orthocoeffs[i*nmodesB+j] *= exp(-alpha*fac);
1654 OrthoExp.
BwdTrans(orthocoeffs,array);
1662 int n_coeffs = inarray.size();
1670 int nmodes0 =
m_base[0]->GetNumModes();
1671 int nmodes1 =
m_base[1]->GetNumModes();
1672 int numMax = nmodes0;
1688 b0, b1, coeff_tmp, bortho0, bortho1, coeff);
1693 for (
int i = 0; i < numMin+1; ++i)
1697 tmp2 = coeff_tmp+cnt,1);
1703 bortho0, bortho1, coeff_tmp, b0, b1, outarray);
1752 int nquad0 =
m_base[0]->GetNumPoints();
1753 int nquad1 =
m_base[1]->GetNumPoints();
1759 for(i = 0; i < nquad1; ++i)
1762 w0.get(),1,outarray.get()+i*nquad0,1);
1765 for(i = 0; i < nquad0; ++i)
1767 Vmath::Vmul(nquad1,outarray.get()+i,nquad0,w1.get(),1,
1768 outarray.get()+i,nquad0);
1776 boost::ignore_unused(standard);
1778 int np1 =
m_base[0]->GetNumPoints();
1779 int np2 =
m_base[1]->GetNumPoints();
1780 int np = max(np1,np2);
1787 for(
int i = 0; i < np-1; ++i)
1790 for(
int j = 0; j < np-1; ++j)
1792 conn[cnt++] = row +j;
1793 conn[cnt++] = row +j+1;
1794 conn[cnt++] = rowp1 +j;
1796 conn[cnt++] = rowp1 +j+1;
1797 conn[cnt++] = rowp1 +j;
1798 conn[cnt++] = row +j+1;
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
#define sign(a, b)
return the sign(b)*a
Describes the specification for a Basis.
Defines a specification for a set of points.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d0, Array< OneD, NekDouble > &outarray_d1)
Calculate the 2D derivative in the local tensor/collapsed coordinate at the physical points.
NekDouble Integral(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &w0, const Array< OneD, const NekDouble > &w1)
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)
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
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)
The base class for all shapes.
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
const LibUtilities::BasisSharedPtr & GetBasis(int dir) const
This function gets the shared point to basis in the dir direction.
void WeakDerivMatrixOp_MatFree(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
int NumBndryCoeffs(void) const
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
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 LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void GetTraceToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space.
int GetTraceNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th trace.
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space.
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLess > m_stdStaticCondMatrixManager
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
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.
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void MassMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
const ConstFactorMap & GetConstFactors() const
MatrixType GetMatrixType() const
NekDouble GetConstFactor(const ConstFactorType &factor) const
bool ConstFactorExists(const ConstFactorType &factor) const
virtual int v_NumBndryCoeffs() const
virtual void v_IProductWRTDerivBase(const int dir, 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.
virtual const LibUtilities::BasisKey v_GetTraceBasisKey(const int i, const int j) const
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_IProductWRTBase_MatOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
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)
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.
virtual void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray)
void v_GeneralMatrixOp_MatOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_GetTraceToElementMap(const int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient=eForwards, int P=-1, int Q=-1)
virtual int v_NumDGBndryCoeffs() const
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
virtual void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
Integrates the specified function over the domain.
virtual void v_ExponentialFilter(Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
virtual void v_GetSimplexEquiSpacedConnectivity(Array< OneD, int > &conn, bool standard=true)
virtual LibUtilities::ShapeType v_DetShapeType() const
virtual void v_FwdTrans_BndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual int v_CalcNumberOfCoefficients(const std::vector< unsigned int > &nummodes, int &modes_offset)
virtual void v_WeakDerivMatrixOp(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey)
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.
virtual bool v_IsBoundaryInteriorExpansion()
virtual void v_LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
virtual int v_GetNverts() const
virtual int v_GetNtraces() const
virtual void v_FillMode(const int mode, Array< OneD, NekDouble > &array)
Fill outarray with mode mode of expansion.
virtual int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false)
virtual NekDouble v_PhysEvaluateBasis(const Array< OneD, const NekDouble > &coords, int mode)
This function evaluates the basis function mode mode at a point coords of the domain.
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)
virtual void v_LocCollapsedToLocCoord(const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi)
virtual void v_GetInteriorMap(Array< OneD, unsigned int > &outarray)
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2)
virtual void v_IProductWRTDerivBase_MatOp(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual int v_GetTraceNcoeffs(const int i) const
virtual int v_GetTraceNumPoints(const int i) const
virtual void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
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)
virtual void v_GetTraceInteriorToElementMap(const int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation edgeOrient=eForwards)
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
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].
static void Dcopy(const int &n, const double *x, const int &incx, double *y, const int &incy)
BLAS level 1: Copy x to y.
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
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 op(A)[m x k], op(B)[k x n], C[m x n] DGEMM perfo...
int getNumberOfCoefficients(int Na, int Nb)
BasisManagerT & BasisManager(void)
std::shared_ptr< Basis > BasisSharedPtr
void InterpCoeff2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
@ eGaussGaussLegendre
1D Gauss-Gauss-Legendre quadrature points
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points .
@ eOrtho_A
Principle Orthogonal Functions .
@ eGLL_Lagrange
Lagrange for SEM basis .
@ eModified_A
Principle Modified Functions .
@ eFourier
Fourier Expansion .
static const NekDouble kNekZeroTol
@ eFactorSVVDGKerDiffCoeff
@ eFactorSVVPowerKerDiffCoeff
const int kSVVDGFiltermodesmin
const int kSVVDGFiltermodesmax
const NekDouble kSVVDGFilter[9][11]
@ ePhysInterpToEquiSpaced
std::map< ConstFactorType, NekDouble > ConstFactorMap
std::shared_ptr< StdSegExp > StdSegExpSharedPtr
The above copyright notice and this permission notice shall be included.
static Array< OneD, NekDouble > NullNekDouble1DArray
std::shared_ptr< DNekMat > DNekMatSharedPtr
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
void Zero(int n, T *x, const int incx)
Zero vector.
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
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.