76 int nquad0 =
m_base[0]->GetNumPoints();
77 int nquad1 =
m_base[1]->GetNumPoints();
79 if (outarray_d0.num_elements() > 0)
82 if(inarray.data() == outarray_d0.data())
85 Vmath::Vcopy(nquad0 * nquad1,inarray.get(),1,wsp.get(),1);
86 Blas::Dgemm(
'N',
'N', nquad0, nquad1, nquad0, 1.0,
87 &(D0->GetPtr())[0], nquad0, &wsp[0], nquad0, 0.0,
88 &outarray_d0[0], nquad0);
92 Blas::Dgemm(
'N',
'N', nquad0, nquad1, nquad0, 1.0,
93 &(D0->GetPtr())[0], nquad0, &inarray[0], nquad0, 0.0,
94 &outarray_d0[0], nquad0);
98 if (outarray_d1.num_elements() > 0)
101 if(inarray.data() == outarray_d1.data())
104 Vmath::Vcopy(nquad0 * nquad1,inarray.get(),1,wsp.get(),1);
105 Blas:: Dgemm(
'N',
'T', nquad0, nquad1, nquad1, 1.0, &wsp[0], nquad0,
106 &(D1->GetPtr())[0], nquad1, 0.0, &outarray_d1[0], nquad0);
110 Blas:: Dgemm(
'N',
'T', nquad0, nquad1, nquad1, 1.0, &inarray[0], nquad0,
111 &(D1->GetPtr())[0], nquad1, 0.0, &outarray_d1[0], nquad0);
128 I[0] =
m_base[0]->GetI(coll);
129 I[1] =
m_base[1]->GetI(coll+1);
140 int nq0 =
m_base[0]->GetNumPoints();
141 int nq1 =
m_base[1]->GetNumPoints();
145 for (i = 0; i < nq1;++i)
147 wsp1[i] =
Blas::Ddot(nq0, &(I[0]->GetPtr())[0], 1,
148 &physvals[i * nq0], 1);
152 val =
Blas::Ddot(nq1, I[1]->GetPtr(), 1, wsp1, 1);
167 int nquad0 =
m_base[0]->GetNumPoints();
168 int nquad1 =
m_base[1]->GetNumPoints();
172 for (i = 0; i < nquad1; ++i)
174 Vmath::Vmul(nquad0, &inarray[0] + i*nquad0, 1, w0.get(),
175 1, &tmp[0] + i*nquad0, 1);
178 for (i = 0; i < nquad0; ++i)
180 Vmath::Vmul(nquad1, &tmp[0]+ i, nquad0, w1.get(), 1,
181 &tmp[0] + i, nquad0);
194 bool doCheckCollDir0,
195 bool doCheckCollDir1)
206 bool doCheckCollDir0,
207 bool doCheckCollDir1)
222 int nquad0 =
m_base[0]->GetNumPoints();
223 int nquad1 =
m_base[1]->GetNumPoints();
224 int nqtot = nquad0*nquad1;
225 int nmodes0 =
m_base[0]->GetNumModes();
226 int nmodes1 =
m_base[1]->GetNumModes();
227 int wspsize = max(max(max(nqtot,
m_ncoeffs),nquad1*nmodes0),nquad0*nmodes1);
236 if(!(
m_base[0]->Collocation() &&
m_base[1]->Collocation()))
253 inarray,outarray,mkey);
267 int nquad0 =
m_base[0]->GetNumPoints();
268 int nquad1 =
m_base[1]->GetNumPoints();
269 int nqtot = nquad0*nquad1;
270 int nmodes0 =
m_base[0]->GetNumModes();
271 int nmodes1 =
m_base[1]->GetNumModes();
272 int wspsize = max(max(max(nqtot,
m_ncoeffs),nquad1*nmodes0),
285 if (!(
m_base[0]->Collocation() &&
m_base[1]->Collocation()))
293 wsp0, wsp2,
true,
true);
309 &wsp1[0], 1, &outarray[0], 1);
314 inarray,outarray,mkey);
void HelmholtzMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
NekDouble GetConstFactor(const ConstFactorType &factor) const
void LaplacianMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
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)
void MultiplyByQuadratureMetric(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)=0
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
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
NekDouble Integral(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &w0, const Array< OneD, const NekDouble > &w1)
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
boost::shared_ptr< DNekMat > DNekMatSharedPtr
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...
bool ConstFactorExists(const ConstFactorType &factor) const
The base class for all shapes.
static const NekDouble kNekZeroTol
virtual NekDouble v_PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
This function evaluates the expansion at a single (arbitrary) point of the domain.
void LaplacianMatrixOp_MatFree_Kernel(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
T Ddot(int n, const Array< OneD, const T > &w, const int incw, const Array< OneD, const T > &x, const int incx, const Array< OneD, const int > &y, const int incy)
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
virtual ~StdExpansion2D()
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)
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Describes the specification for a Basis.
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)=0
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_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)