OpenFPM_pdata  1.1.0
Project that contain the implementation of distributed structures
 All Data Structures Namespaces Functions Variables Typedefs Enumerations Friends Pages
Vector_unit_tests.hpp
1 /*
2  * Vector_unit_tests.hpp
3  *
4  * Created on: Apr 6, 2016
5  * Author: i-bird
6  */
7 
8 #ifndef OPENFPM_NUMERICS_SRC_VECTOR_VECTOR_UNIT_TESTS_HPP_
9 #define OPENFPM_NUMERICS_SRC_VECTOR_VECTOR_UNIT_TESTS_HPP_
10 
11 #include "Vector/Vector.hpp"
12 
13 BOOST_AUTO_TEST_SUITE( vector_test_suite )
14 
15 BOOST_AUTO_TEST_CASE(vector_eigen_parallel)
16 {
17  Vcluster & vcl = create_vcluster();
18 
19  if (vcl.getProcessingUnits() != 3)
20  return;
21 
22  // 3 Processors 9x9 Matrix to invert
23 
24  Vector<double> v(9);
25 
26  if (vcl.getProcessUnitID() == 0)
27  {
28  // v row1
29  v.insert(0,0);
30  v.insert(1,1);
31  v.insert(2,2);
32  }
33  else if (vcl.getProcessUnitID() == 1)
34  {
35  v.insert(3,3);
36  v.insert(4,4);
37  v.insert(5,5);
38  }
39  else if (vcl.getProcessUnitID() == 2)
40  {
41  v.insert(6,6);
42  v.insert(7,7);
43  v.insert(8,8);
44  }
45 
46  Vector<double> v3;
47  v3 = v;
48 
49  // force to sync
50  v.getVec();
51 
52  Vector<double> v2;
53  v2 = v;
54 
55  // Master has the full vector
56 
57  if (vcl.getProcessUnitID() == 0)
58  {
59  BOOST_REQUIRE_EQUAL(v(0),0);
60  BOOST_REQUIRE_EQUAL(v(1),1);
61  BOOST_REQUIRE_EQUAL(v(2),2);
62 
63  BOOST_REQUIRE_EQUAL(v(3),3);
64  BOOST_REQUIRE_EQUAL(v(4),4);
65  BOOST_REQUIRE_EQUAL(v(5),5);
66 
67  BOOST_REQUIRE_EQUAL(v(6),6);
68  BOOST_REQUIRE_EQUAL(v(7),7);
69  BOOST_REQUIRE_EQUAL(v(8),8);
70 
71  // Change the vector on Master
72 
73  v(0) = 8;
74  v(1) = 7;
75  v(2) = 6;
76  v(3) = 5;
77  v(4) = 4;
78  v(5) = 3;
79  v(6) = 2;
80  v(7) = 1;
81  v(8) = 0;
82  }
83 
84  v.scatter();
85  v2.scatter();
86  v3.scatter();
87 
88  if (vcl.getProcessUnitID() == 0)
89  {
90  BOOST_REQUIRE_EQUAL(v(0),8);
91  BOOST_REQUIRE_EQUAL(v(1),7);
92  BOOST_REQUIRE_EQUAL(v(2),6);
93 
94  BOOST_REQUIRE_EQUAL(v2(0),0);
95  BOOST_REQUIRE_EQUAL(v2(1),1);
96  BOOST_REQUIRE_EQUAL(v2(2),2);
97 
98  BOOST_REQUIRE_EQUAL(v3(0),0);
99  BOOST_REQUIRE_EQUAL(v3(1),1);
100  BOOST_REQUIRE_EQUAL(v3(2),2);
101  }
102  else if (vcl.getProcessUnitID() == 1)
103  {
104  BOOST_REQUIRE_EQUAL(v(3),5);
105  BOOST_REQUIRE_EQUAL(v(4),4);
106  BOOST_REQUIRE_EQUAL(v(5),3);
107 
108  BOOST_REQUIRE_EQUAL(v2(3),3);
109  BOOST_REQUIRE_EQUAL(v2(4),4);
110  BOOST_REQUIRE_EQUAL(v2(5),5);
111 
112  BOOST_REQUIRE_EQUAL(v3(3),3);
113  BOOST_REQUIRE_EQUAL(v3(4),4);
114  BOOST_REQUIRE_EQUAL(v3(5),5);
115  }
116  else if (vcl.getProcessUnitID() == 2)
117  {
118  BOOST_REQUIRE_EQUAL(v(6),2);
119  BOOST_REQUIRE_EQUAL(v(7),1);
120  BOOST_REQUIRE_EQUAL(v(8),0);
121 
122  BOOST_REQUIRE_EQUAL(v2(6),6);
123  BOOST_REQUIRE_EQUAL(v2(7),7);
124  BOOST_REQUIRE_EQUAL(v2(8),8);
125 
126  BOOST_REQUIRE_EQUAL(v3(6),6);
127  BOOST_REQUIRE_EQUAL(v3(7),7);
128  BOOST_REQUIRE_EQUAL(v3(8),8);
129  }
130 }
131 
132 #ifdef HAVE_PETSC
133 
134 BOOST_AUTO_TEST_CASE(vector_petsc_parallel)
135 {
136  Vcluster & vcl = create_vcluster();
137 
138  if (vcl.getProcessingUnits() != 3)
139  return;
140 
141  // 3 Processors 9x9 Matrix to invert
142 
144 
145  if (vcl.getProcessUnitID() == 0)
146  {
147  // v row1
148  v.insert(0,0);
149  v.insert(1,1);
150  v.insert(2,2);
151  }
152  else if (vcl.getProcessUnitID() == 1)
153  {
154  v.insert(3,3);
155  v.insert(4,4);
156  v.insert(5,5);
157  }
158  else if (vcl.getProcessUnitID() == 2)
159  {
160  v.insert(6,6);
161  v.insert(7,7);
162  v.insert(8,8);
163  }
164 
166  v3 = v;
167 
168  if (vcl.getProcessUnitID() == 0)
169  {
170  // v row1
171  v(0) = 8;
172  v(1) = 7;
173  v(2) = 6;
174  }
175  else if (vcl.getProcessUnitID() == 1)
176  {
177  v(3) = 5;
178  v(4) = 4;
179  v(5) = 3;
180  }
181  else if (vcl.getProcessUnitID() == 2)
182  {
183  v(6) = 2;
184  v(7) = 1;
185  v(8) = 0;
186  }
187 
188  // Master has the full vector
189 
190  if (vcl.getProcessUnitID() == 0)
191  {
192  BOOST_REQUIRE_EQUAL(v(0),8);
193  BOOST_REQUIRE_EQUAL(v(1),7);
194  BOOST_REQUIRE_EQUAL(v(2),6);
195 
196  BOOST_REQUIRE_EQUAL(v3(0),0);
197  BOOST_REQUIRE_EQUAL(v3(1),1);
198  BOOST_REQUIRE_EQUAL(v3(2),2);
199  }
200  else if (vcl.getProcessUnitID() == 1)
201  {
202  BOOST_REQUIRE_EQUAL(v(3),5);
203  BOOST_REQUIRE_EQUAL(v(4),4);
204  BOOST_REQUIRE_EQUAL(v(5),3);
205 
206  BOOST_REQUIRE_EQUAL(v3(3),3);
207  BOOST_REQUIRE_EQUAL(v3(4),4);
208  BOOST_REQUIRE_EQUAL(v3(5),5);
209  }
210  else if (vcl.getProcessUnitID() == 2)
211  {
212  BOOST_REQUIRE_EQUAL(v(6),2);
213  BOOST_REQUIRE_EQUAL(v(7),1);
214  BOOST_REQUIRE_EQUAL(v(8),0);
215 
216  BOOST_REQUIRE_EQUAL(v3(6),6);
217  BOOST_REQUIRE_EQUAL(v3(7),7);
218  BOOST_REQUIRE_EQUAL(v3(8),8);
219  }
220 
221  // Here we get the petsc vector
222  auto & vp = v.getVec();
223 
224  // We check the correctness of the PETSC vector
225 
226  if (vcl.getProcessUnitID() == 0)
227  {
228  PetscInt ix[] = {0,1,2};
229  PetscScalar y[3];
230 
231  VecGetValues(vp,3,ix,y);
232 
233  BOOST_REQUIRE_EQUAL(y[0],8);
234  BOOST_REQUIRE_EQUAL(y[1],7);
235  BOOST_REQUIRE_EQUAL(y[2],6);
236  }
237  else if (vcl.getProcessUnitID() == 1)
238  {
239  PetscInt ix[] = {3,4,5};
240  PetscScalar y[3];
241 
242  VecGetValues(vp,3,ix,y);
243 
244  BOOST_REQUIRE_EQUAL(y[0],5);
245  BOOST_REQUIRE_EQUAL(y[1],4);
246  BOOST_REQUIRE_EQUAL(y[2],3);
247  }
248  else if (vcl.getProcessUnitID() == 2)
249  {
250  PetscInt ix[] = {6,7,8};
251  PetscScalar y[3];
252 
253  VecGetValues(vp,3,ix,y);
254 
255  BOOST_REQUIRE_EQUAL(y[0],2);
256  BOOST_REQUIRE_EQUAL(y[1],1);
257  BOOST_REQUIRE_EQUAL(y[2],0);
258  }
259 
260 }
261 
262 #endif
263 
264 BOOST_AUTO_TEST_SUITE_END()
265 
266 
267 #endif /* OPENFPM_NUMERICS_SRC_VECTOR_VECTOR_UNIT_TESTS_HPP_ */
size_t getProcessUnitID()
Get the process unit id.
Sparse Matrix implementation stub object when OpenFPM is compiled with no linear algebra support...
Definition: Vector.hpp:39
Implementation of VCluster class.
Definition: VCluster.hpp:36
void scatter()
scatter
Definition: Vector.hpp:144
size_t getProcessingUnits()
Get the total number of processors.