Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_SacadoMPVectorCommTests.cpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Stokhos Package
5 // Copyright (2009) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Eric T. Phipps (etphipp@sandia.gov).
38 //
39 // ***********************************************************************
40 // @HEADER
45 #include "Teuchos_CommHelpers.hpp"
46 #include "Teuchos_DefaultComm.hpp"
47 #include "Teuchos_Array.hpp"
48 #include "Teuchos_Comm.hpp"
49 
50 #include "Sacado.hpp"
51 #include "Stokhos_Sacado.hpp"
53 #include "Sacado_Fad_DFad.hpp"
54 #include "Sacado_mpl_apply.hpp"
55 #include "Sacado_Random.hpp"
56 
57 #include <Kokkos_Core.hpp>
58 
59 //
60 // Currently this doesn't test:
61 // * the device
62 // * threaded storage (needs the device)
63 // * strided storage with non-trivial stride
64 //
65 
66 using Teuchos::RCP;
67 using Teuchos::rcp;
68 
69 // Common setup for unit tests
70 template <typename VecType, typename FadType>
71 struct UnitTestSetup {
72  int sz;
73 
76 
77  typedef typename Sacado::mpl::apply<FadType,VecType>::type FadVecType;
80 
82  sz = 8;
83 
84  // Serializers
86  rcp(new VecSerializerT(
89  }
90 };
91 
92 template <typename VecType>
94  const Teuchos::Array<VecType>& x2,
95  const std::string& tag,
96  Teuchos::FancyOStream& out) {
97 
98  // Check sizes match
99  bool success = (x.size() == x2.size());
100  out << tag << " Vec array size test";
101  if (success)
102  out << " passed";
103  else
104  out << " failed";
105  out << ": \n\tExpected: " << x.size() << ", \n\tGot: " << x2.size()
106  << "." << std::endl;
107 
108  // Check Fads match
109  for (int i=0; i<x.size(); i++) {
110  bool success2 = Sacado::IsEqual<VecType>::eval(x[i], x2[i]);
111  out << tag << " Vec array comparison test " << i;
112  if (success2)
113  out << " passed";
114  else
115  out << " failed";
116  out << ": \n\tExpected: " << x[i] << ", \n\tGot: " << x2[i] << "."
117  << std::endl;
118  success = success && success2;
119  }
120 
121  return success;
122 }
123 
124 template<typename Ordinal>
126  const Teuchos::Comm<Ordinal> &comm,
128  const bool result
129  )
130 {
131  out << "\nChecking that the above test passed in all processes ...";
132  int thisResult = ( result ? 1 : 0 );
133  int sumResult = -1;
134  Teuchos::reduceAll(comm,Teuchos::REDUCE_SUM,Ordinal(1),&thisResult,
135  &sumResult);
136  const bool passed = sumResult==Teuchos::size(comm);
137  if(passed)
138  out << " passed\n";
139  else
140  out << " (sumResult="<<sumResult<<"!=numProcs="<<Teuchos::size(comm)<<") failed\n";
141  return passed;
142 }
143 
144 #define VEC_COMM_TESTS(VecType, FadType, Vec, FAD) \
145 TEUCHOS_UNIT_TEST( Vec##_Comm, Vec_Broadcast ) { \
146  Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
147  comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
148  \
149  int n = 7; \
150  Teuchos::Array<VecType> x(n), x2(n); \
151  for (int i=0; i<n; i++) { \
152  x[i] = VecType(setup.sz, 0.0); \
153  for (int j=0; j<setup.sz; j++) \
154  x[i].fastAccessCoeff(j) = rnd.number(); \
155  } \
156  if (comm->getRank() == 0) \
157  x2 = x; \
158  Teuchos::broadcast(*comm, *setup.vec_serializer, 0, n, &x2[0]); \
159  success = checkVecArrays(x, x2, std::string(#Vec)+" Broadcast", out); \
160  success = checkResultOnAllProcs(*comm, out, success); \
161 } \
162  \
163 TEUCHOS_UNIT_TEST( Vec##_Comm, Vec_GatherAll ) { \
164  Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
165  comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
166  \
167  int n = 7; \
168  int size = comm->getSize(); \
169  int rank = comm->getRank(); \
170  int N = n*size; \
171  Teuchos::Array<VecType> x(n), x2(N), x3(N); \
172  for (int i=0; i<n; i++) { \
173  x[i] = VecType(setup.sz, 0.0); \
174  for (int j=0; j<setup.sz; j++) \
175  x[i].fastAccessCoeff(j) = (rank+1)*(i+1)*(j+1); \
176  } \
177  for (int j=0; j<size; j++) { \
178  for (int i=0; i<n; i++) { \
179  x3[n*j+i] = VecType(setup.sz, 0.0); \
180  for (int k=0; k<setup.sz; k++) \
181  x3[n*j+i].fastAccessCoeff(k) = (j+1)*(i+1)*(k+1); \
182  } \
183  } \
184  Teuchos::gatherAll(*comm, *setup.vec_serializer, \
185  n, &x[0], N, &x2[0]); \
186  success = checkVecArrays(x3, x2, std::string(#Vec)+" Gather All", out); \
187  success = checkResultOnAllProcs(*comm, out, success); \
188 } \
189  \
190 TEUCHOS_UNIT_TEST( Vec##_Comm, Vec_SumAll ) { \
191  Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
192  comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
193  \
194  int n = 7; \
195  int num_proc = comm->getSize(); \
196  \
197  Teuchos::Array<VecType> x(n), sums(n), sums2(n); \
198  for (int i=0; i<n; i++) { \
199  x[i] = VecType(setup.sz, 0.0); \
200  for (int j=0; j<setup.sz; j++) \
201  x[i].fastAccessCoeff(j) = 2.0*(i+1); \
202  } \
203  for (int i=0; i<n; i++) { \
204  sums[i] = VecType(setup.sz, 0.0); \
205  for (int j=0; j<setup.sz; j++) \
206  sums[i].fastAccessCoeff(j) = 2.0*(i+1)*num_proc; \
207  } \
208  Teuchos::reduceAll(*comm, *setup.vec_serializer, \
209  Teuchos::REDUCE_SUM, n, &x[0], &sums2[0]); \
210  success = checkVecArrays(sums, sums2, \
211  std::string(#Vec)+" Sum All", out); \
212  success = checkResultOnAllProcs(*comm, out, success); \
213 } \
214  \
215 TEUCHOS_UNIT_TEST( Vec##_Comm, Vec_MaxAll ) { \
216  Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
217  comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
218  \
219  int n = 7; \
220  int rank = comm->getRank(); \
221  int num_proc = comm->getSize(); \
222  \
223  Teuchos::Array<VecType> x(n), maxs(n), maxs2(n); \
224  for (int i=0; i<n; i++) { \
225  x[i] = VecType(setup.sz, 0.0); \
226  for (int j=0; j<setup.sz; j++) \
227  x[i].fastAccessCoeff(j) = 2.0*(i+1)*(rank+1); \
228  } \
229  for (int i=0; i<n; i++) { \
230  maxs[i] = VecType(setup.sz, 0.0); \
231  for (int j=0; j<setup.sz; j++) \
232  maxs[i].fastAccessCoeff(j) = 2.0*(i+1)*num_proc; \
233  } \
234  Teuchos::reduceAll(*comm, *setup.vec_serializer, \
235  Teuchos::REDUCE_MAX, n, &x[0], &maxs2[0]); \
236  success = checkVecArrays(maxs, maxs2, \
237  std::string(#Vec)+" Max All", out); \
238  success = checkResultOnAllProcs(*comm, out, success); \
239 } \
240  \
241 TEUCHOS_UNIT_TEST( Vec##_Comm, Vec_MinAll ) { \
242  Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
243  comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
244  \
245  int n = 7; \
246  int rank = comm->getRank(); \
247  \
248  Teuchos::Array<VecType> x(n), mins(n), mins2(n); \
249  for (int i=0; i<n; i++) { \
250  x[i] = VecType(setup.sz, 0.0); \
251  for (int j=0; j<setup.sz; j++) \
252  x[i].fastAccessCoeff(j) = 2.0*(i+1)*(rank+1); \
253  } \
254  for (int i=0; i<n; i++) { \
255  mins[i] = VecType(setup.sz, 0.0); \
256  for (int j=0; j<setup.sz; j++) \
257  mins[i].fastAccessCoeff(j) = 2.0*(i+1); \
258  } \
259  Teuchos::reduceAll(*comm, *setup.vec_serializer, \
260  Teuchos::REDUCE_MIN, n, &x[0], &mins2[0]); \
261  success = checkVecArrays(mins, mins2, \
262  std::string(#Vec)+" Min All", out); \
263  success = checkResultOnAllProcs(*comm, out, success); \
264 } \
265  \
266 TEUCHOS_UNIT_TEST( Vec##_Comm, Vec_ScanSum ) { \
267  Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
268  comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
269  \
270  int n = 7; \
271  int rank = comm->getRank(); \
272  \
273  Teuchos::Array<VecType> x(n), sums(n), sums2(n); \
274  for (int i=0; i<n; i++) { \
275  x[i] = VecType(setup.sz, 0.0); \
276  for (int j=0; j<setup.sz; j++) \
277  x[i].fastAccessCoeff(j) = 2.0*(i+1); \
278  } \
279  for (int i=0; i<n; i++) { \
280  sums[i] = VecType(setup.sz, 0.0); \
281  for (int j=0; j<setup.sz; j++) \
282  sums[i].fastAccessCoeff(j) = 2.0*(i+1)*(rank+1); \
283  } \
284  Teuchos::scan(*comm, *setup.vec_serializer, \
285  Teuchos::REDUCE_SUM, n, &x[0], &sums2[0]); \
286  success = checkVecArrays(sums, sums2, \
287  std::string(#Vec)+" Scan Sum", out); \
288  success = checkResultOnAllProcs(*comm, out, success); \
289 } \
290  \
291 TEUCHOS_UNIT_TEST( Vec##_Comm, Vec_ScanMax ) { \
292  Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
293  comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
294  \
295  int n = 7; \
296  int rank = comm->getRank(); \
297  \
298  Teuchos::Array<VecType> x(n), maxs(n), maxs2(n); \
299  for (int i=0; i<n; i++) { \
300  x[i] = VecType(setup.sz, 0.0); \
301  for (int j=0; j<setup.sz; j++) \
302  x[i].fastAccessCoeff(j) = 2.0*(i+1)*(rank+1); \
303  } \
304  for (int i=0; i<n; i++) { \
305  maxs[i] = VecType(setup.sz, 0.0); \
306  for (int j=0; j<setup.sz; j++) \
307  maxs[i].fastAccessCoeff(j) = 2.0*(i+1)*(rank+1); \
308  } \
309  Teuchos::scan(*comm, *setup.vec_serializer, \
310  Teuchos::REDUCE_MAX, n, &x[0], &maxs2[0]); \
311  success = checkVecArrays(maxs, maxs2, \
312  std::string(#Vec)+" Scan Max", out); \
313  success = checkResultOnAllProcs(*comm, out, success); \
314 } \
315  \
316 TEUCHOS_UNIT_TEST( Vec##_Comm, Vec_ScanMin ) { \
317  Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
318  comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
319  \
320  int n = 7; \
321  int rank = comm->getRank(); \
322  \
323  Teuchos::Array<VecType> x(n), mins(n), mins2(n); \
324  for (int i=0; i<n; i++) { \
325  x[i] = VecType(setup.sz, 0.0); \
326  for (int j=0; j<setup.sz; j++) \
327  x[i].fastAccessCoeff(j) = 2.0*(i+1)*(rank+1); \
328  } \
329  for (int i=0; i<n; i++) { \
330  mins[i] = VecType(setup.sz, 0.0); \
331  for (int j=0; j<setup.sz; j++) \
332  mins[i].fastAccessCoeff(j) = 2.0*(i+1); \
333  } \
334  Teuchos::scan(*comm, *setup.vec_serializer, \
335  Teuchos::REDUCE_MIN, n, &x[0], &mins2[0]); \
336  success = checkVecArrays(mins, mins2, \
337  std::string(#Vec)+" Scan Min", out); \
338  success = checkResultOnAllProcs(*comm, out, success); \
339 } \
340  \
341 TEUCHOS_UNIT_TEST( Vec##_Comm, Vec_SendReceive ) { \
342  Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
343  comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
344  \
345  int num_proc = comm->getSize(); \
346  if (num_proc > 1) { \
347  int rank = comm->getRank(); \
348  int n = 7; \
349  Teuchos::Array<VecType> x(n), x2(n); \
350  for (int i=0; i<n; i++) { \
351  x[i] = VecType(setup.sz, 0.0); \
352  for (int j=0; j<setup.sz; j++) \
353  x[i].fastAccessCoeff(j) = 2.0*(i+1)*(j+1); \
354  } \
355  if (rank != 1) \
356  x2 = x; \
357  if (rank == 0) Teuchos::send(*comm, *setup.vec_serializer, \
358  n, &x[0], 1); \
359  if (rank == 1) Teuchos::receive(*comm, *setup.vec_serializer, \
360  0, n, &x2[0]); \
361  success = checkVecArrays(x, x2, \
362  std::string(#Vec)+" Send/Receive", out); \
363  success = checkResultOnAllProcs(*comm, out, success); \
364  } \
365  else \
366  success = true; \
367 } \
368  \
369 TEUCHOS_UNIT_TEST( Vec##_Comm, FadVec_Broadcast ) { \
370  typedef Sacado::mpl::apply<FadType,VecType>::type FadVecType; \
371  Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
372  comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
373  \
374  int n = 7; \
375  int p = 5; \
376  Teuchos::Array<FadVecType> x(n), x2(n); \
377  for (int i=0; i<n; i++) { \
378  VecType f(setup.sz, 0.0); \
379  for (int k=0; k<setup.sz; k++) \
380  f.fastAccessCoeff(k) = rnd.number(); \
381  x[i] = FadVecType(p, f); \
382  for (int j=0; j<p; j++) { \
383  VecType g(setup.sz, 0.0); \
384  for (int k=0; k<setup.sz; k++) \
385  g.fastAccessCoeff(k) = rnd.number(); \
386  x[i].fastAccessDx(j) = g; \
387  } \
388  } \
389  if (comm->getRank() == 0) \
390  x2 = x; \
391  Teuchos::broadcast(*comm, *setup.fad_vec_serializer, 0, n, &x2[0]); \
392  success = checkVecArrays(x, x2, \
393  std::string(#FAD)+"<"+#Vec+"> Broadcast", out); \
394  success = checkResultOnAllProcs(*comm, out, success); \
395 } \
396  \
397 TEUCHOS_UNIT_TEST( Vec##_Comm, FadVec_GatherAll ) { \
398  typedef Sacado::mpl::apply<FadType,VecType>::type FadVecType; \
399  Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
400  comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
401  \
402  int n = 7; \
403  int p = 5; \
404  int size = comm->getSize(); \
405  int rank = comm->getRank(); \
406  int N = n*size; \
407  Teuchos::Array<FadVecType> x(n), x2(N), x3(N); \
408  for (int i=0; i<n; i++) { \
409  VecType f(setup.sz, 0.0); \
410  for (int k=0; k<setup.sz; k++) \
411  f.fastAccessCoeff(k) = (rank+1)*(i+1)*(k+1); \
412  x[i] = FadVecType(p, f); \
413  for (int j=0; j<p; j++) { \
414  x[i].fastAccessDx(j) = f; \
415  } \
416  } \
417  for (int j=0; j<size; j++) { \
418  for (int i=0; i<n; i++) { \
419  VecType f(setup.sz, 0.0); \
420  for (int k=0; k<setup.sz; k++) \
421  f.fastAccessCoeff(k) = (j+1)*(i+1)*(k+1); \
422  x3[n*j+i] = FadVecType(p, f); \
423  for (int k=0; k<p; k++) \
424  x3[n*j+i].fastAccessDx(k) = f; \
425  } \
426  } \
427  Teuchos::gatherAll(*comm, *setup.fad_vec_serializer, \
428  n, &x[0], N, &x2[0]); \
429  success = checkVecArrays(x3, x2, \
430  std::string(#FAD)+"<"+#Vec+"> Gather All", out); \
431  success = checkResultOnAllProcs(*comm, out, success); \
432 } \
433  \
434 TEUCHOS_UNIT_TEST( Vec##_Comm, FadVec_SumAll ) { \
435  typedef Sacado::mpl::apply<FadType,VecType>::type FadVecType; \
436  Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
437  comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
438  \
439  int n = 7; \
440  int p = 5; \
441  int num_proc = comm->getSize(); \
442  \
443  Teuchos::Array<FadVecType> x(n), sums(n), sums2(n); \
444  for (int i=0; i<n; i++) { \
445  VecType f(setup.sz, 0.0); \
446  for (int k=0; k<setup.sz; k++) \
447  f.fastAccessCoeff(k) = 2.0*(i+1); \
448  x[i] = FadVecType(p, f); \
449  for (int j=0; j<p; j++) { \
450  VecType g(setup.sz, 0.0); \
451  for (int k=0; k<setup.sz; k++) \
452  g.fastAccessCoeff(k) = 2.0*(i+1); \
453  x[i].fastAccessDx(j) = g; \
454  } \
455  } \
456  for (int i=0; i<n; i++) { \
457  VecType f(setup.sz, 0.0); \
458  for (int k=0; k<setup.sz; k++) \
459  f.fastAccessCoeff(k) = 2.0*(i+1)*num_proc; \
460  sums[i] = FadVecType(p, f); \
461  for (int j=0; j<p; j++) { \
462  VecType g(setup.sz, 0.0); \
463  for (int k=0; k<setup.sz; k++) \
464  g.fastAccessCoeff(k) = 2.0*(i+1)*num_proc; \
465  sums[i].fastAccessDx(j) = g; \
466  } \
467  } \
468  Teuchos::reduceAll(*comm, *setup.fad_vec_serializer, \
469  Teuchos::REDUCE_SUM, n, &x[0], &sums2[0]); \
470  success = checkVecArrays(sums, sums2, \
471  std::string(#FAD)+"<"+#Vec+"> Sum All", out); \
472  success = checkResultOnAllProcs(*comm, out, success); \
473 } \
474  \
475 TEUCHOS_UNIT_TEST( Vec##_Comm, FadVec_MaxAll ) { \
476  typedef Sacado::mpl::apply<FadType,VecType>::type FadVecType; \
477  Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
478  comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
479  \
480  int n = 8; \
481  int p = 5; \
482  int rank = comm->getRank(); \
483  int num_proc = comm->getSize(); \
484  \
485  Teuchos::Array<FadVecType> x(n), maxs(n), maxs2(n); \
486  for (int i=0; i<n; i++) { \
487  VecType f(setup.sz, 0.0); \
488  for (int k=0; k<setup.sz; k++) \
489  f.fastAccessCoeff(k) = 2.0*(i+1)*(rank+1); \
490  x[i] = FadVecType(p, f); \
491  for (int j=0; j<p; j++) { \
492  x[i].fastAccessDx(j) = f; \
493  } \
494  } \
495  for (int i=0; i<n; i++) { \
496  VecType f(setup.sz, 0.0); \
497  for (int k=0; k<setup.sz; k++) \
498  f.fastAccessCoeff(k) = 2.0*(i+1)*num_proc; \
499  maxs[i] = FadVecType(p, f); \
500  for (int j=0; j<p; j++) \
501  maxs[i].fastAccessDx(j) = f; \
502  } \
503  Teuchos::reduceAll(*comm, *setup.fad_vec_serializer, \
504  Teuchos::REDUCE_MAX, n, &x[0], &maxs2[0]); \
505  success = checkVecArrays(maxs, maxs2, \
506  std::string(#FAD)+"<"+#Vec+"> Max All", out); \
507  success = checkResultOnAllProcs(*comm, out, success); \
508 } \
509  \
510 TEUCHOS_UNIT_TEST( Vec##_Comm, FadVec_MinAll ) { \
511  typedef Sacado::mpl::apply<FadType,VecType>::type FadVecType; \
512  Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
513  comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
514  \
515  int n = 8; \
516  int p = 5; \
517  int rank = comm->getRank(); \
518  \
519  Teuchos::Array<FadVecType> x(n), mins(n), mins2(n); \
520  for (int i=0; i<n; i++) { \
521  VecType f(setup.sz, 0.0); \
522  for (int k=0; k<setup.sz; k++) \
523  f.fastAccessCoeff(k) = 2.0*(i+1)*(rank+1); \
524  x[i] = FadVecType(p, f); \
525  for (int j=0; j<p; j++) { \
526  x[i].fastAccessDx(j) = f; \
527  } \
528  } \
529  for (int i=0; i<n; i++) { \
530  VecType f(setup.sz, 0.0); \
531  for (int k=0; k<setup.sz; k++) \
532  f.fastAccessCoeff(k) = 2.0*(i+1); \
533  mins[i] = FadVecType(p, f); \
534  for (int j=0; j<p; j++) \
535  mins[i].fastAccessDx(j) = f; \
536  } \
537  Teuchos::reduceAll(*comm, *setup.fad_vec_serializer, \
538  Teuchos::REDUCE_MIN, n, &x[0], &mins2[0]); \
539  success = checkVecArrays(mins, mins2, \
540  std::string(#FAD)+"<"+#Vec+"> Min All", out); \
541  success = checkResultOnAllProcs(*comm, out, success); \
542 } \
543  \
544 TEUCHOS_UNIT_TEST( Vec##_Comm, FadVec_ScanSum ) { \
545  typedef Sacado::mpl::apply<FadType,VecType>::type FadVecType; \
546  Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
547  comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
548  \
549  int n = 7; \
550  int p = 5; \
551  int rank = comm->getRank(); \
552  \
553  Teuchos::Array<FadVecType> x(n), sums(n), sums2(n); \
554  for (int i=0; i<n; i++) { \
555  VecType f(setup.sz, 0.0); \
556  for (int k=0; k<setup.sz; k++) \
557  f.fastAccessCoeff(k) = 2.0*(i+1); \
558  x[i] = FadVecType(p, f); \
559  for (int j=0; j<p; j++) { \
560  x[i].fastAccessDx(j) = f; \
561  } \
562  } \
563  for (int i=0; i<n; i++) { \
564  VecType f(setup.sz, 0.0); \
565  for (int k=0; k<setup.sz; k++) \
566  f.fastAccessCoeff(k) = 2.0*(i+1)*(rank+1); \
567  sums[i] = FadVecType(p, f); \
568  for (int j=0; j<p; j++) \
569  sums[i].fastAccessDx(j) = f; \
570  } \
571  Teuchos::scan(*comm, *setup.fad_vec_serializer, \
572  Teuchos::REDUCE_SUM, n, &x[0], &sums2[0]); \
573  success = checkVecArrays(sums, sums2, \
574  std::string(#FAD)+"<"+#Vec+"> Scan Sum", out); \
575  success = checkResultOnAllProcs(*comm, out, success); \
576 } \
577  \
578 TEUCHOS_UNIT_TEST( Vec##_Comm, FadVec_ScanMax ) { \
579  typedef Sacado::mpl::apply<FadType,VecType>::type FadVecType; \
580  Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
581  comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
582  \
583  int n = 7; \
584  int p = 5; \
585  int rank = comm->getRank(); \
586  \
587  Teuchos::Array<FadVecType> x(n), maxs(n), maxs2(n); \
588  for (int i=0; i<n; i++) { \
589  VecType f(setup.sz, 0.0); \
590  for (int k=0; k<setup.sz; k++) \
591  f.fastAccessCoeff(k) = 2.0*(i+1)*(rank+1); \
592  x[i] = FadVecType(p, f); \
593  for (int j=0; j<p; j++) { \
594  x[i].fastAccessDx(j) = f; \
595  } \
596  } \
597  for (int i=0; i<n; i++) { \
598  VecType f(setup.sz, 0.0); \
599  for (int k=0; k<setup.sz; k++) \
600  f.fastAccessCoeff(k) = 2.0*(i+1)*(rank+1); \
601  maxs[i] = FadVecType(p, f); \
602  for (int j=0; j<p; j++) \
603  maxs[i].fastAccessDx(j) = f; \
604  } \
605  Teuchos::scan(*comm, *setup.fad_vec_serializer, \
606  Teuchos::REDUCE_MAX, n, &x[0], &maxs2[0]); \
607  success = checkVecArrays(maxs, maxs2, \
608  std::string(#FAD)+"<"+#Vec+"> Scan Max", out); \
609  success = checkResultOnAllProcs(*comm, out, success); \
610 } \
611  \
612 TEUCHOS_UNIT_TEST( Vec##_Comm, FadVec_ScanMin ) { \
613  typedef Sacado::mpl::apply<FadType,VecType>::type FadVecType; \
614  Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
615  comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
616  \
617  int n = 7; \
618  int p = 5; \
619  int rank = comm->getRank(); \
620  \
621  Teuchos::Array<FadVecType> x(n), mins(n), mins2(n); \
622  for (int i=0; i<n; i++) { \
623  VecType f(setup.sz, 0.0); \
624  for (int k=0; k<setup.sz; k++) \
625  f.fastAccessCoeff(k) = 2.0*(i+1)*(rank+1); \
626  x[i] = FadVecType(p, f); \
627  for (int j=0; j<p; j++) { \
628  x[i].fastAccessDx(j) = f; \
629  } \
630  } \
631  for (int i=0; i<n; i++) { \
632  VecType f(setup.sz, 0.0); \
633  for (int k=0; k<setup.sz; k++) \
634  f.fastAccessCoeff(k) = 2.0*(i+1); \
635  mins[i] = FadVecType(p, f); \
636  for (int j=0; j<p; j++) \
637  mins[i].fastAccessDx(j) = f; \
638  } \
639  Teuchos::scan(*comm, *setup.fad_vec_serializer, \
640  Teuchos::REDUCE_MIN, n, &x[0], &mins2[0]); \
641  success = checkVecArrays(mins, mins2, \
642  std::string(#FAD)+"<"+#Vec+"> Scan Min", out); \
643  success = checkResultOnAllProcs(*comm, out, success); \
644 } \
645  \
646 TEUCHOS_UNIT_TEST( Vec##_Comm, FadVec_SendReceive ) { \
647  typedef Sacado::mpl::apply<FadType,VecType>::type FadVecType; \
648  Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
649  comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
650  \
651  int num_proc = comm->getSize(); \
652  if (num_proc > 1) { \
653  int rank = comm->getRank(); \
654  int n = 7; \
655  int p = 5; \
656  Teuchos::Array<FadVecType> x(n), x2(n); \
657  for (int i=0; i<n; i++) { \
658  VecType f(setup.sz, 0.0); \
659  for (int k=0; k<setup.sz; k++) \
660  f.fastAccessCoeff(k) = 2.0*(i+1)*(k+1); \
661  x[i] = FadVecType(p, f); \
662  for (int j=0; j<p; j++) \
663  x[i].fastAccessDx(j) = f; \
664  } \
665  if (rank != 1) \
666  x2 = x; \
667  if (rank == 0) Teuchos::send(*comm, *setup.fad_vec_serializer, \
668  n, &x[0], 1); \
669  if (rank == 1) Teuchos::receive(*comm, *setup.fad_vec_serializer, \
670  0, n, &x2[0]); \
671  success = checkVecArrays(x, x2, \
672  std::string(#FAD)+"<"+#Vec+"> Send/Receive", out); \
673  success = checkResultOnAllProcs(*comm, out, success); \
674  } \
675  else \
676  success = true; \
677 }
678 
679 namespace DynamicVecTest {
680  Sacado::Random<double> rnd;
681  typedef int Ordinal;
682  typedef Kokkos::DefaultExecutionSpace execution_space;
684  typedef Sacado::Fad::DFad<double> fad_type;
687  VEC_COMM_TESTS(vec_type, fad_type, DynamicVector, DFad)
688 }
689 
690 namespace DynamicStridedVecTest {
691  Sacado::Random<double> rnd;
692  typedef int Ordinal;
693  typedef Kokkos::DefaultExecutionSpace execution_space;
695  typedef Sacado::Fad::DFad<double> fad_type;
698  VEC_COMM_TESTS(vec_type, fad_type, DynamicStridedVector, DFad)
699 }
700 
701 namespace StaticVecTest {
702  Sacado::Random<double> rnd;
703  typedef int Ordinal;
704  typedef Kokkos::DefaultExecutionSpace execution_space;
706  typedef Sacado::Fad::DFad<double> fad_type;
709  VEC_COMM_TESTS(vec_type, fad_type, StaticVector, DFad)
710 }
711 
712 namespace StaticFixedVecTest {
713  Sacado::Random<double> rnd;
714  typedef int Ordinal;
715  typedef Kokkos::DefaultExecutionSpace execution_space;
717  typedef Sacado::Fad::DFad<double> fad_type;
720  VEC_COMM_TESTS(vec_type, fad_type, StaticFixedVector, DFad)
721 }
722 
723 int main( int argc, char* argv[] ) {
724  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
726 }
Sacado::Fad::DFad< double > fad_type
UnitTestSetup< vec_type, fad_type > setup
Stokhos::StaticFixedStorage< int, double, 8, execution_space > storage_type
Kokkos::DefaultExecutionSpace execution_space
Statically allocated storage class.
Sacado::Fad::DFad< double > fad_type
#define VEC_COMM_TESTS(VecType, FadType, Vec, FAD)
Kokkos::DefaultExecutionSpace execution_space
UnitTestSetup< vec_type, fad_type > setup
Teuchos::ValueTypeSerializer< int, VecType > VecSerializerT
Sacado::Fad::DFad< double > fad_type
Stokhos::DynamicStorage< int, double, execution_space > storage_type
UnitTestSetup< vec_type, fad_type > setup
Stokhos::DynamicStridedStorage< int, double, execution_space > storage_type
Stokhos::StaticStorage< int, double, 8, execution_space > storage_type
Sacado::MP::Vector< storage_type > vec_type
static int runUnitTestsFromMain(int argc, char *argv[])
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
RCP< VecSerializerT > vec_serializer
Sacado::Random< double > rnd
Teuchos::ValueTypeSerializer< int, FadVecType > FadVecSerializerT
UnitTestSetup< vec_type, fad_type > setup
Kokkos::DefaultExecutionSpace execution_space
Sacado::Random< double > rnd
int main(int argc, char **argv)
Kokkos::DefaultExecutionSpace execution_space
Sacado::MP::Vector< storage_type > vec_type
Sacado::MP::Vector< storage_type > vec_type
Sacado::mpl::apply< FadType, VecType >::type FadVecType
size_type size() const
bool checkResultOnAllProcs(const Teuchos::Comm< Ordinal > &comm, Teuchos::FancyOStream &out, const bool result)
Sacado::MP::Vector< storage_type > vec_type
Statically allocated storage class.
bool checkVecArrays(const Teuchos::Array< VecType > &x, const Teuchos::Array< VecType > &x2, const std::string &tag, Teuchos::FancyOStream &out)
Sacado::Fad::DFad< double > fad_type
RCP< FadVecSerializerT > fad_vec_serializer