Cyclops Tensor Framework
parallel arithmetic on multidimensional arrays
diag_ctr.cxx
Go to the documentation of this file.
1 /*Copyright (c) 2011, Edgar Solomonik, all rights reserved.*/
2 
9 #include <ctf.hpp>
10 
11 using namespace CTF;
12 
13 int diag_ctr(int n,
14  int m,
15  World & dw){
16  int rank, i, num_pes, pass;
17  int64_t np;
18  double * pairs;
19  int64_t * indices;
20 
21  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
22  MPI_Comm_size(MPI_COMM_WORLD, &num_pes);
23 
24 
25  int shapeN4[] = {NS,NS,NS,NS};
26  int sizeN4[] = {n,m,n,m};
27 
28  //* Creates distributed tensors initialized with zeros
29  Tensor<> A(4, sizeN4, shapeN4, dw);
30 
31  srand48(13*rank);
32 
33  Matrix<> mA(n,m,NS,dw);
34  Matrix<> mB(n,m,NS,dw);
35  A.get_local_data(&np, &indices, &pairs);
36  for (i=0; i<np; i++ ) pairs[i] = drand48()-.5; //(1.E-3)*sin(indices[i]);
37  A.write(np, indices, pairs);
38  delete [] pairs;
39  free(indices);
40  pass = 1;
41  double tr = 0.0;
42  tr += A["aiai"];
43  if (fabs(tr) < 1.E-10){
44  pass = 0;
45  }
46  mA["ai"] = A["aiai"];
47  tr -= mA["ai"];
48 
49  if (fabs(tr) > 1.E-10)
50  pass = 0;
51  if (pass){
52  if (rank == 0)
53  printf("{ sum(ai)A[\"aiai\"]=sum(ai)mA[\"ai\"] } passed \n");
54  } else {
55  if (rank == 0)
56  printf("{ sum(ai)A[\"aiai\"]=sum(ai)mA[\"ai\"] } failed \n");
57  }
58 
59 
60  return pass;
61 }
62 
63 
64 #ifndef TEST_SUITE
65 char* getCmdOption(char ** begin,
66  char ** end,
67  const std::string & option){
68  char ** itr = std::find(begin, end, option);
69  if (itr != end && ++itr != end){
70  return *itr;
71  }
72  return 0;
73 }
74 
75 
76 int main(int argc, char ** argv){
77  int rank, np, n, m;
78  int in_num = argc;
79  char ** input_str = argv;
80 
81  MPI_Init(&argc, &argv);
82  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
83  MPI_Comm_size(MPI_COMM_WORLD, &np);
84 
85  if (getCmdOption(input_str, input_str+in_num, "-n")){
86  n = atoi(getCmdOption(input_str, input_str+in_num, "-n"));
87  if (n < 0) n = 7;
88  } else n = 7;
89 
90  if (getCmdOption(input_str, input_str+in_num, "-m")){
91  m = atoi(getCmdOption(input_str, input_str+in_num, "-m"));
92  if (m < 0) m = 7;
93  } else m = 7;
94 
95  {
96  World dw(argc, argv);
97  diag_ctr(n, m, dw);
98  }
99 
100  MPI_Finalize();
101  return 0;
102 }
108 #endif
int diag_ctr(int n, int m, World &dw)
Definition: diag_ctr.cxx:13
Matrix class which encapsulates a 2D tensor.
Definition: matrix.h:18
def rank(self)
Definition: core.pyx:312
Definition: common.h:37
an instance of the CTF library (world) on a MPI communicator
Definition: world.h:19
string
Definition: core.pyx:456
int main(int argc, char **argv)
Definition: diag_ctr.cxx:76
void get_local_data(int64_t *npair, int64_t **global_idx, dtype **data, bool nonzeros_only=false, bool unpack_sym=false) const
Gives the global indices and values associated with the local data.
Definition: tensor.cxx:159
char * getCmdOption(char **begin, char **end, const std::string &option)
Definition: diag_ctr.cxx:65
Definition: apsp.cxx:17
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
def np(self)
Definition: core.pyx:315