-
Notifications
You must be signed in to change notification settings - Fork 54
/
Copy pathbivar_function.cxx
116 lines (87 loc) · 2.42 KB
/
bivar_function.cxx
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
/*Copyright (c) 2011, Edgar Solomonik, all rights reserved.*/
/** \addtogroup tests
* @{
* \defgroup bivar_function bivar_function
* @{
* \brief tests custom element-wise functions by implementing division elementwise on 4D tensors
*/
#include <ctf.hpp>
using namespace CTF;
double f2(double a, double b){
return a*b+b*a;
}
int bivar_function(int n,
World & dw){
int shapeN4[] = {NS,NS,NS,NS};
int sizeN4[] = {n+1,n,n+2,n+3};
Tensor<> A(4, sizeN4, shapeN4, dw);
Tensor<> B(4, sizeN4, shapeN4, dw);
srand48(dw.rank);
A.fill_random(-.5, .5);
B.fill_random(-.5, .5);
double * all_start_data_A;
int64_t nall_A;
A.read_all(&nall_A, &all_start_data_A);
double * all_start_data_B;
int64_t nall_B;
B.read_all(&nall_B, &all_start_data_B);
CTF::Function<> bfun([](double a, double b){ return a*b + b*a; });
.5*A["ijkl"]+=bfun(A["ijkl"],B["ijkl"]);
double * all_end_data_A;
int64_t nall2_A;
A.read_all(&nall2_A, &all_end_data_A);
int pass = (nall_A == nall2_A);
if (pass){
for (int64_t i=0; i<nall_A; i++){
if (fabs(.5*all_start_data_A[i]+f2(all_start_data_A[i],all_start_data_B[i])-all_end_data_A[i])>=1.E-6) pass =0;
}
}
MPI_Allreduce(MPI_IN_PLACE, &pass, 1, MPI_INT, MPI_MIN, MPI_COMM_WORLD);
if (dw.rank == 0){
if (pass){
printf("{ A[\"ijkl\"] = f2(A[\"ijkl\"], B[\"ijkl\"]) } passed\n");
} else {
printf("{ A[\"ijkl\"] = f2(A[\"ijkl\"], B[\"ijkl\"]) } failed\n");
}
}
delete [] all_start_data_A;
delete [] all_end_data_A;
delete [] all_start_data_B;
return pass;
}
#ifndef TEST_SUITE
char* getCmdOption(char ** begin,
char ** end,
const std::string & option){
char ** itr = std::find(begin, end, option);
if (itr != end && ++itr != end){
return *itr;
}
return 0;
}
int main(int argc, char ** argv){
int rank, np, n;
int const in_num = argc;
char ** input_str = argv;
MPI_Init(&argc, &argv);
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
MPI_Comm_size(MPI_COMM_WORLD, &np);
if (getCmdOption(input_str, input_str+in_num, "-n")){
n = atoi(getCmdOption(input_str, input_str+in_num, "-n"));
if (n < 0) n = 5;
} else n = 5;
{
World dw(MPI_COMM_WORLD, argc, argv);
if (rank == 0){
printf("Computing bivar_function A_ijkl = f(B_ijkl, A_ijkl)\n");
}
bivar_function(n, dw);
}
MPI_Finalize();
return 0;
}
/**
* @}
* @}
*/
#endif