Cyclops Tensor Framework
parallel arithmetic on multidimensional arrays
fast_tensor_ctr.cxx
Go to the documentation of this file.
1 /*Copyright (c) 2014, Edgar Solomonik, all rights reserved.*/
9 #include <ctf.hpp>
10 using namespace CTF;
11 
13  int dim = tsr.order;
14  Tensor<> ptsr(dim, tsr.lens, tsr.sym, *tsr.wrld);
15  char str[dim];
16  char pstr[dim];
17  for (int i=0; i<dim; i++){
18  str[i] = 'a'+i;
19  pstr[i] = 'a'+i;
20  }
21  for (int i=0; i<dim; i++){
22  for (int j=i+1; j<dim; j++){
23  pstr[i] = str[j];
24  pstr[j] = str[i];
25  ptsr[str] += tsr[str];
26  ptsr[str] -= tsr[pstr];
27  if (ptsr.norm2() > 1.E-6) return false;
28  pstr[i] = str[i];
29  pstr[j] = str[j];
30  }
31  }
32  return true;
33 }
34 
35 int64_t fact(int64_t n){
36  int64_t f = 1;
37  for (int64_t i=1; i<=n; i++){
38  f*=i;
39  }
40  return f;
41 }
42 
43 int64_t choose(int64_t n, int64_t k){
44  return fact(n)/(fact(k)*fact(n-k));
45 }
46 
47 int64_t chchoose(int64_t n, int64_t k){
48  return fact(n+k-1)/(fact(k)*fact(n-1));
49 }
50 
51 void chi(char const * idx,
52  int idx_len,
53  int p_len,
54  int q_len,
55  int * npair,
56  char *** idx_p,
57  char *** idx_q){
58  int np;
59 
60 // printf("entering (i<%d>,j<%d>) in chi(k<%d>)\n",p_len,q_len,idx_len);
61 
62  if (p_len+q_len > idx_len){
63  *npair = 0;
64  return;
65  }
66  if (idx_len == 0 || (p_len == 0 && q_len == 0)){
67  *npair = 1;
68  char ** ip = (char**)malloc(sizeof(char*));
69  char ** iq = (char**)malloc(sizeof(char*));
70  *idx_p = ip;
71  *idx_q = iq;
72  return;
73  }
74 
75  np = choose(idx_len, p_len);
76  np *= choose(idx_len-p_len, q_len);
77  *npair = np;
78 // printf("expecting %d pairs\n",np);
79 
80  char ** ip = (char**)malloc(sizeof(char*)*np);
81  char ** iq = (char**)malloc(sizeof(char*)*np);
82 
83  *idx_p = ip;
84  *idx_q = iq;
85 
86  for (int i=0; i<np; i++){
87  ip[i] = (char*)malloc(sizeof(char)*p_len);
88  iq[i] = (char*)malloc(sizeof(char)*q_len);
89  }
90 
91  if (q_len == 0){
92  char ** n1_ip;
93  char ** qnull;
94  int n1_len;
95  chi(idx, idx_len-1, p_len-1, 0, &n1_len, &n1_ip, &qnull);
96 
97  for (int i=0; i<n1_len; i++){
98  memcpy(ip[i], n1_ip[i], sizeof(char)*(p_len-1));
99  ip[i][p_len-1] = idx[idx_len-1];
100  }
101 
102  char ** n2_ip;
103  int n2_len;
104  chi(idx, idx_len-1, p_len, 0, &n2_len, &n2_ip, &qnull);
105  assert(n2_len + n1_len == np);
106 
107  for (int i=0; i<n2_len; i++){
108  memcpy(ip[i+n1_len], n2_ip[i], sizeof(char)*p_len);
109  }
110  } else if (p_len == 0){
111  char ** n1_iq;
112  char ** pnull;
113  int n1_len;
114  chi(idx, idx_len-1, 0, q_len-1, &n1_len, &pnull, &n1_iq);
115 
116  for (int i=0; i<n1_len; i++){
117  memcpy(iq[i], n1_iq[i], sizeof(char)*(q_len-1));
118  iq[i][q_len-1] = idx[idx_len-1];
119  }
120 
121  char ** n2_iq;
122  int n2_len;
123  chi(idx, idx_len-1, 0, q_len, &n2_len, &pnull, &n2_iq);
124  assert(n2_len + n1_len == np);
125 
126  for (int i=0; i<n2_len; i++){
127  memcpy(iq[i+n1_len], n2_iq[i], sizeof(char)*q_len);
128  }
129  } else {
130  char ** n1_ip;
131  char ** n1_iq;
132  int n1_len;
133  chi(idx, idx_len-1, p_len-1, q_len, &n1_len, &n1_ip, &n1_iq);
134 
135  assert(n1_len<=np);
136  for (int i=0; i<n1_len; i++){
137  memcpy(ip[i], n1_ip[i], sizeof(char)*(p_len-1));
138  ip[i][p_len-1] = idx[idx_len-1];
139  memcpy(iq[i], n1_iq[i], sizeof(char)*q_len);
140  }
141 
142  char ** n2_ip;
143  char ** n2_iq;
144  int n2_len;
145  chi(idx, idx_len-1, p_len, q_len-1, &n2_len, &n2_ip, &n2_iq);
146 
147  for (int i=0; i<n2_len; i++){
148  memcpy(ip[i+n1_len], n2_ip[i], sizeof(char)*p_len);
149  memcpy(iq[i+n1_len], n2_iq[i], sizeof(char)*(q_len-1));
150  iq[i+n1_len][q_len-1] = idx[idx_len-1];
151  }
152 
153  char ** n3_ip;
154  char ** n3_iq;
155  int n3_len;
156  chi(idx, idx_len-1, p_len, q_len, &n3_len, &n3_ip, &n3_iq);
157 
158  for (int i=0; i<n3_len; i++){
159  memcpy(ip[i+n1_len+n2_len], n3_ip[i], sizeof(char)*p_len);
160  memcpy(iq[i+n1_len+n2_len], n3_iq[i], sizeof(char)*q_len);
161  }
162 
163  assert(n1_len+n2_len+n3_len==np);
164 
165  }
166  /* printf("exiting (i<%d>,j<%d>) in chi(k<%d>) with npair=%d\n",p_len,q_len,idx_len,*npair);
167  for (int i=0; i<*npair; i++){
168  printf("(");
169  for (int j=0; j<p_len; j++){
170  printf("%c",ip[i][j]);
171  }
172  printf(", ");
173  for (int j=0; j<q_len; j++){
174  printf("%c",iq[i][j]);
175  }
176  printf(")\n");
177  }*/
178 }
179 
180 void chi(char const * idx,
181  int idx_len,
182  int p_len,
183  int * npair,
184  char *** idx_p){
185  char ** idx_q;
186 
187  chi(idx, idx_len, p_len, idx_len-p_len, npair, idx_p, &idx_q);
188 
189 }
190 
191 int fast_tensor_ctr(int n,
192  int s,
193  int t,
194  int v,
195  World &ctf){
196  int rank, i;
197  int64_t * indices, size;
198  double * values;
199 
200  MPI_Comm_rank(ctf.comm, &rank);
201 
202  int len_A[s+v];
203  int len_B[t+v];
204  int len_C[s+t];
205 
206  int sym_A[s+v];
207  int sym_B[t+v];
208  int sym_C[s+t];
209 
210  char idx_A[s+v];
211  char idx_B[t+v];
212  char idx_C[s+t];
213 
214  for (i=0; i<s+v; i++){
215  len_A[i] = n;
216  sym_A[i] = NS;
217  if (i<s)
218  idx_A[i] = 'a'+i;
219  else
220  idx_A[i] = 'a'+(s+t)+(i-s);
221  }
222  for (i=0; i<t+v; i++){
223  len_B[i] = n;
224  sym_B[i] = NS;
225  if (i<t)
226  idx_B[i] = 'a'+s+i;
227  else
228  idx_B[i] = 'a'+(s+t)+(i-t);
229  }
230  for (i=0; i<s+t; i++){
231  len_C[i] = n;
232  sym_C[i] = NS;
233  idx_C[i] = 'a'+i;
234  }
235 
236  Tensor<> A(s+v, len_A, sym_A, ctf, "A", 1);
237  Tensor<> B(t+v, len_B, sym_B, ctf, "B", 1);
238  Tensor<> C(s+t, len_C, sym_C, ctf, "C", 1);
239  Tensor<> C_int(s+t, len_C, sym_C, ctf, "C_psym", 1);
240  Tensor<> C_ans(s+t, len_C, sym_C, ctf, "C_ans", 1);
241 
242  Vector<> vec(n, ctf, "vec", 1);
243  srand48(13);
244  vec.read_local(&size, &indices, &values);
245  for (i=0; i<size; i++){
246  values[i] = drand48();
247  }
248  vec.write(size, indices, values);
249  free(indices);
250  free(values);
251 
252  for (i=0; i<s+v; i++){
253  A[idx_A] += vec[idx_A+i];
254  }
255 
256  vec.read_local(&size, &indices, &values);
257  for (i=0; i<size; i++){
258  values[i] = drand48();
259  }
260  vec.write(size, indices, values);
261  free(indices);
262  free(values);
263 
264  for (i=0; i<t+v; i++){
265  B[idx_B] += vec[idx_B+i];
266  }
267 
268 
269  C_int[idx_C] += A[idx_A]*B[idx_B];
270 
271  int ncperms;
272  char ** idx_As;
273  char ** idx_Bt;
274 
275  chi(idx_C, s+t, s, t, &ncperms, &idx_As, &idx_Bt);
276 
277  for (i=0; i<ncperms; i++){
278  char idx_C_int[s+t];
279  memcpy(idx_C_int, idx_As[i], sizeof(char)*s);
280  memcpy(idx_C_int+s, idx_Bt[i], sizeof(char)*t);
281  C_ans[idx_C] += C_int[idx_C_int];
282  }
283  bool is_C_sym = check_sym(C_ans);
284  if (is_C_sym) printf("C_ans is symmetric\n");
285  else printf("C_ans is NOT symmetric!!\n");
286 
287  int len_Z[s+v+t];
288  int sym_Z[s+v+t];
289  char idx_Z[s+v+t];
290  for (i=0; i<s+v+t; i++){
291  len_Z[i] = n;
292  sym_Z[i] = NS;
293  idx_Z[i] = 'a'+i;
294  }
295 
296  Tensor<> Z_A_ops(s+v+t, len_Z, sym_Z, ctf, "Z_A", 1);
297  Tensor<> Z_B_ops(s+v+t, len_Z, sym_Z, ctf, "Z_B", 1);
298  Tensor<> Z_mults(s+v+t, len_Z, sym_Z, ctf, "Z", 1);
299 
300  int nAperms;
301  char ** idx_Asv;
302 
303  chi(idx_Z, s+t+v, s+v, &nAperms, &idx_Asv);
304 
305  for (i=0; i<nAperms; i++){
306  Z_A_ops[idx_Z] += A[idx_Asv[i]];
307  }
308 
309  int nBperms;
310  char ** idx_Btv;
311 
312  chi(idx_Z, s+t+v, t+v, &nBperms, &idx_Btv);
313 
314  for (i=0; i<nBperms; i++){
315  Z_B_ops[idx_Z] += B[idx_Btv[i]];
316  }
317 
318 
319  Z_mults[idx_Z] = Z_A_ops[idx_Z]*Z_B_ops[idx_Z];
320 
321  bool is_Z_sym = check_sym(Z_mults);
322  if (is_Z_sym) printf("Z_mults is symmetric\n");
323  else printf("Z_mults is NOT symmetric!!\n");
324 
325  memcpy(idx_Z,idx_C,(s+t)*sizeof(char));
326  for (i=s+t; i<s+t+v; i++){
327  idx_Z[i] = idx_Z[s+t-1]+(i-s-t+1);
328  }
329 
330  C[idx_C]+=Z_mults[idx_Z];
331 
332  Tensor<> V(s+t, len_C, sym_C, ctf, "V");
333  for (int r=0; r<v; r++){
334  for (int p=std::max(v-t-r,0); p<=v-r; p++){
335  for (int q=std::max(v-s-r,0); q<=v-p-r; q++){
336  double prefact = (double)(choose(v,r)*choose(v-r,p)*choose(v-p-r,q)*pow(n,v-p-q-r));
337  char idx_kr[r];
338  for (i=0; i<r; i++){
339  idx_kr[i] = 'a'+s+t+i;
340  }
341  char idx_kp[p];
342  for (i=0; i<p; i++){
343  idx_kp[i] = 'a'+s+t+r+i;
344  }
345  char idx_kq[q];
346  for (i=0; i<q; i++){
347  idx_kq[i] = 'a'+s+t+r+p+i;
348  }
349 
350  Tensor<> V_A_ops(s+t+r, len_Z, sym_Z, ctf, "V_A_ops");
351  char idx_VA[s+t+r];
352  memcpy(idx_VA,idx_C,(s+t)*sizeof(char));
353  memcpy(idx_VA+s+t,idx_kr,r*sizeof(char));
354  //printf("r=%d,p=%d,q=%d\n",r,p,q);
355  int nvAperms;
356  char ** idx_VAsvpr;
357  chi(idx_C, s+t, s+v-p-r, &nvAperms, &idx_VAsvpr);
358  for (i=0; i<nvAperms; i++){
359  char idx_VAA[s+v];
360  memcpy(idx_VAA, idx_VAsvpr[i], (s+v-p-r)*sizeof(char));
361  memcpy(idx_VAA+s+v-p-r, idx_kr, r*sizeof(char));
362  memcpy(idx_VAA+s+v-p, idx_kp, p*sizeof(char));
363  V_A_ops[idx_VA] += A[idx_VAA];
364  }
365 
366  Tensor<> V_B_ops(s+t+r, len_Z, sym_Z, ctf, "V_B_ops");
367  char idx_VB[s+t+r];
368  memcpy(idx_VB,idx_C,(s+t)*sizeof(char));
369  memcpy(idx_VB+s+t,idx_kr,r*sizeof(char));
370 
371  int nvBperms;
372  char ** idx_VBtvqr;
373  chi(idx_C, s+t, t+v-q-r, &nvBperms, &idx_VBtvqr);
374  for (i=0; i<nvBperms; i++){
375  char idx_VBB[t+v];
376  memcpy(idx_VBB, idx_VBtvqr[i], (t+v-q-r)*sizeof(char));
377  memcpy(idx_VBB+t+v-q-r, idx_kr, r*sizeof(char));
378  memcpy(idx_VBB+t+v-q, idx_kq, q*sizeof(char));
379  /*for (i=0; i<t+v; i++){
380  printf("index %d of B is %c\n",i, idx_VBB[i]);
381  }
382  for (i=0; i<s+t+r; i++){
383  printf("index %d of V_B is %c\n",i, idx_VB[i]);
384  }*/
385  V_B_ops[idx_VB] += B[idx_VBB];
386  }
387 
388  V[idx_C] += prefact*V_A_ops[idx_VA]*V_B_ops[idx_VB];
389  }
390  }
391  }
392  Tensor<> W(s+t, len_C, sym_C, ctf, "W");
393  for (int r=1; r<=std::min(s,t); r++){
394  char idx_kr[r];
395  for (int i=0; i<r; i++){
396  idx_kr[i] = 'a'+s+t+i;
397  }
398  char idx_kv[v];
399  for (int i=0; i<v; i++){
400  idx_kv[i] = 'a'+s+t+r+i;
401  }
402  Tensor<> U(s+t-r, len_C, sym_C, ctf, "U");
403  char idx_U[s+t-r];
404  char idx_UA[s+v];
405  char idx_UB[t+v];
406  memcpy(idx_U, idx_kr, sizeof(char)*r);
407  memcpy(idx_UA, idx_kr, sizeof(char)*r);
408  memcpy(idx_UB, idx_kr, sizeof(char)*r);
409  memcpy(idx_UA+r, idx_kv, sizeof(char)*v);
410  memcpy(idx_UB+r, idx_kv, sizeof(char)*v);
411  int npermU;
412  char ** idxj, ** idxl;
413  chi(idx_C, s+t-2*r, s-r, t-r, &npermU, &idxj, &idxl);
414  memcpy(idx_U+r,idx_C,s+t-2*r);
415  for (int i=0; i<npermU; i++){
416  memcpy(idx_UA+r+v,idxj[i],s-r);
417  memcpy(idx_UB+r+v,idxl[i],t-r);
418  U[idx_U] += A[idx_UA]*B[idx_UB];
419  }
420  int npermW;
421  char ** idxh, ** idxr;
422  chi(idx_C, s+t, r, s+t-2*r, &npermW, &idxr, &idxh);
423  for (int i=0; i<npermW; i++){
424  memcpy(idx_U,idxr[i],r);
425  memcpy(idx_U+r,idxh[i],s+t-2*r);
426  W[idx_C] += U[idx_U];
427  }
428  }
429 
430  C[idx_C] -= V[idx_C];
431  C[idx_C] -= W[idx_C];
432 
433  C[idx_C] -= C_ans[idx_C];
434 
435  double nrm = C.norm2();
436 
437  printf("error 2-norm is %.4E\n",nrm);
438 
439  int pass = (nrm <=1.E-3);
440 
441  if (rank == 0){
442  if (pass) printf("{ fast symmetric tensor contraction algorithm test } passed\n");
443  else printf("{ fast symmetric tensor contraction algorithm test } failed\n");
444  }
445  return pass;
446 }
447 
448 #ifndef TEST_SUITE
449 char* getCmdOption(char ** begin,
450  char ** end,
451  const std::string & option){
452  char ** itr = std::find(begin, end, option);
453  if (itr != end && ++itr != end){
454  return *itr;
455  }
456  return 0;
457 }
458 
459 
460 int main(int argc, char ** argv){
461  int rank, np, n, s,t,v;
462  int const in_num = argc;
463  char ** input_str = argv;
464 
465  MPI_Init(&argc, &argv);
466  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
467  MPI_Comm_size(MPI_COMM_WORLD, &np);
468 
469  if (getCmdOption(input_str, input_str+in_num, "-n")){
470  n = atoi(getCmdOption(input_str, input_str+in_num, "-n"));
471  if (n < 0) n = 4;
472  } else n = 4;
473 
474  if (getCmdOption(input_str, input_str+in_num, "-s")){
475  s = atoi(getCmdOption(input_str, input_str+in_num, "-s"));
476  if (s < 0) s = 1;
477  } else s = 1;
478 
479  if (getCmdOption(input_str, input_str+in_num, "-t")){
480  t = atoi(getCmdOption(input_str, input_str+in_num, "-t"));
481  if (t < 0) t = 1;
482  } else t = 1;
483 
484  if (getCmdOption(input_str, input_str+in_num, "-v")){
485  v = atoi(getCmdOption(input_str, input_str+in_num, "-v"));
486  if (v < 0) v = 1;
487  } else v = 1;
488 
489  {
490  World dw(MPI_COMM_WORLD, argc, argv);
491  if (rank == 0){
492  printf("Contracting symmetric A of order %d with B of order %d into C of order %d, all with dimension %d\n",s+v,t+v,s+t,n);
493  }
494  int pass = fast_tensor_ctr(n, s, t, v, dw);
495  assert(pass);
496  }
497 
498  MPI_Finalize();
499  return 0;
500 }
506 #endif
507 
int64_t choose(int64_t n, int64_t k)
int main(int argc, char **argv)
bool check_sym(Tensor<> &tsr)
int * sym
symmetries among tensor dimensions
def rank(self)
Definition: core.pyx:312
Vector class which encapsulates a 1D tensor.
Definition: vector.h:14
Definition: common.h:37
int64_t chchoose(int64_t n, int64_t k)
an instance of the CTF library (world) on a MPI communicator
Definition: world.h:19
int order
number of tensor dimensions
void read_local(int64_t *npair, int64_t **global_idx, dtype **data, bool unpack_sym=false) const
Using get_local_data(), which returns an array that must be freed with delete [], is more efficient...
Definition: tensor.cxx:182
string
Definition: core.pyx:456
dtype norm2()
computes the frobenius norm of the tensor (needs sqrt()!)
Definition: tensor.h:811
CTF::World * wrld
distributed processor context on which tensor is defined
int64_t fact(int64_t n)
int * lens
unpadded tensor edge lengths
void chi(char const *idx, int idx_len, int p_len, int q_len, int *npair, char ***idx_p, char ***idx_q)
char * getCmdOption(char **begin, char **end, const std::string &option)
Definition: apsp.cxx:17
int fast_tensor_ctr(int n, int s, int t, int v, World &ctf)
an instance of a tensor within a CTF world
Definition: tensor.h:74
void write(int64_t npair, int64_t const *global_idx, dtype const *data)
writes in values associated with any set of indices The sparse data is defined in coordinate format...
Definition: tensor.cxx:264
MPI_Comm comm
set of processors making up this world
Definition: world.h:22
def np(self)
Definition: core.pyx:315