precomp_mat.txx 9.7 KB

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