precomp_mat.txx 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374
  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_mgr.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 scale_invar_, int max_d): scale_invar(scale_invar_), max_depth(max_d){
  21. if(!scale_invar) 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. for(size_t i=0;i<Type_Count;i++){
  27. perm[i].resize(Perm_Count);
  28. }
  29. perm_r.resize((max_d+PRECOMP_MIN_DEPTH)*Type_Count);
  30. perm_c.resize((max_d+PRECOMP_MIN_DEPTH)*Type_Count);
  31. for(size_t i=0;i<perm_r.size();i++){
  32. perm_r[i].resize(500);
  33. perm_c[i].resize(500);
  34. }
  35. }
  36. template <class T>
  37. Matrix<T>& PrecompMat<T>::Mat(int l, Mat_Type type, size_t indx){
  38. int level=(scale_invar?0:l+PRECOMP_MIN_DEPTH);
  39. assert(level*Type_Count+type<mat.size());
  40. //#pragma omp critical (PrecompMAT)
  41. if(indx>=mat[level*Type_Count+type].size()){
  42. mat[level*Type_Count+type].resize(indx+1);
  43. assert(false); //TODO: this is not thread safe.
  44. }
  45. return mat[level*Type_Count+type][indx];
  46. }
  47. template <class T>
  48. Permutation<T>& PrecompMat<T>::Perm_R(int l, Mat_Type type, size_t indx){
  49. int level=l+PRECOMP_MIN_DEPTH;
  50. assert(level*Type_Count+type<perm_r.size());
  51. //#pragma omp critical (PrecompMAT)
  52. if(indx>=perm_r[level*Type_Count+type].size()){
  53. perm_r[level*Type_Count+type].resize(indx+1);
  54. assert(false); //TODO: this is not thread safe.
  55. }
  56. return perm_r[level*Type_Count+type][indx];
  57. }
  58. template <class T>
  59. Permutation<T>& PrecompMat<T>::Perm_C(int l, Mat_Type type, size_t indx){
  60. int level=l+PRECOMP_MIN_DEPTH;
  61. assert(level*Type_Count+type<perm_c.size());
  62. //#pragma omp critical (PrecompMAT)
  63. if(indx>=perm_c[level*Type_Count+type].size()){
  64. perm_c[level*Type_Count+type].resize(indx+1);
  65. assert(false); //TODO: this is not thread safe.
  66. }
  67. return perm_c[level*Type_Count+type][indx];
  68. }
  69. template <class T>
  70. Permutation<T>& PrecompMat<T>::Perm(Mat_Type type, size_t indx){
  71. assert(indx<Perm_Count);
  72. return perm[type][indx];
  73. }
  74. template <class T>
  75. size_t PrecompMat<T>::CompactData(int level, Mat_Type type, Matrix<char>& comp_data, size_t offset){
  76. struct HeaderData{
  77. size_t total_size;
  78. size_t level;
  79. size_t mat_cnt ;
  80. size_t max_depth;
  81. };
  82. if(comp_data.Dim(0)*comp_data.Dim(1)>offset){
  83. char* indx_ptr=comp_data[0]+offset;
  84. HeaderData& header=*(HeaderData*)indx_ptr; indx_ptr+=sizeof(HeaderData);
  85. if(level==header.level){ // Data already exists.
  86. offset+=header.total_size;
  87. return offset;
  88. }
  89. }
  90. std::vector<Matrix<T> >& mat_=mat[(scale_invar?0:level+PRECOMP_MIN_DEPTH)*Type_Count+type];
  91. size_t mat_cnt=mat_.size();
  92. size_t indx_size=0;
  93. size_t mem_size=0;
  94. int omp_p=omp_get_max_threads();
  95. size_t l0=(scale_invar?0:level);
  96. size_t l1=(scale_invar?max_depth:level+1);
  97. { // Determine memory size.
  98. indx_size+=sizeof(HeaderData); // HeaderData
  99. indx_size+=mat_cnt*(1+(2+2)*(l1-l0))*sizeof(size_t); //Mat, Perm_R, Perm_C.
  100. indx_size=mem::align_ptr(indx_size);
  101. for(size_t j=0;j<mat_cnt;j++){
  102. Matrix <T>& M =Mat (level,type,j);
  103. if(M.Dim(0)>0 && M.Dim(1)>0){
  104. mem_size+=M.Dim(0)*M.Dim(1)*sizeof(T); mem_size=mem::align_ptr(mem_size);
  105. }
  106. for(size_t l=l0;l<l1;l++){
  107. Permutation<T>& Pr=Perm_R(l,type,j);
  108. Permutation<T>& Pc=Perm_C(l,type,j);
  109. if(Pr.Dim()>0){
  110. mem_size+=Pr.Dim()*sizeof(PERM_INT_T); mem_size=mem::align_ptr(mem_size);
  111. mem_size+=Pr.Dim()*sizeof(T); mem_size=mem::align_ptr(mem_size);
  112. }
  113. if(Pc.Dim()>0){
  114. mem_size+=Pc.Dim()*sizeof(PERM_INT_T); mem_size=mem::align_ptr(mem_size);
  115. mem_size+=Pc.Dim()*sizeof(T); mem_size=mem::align_ptr(mem_size);
  116. }
  117. }
  118. }
  119. }
  120. if(comp_data.Dim(0)*comp_data.Dim(1)<offset+indx_size+mem_size){ // Resize if needed.
  121. Matrix<char> old_data;
  122. if(offset>0) old_data=comp_data;
  123. comp_data.Resize(1,offset+indx_size+mem_size);
  124. if(offset>0){
  125. #pragma omp parallel for
  126. for(int tid=0;tid<omp_p;tid++){ // Copy data.
  127. size_t a=(offset*(tid+0))/omp_p;
  128. size_t b=(offset*(tid+1))/omp_p;
  129. mem::memcopy(comp_data[0]+a, old_data[0]+a, b-a);
  130. }
  131. }
  132. }
  133. { // Create indx.
  134. char* indx_ptr=comp_data[0]+offset;
  135. HeaderData& header=*(HeaderData*)indx_ptr; indx_ptr+=sizeof(HeaderData);
  136. Matrix<size_t> offset_indx(mat_cnt,1+(2+2)*(l1-l0), (size_t*)indx_ptr, false);
  137. header.total_size=indx_size+mem_size;
  138. header. level=level ;
  139. header. mat_cnt = mat_cnt ;
  140. header. max_depth=l1-l0 ;
  141. size_t data_offset=offset+indx_size;
  142. for(size_t j=0;j<mat_cnt;j++){
  143. Matrix <T>& M =Mat (level,type,j);
  144. offset_indx[j][0]=data_offset; indx_ptr+=sizeof(size_t);
  145. data_offset+=M.Dim(0)*M.Dim(1)*sizeof(T); mem_size=mem::align_ptr(mem_size);
  146. for(size_t l=l0;l<l1;l++){
  147. Permutation<T>& Pr=Perm_R(l,type,j);
  148. offset_indx[j][1+4*(l-l0)+0]=data_offset;
  149. data_offset+=Pr.Dim()*sizeof(PERM_INT_T); mem_size=mem::align_ptr(mem_size);
  150. offset_indx[j][1+4*(l-l0)+1]=data_offset;
  151. data_offset+=Pr.Dim()*sizeof(T); mem_size=mem::align_ptr(mem_size);
  152. Permutation<T>& Pc=Perm_C(l,type,j);
  153. offset_indx[j][1+4*(l-l0)+2]=data_offset;
  154. data_offset+=Pc.Dim()*sizeof(PERM_INT_T); mem_size=mem::align_ptr(mem_size);
  155. offset_indx[j][1+4*(l-l0)+3]=data_offset;
  156. data_offset+=Pc.Dim()*sizeof(T); mem_size=mem::align_ptr(mem_size);
  157. }
  158. }
  159. }
  160. #pragma omp parallel for
  161. for(int tid=0;tid<omp_p;tid++){ // Copy data.
  162. char* indx_ptr=comp_data[0]+offset;
  163. //HeaderData& header=*(HeaderData*)indx_ptr;
  164. indx_ptr+=sizeof(HeaderData);
  165. Matrix<size_t> offset_indx(mat_cnt,1+(2+2)*(l1-l0), (size_t*)indx_ptr, false);
  166. for(size_t j=0;j<mat_cnt;j++){
  167. Matrix <T>& M =Mat (level,type,j);
  168. if(M.Dim(0)>0 && M.Dim(1)>0){
  169. size_t a=(M.Dim(0)*M.Dim(1)* tid )/omp_p;
  170. size_t b=(M.Dim(0)*M.Dim(1)*(tid+1))/omp_p;
  171. mem::memcopy(comp_data[0]+offset_indx[j][0]+a*sizeof(T), &M[0][a], (b-a)*sizeof(T));
  172. }
  173. for(size_t l=l0;l<l1;l++){
  174. Permutation<T>& Pr=Perm_R(l,type,j);
  175. Permutation<T>& Pc=Perm_C(l,type,j);
  176. if(Pr.Dim()>0){
  177. size_t a=(Pr.Dim()* tid )/omp_p;
  178. size_t b=(Pr.Dim()*(tid+1))/omp_p;
  179. mem::memcopy(comp_data[0]+offset_indx[j][1+4*(l-l0)+0]+a*sizeof(PERM_INT_T), &Pr.perm[a], (b-a)*sizeof(PERM_INT_T));
  180. mem::memcopy(comp_data[0]+offset_indx[j][1+4*(l-l0)+1]+a*sizeof( T), &Pr.scal[a], (b-a)*sizeof( T));
  181. }
  182. if(Pc.Dim()>0){
  183. size_t a=(Pc.Dim()* tid )/omp_p;
  184. size_t b=(Pc.Dim()*(tid+1))/omp_p;
  185. mem::memcopy(comp_data[0]+offset_indx[j][1+4*(l-l0)+2]+a*sizeof(PERM_INT_T), &Pc.perm[a], (b-a)*sizeof(PERM_INT_T));
  186. mem::memcopy(comp_data[0]+offset_indx[j][1+4*(l-l0)+3]+a*sizeof( T), &Pc.scal[a], (b-a)*sizeof( T));
  187. }
  188. }
  189. }
  190. }
  191. return offset+indx_size+mem_size;
  192. }
  193. template <class T>
  194. void PrecompMat<T>::Save2File(const char* fname, bool replace){
  195. FILE* f=fopen(fname,"r");
  196. if(f!=NULL) { //File already exists.
  197. fclose(f);
  198. if(!replace) return;
  199. }
  200. f=fopen(fname,"wb");
  201. if(f==NULL) return;
  202. int tmp;
  203. tmp=sizeof(T);
  204. fwrite(&tmp,sizeof(int),1,f);
  205. tmp=(scale_invar?1:0);
  206. fwrite(&tmp,sizeof(int),1,f);
  207. tmp=max_depth;
  208. fwrite(&tmp,sizeof(int),1,f);
  209. for(size_t i=0;i<mat.size();i++){
  210. int n=mat[i].size();
  211. fwrite(&n,sizeof(int),1,f);
  212. for(int j=0;j<n;j++){
  213. Matrix<T>& M=mat[i][j];
  214. int n1=M.Dim(0);
  215. fwrite(&n1,sizeof(int),1,f);
  216. int n2=M.Dim(1);
  217. fwrite(&n2,sizeof(int),1,f);
  218. if(n1*n2>0)
  219. fwrite(&M[0][0],sizeof(T),n1*n2,f);
  220. }
  221. }
  222. fclose(f);
  223. }
  224. #define MY_FREAD(a,b,c,d) { \
  225. size_t r_cnt=fread(a,b,c,d); \
  226. if(r_cnt!=c){ \
  227. fputs ("Reading error ",stderr); \
  228. exit (-1); \
  229. } }
  230. template <class T>
  231. void PrecompMat<T>::LoadFile(const char* fname, MPI_Comm comm){
  232. Profile::Tic("LoadMatrices",&comm,true,3);
  233. Profile::Tic("ReadFile",&comm,true,4);
  234. size_t f_size=0;
  235. char* f_data=NULL;
  236. int np, myrank;
  237. MPI_Comm_size(comm,&np);
  238. MPI_Comm_rank(comm,&myrank);
  239. if(myrank==0){
  240. FILE* f=fopen(fname,"rb");
  241. if(f==NULL){
  242. f_size=0;
  243. }else{
  244. struct stat fileStat;
  245. if(stat(fname,&fileStat) < 0) f_size=0;
  246. else f_size=fileStat.st_size;
  247. //fseek (f, 0, SEEK_END);
  248. //f_size=ftell (f);
  249. }
  250. if(f_size>0){
  251. f_data=mem::aligned_new<char>(f_size);
  252. fseek (f, 0, SEEK_SET);
  253. MY_FREAD(f_data,sizeof(char),f_size,f);
  254. fclose(f);
  255. }
  256. }
  257. Profile::Toc();
  258. Profile::Tic("Broadcast",&comm,true,4);
  259. MPI_Bcast( &f_size, sizeof(size_t), MPI_BYTE, 0, comm );
  260. if(f_size==0){
  261. Profile::Toc();
  262. Profile::Toc();
  263. return;
  264. }
  265. if(f_data==NULL) f_data=mem::aligned_new<char>(f_size);
  266. char* f_ptr=f_data;
  267. int max_send_size=1000000000;
  268. while(f_size>0){
  269. if(f_size>(size_t)max_send_size){
  270. MPI_Bcast( f_ptr, max_send_size, MPI_BYTE, 0, comm );
  271. f_size-=max_send_size;
  272. f_ptr+=max_send_size;
  273. }else{
  274. MPI_Bcast( f_ptr, f_size, MPI_BYTE, 0, comm );
  275. f_size=0;
  276. }
  277. }
  278. f_ptr=f_data;
  279. {
  280. int tmp;
  281. tmp=*(int*)f_ptr; f_ptr+=sizeof(int);
  282. assert(tmp==sizeof(T));
  283. tmp=*(int*)f_ptr; f_ptr+=sizeof(int);
  284. scale_invar=tmp;
  285. int max_depth_;
  286. max_depth_=*(int*)f_ptr; f_ptr+=sizeof(int);
  287. size_t mat_size=(size_t)Type_Count*(scale_invar?1:max_depth_+PRECOMP_MIN_DEPTH);
  288. if(mat.size()<mat_size){
  289. max_depth=max_depth_;
  290. mat.resize(mat_size);
  291. }
  292. for(size_t i=0;i<mat_size;i++){
  293. int n;
  294. n=*(int*)f_ptr; f_ptr+=sizeof(int);
  295. if(mat[i].size()<(size_t)n)
  296. mat[i].resize(n);
  297. for(int j=0;j<n;j++){
  298. Matrix<T>& M=mat[i][j];
  299. int n1;
  300. n1=*(int*)f_ptr; f_ptr+=sizeof(int);
  301. int n2;
  302. n2=*(int*)f_ptr; f_ptr+=sizeof(int);
  303. if(n1*n2>0){
  304. M.Resize(n1,n2);
  305. mem::memcopy(&M[0][0], f_ptr, sizeof(T)*n1*n2); f_ptr+=sizeof(T)*n1*n2;
  306. }
  307. }
  308. }
  309. // Resize perm_r, perm_c
  310. perm_r.resize((max_depth_+PRECOMP_MIN_DEPTH)*Type_Count);
  311. perm_c.resize((max_depth_+PRECOMP_MIN_DEPTH)*Type_Count);
  312. for(size_t i=0;i<perm_r.size();i++){
  313. perm_r[i].resize(500);
  314. perm_c[i].resize(500);
  315. }
  316. }
  317. mem::aligned_delete<char>(f_data);
  318. Profile::Toc();
  319. Profile::Toc();
  320. }
  321. #undef MY_FREAD
  322. #undef PRECOMP_MIN_DEPTH
  323. template <class T>
  324. std::vector<T>& PrecompMat<T>::RelativeTrgCoord(){
  325. return rel_trg_coord;
  326. }
  327. template <class T>
  328. bool PrecompMat<T>::ScaleInvar(){
  329. return scale_invar;
  330. }
  331. }//end namespace