precomp_mat.txx 9.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313
  1. /**
  2. * \file precomp_mat.txx
  3. * \author Dhairya Malhotra, dhairya.malhotra@gmail.com
  4. * \date 3-07-2011
  5. * \brief This file contains the implementation of the PrecompMat class.
  6. * Handles storage of precomputed translation matrices.
  7. */
  8. #include <sys/stat.h>
  9. #include <stdint.h>
  10. namespace pvfmm{
  11. #define PRECOMP_MIN_DEPTH 40
  12. template <class T>
  13. PrecompMat<T>::PrecompMat(bool homogen, int max_d): homogeneous(homogen), max_depth(max_d){
  14. if(!homogen) mat.resize((max_d+PRECOMP_MIN_DEPTH)*Type_Count); //For each U,V,W,X
  15. else mat.resize(Type_Count);
  16. for(size_t i=0;i<mat.size();i++)
  17. mat[i].resize(500);
  18. perm.resize(Type_Count);
  19. perm_r.resize(Type_Count);
  20. perm_c.resize(Type_Count);
  21. for(size_t i=0;i<Type_Count;i++){
  22. perm[i].resize(Perm_Count);
  23. perm_r[i].resize(500);
  24. perm_c[i].resize(500);
  25. }
  26. }
  27. template <class T>
  28. Matrix<T>& PrecompMat<T>::Mat(int l, Mat_Type type, size_t indx){
  29. int level=(homogeneous?0:l+PRECOMP_MIN_DEPTH);
  30. #pragma omp critical (PrecompMAT)
  31. if(indx>=mat[level*Type_Count+type].size()){
  32. mat[level*Type_Count+type].resize(indx+1);
  33. assert(false); //TODO: this is not thread safe.
  34. }
  35. return mat[level*Type_Count+type][indx];
  36. }
  37. template <class T>
  38. Permutation<T>& PrecompMat<T>::Perm_R(Mat_Type type, size_t indx){
  39. #pragma omp critical (PrecompMAT)
  40. if(indx>=perm_r[type].size()){
  41. perm_r[type].resize(indx+1);
  42. assert(false); //TODO: this is not thread safe.
  43. }
  44. return perm_r[type][indx];
  45. }
  46. template <class T>
  47. Permutation<T>& PrecompMat<T>::Perm_C(Mat_Type type, size_t indx){
  48. #pragma omp critical (PrecompMAT)
  49. if(indx>=perm_c[type].size()){
  50. perm_c[type].resize(indx+1);
  51. assert(false); //TODO: this is not thread safe.
  52. }
  53. return perm_c[type][indx];
  54. }
  55. template <class T>
  56. Permutation<T>& PrecompMat<T>::Perm(Mat_Type type, size_t indx){
  57. assert(indx<Perm_Count);
  58. return perm[type][indx];
  59. }
  60. template <class T>
  61. size_t PrecompMat<T>::CompactData(int l, Mat_Type type, Matrix<char>& comp_data, size_t offset){
  62. std::vector<Matrix<T> >& mat_=mat[(homogeneous?0:l+PRECOMP_MIN_DEPTH)*Type_Count+type];
  63. size_t mat_cnt=mat_.size();
  64. size_t indx_size=0;
  65. size_t mem_size=0;
  66. int omp_p=omp_get_max_threads();
  67. { // Determine memory size.
  68. indx_size+=2*sizeof(size_t); //total_size, mat_cnt
  69. indx_size+=mat_cnt*(1+2+2)*sizeof(size_t); //Mat, Perm_R, Perm_C.
  70. indx_size=((uintptr_t)indx_size+(uintptr_t)(MEM_ALIGN-1)) & ~ (uintptr_t)(MEM_ALIGN-1);
  71. for(size_t j=0;j<mat_cnt;j++){
  72. Matrix <T>& M =Mat( l,type,j);
  73. Permutation<T>& Pr=Perm_R(type,j);
  74. Permutation<T>& Pc=Perm_C(type,j);
  75. if(M.Dim(0)>0 && M.Dim(1)>0){
  76. mem_size+=M.Dim(0)*M.Dim(1)*sizeof(T); mem_size=((uintptr_t)mem_size+(uintptr_t)(MEM_ALIGN-1)) & ~ (uintptr_t)(MEM_ALIGN-1);
  77. }
  78. if(Pr.Dim()>0){
  79. mem_size+=Pr.Dim()*sizeof(PERM_INT_T); mem_size=((uintptr_t)mem_size+(uintptr_t)(MEM_ALIGN-1)) & ~ (uintptr_t)(MEM_ALIGN-1);
  80. mem_size+=Pr.Dim()*sizeof(T); mem_size=((uintptr_t)mem_size+(uintptr_t)(MEM_ALIGN-1)) & ~ (uintptr_t)(MEM_ALIGN-1);
  81. }
  82. if(Pc.Dim()>0){
  83. mem_size+=Pc.Dim()*sizeof(PERM_INT_T); mem_size=((uintptr_t)mem_size+(uintptr_t)(MEM_ALIGN-1)) & ~ (uintptr_t)(MEM_ALIGN-1);
  84. mem_size+=Pc.Dim()*sizeof(T); mem_size=((uintptr_t)mem_size+(uintptr_t)(MEM_ALIGN-1)) & ~ (uintptr_t)(MEM_ALIGN-1);
  85. }
  86. }
  87. }
  88. if(comp_data.Dim(0)*comp_data.Dim(1)<offset+indx_size+mem_size){ // Resize if needed.
  89. Matrix<char> old_data;
  90. if(offset>0) old_data=comp_data;
  91. comp_data.Resize(1,offset+indx_size+mem_size);
  92. if(offset>0) mem::memcopy(comp_data[0], old_data[0], offset); //TODO: This will affect NUMA.
  93. }
  94. { // Create indx.
  95. char* indx_ptr=comp_data[0]+offset;
  96. size_t data_offset=offset+indx_size;
  97. ((size_t*)indx_ptr)[0]=indx_size+mem_size; indx_ptr+=sizeof(size_t);
  98. ((size_t*)indx_ptr)[0]= mat_cnt ; indx_ptr+=sizeof(size_t);
  99. for(size_t j=0;j<mat_cnt;j++){
  100. Matrix <T>& M =Mat( l,type,j);
  101. Permutation<T>& Pr=Perm_R(type,j);
  102. Permutation<T>& Pc=Perm_C(type,j);
  103. ((size_t*)indx_ptr)[0]=data_offset; indx_ptr+=sizeof(size_t);
  104. data_offset+=M.Dim(0)*M.Dim(1)*sizeof(T); data_offset=((uintptr_t)data_offset+(uintptr_t)(MEM_ALIGN-1)) & ~ (uintptr_t)(MEM_ALIGN-1);
  105. ((size_t*)indx_ptr)[0]=data_offset; indx_ptr+=sizeof(size_t);
  106. data_offset+=Pr.Dim()*sizeof(PERM_INT_T); data_offset=((uintptr_t)data_offset+(uintptr_t)(MEM_ALIGN-1)) & ~ (uintptr_t)(MEM_ALIGN-1);
  107. ((size_t*)indx_ptr)[0]=data_offset; indx_ptr+=sizeof(size_t);
  108. data_offset+=Pr.Dim()*sizeof(T); data_offset=((uintptr_t)data_offset+(uintptr_t)(MEM_ALIGN-1)) & ~ (uintptr_t)(MEM_ALIGN-1);
  109. ((size_t*)indx_ptr)[0]=data_offset; indx_ptr+=sizeof(size_t);
  110. data_offset+=Pc.Dim()*sizeof(PERM_INT_T); data_offset=((uintptr_t)data_offset+(uintptr_t)(MEM_ALIGN-1)) & ~ (uintptr_t)(MEM_ALIGN-1);
  111. ((size_t*)indx_ptr)[0]=data_offset; indx_ptr+=sizeof(size_t);
  112. data_offset+=Pc.Dim()*sizeof(T); data_offset=((uintptr_t)data_offset+(uintptr_t)(MEM_ALIGN-1)) & ~ (uintptr_t)(MEM_ALIGN-1);
  113. }
  114. }
  115. { // Copy data.
  116. char* indx_ptr=comp_data[0]+offset;
  117. size_t& total_size=((size_t*)indx_ptr)[0]; indx_ptr+=sizeof(size_t);
  118. size_t& mat_cnt =((size_t*)indx_ptr)[0]; indx_ptr+=sizeof(size_t);
  119. Vector<size_t> data_offset((1+2+2)*mat_cnt, (size_t*)indx_ptr, false);
  120. offset+=total_size;
  121. for(size_t j=0;j<mat_cnt;j++){
  122. Matrix <T>& M =Mat( l,type,j);
  123. Permutation<T>& Pr=Perm_R(type,j);
  124. Permutation<T>& Pc=Perm_C(type,j);
  125. if(M.Dim(0)>0 && M.Dim(1)>0){
  126. #pragma omp parallel for
  127. for(int tid=0;tid<omp_p;tid++){
  128. size_t a=(M.Dim(0)*M.Dim(1)* tid )/omp_p;
  129. size_t b=(M.Dim(0)*M.Dim(1)*(tid+1))/omp_p;
  130. mem::memcopy(comp_data[0]+data_offset[5*j+0]+a*sizeof(T), &M[0][a], (b-a)*sizeof(T));
  131. }
  132. }
  133. if(Pr.Dim()>0){
  134. mem::memcopy(comp_data[0]+data_offset[5*j+1], &Pr.perm[0], Pr.Dim()*sizeof(PERM_INT_T));
  135. mem::memcopy(comp_data[0]+data_offset[5*j+2], &Pr.scal[0], Pr.Dim()*sizeof( T));
  136. }
  137. if(Pc.Dim()>0){
  138. mem::memcopy(comp_data[0]+data_offset[5*j+3], &Pc.perm[0], Pc.Dim()*sizeof(PERM_INT_T));
  139. mem::memcopy(comp_data[0]+data_offset[5*j+4], &Pc.scal[0], Pc.Dim()*sizeof( T));
  140. }
  141. }
  142. }
  143. return offset;
  144. }
  145. template <class T>
  146. void PrecompMat<T>::Save2File(const char* fname, bool replace){
  147. FILE* f=fopen(fname,"r");
  148. if(f!=NULL) { //File already exists.
  149. fclose(f);
  150. if(!replace) return;
  151. }
  152. f=fopen(fname,"wb");
  153. if(f==NULL) return;
  154. int tmp;
  155. tmp=sizeof(T);
  156. fwrite(&tmp,sizeof(int),1,f);
  157. tmp=(homogeneous?1:0);
  158. fwrite(&tmp,sizeof(int),1,f);
  159. tmp=max_depth;
  160. fwrite(&tmp,sizeof(int),1,f);
  161. for(size_t i=0;i<mat.size();i++){
  162. int n=mat[i].size();
  163. fwrite(&n,sizeof(int),1,f);
  164. for(int j=0;j<n;j++){
  165. Matrix<T>& M=mat[i][j];
  166. int n1=M.Dim(0);
  167. fwrite(&n1,sizeof(int),1,f);
  168. int n2=M.Dim(1);
  169. fwrite(&n2,sizeof(int),1,f);
  170. if(n1*n2>0)
  171. fwrite(&M[0][0],sizeof(T),n1*n2,f);
  172. }
  173. }
  174. fclose(f);
  175. }
  176. #define MY_FREAD(a,b,c,d) { \
  177. size_t r_cnt=fread(a,b,c,d); \
  178. if(r_cnt!=c){ \
  179. fputs ("Reading error ",stderr); \
  180. exit (-1); \
  181. } }
  182. template <class T>
  183. void PrecompMat<T>::LoadFile(const char* fname, MPI_Comm comm){
  184. Profile::Tic("LoadMatrices",&comm,true,3);
  185. Profile::Tic("ReadFile",&comm,true,4);
  186. size_t f_size=0;
  187. char* f_data=NULL;
  188. int np, myrank;
  189. MPI_Comm_size(comm,&np);
  190. MPI_Comm_rank(comm,&myrank);
  191. if(myrank==0){
  192. FILE* f=fopen(fname,"rb");
  193. if(f==NULL){
  194. f_size=0;
  195. }else{
  196. struct stat fileStat;
  197. if(stat(fname,&fileStat) < 0) f_size=0;
  198. else f_size=fileStat.st_size;
  199. //fseek (f, 0, SEEK_END);
  200. //f_size=ftell (f);
  201. }
  202. if(f_size>0){
  203. f_data=new char[f_size];
  204. fseek (f, 0, SEEK_SET);
  205. MY_FREAD(f_data,sizeof(char),f_size,f);
  206. fclose(f);
  207. }
  208. }
  209. Profile::Toc();
  210. Profile::Tic("Broadcast",&comm,true,4);
  211. MPI_Bcast( &f_size, sizeof(size_t), MPI_BYTE, 0, comm );
  212. if(f_size==0){
  213. Profile::Toc();
  214. Profile::Toc();
  215. return;
  216. }
  217. if(f_data==NULL) f_data=new char[f_size];
  218. char* f_ptr=f_data;
  219. int max_send_size=1000000000;
  220. while(f_size>0){
  221. if(f_size>(size_t)max_send_size){
  222. MPI_Bcast( f_ptr, max_send_size, MPI_BYTE, 0, comm );
  223. f_size-=max_send_size;
  224. f_ptr+=max_send_size;
  225. }else{
  226. MPI_Bcast( f_ptr, f_size, MPI_BYTE, 0, comm );
  227. f_size=0;
  228. }
  229. }
  230. f_ptr=f_data;
  231. {
  232. int tmp;
  233. tmp=*(int*)f_ptr; f_ptr+=sizeof(int);
  234. assert(tmp==sizeof(T));
  235. tmp=*(int*)f_ptr; f_ptr+=sizeof(int);
  236. homogeneous=tmp;
  237. int max_depth_;
  238. max_depth_=*(int*)f_ptr; f_ptr+=sizeof(int);
  239. size_t mat_size=(size_t)Type_Count*(homogeneous?1:max_depth_+PRECOMP_MIN_DEPTH);
  240. if(mat.size()<mat_size){
  241. max_depth=max_depth_;
  242. mat.resize(mat_size);
  243. }
  244. for(size_t i=0;i<mat_size;i++){
  245. int n;
  246. n=*(int*)f_ptr; f_ptr+=sizeof(int);
  247. if(mat[i].size()<(size_t)n)
  248. mat[i].resize(n);
  249. for(int j=0;j<n;j++){
  250. Matrix<T>& M=mat[i][j];
  251. int n1;
  252. n1=*(int*)f_ptr; f_ptr+=sizeof(int);
  253. int n2;
  254. n2=*(int*)f_ptr; f_ptr+=sizeof(int);
  255. if(n1*n2>0){
  256. M.Resize(n1,n2);
  257. mem::memcopy(&M[0][0], f_ptr, sizeof(T)*n1*n2); f_ptr+=sizeof(T)*n1*n2;
  258. }
  259. }
  260. }
  261. }
  262. if(f_data!=NULL) delete[] f_data;
  263. Profile::Toc();
  264. Profile::Toc();
  265. }
  266. #undef MY_FREAD
  267. #undef PRECOMP_MIN_DEPTH
  268. template <class T>
  269. std::vector<T>& PrecompMat<T>::RelativeTrgCoord(){
  270. return rel_trg_coord;
  271. }
  272. template <class T>
  273. bool PrecompMat<T>::Homogen(){
  274. return homogeneous;
  275. }
  276. }//end namespace