Cyclops Tensor Framework
parallel arithmetic on multidimensional arrays
diag_sym.cxx
Go to the documentation of this file.
1 /*Copyright (c) 2011, Edgar Solomonik, all rights reserved.*/
2 
10 #include <ctf.hpp>
11 
12 using namespace CTF;
13 
14 int diag_sym(int n,
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[] = {SY,NS,SY,NS};
26  int sizeN4[] = {n,n,n,n};
27 
28  //* Creates distributed tensors initialized with zeros
29  Tensor<> A(4, sizeN4, shapeN4, dw);
30  Tensor<> B(4, sizeN4, shapeN4, dw);
31  Tensor<> C(4, sizeN4, shapeN4, dw);
32 
33  srand48(13*rank);
34 
35  Matrix<> mA(n,n,NS,dw);
36  Matrix<> mB(n,n,NS,dw);
37  mA.get_local_data(&np, &indices, &pairs);
38  for (i=0; i<np; i++ ) pairs[i] = drand48()-.5; //(1.E-3)*sin(indices[i]);
39  mA.write(np, indices, pairs);
40  delete [] pairs;
41  free(indices);
42  mB.get_local_data(&np, &indices, &pairs);
43  for (i=0; i<np; i++ ) pairs[i] = drand48()-.5; //(1.E-3)*sin(indices[i]);
44  mB.write(np, indices, pairs);
45  delete [] pairs;
46  free(indices);
47 
48  A["abij"] = mA["ii"];
49  B["abij"] = mA["jj"];
50  A["abij"] -= mB["aa"];
51  B["abij"] -= mB["bb"];
52  C["abij"] = A["abij"]-B["abij"];
53 
54  double norm = C.norm2();
55 
56  if (norm < 1.E-10){
57  pass = 1;
58  if (rank == 0)
59  printf("{ (A[\"(ab)(ij)\"]=mA[\"ii\"]-mB[\"aa\"]=mA[\"jj\"]-mB[\"bb\"] } passed \n");
60  } else {
61  pass = 0;
62  if (rank == 0)
63  printf("{ (A[\"(ab)(ij)\"]=mA[\"ii\"]-mB[\"aa\"]=mA[\"jj\"]-mB[\"bb\"] } failed \n");
64  }
65  return pass;
66 }
67 
68 
69 #ifndef TEST_SUITE
70 char* getCmdOption(char ** begin,
71  char ** end,
72  const std::string & option){
73  char ** itr = std::find(begin, end, option);
74  if (itr != end && ++itr != end){
75  return *itr;
76  }
77  return 0;
78 }
79 
80 
81 int main(int argc, char ** argv){
82  int rank, np, n;
83  int in_num = argc;
84  char ** input_str = argv;
85 
86  MPI_Init(&argc, &argv);
87  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
88  MPI_Comm_size(MPI_COMM_WORLD, &np);
89 
90  if (getCmdOption(input_str, input_str+in_num, "-n")){
91  n = atoi(getCmdOption(input_str, input_str+in_num, "-n"));
92  if (n < 0) n = 7;
93  } else n = 7;
94 
95 
96  {
97  World dw(argc, argv);
98  diag_sym(n, dw);
99  }
100 
101  MPI_Finalize();
102  return 0;
103 }
109 #endif
char * getCmdOption(char **begin, char **end, const std::string &option)
Definition: diag_sym.cxx:70
int main(int argc, char **argv)
Definition: diag_sym.cxx:81
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
dtype norm2()
computes the frobenius norm of the tensor (needs sqrt()!)
Definition: tensor.h:811
int diag_sym(int n, World &dw)
Definition: diag_sym.cxx:14
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
Definition: apsp.cxx:17
an instance of a tensor within a CTF world
Definition: tensor.h:74
Definition: common.h:37
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