46 namespace MultiRegions
49 ExpListHomogeneous1D::ExpListHomogeneous1D():
62 const bool dealiasing):
67 m_dealiasing(dealiasing)
70 "Homogeneous Basis is a null basis");
75 m_comm->GetColumnComm()->GetColumnComm() :
80 "HomModesZ should be a multiple of npz.");
90 "HomModesZ/npz should be an even integer.");
119 ASSERTL0(
false,
"Dealiasing available just in combination "
131 m_transposition(In.m_transposition),
132 m_StripZcomm(In.m_StripZcomm),
133 m_useFFT(In.m_useFFT),
136 m_tmpOUT(In.m_tmpOUT),
137 m_homogeneousBasis(In.m_homogeneousBasis),
139 m_homogeneous1DBlockMat(In.m_homogeneous1DBlockMat),
140 m_dealiasing(In.m_dealiasing),
141 m_padsize(In.m_padsize)
147 const std::vector<unsigned int> &eIDs):
149 m_transposition(In.m_transposition),
150 m_useFFT(In.m_useFFT),
153 m_tmpOUT(In.m_tmpOUT),
154 m_homogeneousBasis(In.m_homogeneousBasis),
157 m_dealiasing(In.m_dealiasing),
158 m_padsize(In.m_padsize)
202 int num_dofs = inarray1.num_elements();
220 int num_points_per_plane = num_dofs/
m_planes.num_elements();
222 if(!
m_session->DefinesSolverInfo(
"HomoStrip"))
224 num_proc =
m_comm->GetColumnComm()->GetSize();
230 int num_dfts_per_proc = num_points_per_plane / num_proc
231 + (num_points_per_plane % num_proc > 0);
249 for(
int i = 0 ; i < num_dfts_per_proc ; i++)
253 Vmath::Vcopy(N, &(ShufV1[i*N]), 1, &(ShufV1_PAD_coef[0]), 1);
254 Vmath::Vcopy(N, &(ShufV2[i*N]), 1, &(ShufV2_PAD_coef[0]), 1);
257 m_FFT_deal->FFTBwdTrans(ShufV1_PAD_coef, ShufV1_PAD_phys);
258 m_FFT_deal->FFTBwdTrans(ShufV2_PAD_coef, ShufV2_PAD_phys);
264 ShufV1V2_PAD_phys, 1);
268 m_FFT_deal->FFTFwdTrans(ShufV1V2_PAD_phys, ShufV1V2_PAD_coef);
273 &(ShufV1V2[i*N]), 1);
305 int ndim = inarray1.num_elements();
306 ASSERTL1( inarray2.num_elements() % ndim == 0,
307 "Wrong dimensions for DealiasedDotProd.");
308 int nvec = inarray2.num_elements() / ndim;
310 int num_dofs = inarray1[0].num_elements();
313 int num_points_per_plane = num_dofs/
m_planes.num_elements();
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.num_elements(); ++n)
471 m_planes[n]->FwdTrans(inarray+cnt, tmparray = outarray + cnt1,
490 int cnt = 0, cnt1 = 0;
494 for(
int n = 0; n <
m_planes.num_elements(); ++n)
496 m_planes[n]->FwdTrans_IterPerExp(inarray+cnt, tmparray = outarray + cnt1);
513 int cnt = 0, cnt1 = 0;
517 for(
int n = 0; n <
m_planes.num_elements(); ++n)
519 m_planes[n]->FwdTrans_BndConstrained(inarray+cnt, tmparray = outarray + cnt1);
534 int cnt = 0, cnt1 = 0;
537 for(
int n = 0; n <
m_planes.num_elements(); ++n)
539 m_planes[n]->BwdTrans(inarray+cnt, tmparray = outarray + cnt1,
542 cnt1 +=
m_planes[n]->GetTotPoints();
555 int cnt = 0, cnt1 = 0;
558 for(
int n = 0; n <
m_planes.num_elements(); ++n)
560 m_planes[n]->BwdTrans_IterPerExp(inarray+cnt, tmparray = outarray + cnt1);
563 cnt1 +=
m_planes[n]->GetTotPoints();
576 int cnt = 0, cnt1 = 0;
589 for(
int n = 0; n <
m_planes.num_elements(); ++n)
591 m_planes[n]->IProductWRTBase(tmpIn+cnt, tmparray = outarray + cnt1,coeffstate);
603 int cnt = 0, cnt1 = 0;
616 for(
int n = 0; n <
m_planes.num_elements(); ++n)
618 m_planes[n]->IProductWRTBase_IterPerExp(tmpIn+cnt, tmparray = outarray + cnt1);
638 num_dofs = inarray.num_elements();
642 num_dofs = outarray.num_elements();
647 int num_points_per_plane = num_dofs/
m_planes.num_elements();
648 int num_dfts_per_proc;
649 if(!
m_session->DefinesSolverInfo(
"HomoStrip"))
651 int nP =
m_comm->GetColumnComm()->GetSize();
652 num_dfts_per_proc = num_points_per_plane / nP
653 + (num_points_per_plane % nP > 0);
658 num_dfts_per_proc = num_points_per_plane / nP
659 + (num_points_per_plane % nP > 0);
677 for(
int i = 0 ; i < num_dfts_per_proc ; i++)
684 for(
int i = 0 ; i < num_dfts_per_proc ; i++)
697 fft_out,1,outarray,1);
727 int nrows = blkmat->GetRows();
728 int ncols = blkmat->GetColumns();
772 return matrixIter->second;
779 boost::ignore_unused(coeffstate);
784 int num_trans_per_proc = 0;
793 n_exp =
m_planes[0]->GetTotPoints();
796 num_trans_per_proc = n_exp/
m_comm->GetColumnComm()->GetSize() + (n_exp%
m_comm->GetColumnComm()->GetSize() > 0);
863 for(
int i = 0; i < num_trans_per_proc; ++i)
865 BlkMatrix->SetBlock(i,i,loc_mat);
873 std::vector<LibUtilities::FieldDefinitionsSharedPtr> returnval;
878 std::vector<NekDouble> HomoLen;
879 HomoLen.push_back(
m_lhom);
881 std::vector<unsigned int> StripsIDs;
884 m_session->MatchSolverInfo(
"HomoStrip",
"True",strips,
false);
890 std::vector<unsigned int> PlanesIDs;
900 for(
int i = 0; i <
m_planes.num_elements(); i++)
905 m_planes[0]->GeneralGetFieldDefinitions(returnval, 1, HomoBasis,
906 HomoLen, strips, StripsIDs, PlanesIDs);
915 std::vector<NekDouble> HomoLen;
916 HomoLen.push_back(
m_lhom);
918 std::vector<unsigned int> StripsIDs;
921 m_session->MatchSolverInfo(
"HomoStrip",
"True",strips,
false);
927 std::vector<unsigned int> PlanesIDs;
935 for(
int i = 0; i <
m_planes.num_elements(); i++)
941 m_planes[0]->GeneralGetFieldDefinitions(fielddef, 1, HomoBasis,
942 HomoLen, strips, StripsIDs, PlanesIDs);
955 int ncoeffs_per_plane =
m_planes[0]->GetNcoeffs();
961 for(i = 0; i <
m_planes[0]->GetExpSize(); ++i)
967 for(i = 0; i < fielddef->m_elementIDs.size(); ++i)
970 int datalen = (*m_exp)[eid]->GetNcoeffs();
972 for(n = 0; n <
m_planes.num_elements(); ++n)
974 fielddata.insert(fielddata.end(),&coeffs[
m_coeff_offset[eid]+n*ncoeffs_per_plane],&coeffs[
m_coeff_offset[eid]+n*ncoeffs_per_plane]+datalen);
987 std::vector<NekDouble> &fielddata,
994 int datalen = fielddata.size()/fielddef->m_fields.size();
995 std::vector<unsigned int> fieldDefHomoZids;
999 for(i = 0; i < fielddef->m_fields.size(); ++i)
1001 if(fielddef->m_fields[i] == field)
1008 if(i == fielddef->m_fields.size())
1010 cout <<
"Field "<< field<<
"not found in data file. " << endl;
1015 int modes_offset = 0;
1016 int planes_offset = 0;
1023 for (i = 0; i <
m_planes.num_elements(); ++i)
1028 for (i = 0; i <
m_planes[0]->GetExpSize(); ++i)
1034 if(fielddef->m_numHomogeneousDir)
1036 nzmodes = fielddef->m_homogeneousZIDs.size();
1037 fieldDefHomoZids = fielddef->m_homogeneousZIDs;
1042 fieldDefHomoZids.push_back(0);
1046 int ncoeffs_per_plane =
m_planes[0]->GetNcoeffs();
1048 for(i = 0; i < fielddef->m_elementIDs.size(); ++i)
1050 if(fielddef->m_uniOrder ==
true)
1056 fielddef->m_numModes,
1065 offset += datalen*nzmodes;
1066 modes_offset += (*m_exp)[0]->GetNumBases() +
1067 fielddef->m_numHomogeneousDir;
1071 int eid = it->second;
1072 bool sameBasis =
true;
1073 for (
int j = 0; j < fielddef->m_basis.size()-1; ++j)
1075 if (fielddef->m_basis[j] != (*
m_exp)[eid]->GetBasisType(j))
1082 for(n = 0; n < nzmodes; ++n, offset += datalen)
1093 planes_offset = it->second;
1100 (*m_exp)[eid]->ExtractDataToCoeffs(&fielddata[offset], fielddef->m_numModes,modes_offset,&coeffs[
m_coeff_offset[eid] + planes_offset*ncoeffs_per_plane], fielddef->m_basis);
1103 modes_offset += (*m_exp)[0]->GetNumBases() + fielddef->m_numHomogeneousDir;
1113 int fromNcoeffs_per_plane = fromExpList->GetPlane(0)->GetNcoeffs();
1114 int toNcoeffs_per_plane =
m_planes[0]->GetNcoeffs();
1117 for(i = 0; i <
m_planes.num_elements(); ++i)
1119 m_planes[i]->ExtractCoeffsToCoeffs(fromExpList->GetPlane(i),fromcoeffs_tmp = fromCoeffs + fromNcoeffs_per_plane*i, tocoeffs_tmp = toCoeffs + toNcoeffs_per_plane*i);
1129 m_planes[0]->WriteVtkPieceData(outfile, expansion, var);
1134 int nq = (*m_exp)[expansion]->GetTotPoints();
1135 int npoints_per_plane =
m_planes[0]->GetTotPoints();
1138 int outputExtraPlane = 0;
1144 outputExtraPlane = 1;
1155 int fromRank = (rank+1) % size;
1156 int toRank = (rank == 0) ? size-1 : rank-1;
1162 fromRank, extraPlane);
1167 outfile <<
" <DataArray type=\"Float64\" Name=\""
1168 << var <<
"\">" << endl;
1170 for (
int n = 0; n <
m_planes.num_elements(); ++n)
1173 for(i = 0; i < nq; ++i)
1178 if (outputExtraPlane)
1180 for(i = 0; i < nq; ++i)
1183 0 : extraPlane[i]) <<
" ";
1187 outfile <<
" </DataArray>" << endl;
1195 cnt1 =
m_planes[0]->Get1DScaledTotPoints(scale);
1197 ASSERTL1(
m_planes.num_elements()*cnt1 <= outarray.num_elements(),
"size of outarray does not match internal estimage");
1200 for(
int i = 0; i <
m_planes.num_elements(); i++)
1203 m_planes[i]->PhysInterp1DScaled(scale,inarray+i*cnt,
1204 tmparray = outarray+i*cnt1);
1213 cnt =
m_planes[0]->Get1DScaledTotPoints(scale);
1214 cnt1 =
m_planes[0]->GetTotPoints();
1216 ASSERTL1(
m_planes.num_elements()*cnt <= inarray.num_elements(),
"size of outarray does not match internal estimage");
1219 for(
int i = 0; i <
m_planes.num_elements(); i++)
1221 m_planes[i]->PhysGalerkinProjection1DScaled(scale,inarray+i*cnt,
1222 tmparray = outarray+i*cnt1);
1231 int nT_pts = inarray.num_elements();
1232 int nP_pts = nT_pts/
m_planes.num_elements();
1240 for(
int i = 0; i <
m_planes.num_elements(); i++)
1242 m_planes[i]->PhysDeriv(inarray + i*nP_pts ,tmp2 = out_d0 + i*nP_pts , tmp3 = out_d1 + i*nP_pts );
1252 temparray = inarray;
1279 for(
int i = 0; i <
m_planes.num_elements(); i++)
1283 Vmath::Smul(nP_pts,beta,tmp1 = temparray + i*nP_pts,1,tmp2 = outarray + (i-
int(
sign))*nP_pts,1);
1300 if(!
m_session->DefinesSolverInfo(
"HomoStrip"))
1303 "Parallelisation in the homogeneous direction "
1304 "implemented just for Fourier basis");
1309 "Parallelisation in the homogeneous direction "
1310 "implemented just for Fourier basis");
1315 ASSERTL0(
false,
"Semi-phyisical time-stepping not "
1316 "implemented yet for non-Fourier "
1325 for(
int i = 0; i < nP_pts; i++)
1342 int nT_pts = inarray.num_elements();
1343 int nP_pts = nT_pts/
m_planes.num_elements();
1354 for(
int i=0; i <
m_planes.num_elements(); i++)
1356 m_planes[i]->PhysDeriv(edir, inarray + i*nP_pts ,tmp2 = out_d + i*nP_pts);
1366 temparray = inarray;
1392 for(
int i = 0; i <
m_planes.num_elements(); i++)
1396 Vmath::Smul(nP_pts,beta,tmp1 = temparray + i*nP_pts,1,tmp2 = outarray + (i-
int(
sign))*nP_pts,1);
1412 if(!
m_session->DefinesSolverInfo(
"HomoStrip"))
1415 "Parallelisation in the homogeneous direction "
1416 "implemented just for Fourier basis");
1421 "Parallelisation in the homogeneous direction "
1422 "implemented just for Fourier basis");
1427 ASSERTL0(
false,
"Semi-phyisical time-stepping not implemented yet for non-Fourier basis");
1435 for(
int i = 0; i < nP_pts; 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...
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
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.
void Homogeneous1DTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool IsForwards, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
LibUtilities::TranspositionSharedPtr m_transposition
void HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
virtual void v_ExtractCoeffsToCoeffs(const std::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
virtual NekDouble v_GetHomoLen(void)
Array< OneD, ExpListSharedPtr > m_planes
std::unordered_map< int, int > m_zIdToPlane
LibUtilities::CommSharedPtr m_StripZcomm
void HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
Homo1DBlockMatrixMapShPtr m_homogeneous1DBlockMat
ExpListHomogeneous1D()
Default constructor.
LibUtilities::BasisSharedPtr m_homogeneousBasis
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
DNekBlkMatSharedPtr GenHomogeneous1DBlockMatrix(Homogeneous1DMatType mattype, CoeffState coeffstate=eLocal) const
virtual void v_DealiasedDotProd(const Array< OneD, Array< OneD, NekDouble > > &inarray1, const Array< OneD, Array< OneD, NekDouble > > &inarray2, Array< OneD, Array< OneD, NekDouble > > &outarray, CoeffState coeffstate=eLocal)
virtual void v_FwdTrans_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
DNekBlkMatSharedPtr GetHomogeneous1DBlockMatrix(Homogeneous1DMatType mattype, CoeffState coeffstate=eLocal) const
LibUtilities::NektarFFTSharedPtr m_FFT_deal
virtual Array< OneD, const unsigned int > v_GetZIDs(void)
virtual void v_HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, 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 void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
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)
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, CoeffState coeffstate=eLocal, 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_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_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
bool m_useFFT
FFT variables.
virtual LibUtilities::TranspositionSharedPtr v_GetTransposition(void)
virtual void v_DealiasedProd(const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
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
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*y.
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)