42 namespace LocalRegions
50 StdNodalTriExp(Ba,Bb,Ntype),
54 std::bind(&
NodalTriExp::CreateMatrix, this, std::placeholders::_1),
55 std::string(
"NodalTriExpMatrix")),
56 m_staticCondMatrixManager(
57 std::bind(&
NodalTriExp::CreateStaticCondMatrix, this, std::placeholders::_1),
58 std::string(
"NodalTriExpStaticCondMatrix"))
65 StdRegions::StdNodalTriExp(T),
68 m_matrixManager(T.m_matrixManager),
69 m_staticCondMatrixManager(T.m_staticCondMatrixManager)
100 int nquad0 =
m_base[0]->GetNumPoints();
101 int nquad1 =
m_base[1]->GetNumPoints();
109 Vmath::Vmul(nquad0*nquad1, jac, 1, inarray, 1,tmp, 1);
113 Vmath::Smul(nquad0*nquad1, jac[0], inarray, 1, tmp, 1);
117 ival = StdNodalTriExp::v_Integral(tmp);
124 bool multiplybyweights)
126 int nquad0 =
m_base[0]->GetNumPoints();
127 int nquad1 =
m_base[1]->GetNumPoints();
128 int order1 =
m_base[1]->GetNumModes();
130 if(multiplybyweights)
136 StdTriExp::IProductWRTBase_SumFacKernel(
m_base[0]->GetBdata(),
m_base[1]->GetBdata(),tmp,outarray,wsp);
143 StdTriExp::IProductWRTBase_SumFacKernel(
m_base[0]->GetBdata(),
m_base[1]->GetBdata(),inarray,outarray,wsp);
156 m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
163 int nquad0 =
m_base[0]->GetNumPoints();
164 int nquad1 =
m_base[1]->GetNumPoints();
165 int nqtot = nquad0*nquad1;
194 ASSERTL1((dir==0)||(dir==1)||(dir==2),
"Invalid direction.");
196 "Invalid direction.");
198 int nquad0 =
m_base[0]->GetNumPoints();
199 int nquad1 =
m_base[1]->GetNumPoints();
200 int nqtot = nquad0*nquad1;
218 for (
int i = 0; i < nquad1; ++i)
220 gfac0[i] = 2.0/(1-z1[i]);
222 for (
int i = 0; i < nquad0; ++i)
224 gfac1[i] = 0.5*(1+z0[i]);
227 for (
int i = 0; i < nquad1; ++i)
229 Vmath::Smul(nquad0, gfac0[i], &inarray[0]+i*nquad0, 1,
230 &tmp0[0]+i*nquad0, 1);
233 for (
int i = 0; i < nquad1; ++i)
235 Vmath::Vmul(nquad0, &gfac1[0], 1, &tmp0[0]+i*nquad0, 1,
236 &tmp1[0]+i*nquad0, 1);
241 Vmath::Vmul(nqtot,&df[2*dir][0], 1,&tmp0[0], 1,&tmp0[0],1);
242 Vmath::Vmul(nqtot,&df[2*dir+1][0],1,&tmp1[0], 1,&tmp1[0],1);
243 Vmath::Vmul(nqtot,&df[2*dir+1][0],1,&inarray[0],1,&tmp2[0],1);
247 Vmath::Smul(nqtot, df[2*dir][0], tmp0, 1, tmp0, 1);
248 Vmath::Smul(nqtot, df[2*dir+1][0], tmp1, 1, tmp1, 1);
249 Vmath::Smul(nqtot, df[2*dir+1][0], inarray, 1, tmp2, 1);
280 ASSERTL1(
false,
"input dir is out of range");
289 m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
305 int nquad0 =
m_base[0]->GetNumPoints();
306 int nquad1 =
m_base[1]->GetNumPoints();
307 int nqtot = nquad0*nquad1;
314 StdNodalTriExp::v_PhysDeriv(inarray, diff0, diff1);
321 Vmath::Vvtvp (nqtot,df[1],1,diff1,1, out_d0, 1, out_d0,1);
327 Vmath::Vvtvp (nqtot,df[3],1,diff1,1, out_d1, 1, out_d1,1);
333 Vmath::Vvtvp (nqtot,df[5],1,diff1,1, out_d2, 1, out_d2,1);
340 Vmath::Smul (nqtot, df[0][0], diff0, 1, out_d0, 1);
341 Blas::Daxpy (nqtot, df[1][0], diff1, 1, out_d0, 1);
346 Vmath::Smul (nqtot, df[2][0], diff0, 1, out_d1, 1);
347 Blas::Daxpy (nqtot, df[3][0], diff1, 1, out_d1, 1);
352 Vmath::Smul (nqtot, df[4][0], diff0, 1, out_d2, 1);
353 Blas::Daxpy (nqtot, df[5][0], diff1, 1, out_d2, 1);
393 if(inarray.get() == outarray.get())
399 m_ncoeffs, tmp.get(), 1, 0.0, outarray.get(), 1);
404 m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
421 ASSERTL1(Lcoords[0] >= -1.0 && Lcoords[1] <= 1.0 &&
422 Lcoords[1] >= -1.0 && Lcoords[1] <=1.0,
423 "Local coordinates are not in region [-1,1]");
427 for(i = 0; i <
m_geom->GetCoordim(); ++i)
429 coords[i] =
m_geom->GetCoord(i,Lcoords);
441 return tmp->GetStdMatrix(mkey);
452 m_geom->GetLocCoords(coord,Lcoord);
454 return StdNodalTriExp::v_PhysEvaluate(Lcoord, physvals);
515 ASSERTL1(
m_geom->GetCoordim() == 2,
"Standard Region Laplacian is only set up for Quads in two-dimensional");
531 int rows = lap00.GetRows();
532 int cols = lap00.GetColumns();
536 (*lap) = gmat[0][0] * lap00 +
537 gmat[1][0] * (lap01 +
Transpose(lap01)) +
554 int rows = LapMat.GetRows();
555 int cols = LapMat.GetColumns();
560 (*helm) = LapMat + factor*MassMat;
581 unsigned int nint = (
unsigned int)(
m_ncoeffs - nbdry);
582 unsigned int exp_size[] = {nbdry,nint};
583 unsigned int nblks = 2;
594 goto UseLocRegionsMatrix;
600 goto UseLocRegionsMatrix;
605 factor = mat->Scale();
606 goto UseStdRegionsMatrix;
639 for(i = 0; i < nbdry; ++i)
641 for(j = 0; j < nbdry; ++j)
643 (*A)(i,j) = mat(bmap[i],bmap[j]);
646 for(j = 0; j < nint; ++j)
648 (*B)(i,j) = mat(bmap[i],imap[j]);
652 for(i = 0; i < nint; ++i)
654 for(j = 0; j < nbdry; ++j)
656 (*C)(i,j) = mat(imap[i],bmap[j]);
659 for(j = 0; j < nint; ++j)
661 (*D)(i,j) = mat(imap[i],imap[j]);
670 (*A) = (*A) - (*B)*(*C);
699 2,
m_base[0]->GetPointsKey());
701 2,
m_base[1]->GetPointsKey());
722 returnval = StdNodalTriExp::v_GenMatrix(mkey);
737 int nqe =
m_base[0]->GetNumPoints();
742 for (i = 0; i < dim; ++i)
768 Vmath::Fill(nqe,df[2*i+1][0] + df[2*i][0],normal[i],1);
778 ASSERTL0(
false,
"Edge is out of range (edge < 3)");
785 fac += normal[i][0]*normal[i][0];
800 int nquad0 = ptsKeys[0].GetNumPoints();
801 int nquad1 = ptsKeys[1].GetNumPoints();
814 for(j = 0; j < nquad0; ++j)
819 normals[i*nquad0+j] = -df[2*i+1][j]*edgejac[j];
822 from_key = ptsKeys[0];
825 for(j = 0; j < nquad1; ++j)
827 edgejac[j] = jac[nquad0*j+nquad0-1];
830 normals[i*nquad1+j] = (df[2*i][nquad0*j + nquad0-1] + df[2*i+1][nquad0*j + nquad0-1])*edgejac[j];
833 from_key = ptsKeys[1];
836 for(j = 0; j < nquad1; ++j)
838 edgejac[j] = jac[nquad0*j];
841 normals[i*nquad1+j] = -df[2*i][nquad0*j]*edgejac[j];
844 from_key = ptsKeys[1];
847 ASSERTL0(
false,
"edge is out of range (edge < 3)");
869 Vmath::Vvtvp(nqe,normal[i],1, normal[i],1,work,1,work,1);
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Describes the specification for a Basis.
Defines a specification for a set of points.
PointsType GetPointsType() const
unsigned int GetNumPoints() const
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey)
std::map< int, NormalVector > m_edgeNormals
std::map< int, Array< OneD, NekDouble > > m_elmtBndNormDirElmtLen
the element length in each element boundary(Vertex, edge or face) normal direction calculated based o...
SpatialDomains::GeometrySharedPtr GetGeom() const
SpatialDomains::GeometrySharedPtr m_geom
void AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3)
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
~NodalTriExp()
Destructor.
DNekScalMatSharedPtr CreateMatrix(const MatrixKey &mkey)
void GetCoord(const Array< OneD, const NekDouble > &Lcoords, Array< OneD, NekDouble > &coords)
void PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Differentiation Methods.
DNekMatSharedPtr CreateStdMatrix(const StdRegions::StdMatrixKey &mkey)
void GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
void GeneralMatrixOp_MatOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Inner product of inarray over region with respect to the expansion basis (this)->_Base[0] and return ...
void IProductWRTDerivBase_MatOp(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Forward transform from physical quadrature space stored in inarray and evaluate the expansion coeffic...
DNekScalBlkMatSharedPtr CreateStaticCondMatrix(const MatrixKey &mkey)
void IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
NekDouble PhysEvaluate(const Array< OneD, const NekDouble > &coord, const Array< OneD, const NekDouble > &physvals)
NekDouble Integral(const Array< OneD, const NekDouble > &inarray)
Integrate the physical point list inarray over region.
void IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
virtual void v_AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
void v_ComputeTraceNormal(const int edge)
void IProductWRTBase_MatOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual StdRegions::StdExpansionSharedPtr v_GetStdExp(void) const
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey)
NodalTriExp(const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::PointsType Ntype, const SpatialDomains::TriGeomSharedPtr &geom)
Constructor using BasisKey class for quadrature points and order definition.
virtual StdRegions::StdExpansionSharedPtr v_GetLinStdExp(void) const
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
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)
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
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)
const LibUtilities::PointsKeyVector GetPointsKeys() const
DNekBlkMatSharedPtr GetStdStaticCondMatrix(const StdMatrixKey &mkey)
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Array< OneD, LibUtilities::BasisSharedPtr > m_base
LibUtilities::ShapeType GetShapeType() const
MatrixType GetMatrixType() const
NekDouble GetConstFactor(const ConstFactorType &factor) const
void NodalToModalTranspose(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
LibUtilities::PointsKey m_nodalPointsKey
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 Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
int getNumberOfCoefficients(int Na)
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 ...
std::vector< PointsKey > PointsKeyVector
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
GeomType
Indicates the type of element geometry.
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eNoGeomType
No type defined.
@ eMovingRegular
Currently unused.
@ eDeformed
Geometry is curved or has non-constant factors.
std::shared_ptr< TriGeom > TriGeomSharedPtr
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
std::shared_ptr< StdNodalTriExp > StdNodalTriExpSharedPtr
static ConstFactorMap NullConstFactorMap
static VarCoeffMap NullVarCoeffMap
The above copyright notice and this permission notice shall be included.
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
std::shared_ptr< DNekMat > DNekMatSharedPtr
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
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 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
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add 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 Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/y.
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 Reverse(int n, const T *x, const int incx, T *y, const int incy)
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
scalarT< T > sqrt(scalarT< T > in)