Cyclops Tensor Framework
parallel arithmetic on multidimensional arrays
fast_as_as_sy_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 
12 //gets parity of (a,b) in chi(c)
13 int parity(char const * a, char const * b, char const * c, int len_A, int len_B){
14  char ab[len_A+len_B];
15  for (int i=0; i<len_A; i++){
16  ab[i] = a[i];
17  }
18  for (int i=0; i<len_B; i++){
19  ab[len_A+i] = b[i];
20  }
21  /*for (int i=0; i<len_A+len_B; i++){
22  printf("%c", ab[i]);
23  }
24  printf("\n");
25  for (int i=0; i<len_A+len_B; i++){
26  printf("%c", c[i]);
27  }
28  printf("\n");*/
29  int par = 0;
30  for (int i=0; i<len_A+len_B; i++){
31  if (ab[i] != c[i]){
32  int j;
33  for (j=i+1; j<len_A+len_B; j++){
34  if (ab[j] == c[i]){
35  ab[j] = ab[i];
36  ab[i] = c[i];
37  par++;
38  break;
39  }
40  }
41  assert(j<len_A+len_B);
42  }
43  }
44  //printf("parity = %d\n", par);
45  return par;
46 }
47 
48 //gets parity of a in chi(c)
49 int parity(char const * a, char const * c, int len_A, int len_C){
50  if (len_A == len_C) return parity(a, NULL, c, len_A, 0);
51  else {
52  char b[len_C-len_A];
53  int ib = 0;
54  for (int i=0; i<len_C; i++){
55  int j;
56  for (j=0; j<len_A; j++){
57  if (c[i] == a[j]) break;
58  }
59  if (j==len_A){
60  b[ib] = c[i];
61  ib++;
62  }
63  }
64  assert(ib == len_C-len_A);
65  return parity(a, b, c, len_A, len_C-len_A);
66  }
67 }
68 
69 double sign(int par){
70  //printf("par = %d, %d\n",par, par%2);
71  if (par%2 == 1) return -1.0;
72  else return 1.0;
73 }
74 
75 
77  int64_t * indices, size;
78  double * values;
79  int dim = tsr.order;
80  if (dim == 1){
81  srand48(seed);
82  tsr.read_local(&size, &indices, &values);
83  for (int64_t i=0; i<size; i++){
84  values[i] = drand48();
85  }
86  tsr.write(size, indices, values);
87  free(indices);
88  free(values);
89  } else {
90  Tensor<> tsr_n1(dim-1, tsr.lens, tsr.sym, *tsr.wrld);
91  get_rand_as_tsr(tsr_n1, seed+1);
92  Vector<> v(tsr.lens[0], *tsr.wrld);
93  srand48(2*seed);
94  v.read_local(&size, &indices, &values);
95  for (int64_t i=0; i<size; i++){
96  values[i] = drand48();
97  }
98  v.write(size, indices, values);
99  free(indices);
100  free(values);
101  char str[dim];
102  char str_n1[dim-1];
103  for (int i=0; i<dim; i++){
104  str[i] = 'a'+i;
105  }
106  double sgn = 1.0;
107  for (int i=0; i<dim; i++){
108  for (int j=0; j<dim-1; j++){
109  if (j>=i) str_n1[j] = str[j+1];
110  else str_n1[j] =str[j];
111  }
112  //printf("sgn = %lf %lf\n",sgn,sign(parity(str+i, str_n1, str, 1, dim-1)));
113  assert(sign(parity(str+i, str_n1, str, 1, dim-1)) == sgn);
114  tsr[str] += sgn*v[str+i]*tsr_n1[str_n1];
115  sgn *= -1.0;
116  }
117  }
118 }
119 
121  int dim = tsr.order;
122  Tensor<> ptsr(dim, tsr.lens, tsr.sym, *tsr.wrld);
123  char str[dim];
124  char pstr[dim];
125  for (int i=0; i<dim; i++){
126  str[i] = 'a'+i;
127  pstr[i] = 'a'+i;
128  }
129  for (int i=0; i<dim; i++){
130  for (int j=i+1; j<dim; j++){
131  pstr[i] = str[j];
132  pstr[j] = str[i];
133  ptsr[str] += tsr[str];
134  ptsr[str] += tsr[pstr];
135  if (ptsr.norm2() > 1.E-6) return false;
136  pstr[i] = str[i];
137  pstr[j] = str[j];
138  }
139  }
140  return true;
141 }
142 
144  int dim = tsr.order;
145  Tensor<> ptsr(dim, tsr.lens, tsr.sym, *tsr.wrld);
146  char str[dim];
147  char pstr[dim];
148  for (int i=0; i<dim; i++){
149  str[i] = 'a'+i;
150  pstr[i] = 'a'+i;
151  }
152  for (int i=0; i<dim; i++){
153  for (int j=i+1; j<dim; j++){
154  pstr[i] = str[j];
155  pstr[j] = str[i];
156  ptsr[str] += tsr[str];
157  ptsr[str] -= tsr[pstr];
158  if (ptsr.norm2() > 1.E-6) return false;
159  pstr[i] = str[i];
160  pstr[j] = str[j];
161  }
162  }
163  return true;
164 }
165 
166 int64_t fact(int64_t n){
167  int64_t f = 1;
168  for (int64_t i=1; i<=n; i++){
169  f*=i;
170  }
171  return f;
172 }
173 
174 int64_t choose(int64_t n, int64_t k){
175  return fact(n)/(fact(k)*fact(n-k));
176 }
177 
178 int64_t chchoose(int64_t n, int64_t k){
179  return fact(n+k-1)/(fact(k)*fact(n-1));
180 }
181 
182 void chi(char const * idx,
183  int idx_len,
184  int p_len,
185  int q_len,
186  int * npair,
187  char *** idx_p,
188  char *** idx_q){
189  int np;
190 
191 // printf("entering (i<%d>,j<%d>) in chi(k<%d>)\n",p_len,q_len,idx_len);
192 
193  if (p_len+q_len > idx_len){
194  *npair = 0;
195  return;
196  }
197  if (idx_len == 0 || (p_len == 0 && q_len == 0)){
198  *npair = 1;
199  char ** ip = (char**)malloc(sizeof(char*));
200  char ** iq = (char**)malloc(sizeof(char*));
201  *idx_p = ip;
202  *idx_q = iq;
203  return;
204  }
205 
206  np = choose(idx_len, p_len);
207  np *= choose(idx_len-p_len, q_len);
208  *npair = np;
209 // printf("expecting %d pairs\n",np);
210 
211  char ** ip = (char**)malloc(sizeof(char*)*np);
212  char ** iq = (char**)malloc(sizeof(char*)*np);
213 
214  *idx_p = ip;
215  *idx_q = iq;
216 
217  for (int i=0; i<np; i++){
218  ip[i] = (char*)malloc(sizeof(char)*p_len);
219  iq[i] = (char*)malloc(sizeof(char)*q_len);
220  }
221 
222  if (q_len == 0){
223  char ** n1_ip;
224  char ** qnull;
225  int n1_len;
226  chi(idx, idx_len-1, p_len-1, 0, &n1_len, &n1_ip, &qnull);
227 
228  for (int i=0; i<n1_len; i++){
229  memcpy(ip[i], n1_ip[i], sizeof(char)*(p_len-1));
230  ip[i][p_len-1] = idx[idx_len-1];
231  }
232 
233  char ** n2_ip;
234  int n2_len;
235  chi(idx, idx_len-1, p_len, 0, &n2_len, &n2_ip, &qnull);
236  assert(n2_len + n1_len == np);
237 
238  for (int i=0; i<n2_len; i++){
239  memcpy(ip[i+n1_len], n2_ip[i], sizeof(char)*p_len);
240  }
241  } else if (p_len == 0){
242  char ** n1_iq;
243  char ** pnull;
244  int n1_len;
245  chi(idx, idx_len-1, 0, q_len-1, &n1_len, &pnull, &n1_iq);
246 
247  for (int i=0; i<n1_len; i++){
248  memcpy(iq[i], n1_iq[i], sizeof(char)*(q_len-1));
249  iq[i][q_len-1] = idx[idx_len-1];
250  }
251 
252  char ** n2_iq;
253  int n2_len;
254  chi(idx, idx_len-1, 0, q_len, &n2_len, &pnull, &n2_iq);
255  assert(n2_len + n1_len == np);
256 
257  for (int i=0; i<n2_len; i++){
258  memcpy(iq[i+n1_len], n2_iq[i], sizeof(char)*q_len);
259  }
260  } else {
261  char ** n1_ip;
262  char ** n1_iq;
263  int n1_len;
264  chi(idx, idx_len-1, p_len-1, q_len, &n1_len, &n1_ip, &n1_iq);
265 
266  assert(n1_len<=np);
267  for (int i=0; i<n1_len; i++){
268  memcpy(ip[i], n1_ip[i], sizeof(char)*(p_len-1));
269  ip[i][p_len-1] = idx[idx_len-1];
270  memcpy(iq[i], n1_iq[i], sizeof(char)*q_len);
271  }
272 
273  char ** n2_ip;
274  char ** n2_iq;
275  int n2_len;
276  chi(idx, idx_len-1, p_len, q_len-1, &n2_len, &n2_ip, &n2_iq);
277 
278  for (int i=0; i<n2_len; i++){
279  memcpy(ip[i+n1_len], n2_ip[i], sizeof(char)*p_len);
280  memcpy(iq[i+n1_len], n2_iq[i], sizeof(char)*(q_len-1));
281  iq[i+n1_len][q_len-1] = idx[idx_len-1];
282  }
283 
284  char ** n3_ip;
285  char ** n3_iq;
286  int n3_len;
287  chi(idx, idx_len-1, p_len, q_len, &n3_len, &n3_ip, &n3_iq);
288 
289  for (int i=0; i<n3_len; i++){
290  memcpy(ip[i+n1_len+n2_len], n3_ip[i], sizeof(char)*p_len);
291  memcpy(iq[i+n1_len+n2_len], n3_iq[i], sizeof(char)*q_len);
292  }
293 
294  assert(n1_len+n2_len+n3_len==np);
295 
296  }
297  /* printf("exiting (i<%d>,j<%d>) in chi(k<%d>) with npair=%d\n",p_len,q_len,idx_len,*npair);
298  for (int i=0; i<*npair; i++){
299  printf("(");
300  for (int j=0; j<p_len; j++){
301  printf("%c",ip[i][j]);
302  }
303  printf(", ");
304  for (int j=0; j<q_len; j++){
305  printf("%c",iq[i][j]);
306  }
307  printf(")\n");
308  }*/
309 }
310 
311 void chi(char const * idx,
312  int idx_len,
313  int p_len,
314  int * npair,
315  char *** idx_p){
316  char ** idx_q;
317 
318  chi(idx, idx_len, p_len, idx_len-p_len, npair, idx_p, &idx_q);
319 
320 }
321 
322 int fast_tensor_ctr(int n,
323  int s,
324  int t,
325  int v,
326  World &ctf){
327  int rank, i;
328  int64_t * indices, size;
329  double * values;
330 
331  MPI_Comm_rank(ctf.comm, &rank);
332  //MPI_Comm_size(MPI_COMM_WORLD, &num_pes);
333 
334  int len_A[s+v];
335  int len_B[t+v];
336  int len_C[s+t];
337 
338  int sym_A[s+v];
339  int sym_B[t+v];
340  int sym_C[s+t];
341 
342  char idx_A[s+v];
343  char idx_B[t+v];
344  char idx_C[s+t];
345 
346  for (i=0; i<s+v; i++){
347  len_A[i] = n;
348  sym_A[i] = NS;
349  if (i<s)
350  idx_A[i] = 'a'+i;
351  else
352  idx_A[i] = 'a'+(s+t)+(i-s);
353  }
354  for (i=0; i<t+v; i++){
355  len_B[i] = n;
356  sym_B[i] = NS;
357  /*if (i<t)
358  idx_B[i] = 'a'+s+i;
359  else
360  idx_B[i] = 'a'+(s+t)+(i-t);*/
361  if (i<v)
362  idx_B[i] = 'a'+(s+t)+i;
363  else
364  idx_B[i] = 'a'+s+(i-v);
365  }
366  for (i=0; i<s+t; i++){
367  len_C[i] = n;
368  sym_C[i] = NS;
369  idx_C[i] = 'a'+i;
370  }
371 
372  Tensor<> A(s+v, len_A, sym_A, ctf, "A", 1);
373  Tensor<> B(t+v, len_B, sym_B, ctf, "B", 1);
374  Tensor<> C(s+t, len_C, sym_C, ctf, "C", 1);
375  Tensor<> C_int(s+t, len_C, sym_C, ctf, "C_psym", 1);
376  Tensor<> C_ans(s+t, len_C, sym_C, ctf, "C_ans", 1);
377 
378  /*
379  srand48(13);
380  vec.read_local(&size, &indices, &values);
381  for (i=0; i<size; i++){
382  values[i] = drand48();
383  }
384  vec.write(size, indices, values);
385  free(indices);
386  free(values);
387 
388  for (i=0; i<s+v; i++){
389  A[idx_A] += vec[idx_A+i];
390  }*/
391  get_rand_as_tsr(A, 13);
392  //A.print();
393  assert(check_asym(A));
394  if (A.order > 1)
395  assert(!check_sym(A));
396 
397  Vector<> vec(n, ctf, "vec", 1);
398  vec.read_local(&size, &indices, &values);
399  for (i=0; i<size; i++){
400  values[i] = drand48();
401  }
402  vec.write(size, indices, values);
403  free(indices);
404  free(values);
405 
406  for (i=0; i<t+v; i++){
407  B[idx_B] += vec[idx_B+i];
408  }
409  //get_rand_as_tsr(B, 21);
410  assert(check_sym(B));
411  if (B.order > 1)
412  assert(!check_asym(B));
413 
414  C_int[idx_C] += A[idx_A]*B[idx_B];
415 
416  int ncperms;
417  char ** idx_As;
418  char ** idx_Bt;
419 
420  chi(idx_C, s+t, s, t, &ncperms, &idx_As, &idx_Bt);
421 
422  for (i=0; i<ncperms; i++){
423  char idx_C_int[s+t];
424  memcpy(idx_C_int, idx_As[i], sizeof(char)*s);
425  memcpy(idx_C_int+s, idx_Bt[i], sizeof(char)*t);
426  C_ans[idx_C] += sign(parity(idx_As[i], idx_Bt[i], idx_C, s, t))*C_int[idx_C_int];
427  }
428  bool is_C_asym = check_asym(C_ans);
429  if (is_C_asym) printf("C_ans is antisymmetric\n");
430  else printf("C_ans is NOT antisymmetric!!\n");
431 
432  int len_Z[s+v+t];
433  int sym_Z[s+v+t];
434  char idx_Z[s+v+t];
435  for (i=0; i<s+v+t; i++){
436  len_Z[i] = n;
437  sym_Z[i] = NS;
438  idx_Z[i] = 'a'+i;
439  }
440 
441  Tensor<> Z_A_ops(s+v+t, len_Z, sym_Z, ctf, "Z_A", 1);
442  Tensor<> Z_B_ops(s+v+t, len_Z, sym_Z, ctf, "Z_B", 1);
443  Tensor<> Z_mults(s+v+t, len_Z, sym_Z, ctf, "Z", 1);
444 
445  int nAperms;
446  char ** idx_Asv;
447 
448  chi(idx_Z, s+t+v, s+v, &nAperms, &idx_Asv);
449 
450  for (i=0; i<nAperms; i++){
451  Z_A_ops[idx_Z] += sign(parity(idx_Asv[i], idx_Z, s+v, s+t+v))*A[idx_Asv[i]];
452  }
453  bool is_A_asym = check_asym(Z_A_ops);
454  if (is_A_asym) printf("Z_A_ops is antisymmetric\n");
455  else printf("Z_A_ops is NOT antisymmetric!!\n");
456 
457  int nBperms;
458  char ** idx_Btv;
459 
460  chi(idx_Z, s+t+v, t+v, &nBperms, &idx_Btv);
461 
462  for (i=0; i<nBperms; i++){
463  Z_B_ops[idx_Z] += B[idx_Btv[i]];
464  }
465  bool is_B_sym = check_sym(Z_B_ops);
466  if (is_B_sym) printf("Z_B_ops is symmetric\n");
467  else printf("Z_B_ops is NOT symmetric!!\n");
468 
469 
470  Z_mults[idx_Z] = Z_A_ops[idx_Z]*Z_B_ops[idx_Z];
471 
472  bool is_Z_asym = check_asym(Z_mults);
473  if (is_Z_asym) printf("Z_mults is antisymmetric\n");
474  else printf("Z_mults is NOT antisymmetric!!\n");
475 
476  memcpy(idx_Z,idx_C,(s+t)*sizeof(char));
477  for (i=s+t; i<s+t+v; i++){
478  idx_Z[i] = idx_Z[s+t-1]+(i-s-t+1);
479  }
480 
481  C[idx_C]+=sign(t*v)*Z_mults[idx_Z];
482 
483  Tensor<> V(s+t, len_C, sym_C, ctf, "V");
484  for (int r=0; r<v; r++){
485  for (int p=std::max(v-t-r,0); p<=v-r; p++){
486  for (int q=std::max(v-s-r,0); q<=v-p-r; q++){
487  double prefact = (double)(choose(v,r)*choose(v-r,p)*choose(v-p-r,q)*pow(n,v-p-q-r));
488  double sgn_V = sign(p+1);
489 // printf("sgn_V = %lf\n",sgn_V);
490  prefact*= sgn_V;
491 
492  char idx_kr[r];
493  for (i=0; i<r; i++){
494  idx_kr[i] = 'a'+s+t+i;
495  }
496  char idx_kp[p];
497  for (i=0; i<p; i++){
498  idx_kp[i] = 'a'+s+t+r+i;
499  }
500  char idx_kq[q];
501  for (i=0; i<q; i++){
502  idx_kq[i] = 'a'+s+t+r+p+i;
503  }
504 
505  Tensor<> V_A_ops(s+t+r, len_Z, sym_Z, ctf, "V_A_ops");
506  char idx_VA[s+t+r];
507  memcpy(idx_VA,idx_C,(s+t)*sizeof(char));
508  memcpy(idx_VA+s+t,idx_kr,r*sizeof(char));
509  //printf("r=%d,p=%d,q=%d\n",r,p,q);
510  int nvAperms;
511  char ** idx_VAsvpr;
512  chi(idx_C, s+t, s+v-p-r, &nvAperms, &idx_VAsvpr);
513  for (i=0; i<nvAperms; i++){
514  char idx_VAA[s+v];
515  memcpy(idx_VAA, idx_VAsvpr[i], (s+v-p-r)*sizeof(char));
516  memcpy(idx_VAA+s+v-p-r, idx_kr, r*sizeof(char));
517  memcpy(idx_VAA+s+v-p, idx_kp, p*sizeof(char));
518  double sgn_VA = sign(parity(idx_VAsvpr[i], idx_C, s+v-p-r, s+t));
519  V_A_ops[idx_VA] += sgn_VA*A[idx_VAA];
520  }
521 
522  Tensor<> V_B_ops(s+t+r, len_Z, sym_Z, ctf, "V_B_ops");
523  char idx_VB[s+t+r];
524  memcpy(idx_VB,idx_C,(s+t)*sizeof(char));
525  memcpy(idx_VB+s+t,idx_kr,r*sizeof(char));
526 
527 
528  int nvBperms;
529  char ** idx_VBtvqr;
530  chi(idx_C, s+t, t+v-q-r, &nvBperms, &idx_VBtvqr);
531  for (i=0; i<nvBperms; i++){
532  char idx_VBB[t+v];
533  memcpy(idx_VBB, idx_VBtvqr[i], (t+v-q-r)*sizeof(char));
534  memcpy(idx_VBB+t+v-q-r, idx_kr, r*sizeof(char));
535  memcpy(idx_VBB+t+v-q, idx_kq, q*sizeof(char));
536  /*for (i=0; i<t+v; i++){
537  printf("index %d of B is %c\n",i, idx_VBB[i]);
538  }
539  for (i=0; i<s+t+r; i++){
540  printf("index %d of V_B is %c\n",i, idx_VB[i]);
541  }*/
542  V_B_ops[idx_VB] += B[idx_VBB];
543  }
544 
545  V[idx_C] += prefact*V_A_ops[idx_VA]*V_B_ops[idx_VB];
546  }
547  }
548  }
549  Tensor<> W(s+t, len_C, sym_C, ctf, "W");
550  for (int r=1; r<=std::min(s,t); r++){
551  char idx_kr[r];
552  for (int i=0; i<r; i++){
553  idx_kr[i] = 'a'+s+t+i;
554  }
555  char idx_kv[v];
556  for (int i=0; i<v; i++){
557  idx_kv[i] = 'a'+s+t+r+i;
558  }
559  Tensor<> U(s+t-r, len_C, sym_C, ctf, "U");
560  char idx_U[s+t-r];
561  char idx_UA[s+v];
562  char idx_UB[t+v];
563  memcpy(idx_U, idx_kr, sizeof(char)*r);
564  memcpy(idx_UA, idx_kr, sizeof(char)*r);
565  //memcpy(idx_UB, idx_kr, sizeof(char)*r);
566  memcpy(idx_UB+t+v-r, idx_kr, sizeof(char)*r);
567  int npermU;
568  char ** idxj, ** idxl;
569  chi(idx_C, s+t-2*r, s-r, t-r, &npermU, &idxj, &idxl);
570  memcpy(idx_U+r,idx_C,s+t-2*r);
571  for (int i=0; i<npermU; i++){
572  memcpy(idx_UA+r,idxj[i],s-r);
573  //memcpy(idx_UB+r,idxl[i],t-r);
574  memcpy(idx_UB+v,idxl[i],t-r);
575  memcpy(idx_UA+s, idx_kv, sizeof(char)*v);
576  //memcpy(idx_UB+t, idx_kv, sizeof(char)*v);
577  memcpy(idx_UB, idx_kv, sizeof(char)*v);
578 // double sgnU = sign(parity(idxj[i], idxl[i], idx_C, s-r, t-r));
579  U[idx_U] += A[idx_UA]*B[idx_UB];
580  }
581  int npermW1;
582  char ** idxh1;
583  chi(idx_C, s+t, s+t-r, &npermW1,&idxh1);
584  for (int j=0; j<npermW1; j++){
585  double sgn1 = sign(parity(idxh1[j], idx_C, s+t-r, s+t));
586  int npermW;
587  char ** idxh, ** idxr;
588  chi(idxh1[j], s+t-r, r, s+t-2*r, &npermW, &idxr, &idxh);
589  printf("npermW = %d\n",npermW);
590  for (int i=0; i<npermW; i++){
591  memcpy(idx_U,idxr[i],r);
592  memcpy(idx_U+r,idxh[i],s+t-2*r);
593  W[idx_C] += sgn1*sign(parity(idxr[i], idxh[i],idxh1[j], r,s+t-2*r))*U[idx_U];
594  }
595  }
596  }
597  assert(check_asym(C));
598  assert(check_asym(V));
599  assert(check_asym(W));
600 
601  C[idx_C] -= V[idx_C];
602  C[idx_C] -= W[idx_C];
603 
604  C[idx_C] -= C_ans[idx_C];
605 
606  double nrm = C.norm2();
607 
608  printf("error 2-norm is %.4E\n",nrm);
609 
610  int pass = (nrm <=1.E-3);
611 
612  if (rank == 0){
613  if (pass) printf("{ fast symmetric tensor contraction algorithm test } passed\n");
614  else printf("{ fast symmetric tensor contraction algorithm test } failed\n");
615  }
616  return pass;
617 }
618 
619 #ifndef TEST_SUITE
620 char* getCmdOption(char ** begin,
621  char ** end,
622  const std::string & option){
623  char ** itr = std::find(begin, end, option);
624  if (itr != end && ++itr != end){
625  return *itr;
626  }
627  return 0;
628 }
629 
630 
631 int main(int argc, char ** argv){
632  int rank, np, n, s,t,v;
633  int const in_num = argc;
634  char ** input_str = argv;
635 
636  MPI_Init(&argc, &argv);
637  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
638  MPI_Comm_size(MPI_COMM_WORLD, &np);
639 
640  if (getCmdOption(input_str, input_str+in_num, "-n")){
641  n = atoi(getCmdOption(input_str, input_str+in_num, "-n"));
642  if (n < 0) n = 4;
643  } else n = 4;
644 
645  if (getCmdOption(input_str, input_str+in_num, "-s")){
646  s = atoi(getCmdOption(input_str, input_str+in_num, "-s"));
647  if (s < 0) s = 1;
648  } else s = 1;
649 
650  if (getCmdOption(input_str, input_str+in_num, "-t")){
651  t = atoi(getCmdOption(input_str, input_str+in_num, "-t"));
652  if (t < 0) t = 1;
653  } else t = 1;
654 
655  if (getCmdOption(input_str, input_str+in_num, "-v")){
656  v = atoi(getCmdOption(input_str, input_str+in_num, "-v"));
657  if (v < 0) v = 1;
658  } else v = 1;
659 
660  {
661  World dw(MPI_COMM_WORLD, argc, argv);
662  if (rank == 0){
663  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);
664  }
665  int pass = fast_tensor_ctr(n, s, t, v, dw);
666  assert(pass);
667  }
668 
669  MPI_Finalize();
670  return 0;
671 }
677 #endif
678 
int * sym
symmetries among tensor dimensions
int64_t chchoose(int64_t n, int64_t k)
def rank(self)
Definition: core.pyx:312
int64_t fact(int64_t n)
double sign(int par)
Vector class which encapsulates a 1D tensor.
Definition: vector.h:14
Definition: common.h:37
an instance of the CTF library (world) on a MPI communicator
Definition: world.h:19
bool check_asym(Tensor<> &tsr)
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
int * lens
unpadded tensor edge lengths
def seed(seed)
Definition: random.pyx:15
bool check_sym(Tensor<> tsr)
void chi(char const *idx, int idx_len, int p_len, int q_len, int *npair, char ***idx_p, char ***idx_q)
int fast_tensor_ctr(int n, int s, int t, int v, World &ctf)
void get_rand_as_tsr(Tensor<> &tsr, int seed)
int parity(char const *a, char const *b, char const *c, int len_A, int len_B)
int main(int argc, char **argv)
char * getCmdOption(char **begin, char **end, const std::string &option)
Definition: apsp.cxx:17
int64_t choose(int64_t n, int64_t k)
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