/** * \file matrix.txx * \author Dhairya Malhotra, dhairya.malhotra@gmail.com * \date 2-11-2011 * \brief This file contains inplementation of the class Matrix. */ #include #include #include #include #include #include #include #include #include #include namespace pvfmm{ template std::ostream& operator<<(std::ostream& output, const Matrix& M){ std::ios::fmtflags f(std::cout.flags()); output< Matrix::Matrix(){ dim[0]=0; dim[1]=0; own_data=true; data_ptr=NULL; dev.dev_ptr=(uintptr_t)NULL; } template Matrix::Matrix(size_t dim1, size_t dim2, T* data_, bool own_data_){ dim[0]=dim1; dim[1]=dim2; own_data=own_data_; if(own_data){ if(dim[0]*dim[1]>0){ data_ptr=mem::aligned_new(dim[0]*dim[1]); #if !defined(__MIC__) || !defined(__INTEL_OFFLOAD) Profile::Add_MEM(dim[0]*dim[1]*sizeof(T)); #endif if(data_!=NULL) mem::memcopy(data_ptr,data_,dim[0]*dim[1]*sizeof(T)); }else data_ptr=NULL; }else data_ptr=data_; dev.dev_ptr=(uintptr_t)NULL; } template Matrix::Matrix(const Matrix& M){ dim[0]=M.dim[0]; dim[1]=M.dim[1]; own_data=true; if(dim[0]*dim[1]>0){ data_ptr=mem::aligned_new(dim[0]*dim[1]); #if !defined(__MIC__) || !defined(__INTEL_OFFLOAD) Profile::Add_MEM(dim[0]*dim[1]*sizeof(T)); #endif mem::memcopy(data_ptr,M.data_ptr,dim[0]*dim[1]*sizeof(T)); }else data_ptr=NULL; dev.dev_ptr=(uintptr_t)NULL; } template Matrix::~Matrix(){ FreeDevice(false); if(own_data){ if(data_ptr!=NULL){ mem::aligned_delete(data_ptr); #if !defined(__MIC__) || !defined(__INTEL_OFFLOAD) Profile::Add_MEM(-dim[0]*dim[1]*sizeof(T)); #endif } } data_ptr=NULL; dim[0]=0; dim[1]=0; } template void Matrix::Swap(Matrix& M){ size_t dim_[2]={dim[0],dim[1]}; T* data_ptr_=data_ptr; bool own_data_=own_data; Device dev_=dev; Vector dev_sig_=dev_sig; dim[0]=M.dim[0]; dim[1]=M.dim[1]; data_ptr=M.data_ptr; own_data=M.own_data; dev=M.dev; dev_sig=M.dev_sig; M.dim[0]=dim_[0]; M.dim[1]=dim_[1]; M.data_ptr=data_ptr_; M.own_data=own_data_; M.dev=dev_; M.dev_sig=dev_sig_; } template void Matrix::ReInit(size_t dim1, size_t dim2, T* data_, bool own_data_){ Matrix tmp(dim1,dim2,data_,own_data_); this->Swap(tmp); } template typename Matrix::Device& Matrix::AllocDevice(bool copy){ size_t len=dim[0]*dim[1]; if(dev.dev_ptr==(uintptr_t)NULL && len>0) // Allocate data on device. dev.dev_ptr=DeviceWrapper::alloc_device((char*)data_ptr, len*sizeof(T)); if(dev.dev_ptr!=(uintptr_t)NULL && copy) // Copy data to device dev.lock_idx=DeviceWrapper::host2device((char*)data_ptr,(char*)data_ptr,dev.dev_ptr,len*sizeof(T)); dev.dim[0]=dim[0]; dev.dim[1]=dim[1]; return dev; } template void Matrix::Device2Host(T* host_ptr){ dev.lock_idx=DeviceWrapper::device2host((char*)data_ptr,dev.dev_ptr,(char*)(host_ptr==NULL?data_ptr:host_ptr),dim[0]*dim[1]*sizeof(T)); //#if defined(PVFMM_HAVE_CUDA) // cudaEventCreate(&lock); // cudaEventRecord(lock, 0); //#endif } template void Matrix::Device2HostWait(){ //#if defined(PVFMM_HAVE_CUDA) // cudaEventSynchronize(lock); // cudaEventDestroy(lock); //#endif DeviceWrapper::wait(dev.lock_idx); dev.lock_idx=-1; } template void Matrix::FreeDevice(bool copy){ if(dev.dev_ptr==(uintptr_t)NULL) return; if(copy) DeviceWrapper::device2host((char*)data_ptr,dev.dev_ptr,(char*)data_ptr,dim[0]*dim[1]*sizeof(T)); DeviceWrapper::free_device((char*)data_ptr, dev.dev_ptr); dev.dev_ptr=(uintptr_t)NULL; dev.dim[0]=0; dev.dim[1]=0; } template void Matrix::Write(const char* fname){ FILE* f1=fopen(fname,"wb+"); if(f1==NULL){ std::cout<<"Unable to open file for writing:"< size_t Matrix::Dim(size_t i) const{ return dim[i]; } template void Matrix::Resize(size_t i, size_t j){ if(dim[0]==i && dim[1]==j) return; FreeDevice(false); if(own_data){ if(data_ptr!=NULL){ mem::aligned_delete(data_ptr); #if !defined(__MIC__) || !defined(__INTEL_OFFLOAD) Profile::Add_MEM(-dim[0]*dim[1]*sizeof(T)); #endif } } dim[0]=i; dim[1]=j; if(own_data){ if(dim[0]*dim[1]>0){ data_ptr=mem::aligned_new(dim[0]*dim[1]); #if !defined(__MIC__) || !defined(__INTEL_OFFLOAD) Profile::Add_MEM(dim[0]*dim[1]*sizeof(T)); #endif }else data_ptr=NULL; } } template void Matrix::SetZero(){ if(dim[0]*dim[1]>0) memset(data_ptr,0,dim[0]*dim[1]*sizeof(T)); } template Matrix& Matrix::operator=(const Matrix& M){ if(this!=&M){ FreeDevice(false); if(own_data && dim[0]*dim[1]!=M.dim[0]*M.dim[1]){ if(data_ptr!=NULL){ mem::aligned_delete(data_ptr); data_ptr=NULL; #if !defined(__MIC__) || !defined(__INTEL_OFFLOAD) Profile::Add_MEM(-dim[0]*dim[1]*sizeof(T)); #endif } if(M.dim[0]*M.dim[1]>0){ data_ptr=mem::aligned_new(M.dim[0]*M.dim[1]); #if !defined(__MIC__) || !defined(__INTEL_OFFLOAD) Profile::Add_MEM(M.dim[0]*M.dim[1]*sizeof(T)); #endif } } dim[0]=M.dim[0]; dim[1]=M.dim[1]; mem::memcopy(data_ptr,M.data_ptr,dim[0]*dim[1]*sizeof(T)); } return *this; } template Matrix& Matrix::operator+=(const Matrix& M){ assert(M.Dim(0)==Dim(0) && M.Dim(1)==Dim(1)); Profile::Add_FLOP(dim[0]*dim[1]); for(size_t i=0;i Matrix& Matrix::operator-=(const Matrix& M){ assert(M.Dim(0)==Dim(0) && M.Dim(1)==Dim(1)); Profile::Add_FLOP(dim[0]*dim[1]); for(size_t i=0;i Matrix Matrix::operator+(const Matrix& M2){ Matrix& M1=*this; assert(M2.Dim(0)==M1.Dim(0) && M2.Dim(1)==M1.Dim(1)); Profile::Add_FLOP(dim[0]*dim[1]); Matrix M_r(M1.Dim(0),M1.Dim(1),NULL); for(size_t i=0;i Matrix Matrix::operator-(const Matrix& M2){ Matrix& M1=*this; assert(M2.Dim(0)==M1.Dim(0) && M2.Dim(1)==M1.Dim(1)); Profile::Add_FLOP(dim[0]*dim[1]); Matrix M_r(M1.Dim(0),M1.Dim(1),NULL); for(size_t i=0;i inline T& Matrix::operator()(size_t i,size_t j) const{ assert(i inline T* Matrix::operator[](size_t i) const{ assert(i Matrix Matrix::operator*(const Matrix& M){ assert(dim[1]==M.dim[0]); Profile::Add_FLOP(2*(((long long)dim[0])*dim[1])*M.dim[1]); Matrix M_r(dim[0],M.dim[1],NULL); if(M.Dim(0)*M.Dim(1)==0 || this->Dim(0)*this->Dim(1)==0) return M_r; mat::gemm('N','N',M.dim[1],dim[0],dim[1], 1.0,M.data_ptr,M.dim[1],data_ptr,dim[1],0.0,M_r.data_ptr,M_r.dim[1]); return M_r; } template void Matrix::GEMM(Matrix& M_r, const Matrix& A, const Matrix& B, T beta){ if(A.Dim(0)*A.Dim(1)==0 || B.Dim(0)*B.Dim(1)==0) return; assert(A.dim[1]==B.dim[0]); assert(M_r.dim[0]==A.dim[0]); assert(M_r.dim[1]==B.dim[1]); #if !defined(__MIC__) || !defined(__INTEL_OFFLOAD) Profile::Add_FLOP(2*(((long long)A.dim[0])*A.dim[1])*B.dim[1]); #endif mat::gemm('N','N',B.dim[1],A.dim[0],A.dim[1], 1.0,B.data_ptr,B.dim[1],A.data_ptr,A.dim[1],beta,M_r.data_ptr,M_r.dim[1]); } // cublasgemm wrapper #if defined(PVFMM_HAVE_CUDA) template void Matrix::CUBLASGEMM(Matrix& M_r, const Matrix& A, const Matrix& B, T beta){ if(A.Dim(0)*A.Dim(1)==0 || B.Dim(0)*B.Dim(1)==0) return; assert(A.dim[1]==B.dim[0]); assert(M_r.dim[0]==A.dim[0]); assert(M_r.dim[1]==B.dim[1]); Profile::Add_FLOP(2*(((long long)A.dim[0])*A.dim[1])*B.dim[1]); mat::cublasgemm('N', 'N', B.dim[1], A.dim[0], A.dim[1], 1.0, B.data_ptr, B.dim[1], A.data_ptr, A.dim[1], beta, M_r.data_ptr, M_r.dim[1]); } #endif #define myswap(t,a,b) {t c=a;a=b;b=c;} template void Matrix::RowPerm(const Permutation& P){ Matrix& M=*this; if(P.Dim()==0) return; assert(M.Dim(0)==P.Dim()); size_t d0=M.Dim(0); size_t d1=M.Dim(1); #pragma omp parallel for for(size_t i=0;i P_=P; for(size_t i=0;i void Matrix::ColPerm(const Permutation& P){ Matrix& M=*this; if(P.Dim()==0) return; assert(M.Dim(1)==P.Dim()); size_t d0=M.Dim(0); size_t d1=M.Dim(1); int omp_p=omp_get_max_threads(); Matrix M_buff(omp_p,d1); const size_t* perm_=&(P.perm[0]); const T* scal_=&(P.scal[0]); #pragma omp parallel for for(size_t i=0;i Matrix Matrix::Transpose(){ Matrix& M=*this; size_t d0=M.dim[0]; size_t d1=M.dim[1]; Matrix M_r(d1,d0,NULL); const size_t blk0=((d0+B1-1)/B1); const size_t blk1=((d1+B1-1)/B1); const size_t blks=blk0*blk1; // #pragma omp parallel for for(size_t k=0;k=d0) d0_=d0; size_t d1_=j+B1; if(d1_>=d1) d1_=d1; for(size_t ii=i;ii=d0) d0__=d0; size_t d1__=jj+B2; if(d1__>=d1) d1__=d1; for(size_t iii=ii;iii void Matrix::Transpose(Matrix& M_r, const Matrix& M){ size_t d0=M.dim[0]; size_t d1=M.dim[1]; M_r.Resize(d1, d0); const size_t blk0=((d0+B1-1)/B1); const size_t blk1=((d1+B1-1)/B1); const size_t blks=blk0*blk1; #pragma omp parallel for for(size_t k=0;k=d0) d0_=d0; size_t d1_=j+B1; if(d1_>=d1) d1_=d1; for(size_t ii=i;ii=d0) d0__=d0; size_t d1__=jj+B2; if(d1__>=d1) d1__=d1; for(size_t iii=ii;iii void Matrix::SVD(Matrix& tU, Matrix& tS, Matrix& tVT){ pvfmm::Matrix& M=*this; pvfmm::Matrix M_=M; int n=M.Dim(0); int m=M.Dim(1); int k = (mn?m:n); int wssize1 = 5*(mwssize1?wssize:wssize1); T* wsbuf = mem::aligned_new(wssize); pvfmm::mat::svd(&JOBU, &JOBVT, &m, &n, &M[0][0], &m, &tS[0][0], &tVT[0][0], &m, &tU[0][0], &k, wsbuf, &wssize, &INFO); mem::aligned_delete(wsbuf); if(INFO!=0) std::cout< Matrix Matrix::pinv(T eps){ if(eps<0){ eps=1.0; while(eps+(T)1.0>1.0) eps*=0.5; eps=sqrt(eps); } Matrix M_r(dim[1],dim[0]); mat::pinv(data_ptr,dim[0],dim[1],eps,M_r.data_ptr); this->Resize(0,0); return M_r; } template std::ostream& operator<<(std::ostream& output, const Permutation& P){ output< Permutation::Permutation(size_t size){ perm.Resize(size); scal.Resize(size); for(size_t i=0;i Permutation Permutation::RandPerm(size_t size){ Permutation P(size); for(size_t i=0;i Matrix Permutation::GetMatrix() const{ size_t size=perm.Dim(); Matrix M_r(size,size,NULL); for(size_t i=0;i size_t Permutation::Dim() const{ return perm.Dim(); } template Permutation Permutation::Transpose(){ size_t size=perm.Dim(); Permutation P_r(size); Vector& perm_r=P_r.perm; Vector& scal_r=P_r.scal; for(size_t i=0;i Permutation Permutation::operator*(const Permutation& P){ size_t size=perm.Dim(); assert(P.Dim()==size); Permutation P_r(size); Vector& perm_r=P_r.perm; Vector& scal_r=P_r.scal; for(size_t i=0;i Matrix Permutation::operator*(const Matrix& M){ if(Dim()==0) return M; assert(M.Dim(0)==Dim()); size_t d0=M.Dim(0); size_t d1=M.Dim(1); Matrix M_r(d0,d1,NULL); for(size_t i=0;i Matrix operator*(const Matrix& M, const Permutation& P){ if(P.Dim()==0) return M; assert(M.Dim(1)==P.Dim()); size_t d0=M.Dim(0); size_t d1=M.Dim(1); Matrix M_r(d0,d1,NULL); for(size_t i=0;i