46 namespace MultiRegions
63 const bool dealiasing):
68 m_dealiasing(dealiasing)
71 m_comm = pSession->GetComm();
74 "Homogeneous Basis is a null basis");
79 m_comm->GetColumnComm()->GetColumnComm() :
84 "HomModesZ should be a multiple of npz.");
94 "HomModesZ/npz should be an even integer.");
123 ASSERTL0(
false,
"Dealiasing available just in combination "
135 m_transposition(In.m_transposition),
136 m_StripZcomm(In.m_StripZcomm),
137 m_useFFT(In.m_useFFT),
140 m_tmpOUT(In.m_tmpOUT),
141 m_homogeneousBasis(In.m_homogeneousBasis),
143 m_homogeneous1DBlockMat(In.m_homogeneous1DBlockMat),
144 m_dealiasing(In.m_dealiasing),
145 m_padsize(In.m_padsize)
151 const std::vector<unsigned int> &eIDs):
153 m_transposition(In.m_transposition),
154 m_useFFT(In.m_useFFT),
157 m_tmpOUT(In.m_tmpOUT),
158 m_homogeneousBasis(In.m_homogeneousBasis),
161 m_dealiasing(In.m_dealiasing),
162 m_padsize(In.m_padsize)
203 int num_dofs = inarray1.size();
221 int num_points_per_plane = num_dofs/
m_planes.size();
223 if(!
m_session->DefinesSolverInfo(
"HomoStrip"))
225 num_proc =
m_comm->GetColumnComm()->GetSize();
231 int num_dfts_per_proc = num_points_per_plane / num_proc
232 + (num_points_per_plane % num_proc > 0);
250 for(
int i = 0 ; i < num_dfts_per_proc ; i++)
254 Vmath::Vcopy(N, &(ShufV1[i*N]), 1, &(ShufV1_PAD_coef[0]), 1);
255 Vmath::Vcopy(N, &(ShufV2[i*N]), 1, &(ShufV2_PAD_coef[0]), 1);
258 m_FFT_deal->FFTBwdTrans(ShufV1_PAD_coef, ShufV1_PAD_phys);
259 m_FFT_deal->FFTBwdTrans(ShufV2_PAD_coef, ShufV2_PAD_phys);
265 ShufV1V2_PAD_phys, 1);
269 m_FFT_deal->FFTFwdTrans(ShufV1V2_PAD_phys, ShufV1V2_PAD_coef);
274 &(ShufV1V2[i*N]), 1);
305 int ndim = inarray1.size();
306 ASSERTL1( inarray2.size() % ndim == 0,
307 "Wrong dimensions for DealiasedDotProd.");
308 int nvec = inarray2.size() / ndim;
310 int num_dofs = inarray1[0].size();
313 int num_points_per_plane = num_dofs/
m_planes.size();
315 if(!
m_session->DefinesSolverInfo(
"HomoStrip"))
317 num_proc =
m_comm->GetColumnComm()->GetSize();
323 int num_dfts_per_proc = num_points_per_plane / num_proc
324 + (num_points_per_plane % num_proc > 0);
331 for (
int i = 0; i < ndim; i++)
335 for (
int i = 0; i < ndim*nvec; i++)
342 for (
int i = 0; i < ndim; i++)
347 for (
int i = 0; i < ndim*nvec; i++)
358 for (
int i = 0; i < ndim; i++)
361 (num_dfts_per_proc*N,0.0);
369 for (
int i = 0; i < ndim*nvec; i++)
372 (num_dfts_per_proc*N,0.0);
380 for (
int i = 0; i < nvec; i++)
383 (num_dfts_per_proc*N,0.0);
387 for (
int i = 0; i < ndim; i++)
392 for (
int i = 0; i < ndim*nvec; i++)
399 for(
int i = 0 ; i < num_dfts_per_proc ; i++)
401 for (
int j = 0; j < ndim; j++)
406 &(ShufV1_PAD_coef[0]), 1);
408 m_FFT_deal->FFTBwdTrans(ShufV1_PAD_coef, ShufV1_PAD_phys[j]);
410 for (
int j = 0; j < ndim*nvec; j++)
413 &(ShufV2_PAD_coef[0]), 1);
414 m_FFT_deal->FFTBwdTrans(ShufV2_PAD_coef, ShufV2_PAD_phys[j]);
419 for (
int j = 0; j < nvec; j++)
422 for (
int k = 0; k < ndim; k++)
425 ShufV2_PAD_phys[j*ndim+k], 1,
426 ShufV1V2_PAD_phys, 1,
427 ShufV1V2_PAD_phys, 1);
431 m_FFT_deal->FFTFwdTrans(ShufV1V2_PAD_phys, ShufV1V2_PAD_coef);
435 &(ShufV1V2[j][i*N]), 1);
442 for (
int j = 0; j < nvec; j++)
452 for (
int j = 0; j < nvec; j++)
466 int cnt = 0, cnt1 = 0;
469 for(
int n = 0; n <
m_planes.size(); ++n)
471 m_planes[n]->FwdTrans(inarray+cnt, tmparray = outarray + cnt1);
489 int cnt = 0, cnt1 = 0;
493 for(
int n = 0; n <
m_planes.size(); ++n)
495 m_planes[n]->FwdTrans_IterPerExp(inarray+cnt, tmparray = outarray + cnt1);
512 int cnt = 0, cnt1 = 0;
516 for(
int n = 0; n <
m_planes.size(); ++n)
518 m_planes[n]->FwdTrans_BndConstrained(inarray+cnt, tmparray = outarray + cnt1);
533 int cnt = 0, cnt1 = 0;
536 for(
int n = 0; n <
m_planes.size(); ++n)
538 m_planes[n]->BwdTrans(inarray+cnt, tmparray = outarray + cnt1);
540 cnt1 +=
m_planes[n]->GetTotPoints();
553 int cnt = 0, cnt1 = 0;
556 for(
int n = 0; n <
m_planes.size(); ++n)
558 m_planes[n]->BwdTrans_IterPerExp(inarray+cnt, tmparray = outarray + cnt1);
561 cnt1 +=
m_planes[n]->GetTotPoints();
574 int cnt = 0, cnt1 = 0;
587 for(
int n = 0; n <
m_planes.size(); ++n)
589 m_planes[n]->IProductWRTBase(tmpIn+cnt, tmparray = outarray + cnt1);
601 int cnt = 0, cnt1 = 0;
614 for(
int n = 0; n <
m_planes.size(); ++n)
616 m_planes[n]->IProductWRTBase_IterPerExp(tmpIn+cnt, tmparray = outarray + cnt1);
635 num_dofs = inarray.size();
639 num_dofs = outarray.size();
644 int num_points_per_plane = num_dofs/
m_planes.size();
645 int num_dfts_per_proc;
646 if(!
m_session->DefinesSolverInfo(
"HomoStrip"))
648 int nP =
m_comm->GetColumnComm()->GetSize();
649 num_dfts_per_proc = num_points_per_plane / nP
650 + (num_points_per_plane % nP > 0);
655 num_dfts_per_proc = num_points_per_plane / nP
656 + (num_points_per_plane % nP > 0);
674 for(
int i = 0 ; i < num_dfts_per_proc ; i++)
681 for(
int i = 0 ; i < num_dfts_per_proc ; i++)
694 fft_out,1,outarray,1);
724 int nrows = blkmat->GetRows();
725 int ncols = blkmat->GetColumns();
769 return matrixIter->second;
779 int num_trans_per_proc = 0;
788 n_exp =
m_planes[0]->GetTotPoints();
791 num_trans_per_proc = n_exp/
m_comm->GetColumnComm()->GetSize() + (n_exp%
m_comm->GetColumnComm()->GetSize() > 0);
858 for(
int i = 0; i < num_trans_per_proc; ++i)
860 BlkMatrix->SetBlock(i,i,loc_mat);
868 std::vector<LibUtilities::FieldDefinitionsSharedPtr> returnval;
873 std::vector<NekDouble> HomoLen;
874 HomoLen.push_back(
m_lhom);
876 std::vector<unsigned int> StripsIDs;
879 m_session->MatchSolverInfo(
"HomoStrip",
"True",strips,
false);
885 std::vector<unsigned int> PlanesIDs;
895 for(
int i = 0; i <
m_planes.size(); i++)
900 m_planes[0]->GeneralGetFieldDefinitions(returnval, 1, HomoBasis,
901 HomoLen, strips, StripsIDs, PlanesIDs);
910 std::vector<NekDouble> HomoLen;
911 HomoLen.push_back(
m_lhom);
913 std::vector<unsigned int> StripsIDs;
916 m_session->MatchSolverInfo(
"HomoStrip",
"True",strips,
false);
922 std::vector<unsigned int> PlanesIDs;
930 for(
int i = 0; i <
m_planes.size(); i++)
936 m_planes[0]->GeneralGetFieldDefinitions(fielddef, 1, HomoBasis,
937 HomoLen, strips, StripsIDs, PlanesIDs);
950 int ncoeffs_per_plane =
m_planes[0]->GetNcoeffs();
956 for(i = 0; i <
m_planes[0]->GetExpSize(); ++i)
962 for(i = 0; i < fielddef->m_elementIDs.size(); ++i)
965 int datalen = (*m_exp)[eid]->GetNcoeffs();
967 for(n = 0; n <
m_planes.size(); ++n)
969 fielddata.insert(fielddata.end(),&coeffs[
m_coeff_offset[eid]+n*ncoeffs_per_plane],&coeffs[
m_coeff_offset[eid]+n*ncoeffs_per_plane]+datalen);
982 std::vector<NekDouble> &fielddata,
989 int datalen = fielddata.size()/fielddef->m_fields.size();
990 std::vector<unsigned int> fieldDefHomoZids;
994 for(i = 0; i < fielddef->m_fields.size(); ++i)
996 if(fielddef->m_fields[i] == field)
1003 if(i == fielddef->m_fields.size())
1005 cout <<
"Field "<< field<<
"not found in data file. " << endl;
1010 int modes_offset = 0;
1011 int planes_offset = 0;
1018 for (i = 0; i <
m_planes.size(); ++i)
1023 for (i = 0; i <
m_planes[0]->GetExpSize(); ++i)
1029 if(fielddef->m_numHomogeneousDir)
1031 nzmodes = fielddef->m_homogeneousZIDs.size();
1032 fieldDefHomoZids = fielddef->m_homogeneousZIDs;
1037 fieldDefHomoZids.push_back(0);
1041 int ncoeffs_per_plane =
m_planes[0]->GetNcoeffs();
1043 for(i = 0; i < fielddef->m_elementIDs.size(); ++i)
1045 if(fielddef->m_uniOrder ==
true)
1051 fielddef->m_numModes,
1060 offset += datalen*nzmodes;
1061 modes_offset += (*m_exp)[0]->GetNumBases() +
1062 fielddef->m_numHomogeneousDir;
1066 int eid = it->second;
1067 bool sameBasis =
true;
1068 for (
int j = 0; j < fielddef->m_basis.size()-1; ++j)
1070 if (fielddef->m_basis[j] != (*
m_exp)[eid]->GetBasisType(j))
1077 for(n = 0; n < nzmodes; ++n, offset += datalen)
1088 planes_offset = it->second;
1095 (*m_exp)[eid]->ExtractDataToCoeffs(&fielddata[offset], fielddef->m_numModes,modes_offset,&coeffs[
m_coeff_offset[eid] + planes_offset*ncoeffs_per_plane], fielddef->m_basis);
1098 modes_offset += (*m_exp)[0]->GetNumBases() + fielddef->m_numHomogeneousDir;
1108 int fromNcoeffs_per_plane = fromExpList->GetPlane(0)->GetNcoeffs();
1109 int toNcoeffs_per_plane =
m_planes[0]->GetNcoeffs();
1112 for(i = 0; i <
m_planes.size(); ++i)
1114 m_planes[i]->ExtractCoeffsToCoeffs(fromExpList->GetPlane(i),fromcoeffs_tmp = fromCoeffs + fromNcoeffs_per_plane*i, tocoeffs_tmp = toCoeffs + toNcoeffs_per_plane*i);
1124 m_planes[0]->WriteVtkPieceData(outfile, expansion, var);
1129 int nq = (*m_exp)[expansion]->GetTotPoints();
1130 int npoints_per_plane =
m_planes[0]->GetTotPoints();
1133 int outputExtraPlane = 0;
1139 outputExtraPlane = 1;
1150 int fromRank = (rank+1) % size;
1151 int toRank = (rank == 0) ? size-1 : rank-1;
1157 fromRank, extraPlane);
1162 outfile <<
" <DataArray type=\"Float64\" Name=\""
1163 << var <<
"\">" << endl;
1165 for (
int n = 0; n <
m_planes.size(); ++n)
1168 for(i = 0; i < nq; ++i)
1173 if (outputExtraPlane)
1175 for(i = 0; i < nq; ++i)
1178 0 : extraPlane[i]) <<
" ";
1182 outfile <<
" </DataArray>" << endl;
1190 cnt1 =
m_planes[0]->Get1DScaledTotPoints(scale);
1192 ASSERTL1(
m_planes.size()*cnt1 <= outarray.size(),
"size of outarray does not match internal estimage");
1195 for(
int i = 0; i <
m_planes.size(); i++)
1198 m_planes[i]->PhysInterp1DScaled(scale,inarray+i*cnt,
1199 tmparray = outarray+i*cnt1);
1208 cnt =
m_planes[0]->Get1DScaledTotPoints(scale);
1209 cnt1 =
m_planes[0]->GetTotPoints();
1211 ASSERTL1(
m_planes.size()*cnt <= inarray.size(),
"size of outarray does not match internal estimage");
1214 for(
int i = 0; i <
m_planes.size(); i++)
1216 m_planes[i]->PhysGalerkinProjection1DScaled(scale,inarray+i*cnt,
1217 tmparray = outarray+i*cnt1);
1226 int nT_pts = inarray.size();
1227 int nP_pts = nT_pts/
m_planes.size();
1235 for(
int i = 0; i <
m_planes.size(); i++)
1237 m_planes[i]->PhysDeriv(inarray + i*nP_pts ,tmp2 = out_d0 + i*nP_pts , tmp3 = out_d1 + i*nP_pts );
1247 temparray = inarray;
1274 for(
int i = 0; i <
m_planes.size(); i++)
1278 Vmath::Smul(nP_pts,beta,tmp1 = temparray + i*nP_pts,1,tmp2 = outarray + (i-
int(
sign))*nP_pts,1);
1295 if(!
m_session->DefinesSolverInfo(
"HomoStrip"))
1298 "Parallelisation in the homogeneous direction "
1299 "implemented just for Fourier basis");
1304 "Parallelisation in the homogeneous direction "
1305 "implemented just for Fourier basis");
1310 ASSERTL0(
false,
"Semi-phyisical time-stepping not "
1311 "implemented yet for non-Fourier "
1320 for(
int i = 0; i < nP_pts; i++)
1337 int nT_pts = inarray.size();
1338 int nP_pts = nT_pts/
m_planes.size();
1349 for(
int i=0; i <
m_planes.size(); i++)
1351 m_planes[i]->PhysDeriv(edir, inarray + i*nP_pts ,tmp2 = out_d + i*nP_pts);
1361 temparray = inarray;
1387 for(
int i = 0; i <
m_planes.size(); i++)
1391 Vmath::Smul(nP_pts,beta,tmp1 = temparray + i*nP_pts,1,tmp2 = outarray + (i-
int(
sign))*nP_pts,1);
1407 if(!
m_session->DefinesSolverInfo(
"HomoStrip"))
1410 "Parallelisation in the homogeneous direction "
1411 "implemented just for Fourier basis");
1416 "Parallelisation in the homogeneous direction "
1417 "implemented just for Fourier basis");
1422 ASSERTL0(
false,
"Semi-phyisical time-stepping not implemented yet for non-Fourier basis");
1430 for(
int i = 0; i < nP_pts; i++)
1485 for (i = 0; i < (*m_exp).size(); ++i)
#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.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
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.
Abstraction of a two-dimensional multi-elemental expansion which is merely a collection of local expa...
void PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
virtual void v_FwdTrans_BndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
virtual void v_PhysInterp1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
NekDouble m_lhom
Width of homogeneous direction.
LibUtilities::TranspositionSharedPtr m_transposition
virtual void v_ExtractCoeffsToCoeffs(const std::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
virtual void v_DealiasedProd(const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray)
virtual NekDouble v_GetHomoLen(void)
void HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_DealiasedDotProd(const Array< OneD, Array< OneD, NekDouble > > &inarray1, const Array< OneD, Array< OneD, NekDouble > > &inarray2, Array< OneD, Array< OneD, NekDouble > > &outarray)
Array< OneD, ExpListSharedPtr > m_planes
std::unordered_map< int, int > m_zIdToPlane
LibUtilities::CommSharedPtr m_StripZcomm
Homo1DBlockMatrixMapShPtr m_homogeneous1DBlockMat
NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
LibUtilities::BasisSharedPtr m_homogeneousBasis
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
void HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
virtual void v_FwdTrans_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
LibUtilities::NektarFFTSharedPtr m_FFT_deal
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual Array< OneD, const unsigned int > v_GetZIDs(void)
void Homogeneous1DTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool IsForwards, bool Shuff=true, bool UnShuff=true)
Array< OneD, NekDouble > m_tmpIN
virtual void v_SetHomoLen(const NekDouble lhom)
LibUtilities::NektarFFTSharedPtr m_FFT
virtual void v_IProductWRTBase_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual ~ExpListHomogeneous1D()
Destructor.
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)
virtual std::vector< LibUtilities::FieldDefinitionsSharedPtr > v_GetFieldDefinitions(void)
DNekBlkMatSharedPtr GenHomogeneous1DBlockMatrix(Homogeneous1DMatType mattype) const
virtual void v_WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var)
virtual void v_HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
virtual void v_ExtractDataToCoeffs(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
Extract data from raw field data into expansion list.
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
DNekBlkMatSharedPtr GetHomogeneous1DBlockMatrix(Homogeneous1DMatType mattype) const
virtual void v_PhysGalerkinProjection1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_BwdTrans_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Array< OneD, NekDouble > m_tmpOUT
virtual void v_HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
bool m_useFFT
FFT variables.
virtual LibUtilities::TranspositionSharedPtr v_GetTransposition(void)
ExpListHomogeneous1D(const ExpansionType type)
Default constructor.
Base class for all multi-elemental spectral/hp expansions.
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
LibUtilities::CommSharedPtr m_comm
Communicator.
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
std::unordered_map< int, int > m_elmtToExpId
Mapping from geometry ID of element to index inside m_exp.
LibUtilities::SessionReaderSharedPtr m_session
Session.
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
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)
Class representing a segment element in reference space.
static BasisSharedPtr NullBasisSharedPtr
BasisManagerT & BasisManager(void)
NektarFFTFactory & GetNektarFFTFactory()
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< FieldDefinitions > FieldDefinitionsSharedPtr
std::shared_ptr< Transposition > TranspositionSharedPtr
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
int GetNumberOfCoefficients(ShapeType shape, std::vector< unsigned int > &modes, int offset)
@ eFourierSingleMode
Fourier ModifiedExpansion with just the first mode .
@ eFourierHalfModeIm
Fourier Modified expansions with just the imaginary part of the first mode
@ eFourierHalfModeRe
Fourier Modified expansions with just the real part of the first mode
@ eFourier
Fourier Expansion .
std::map< Homogeneous1DMatType, DNekBlkMatSharedPtr > Homo1DBlockMatrixMap
A map between homo matrix keys and their associated block matrices.
static const NekDouble kNekZeroTol
The above copyright notice and this permission notice shall be included.
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
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 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 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 Vcopy(int n, const T *x, const int incx, T *y, const int incy)