23 template <
typename ArrayType>
26 const std::string& tag,
28 typedef typename ArrayType::value_type
FadType;
31 bool success = (x.size() == x2.size());
32 out << tag <<
" Fad array size test";
37 out <<
": \n\tExpected: " << x.size() <<
", \n\tGot: " << x2.size()
41 const int sz = x.size();
42 for (
int i=0;
i<sz;
i++) {
44 out << tag <<
" Fad array comparison test " <<
i;
49 out <<
": \n\tExpected: " << x[
i] <<
", \n\tGot: " << x2[
i] <<
"."
51 success = success && success2;
57 template<
typename Ordinal>
64 out <<
"\nChecking that the above test passed in all processes ...";
65 int thisResult = ( result ? 1 : 0 );
69 const bool passed = sumResult==Teuchos::size(comm);
73 out <<
" (sumResult="<<sumResult<<
"!=numProcs="<<Teuchos::size(comm)<<
") failed\n";
77 #define FAD_BASE_COMM_TESTS(FadType, FAD) \
78 TEUCHOS_UNIT_TEST( FAD##_Comm, Fad_Broadcast ) { \
79 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
80 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
84 ValueTypeSerializer<int,FadType> fts( \
85 rcp(new ValueTypeSerializer<int,double>), p); \
87 Teuchos::Array<FadType> x(n), x2(n), x3(n); \
88 for (int i=0; i<n; i++) { \
89 x[i] = FadType(p, rnd.number()); \
90 for (int j=0; j<p; j++) \
91 x[i].fastAccessDx(j) = rnd.number(); \
93 for (int i=0; i<n; i++) { \
94 x2[i] = FadType(p, 0.0); \
96 if (comm->getRank() == 0) { \
101 Teuchos::broadcast(*comm, 0, n, &x2[0]); \
102 bool success1 = checkFadArrays( \
103 x, x2, std::string(#FAD)+" Broadcast", out); \
104 success1 = checkResultOnAllProcs(*comm, out, success1); \
106 Teuchos::broadcast(*comm, fts, 0, n, &x3[0]); \
107 bool success2 = checkFadArrays( \
108 x, x3, std::string(#FAD)+" Broadcast FTS", out); \
109 success2 = checkResultOnAllProcs(*comm, out, success2); \
111 success = success1 && success2; \
114 TEUCHOS_UNIT_TEST( FAD##_Comm, Fad_GatherAll ) { \
115 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
116 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
120 int size = comm->getSize(); \
121 int rank = comm->getRank(); \
123 ValueTypeSerializer<int,FadType> fts( \
124 rcp(new ValueTypeSerializer<int,double>), p); \
126 Teuchos::Array<FadType> x(n), x2(N), x3(N), x4(N); \
127 for (int i=0; i<n; i++) { \
128 x[i] = FadType(p, (rank+1)*(i+1)); \
129 for (int j=0; j<p; j++) \
130 x[i].fastAccessDx(j) = (rank+1)*(i+1)*(j+1); \
132 for (int i=0; i<N; i++) { \
133 x2[i] = FadType(p, 0.0); \
135 for (int j=0; j<size; j++) { \
136 for (int i=0; i<n; i++) { \
137 x3[n*j+i] = FadType(p, (j+1)*(i+1)); \
138 for (int k=0; k<p; k++) \
139 x3[n*j+i].fastAccessDx(k) = (j+1)*(i+1)*(k+1); \
143 Teuchos::gatherAll(*comm, n, &x[0], N, &x2[0]); \
144 bool success1 = checkFadArrays( \
145 x3, x2, std::string(#FAD)+" Gather All", out); \
146 success1 = checkResultOnAllProcs(*comm, out, success1); \
148 Teuchos::gatherAll(*comm, fts, n, &x[0], N, &x4[0]); \
149 bool success2 = checkFadArrays( \
150 x3, x4, std::string(#FAD)+" Gather All FTS", out); \
151 success2 = checkResultOnAllProcs(*comm, out, success2); \
153 success = success1 && success2; \
156 TEUCHOS_UNIT_TEST( FAD##_Comm, Fad_SumAll ) { \
157 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
158 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
162 int num_proc = comm->getSize(); \
163 ValueTypeSerializer<int,FadType> fts( \
164 rcp(new ValueTypeSerializer<int,double>), p); \
166 Teuchos::Array<FadType> x(n), sums(n), sums2(n), sums3(n); \
167 for (int i=0; i<n; i++) { \
168 x[i] = FadType(p, 1.0*(i+1)); \
169 for (int j=0; j<p; j++) \
170 x[i].fastAccessDx(j) = 2.0*(i+1); \
172 for (int i=0; i<n; i++) { \
173 sums[i] = FadType(p, 1.0*(i+1)*num_proc); \
174 for (int j=0; j<p; j++) \
175 sums[i].fastAccessDx(j) = 2.0*(i+1)*num_proc; \
177 for (int i=0; i<n; i++) { \
178 sums2[i] = FadType(p, 0.0); \
181 Teuchos::reduceAll(*comm, Teuchos::REDUCE_SUM, n, &x[0], &sums2[0]); \
182 bool success1 = checkFadArrays( \
183 sums, sums2, std::string(#FAD)+" Sum All", out); \
184 success1 = checkResultOnAllProcs(*comm, out, success1); \
186 Teuchos::reduceAll(*comm, fts, Teuchos::REDUCE_SUM, n, &x[0], &sums3[0]); \
187 bool success2 = checkFadArrays( \
188 sums, sums3, std::string(#FAD)+" Sum All FTS", out); \
189 success2 = checkResultOnAllProcs(*comm, out, success2); \
191 success = success1 && success2; \
194 TEUCHOS_UNIT_TEST( FAD##_Comm, Fad_MaxAll ) { \
195 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
196 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
200 int rank = comm->getRank(); \
201 int num_proc = comm->getSize(); \
202 ValueTypeSerializer<int,FadType> fts( \
203 rcp(new ValueTypeSerializer<int,double>), p); \
205 Teuchos::Array<FadType> x(n), maxs(n), maxs2(n), maxs3(n); \
206 for (int i=0; i<n; i++) { \
207 x[i] = FadType(p, 1.0*(i+1)*(rank+1)); \
208 for (int j=0; j<p; j++) \
209 x[i].fastAccessDx(j) = 2.0*(i+1)*(rank+1); \
211 for (int i=0; i<n; i++) { \
212 maxs[i] = FadType(p, 1.0*(i+1)*num_proc); \
213 for (int j=0; j<p; j++) \
214 maxs[i].fastAccessDx(j) = 2.0*(i+1)*num_proc; \
216 for (int i=0; i<n; i++) { \
217 maxs2[i] = FadType(p, 0.0); \
220 Teuchos::reduceAll(*comm, Teuchos::REDUCE_MAX, n, &x[0], &maxs2[0]); \
221 bool success1 = checkFadArrays( \
222 maxs, maxs2, std::string(#FAD)+" Max All", out); \
223 success1 = checkResultOnAllProcs(*comm, out, success1); \
225 Teuchos::reduceAll(*comm, fts, Teuchos::REDUCE_MAX, n, &x[0], &maxs3[0]); \
226 bool success2 = checkFadArrays( \
227 maxs, maxs3, std::string(#FAD)+" Max All FTS", out); \
228 success2 = checkResultOnAllProcs(*comm, out, success2); \
230 success = success1 && success2; \
233 TEUCHOS_UNIT_TEST( FAD##_Comm, Fad_MinAll ) { \
234 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
235 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
239 int rank = comm->getRank(); \
240 ValueTypeSerializer<int,FadType> fts( \
241 rcp(new ValueTypeSerializer<int,double>), p); \
243 Teuchos::Array<FadType> x(n), mins(n), mins2(n), mins3(n); \
244 for (int i=0; i<n; i++) { \
245 x[i] = FadType(p, 1.0*(i+1)*(rank+1)); \
246 for (int j=0; j<p; j++) \
247 x[i].fastAccessDx(j) = 2.0*(i+1)*(rank+1); \
249 for (int i=0; i<n; i++) { \
250 mins[i] = FadType(p, 1.0*(i+1)); \
251 for (int j=0; j<p; j++) \
252 mins[i].fastAccessDx(j) = 2.0*(i+1); \
254 for (int i=0; i<n; i++) { \
255 mins2[i] = FadType(p, 0.0); \
258 Teuchos::reduceAll(*comm, Teuchos::REDUCE_MIN, n, &x[0], &mins2[0]); \
259 bool success1 = checkFadArrays( \
260 mins, mins2, std::string(#FAD)+" Min All", out); \
261 success1 = checkResultOnAllProcs(*comm, out, success1); \
263 Teuchos::reduceAll(*comm, fts, Teuchos::REDUCE_MIN, n, &x[0], &mins3[0]); \
264 bool success2 = checkFadArrays( \
265 mins, mins3, std::string(#FAD)+" Min All FTS", out); \
266 success2 = checkResultOnAllProcs(*comm, out, success2); \
268 success = success1 && success2; \
271 TEUCHOS_UNIT_TEST( FAD##_Comm, Fad_ScanSum ) { \
272 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
273 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
277 int rank = comm->getRank(); \
278 ValueTypeSerializer<int,FadType> fts( \
279 rcp(new ValueTypeSerializer<int,double>), p); \
281 Teuchos::Array<FadType> x(n), sums(n), sums2(n), sums3(n); \
282 for (int i=0; i<n; i++) { \
283 x[i] = FadType(p, 1.0*(i+1)); \
284 for (int j=0; j<p; j++) \
285 x[i].fastAccessDx(j) = 2.0*(i+1); \
287 for (int i=0; i<n; i++) { \
288 sums[i] = FadType(p, 1.0*(i+1)*(rank+1)); \
289 for (int j=0; j<p; j++) \
290 sums[i].fastAccessDx(j) = 2.0*(i+1)*(rank+1); \
292 for (int i=0; i<n; i++) { \
293 sums2[i] = FadType(p, 0.0); \
296 Teuchos::scan(*comm, Teuchos::REDUCE_SUM, n, &x[0], &sums2[0]); \
297 bool success1 = checkFadArrays( \
298 sums, sums2, std::string(#FAD)+" Scan Sum", out); \
299 success1 = checkResultOnAllProcs(*comm, out, success1); \
301 Teuchos::scan(*comm, fts, Teuchos::REDUCE_SUM, n, &x[0], &sums3[0]); \
302 bool success2 = checkFadArrays( \
303 sums, sums3, std::string(#FAD)+" Scan Sum FTS", out); \
304 success2 = checkResultOnAllProcs(*comm, out, success2); \
306 success = success1 && success2; \
309 TEUCHOS_UNIT_TEST( FAD##_Comm, Fad_ScanMax ) { \
310 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
311 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
315 int rank = comm->getRank(); \
316 ValueTypeSerializer<int,FadType> fts( \
317 rcp(new ValueTypeSerializer<int,double>), p); \
319 Teuchos::Array<FadType> x(n), maxs(n), maxs2(n), maxs3(n); \
320 for (int i=0; i<n; i++) { \
321 x[i] = FadType(p, 1.0*(i+1)*(rank+1)); \
322 for (int j=0; j<p; j++) \
323 x[i].fastAccessDx(j) = 2.0*(i+1)*(rank+1); \
325 for (int i=0; i<n; i++) { \
326 maxs[i] = FadType(p, 1.0*(i+1)*(rank+1)); \
327 for (int j=0; j<p; j++) \
328 maxs[i].fastAccessDx(j) = 2.0*(i+1)*(rank+1); \
330 for (int i=0; i<n; i++) { \
331 maxs2[i] = FadType(p, 0.0); \
334 Teuchos::scan(*comm, Teuchos::REDUCE_MAX, n, &x[0], &maxs2[0]); \
335 bool success1 = checkFadArrays( \
336 maxs, maxs2, std::string(#FAD)+" Scan Max", out); \
337 success1 = checkResultOnAllProcs(*comm, out, success1); \
339 Teuchos::scan(*comm, fts, Teuchos::REDUCE_MAX, n, &x[0], &maxs3[0]); \
340 bool success2 = checkFadArrays( \
341 maxs, maxs3, std::string(#FAD)+" Scan Max FTS", out); \
342 success2 = checkResultOnAllProcs(*comm, out, success2); \
344 success = success1 && success2; \
347 TEUCHOS_UNIT_TEST( FAD##_Comm, Fad_ScanMin ) { \
348 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
349 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
353 int rank = comm->getRank(); \
354 ValueTypeSerializer<int,FadType> fts( \
355 rcp(new ValueTypeSerializer<int,double>), p); \
357 Teuchos::Array<FadType> x(n), mins(n), mins2(n), mins3(n); \
358 for (int i=0; i<n; i++) { \
359 x[i] = FadType(p, 1.0*(i+1)*(rank+1)); \
360 for (int j=0; j<p; j++) \
361 x[i].fastAccessDx(j) = 2.0*(i+1)*(rank+1); \
363 for (int i=0; i<n; i++) { \
364 mins[i] = FadType(p, 1.0*(i+1)); \
365 for (int j=0; j<p; j++) \
366 mins[i].fastAccessDx(j) = 2.0*(i+1); \
368 for (int i=0; i<n; i++) { \
369 mins2[i] = FadType(p, 0.0); \
372 Teuchos::scan(*comm, Teuchos::REDUCE_MIN, n, &x[0], &mins2[0]); \
373 bool success1 = checkFadArrays( \
374 mins, mins2, std::string(#FAD)+" Scan Min", out); \
375 success1 = checkResultOnAllProcs(*comm, out, success1); \
377 Teuchos::scan(*comm, fts, Teuchos::REDUCE_MIN, n, &x[0], &mins3[0]); \
378 bool success2 = checkFadArrays( \
379 mins, mins3, std::string(#FAD)+" Scan Min FTS", out); \
380 success2 = checkResultOnAllProcs(*comm, out, success2); \
382 success = success1 && success2; \
385 TEUCHOS_UNIT_TEST( FAD##_Comm, Fad_SendReceive ) { \
386 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
387 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
389 int num_proc = comm->getSize(); \
390 if (num_proc > 1) { \
391 int rank = comm->getRank(); \
394 ValueTypeSerializer<int,FadType> fts( \
395 rcp(new ValueTypeSerializer<int,double>), p); \
397 Teuchos::Array<FadType> x(n), x2(n), x3(n); \
398 for (int i=0; i<n; i++) { \
399 x[i] = FadType(p, 1.0*(i+1)); \
400 for (int j=0; j<p; j++) \
401 x[i].fastAccessDx(j) = 2.0*(i+1)*(j+1); \
403 for (int i=0; i<n; i++) { \
404 x2[i] = FadType(p, 0.0); \
411 if (rank == 0) Teuchos::send(*comm, n, &x[0], 1); \
412 if (rank == 1) Teuchos::receive(*comm, 0, n, &x2[0]); \
413 bool success1 = checkFadArrays( \
414 x, x2, std::string(#FAD)+" Send/Receive", out); \
415 success1 = checkResultOnAllProcs(*comm, out, success1); \
417 if (rank == 0) Teuchos::send(*comm, fts, n, &x[0], 1); \
418 if (rank == 1) Teuchos::receive(*comm, fts, 0, n, &x3[0]); \
419 bool success2 = checkFadArrays( \
420 x, x3, std::string(#FAD)+" Send/Receive FTS", out); \
421 success2 = checkResultOnAllProcs(*comm, out, success2); \
423 success = success1 && success2; \
429 TEUCHOS_UNIT_TEST( FAD##_Comm, FadFad_Broadcast ) { \
430 typedef Sacado::mpl::apply<FadType,FadType>::type FadFadType; \
431 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
432 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
437 RCP< ValueTypeSerializer<int,FadType> > fts = \
438 rcp(new ValueTypeSerializer<int,FadType>( \
439 rcp(new ValueTypeSerializer<int,double>), p1)); \
440 ValueTypeSerializer<int,FadFadType> ffts(fts, p2); \
442 Teuchos::Array<FadFadType> x(n), x2(n), x3(n); \
443 for (int i=0; i<n; i++) { \
444 FadType f(p1, rnd.number()); \
445 for (int k=0; k<p1; k++) \
446 f.fastAccessDx(k) = rnd.number(); \
447 x[i] = FadFadType(p2, f); \
448 for (int j=0; j<p2; j++) { \
449 FadType g(p1, rnd.number()); \
450 for (int k=0; k<p1; k++) \
451 g.fastAccessDx(k) = rnd.number(); \
452 x[i].fastAccessDx(j) = g; \
455 for (int i=0; i<n; i++) { \
456 x2[i] = FadFadType(p2, FadType(p1, 0.0)); \
457 for (int j=0; j<p2; j++) \
458 x2[i].fastAccessDx(j) = FadType(p1, 0.0); \
460 if (comm->getRank() == 0) { \
465 Teuchos::broadcast(*comm, 0, n, &x2[0]); \
466 bool success1 = checkFadArrays( \
467 x, x2, std::string(#FAD)+"<"+#FAD+"> Broadcast", out); \
468 success1 = checkResultOnAllProcs(*comm, out, success1); \
470 Teuchos::broadcast(*comm, ffts, 0, n, &x3[0]); \
471 bool success2 = checkFadArrays( \
472 x, x3, std::string(#FAD)+"<"+#FAD+"> Broadcast FTS", out); \
473 success2 = checkResultOnAllProcs(*comm, out, success2); \
475 success = success1 && success2; \
478 TEUCHOS_UNIT_TEST( FAD##_Comm, FadFad_GatherAll ) { \
479 typedef Sacado::mpl::apply<FadType,FadType>::type FadFadType; \
480 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
481 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
486 int size = comm->getSize(); \
487 int rank = comm->getRank(); \
489 RCP< ValueTypeSerializer<int,FadType> > fts = \
490 rcp(new ValueTypeSerializer<int,FadType>( \
491 rcp(new ValueTypeSerializer<int,double>), p1)); \
492 ValueTypeSerializer<int,FadFadType> ffts(fts, p2); \
494 Teuchos::Array<FadFadType> x(n), x2(N), x3(N), x4(N); \
495 for (int i=0; i<n; i++) { \
496 FadType f(p1, (rank+1)*(i+1)); \
497 for (int k=0; k<p1; k++) \
498 f.fastAccessDx(k) = (rank+1)*(i+1)*(k+1); \
499 x[i] = FadFadType(p2, f); \
500 for (int j=0; j<p2; j++) { \
501 x[i].fastAccessDx(j) = f; \
504 for (int i=0; i<N; i++) { \
505 x2[i] = FadFadType(p2, FadType(p1, 0.0)); \
506 for (int j=0; j<p2; j++) \
507 x2[i].fastAccessDx(j) = FadType(p1, 0.0); \
509 for (int j=0; j<size; j++) { \
510 for (int i=0; i<n; i++) { \
511 FadType f(p1, (j+1)*(i+1)); \
512 for (int k=0; k<p1; k++) \
513 f.fastAccessDx(k) = (j+1)*(i+1)*(k+1); \
514 x3[n*j+i] = FadFadType(p2, f); \
515 for (int k=0; k<p2; k++) \
516 x3[n*j+i].fastAccessDx(k) = f; \
520 Teuchos::gatherAll(*comm, n, &x[0], N, &x2[0]); \
521 bool success1 = checkFadArrays( \
522 x3, x2, std::string(#FAD)+"<"+#FAD+"> Gather All", out); \
523 success1 = checkResultOnAllProcs(*comm, out, success1); \
525 Teuchos::gatherAll(*comm, ffts, n, &x[0], N, &x4[0]); \
526 bool success2 = checkFadArrays( \
527 x3, x4, std::string(#FAD)+"<"+#FAD+"> Gather All FTS", out); \
528 success2 = checkResultOnAllProcs(*comm, out, success2); \
530 success = success1 && success2; \
533 TEUCHOS_UNIT_TEST( FAD##_Comm, FadFad_SumAll ) { \
534 typedef Sacado::mpl::apply<FadType,FadType>::type FadFadType; \
535 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
536 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
541 int num_proc = comm->getSize(); \
542 RCP< ValueTypeSerializer<int,FadType> > fts = \
543 rcp(new ValueTypeSerializer<int,FadType>( \
544 rcp(new ValueTypeSerializer<int,double>), p1)); \
545 ValueTypeSerializer<int,FadFadType> ffts(fts, p2); \
547 Teuchos::Array<FadFadType> x(n), sums(n), sums2(n), sums3(n); \
548 for (int i=0; i<n; i++) { \
549 FadType f(p1, 1.0*(i+1)); \
550 for (int k=0; k<p1; k++) \
551 f.fastAccessDx(k) = 2.0*(i+1); \
552 x[i] = FadFadType(p2, f); \
553 for (int j=0; j<p2; j++) { \
554 x[i].fastAccessDx(j) = f; \
557 for (int i=0; i<n; i++) { \
558 FadType f(p1, 1.0*(i+1)*num_proc); \
559 for (int k=0; k<p1; k++) \
560 f.fastAccessDx(k) = 2.0*(i+1)*num_proc; \
561 sums[i] = FadFadType(p2, f); \
562 for (int j=0; j<p2; j++) \
563 sums[i].fastAccessDx(j) = f; \
565 for (int i=0; i<n; i++) { \
566 sums2[i] = FadFadType(p2, FadType(p1, 0.0)); \
567 for (int j=0; j<p2; j++) \
568 sums2[i].fastAccessDx(j) = FadType(p1, 0.0); \
571 Teuchos::reduceAll(*comm, Teuchos::REDUCE_SUM, n, &x[0], &sums2[0]); \
572 bool success1 = checkFadArrays( \
573 sums, sums2, std::string(#FAD)+"<"+#FAD+"> Sum All", out); \
574 success1 = checkResultOnAllProcs(*comm, out, success1); \
576 Teuchos::reduceAll(*comm, ffts, Teuchos::REDUCE_SUM, n, &x[0], &sums3[0]); \
577 bool success2 = checkFadArrays( \
578 sums, sums3, std::string(#FAD)+"<"+#FAD+"> Sum All", out); \
579 success2 = checkResultOnAllProcs(*comm, out, success2); \
581 success = success1 && success2; \
584 TEUCHOS_UNIT_TEST( FAD##_Comm, FadFad_MaxAll ) { \
585 typedef Sacado::mpl::apply<FadType,FadType>::type FadFadType; \
586 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
587 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
592 int rank = comm->getRank(); \
593 int num_proc = comm->getSize(); \
594 RCP< ValueTypeSerializer<int,FadType> > fts = \
595 rcp(new ValueTypeSerializer<int,FadType>( \
596 rcp(new ValueTypeSerializer<int,double>), p1)); \
597 ValueTypeSerializer<int,FadFadType> ffts(fts, p2); \
599 Teuchos::Array<FadFadType> x(n), maxs(n), maxs2(n), maxs3(n); \
600 for (int i=0; i<n; i++) { \
601 FadType f(p1, 1.0*(i+1)*(rank+1)); \
602 for (int k=0; k<p1; k++) \
603 f.fastAccessDx(k) = 2.0*(i+1)*(rank+1); \
604 x[i] = FadFadType(p2, f); \
605 for (int j=0; j<p2; j++) { \
606 x[i].fastAccessDx(j) = f; \
609 for (int i=0; i<n; i++) { \
610 FadType f(p1, 1.0*(i+1)*num_proc); \
611 for (int k=0; k<p1; k++) \
612 f.fastAccessDx(k) = 2.0*(i+1)*num_proc; \
613 maxs[i] = FadFadType(p2, f); \
614 for (int j=0; j<p2; j++) \
615 maxs[i].fastAccessDx(j) = f; \
617 for (int i=0; i<n; i++) { \
618 maxs2[i] = FadFadType(p2, FadType(p1, 0.0)); \
619 for (int j=0; j<p2; j++) \
620 maxs2[i].fastAccessDx(j) = FadType(p1, 0.0); \
623 Teuchos::reduceAll(*comm, Teuchos::REDUCE_MAX, n, &x[0], &maxs2[0]); \
624 bool success1 = checkFadArrays( \
625 maxs, maxs2, std::string(#FAD)+"<"+#FAD+"> Max All", out); \
626 success1 = checkResultOnAllProcs(*comm, out, success1); \
628 Teuchos::reduceAll(*comm, ffts, Teuchos::REDUCE_MAX, n, &x[0], &maxs3[0]); \
629 bool success2 = checkFadArrays( \
630 maxs, maxs3, std::string(#FAD)+"<"+#FAD+"> Max All FTS", out); \
631 success2 = checkResultOnAllProcs(*comm, out, success2); \
633 success = success1 && success2; \
636 TEUCHOS_UNIT_TEST( FAD##_Comm, FadFad_MinAll ) { \
637 typedef Sacado::mpl::apply<FadType,FadType>::type FadFadType; \
638 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
639 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
644 int rank = comm->getRank(); \
645 RCP< ValueTypeSerializer<int,FadType> > fts = \
646 rcp(new ValueTypeSerializer<int,FadType>( \
647 rcp(new ValueTypeSerializer<int,double>), p1)); \
648 ValueTypeSerializer<int,FadFadType> ffts(fts, p2); \
650 Teuchos::Array<FadFadType> x(n), mins(n), mins2(n), mins3(n); \
651 for (int i=0; i<n; i++) { \
652 FadType f(p1, 1.0*(i+1)*(rank+1)); \
653 for (int k=0; k<p1; k++) \
654 f.fastAccessDx(k) = 2.0*(i+1)*(rank+1); \
655 x[i] = FadFadType(p2, f); \
656 for (int j=0; j<p2; j++) { \
657 x[i].fastAccessDx(j) = f; \
660 for (int i=0; i<n; i++) { \
661 FadType f(p1, 1.0*(i+1)); \
662 for (int k=0; k<p1; k++) \
663 f.fastAccessDx(k) = 2.0*(i+1); \
664 mins[i] = FadFadType(p2, f); \
665 for (int j=0; j<p2; j++) \
666 mins[i].fastAccessDx(j) = f; \
668 for (int i=0; i<n; i++) { \
669 mins2[i] = FadFadType(p2, FadType(p1, 0.0)); \
670 for (int j=0; j<p2; j++) \
671 mins2[i].fastAccessDx(j) = FadType(p1, 0.0); \
674 Teuchos::reduceAll(*comm, Teuchos::REDUCE_MIN, n, &x[0], &mins2[0]); \
675 bool success1 = checkFadArrays( \
676 mins, mins2, std::string(#FAD)+"<"+#FAD+"> Min All", out); \
677 success1 = checkResultOnAllProcs(*comm, out, success1); \
679 Teuchos::reduceAll(*comm, ffts, Teuchos::REDUCE_MIN, n, &x[0], &mins3[0]); \
680 bool success2 = checkFadArrays( \
681 mins, mins3, std::string(#FAD)+"<"+#FAD+"> Min All FTS", out); \
682 success2 = checkResultOnAllProcs(*comm, out, success2); \
684 success = success1 && success2; \
687 TEUCHOS_UNIT_TEST( FAD##_Comm, FadFad_ScanSum ) { \
688 typedef Sacado::mpl::apply<FadType,FadType>::type FadFadType; \
689 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
690 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
695 int rank = comm->getRank(); \
696 RCP< ValueTypeSerializer<int,FadType> > fts = \
697 rcp(new ValueTypeSerializer<int,FadType>( \
698 rcp(new ValueTypeSerializer<int,double>), p1)); \
699 ValueTypeSerializer<int,FadFadType> ffts(fts, p2); \
701 Teuchos::Array<FadFadType> x(n), sums(n), sums2(n), sums3(n); \
702 for (int i=0; i<n; i++) { \
703 FadType f(p1, 1.0*(i+1)); \
704 for (int k=0; k<p1; k++) \
705 f.fastAccessDx(k) = 2.0*(i+1); \
706 x[i] = FadFadType(p2, f); \
707 for (int j=0; j<p2; j++) { \
708 x[i].fastAccessDx(j) = f; \
711 for (int i=0; i<n; i++) { \
712 FadType f(p1, 1.0*(i+1)*(rank+1)); \
713 for (int k=0; k<p1; k++) \
714 f.fastAccessDx(k) = 2.0*(i+1)*(rank+1); \
715 sums[i] = FadFadType(p2, f); \
716 for (int j=0; j<p2; j++) \
717 sums[i].fastAccessDx(j) = f; \
719 for (int i=0; i<n; i++) { \
720 sums2[i] = FadFadType(p2, FadType(p1, 0.0)); \
721 for (int j=0; j<p2; j++) \
722 sums2[i].fastAccessDx(j) = FadType(p1, 0.0); \
725 Teuchos::scan(*comm, Teuchos::REDUCE_SUM, n, &x[0], &sums2[0]); \
726 bool success1 = checkFadArrays( \
727 sums, sums2, std::string(#FAD)+"<"+#FAD+"> Scan Sum", out); \
728 success1 = checkResultOnAllProcs(*comm, out, success1); \
730 Teuchos::scan(*comm, ffts, Teuchos::REDUCE_SUM, n, &x[0], &sums3[0]); \
731 bool success2 = checkFadArrays( \
732 sums, sums3, std::string(#FAD)+"<"+#FAD+"> Scan Sum FTS", out); \
733 success2 = checkResultOnAllProcs(*comm, out, success2); \
735 success = success1 && success2; \
738 TEUCHOS_UNIT_TEST( FAD##_Comm, FadFad_ScanMax ) { \
739 typedef Sacado::mpl::apply<FadType,FadType>::type FadFadType; \
740 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
741 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
746 int rank = comm->getRank(); \
747 RCP< ValueTypeSerializer<int,FadType> > fts = \
748 rcp(new ValueTypeSerializer<int,FadType>( \
749 rcp(new ValueTypeSerializer<int,double>), p1)); \
750 ValueTypeSerializer<int,FadFadType> ffts(fts, p2); \
752 Teuchos::Array<FadFadType> x(n), maxs(n), maxs2(n), maxs3(n); \
753 for (int i=0; i<n; i++) { \
754 FadType f(p1, 1.0*(i+1)*(rank+1)); \
755 for (int k=0; k<p1; k++) \
756 f.fastAccessDx(k) = 2.0*(i+1)*(rank+1); \
757 x[i] = FadFadType(p2, f); \
758 for (int j=0; j<p2; j++) { \
759 x[i].fastAccessDx(j) = f; \
762 for (int i=0; i<n; i++) { \
763 FadType f(p1, 1.0*(i+1)*(rank+1)); \
764 for (int k=0; k<p1; k++) \
765 f.fastAccessDx(k) = 2.0*(i+1)*(rank+1); \
766 maxs[i] = FadFadType(p2, f); \
767 for (int j=0; j<p2; j++) \
768 maxs[i].fastAccessDx(j) = f; \
770 for (int i=0; i<n; i++) { \
771 maxs2[i] = FadFadType(p2, FadType(p1, 0.0)); \
772 for (int j=0; j<p2; j++) \
773 maxs2[i].fastAccessDx(j) = FadType(p1, 0.0); \
776 Teuchos::scan(*comm, Teuchos::REDUCE_MAX, n, &x[0], &maxs2[0]); \
777 bool success1 = checkFadArrays( \
778 maxs, maxs2, std::string(#FAD)+"<"+#FAD+"> Scan Max", out); \
779 success1 = checkResultOnAllProcs(*comm, out, success1); \
781 Teuchos::scan(*comm, ffts, Teuchos::REDUCE_MAX, n, &x[0], &maxs3[0]); \
782 bool success2 = checkFadArrays( \
783 maxs, maxs3, std::string(#FAD)+"<"+#FAD+"> Scan Max FTS", out); \
784 success2 = checkResultOnAllProcs(*comm, out, success2); \
786 success = success1 && success2; \
789 TEUCHOS_UNIT_TEST( FAD##_Comm, FadFad_ScanMin ) { \
790 typedef Sacado::mpl::apply<FadType,FadType>::type FadFadType; \
791 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
792 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
797 int rank = comm->getRank(); \
798 RCP< ValueTypeSerializer<int,FadType> > fts = \
799 rcp(new ValueTypeSerializer<int,FadType>( \
800 rcp(new ValueTypeSerializer<int,double>), p1)); \
801 ValueTypeSerializer<int,FadFadType> ffts(fts, p2); \
803 Teuchos::Array<FadFadType> x(n), mins(n), mins2(n), mins3(n); \
804 for (int i=0; i<n; i++) { \
805 FadType f(p1, 1.0*(i+1)*(rank+1)); \
806 for (int k=0; k<p1; k++) \
807 f.fastAccessDx(k) = 2.0*(i+1)*(rank+1); \
808 x[i] = FadFadType(p2, f); \
809 for (int j=0; j<p2; j++) { \
810 x[i].fastAccessDx(j) = f; \
813 for (int i=0; i<n; i++) { \
814 FadType f(p1, 1.0*(i+1)); \
815 for (int k=0; k<p1; k++) \
816 f.fastAccessDx(k) = 2.0*(i+1); \
817 mins[i] = FadFadType(p2, f); \
818 for (int j=0; j<p2; j++) \
819 mins[i].fastAccessDx(j) = f; \
821 for (int i=0; i<n; i++) { \
822 mins2[i] = FadFadType(p2, FadType(p1, 0.0)); \
823 for (int j=0; j<p2; j++) \
824 mins2[i].fastAccessDx(j) = FadType(p1, 0.0); \
827 Teuchos::scan(*comm, Teuchos::REDUCE_MIN, n, &x[0], &mins2[0]); \
828 bool success1 = checkFadArrays( \
829 mins, mins2, std::string(#FAD)+"<"+#FAD+"> Scan Min", out); \
830 success1 = checkResultOnAllProcs(*comm, out, success1); \
832 Teuchos::scan(*comm, ffts, Teuchos::REDUCE_MIN, n, &x[0], &mins3[0]); \
833 bool success2 = checkFadArrays( \
834 mins, mins3, std::string(#FAD)+"<"+#FAD+"> Scan Min FTS", out); \
835 success2 = checkResultOnAllProcs(*comm, out, success2); \
837 success = success1 && success2; \
840 TEUCHOS_UNIT_TEST( FAD##_Comm, FadFad_SendReceive ) { \
841 typedef Sacado::mpl::apply<FadType,FadType>::type FadFadType; \
842 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
843 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
845 int num_proc = comm->getSize(); \
846 if (num_proc > 1) { \
847 int rank = comm->getRank(); \
851 RCP< ValueTypeSerializer<int,FadType> > fts = \
852 rcp(new ValueTypeSerializer<int,FadType>( \
853 rcp(new ValueTypeSerializer<int,double>), p1)); \
854 ValueTypeSerializer<int,FadFadType> ffts(fts, p2); \
856 Teuchos::Array<FadFadType> x(n), x2(n), x3(n); \
857 for (int i=0; i<n; i++) { \
858 FadType f(p1, 1.0*(i+1)); \
859 for (int k=0; k<p1; k++) \
860 f.fastAccessDx(k) = 2.0*(i+1)*(k+1); \
861 x[i] = FadFadType(p2, f); \
862 for (int j=0; j<p2; j++) \
863 x[i].fastAccessDx(j) = f; \
865 for (int i=0; i<n; i++) { \
866 x2[i] = FadFadType(p2, FadType(p1, 0.0)); \
867 for (int j=0; j<p2; j++) \
868 x2[i].fastAccessDx(j) = FadType(p1, 0.0); \
875 if (rank == 0) Teuchos::send(*comm, n, &x[0], 1); \
876 if (rank == 1) Teuchos::receive(*comm, 0, n, &x2[0]); \
877 bool success1 = checkFadArrays( \
878 x, x2, std::string(#FAD)+"<"+#FAD+"> Send/Receive", out); \
879 success1 = checkResultOnAllProcs(*comm, out, success1); \
881 if (rank == 0) Teuchos::send(*comm, ffts, n, &x[0], 1); \
882 if (rank == 1) Teuchos::receive(*comm, ffts, 0, n, &x3[0]); \
883 bool success2 = checkFadArrays( \
884 x, x3, std::string(#FAD)+"<"+#FAD+"> Send/Receive FTS", out); \
885 success2 = checkResultOnAllProcs(*comm, out, success2); \
887 success = success1 && success2; \
893 #if defined(HAVE_SACADO_KOKKOS) && defined(HAVE_SACADO_TEUCHOSKOKKOSCOMM)
895 #include "Kokkos_Core.hpp"
897 #define FAD_KOKKOS_COMM_TESTS_DEV(FadType, FAD, Device) \
898 TEUCHOS_UNIT_TEST( FAD##_Comm_Kokkos_##Device, Fad_Broadcast ) { \
899 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
900 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
905 ValueTypeSerializer<int,FadType> fts( \
906 rcp(new ValueTypeSerializer<int,double>), p); \
908 typedef Kokkos::View<FadType*,Device> ViewType; \
909 typedef ViewType::HostMirror HostViewType; \
910 ViewType x("x",n,p+1), x2("x2",n,p+1), x3("x3",n,p+1); \
911 HostViewType h_x = Kokkos::create_mirror_view(x); \
912 HostViewType h_x2 = Kokkos::create_mirror_view(x2); \
913 HostViewType h_x3 = Kokkos::create_mirror_view(x3); \
914 for (int i=0; i<n; i++) { \
915 h_x[i] = FadType(p, rnd.number()); \
916 for (int j=0; j<p; j++) \
917 h_x[i].fastAccessDx(j) = rnd.number(); \
919 for (int i=0; i<n; i++) { \
920 h_x2[i] = FadType(p, 0.0); \
922 Kokkos::deep_copy(x, h_x); \
923 Kokkos::deep_copy(x2, h_x2); \
924 if (comm->getRank() == 0) { \
934 const bool accessible = \
935 Kokkos::Impl::MemorySpaceAccess< \
937 typename Device::memory_space >::accessible; \
939 Teuchos::broadcast(*comm, 0, n, x2); \
940 Kokkos::deep_copy(h_x2, x2); \
943 Teuchos::broadcast(*comm, 0, n, h_x2); \
944 bool success1 = checkFadArrays( \
945 h_x, h_x2, std::string(#FAD)+" Broadcast", out); \
946 success1 = checkResultOnAllProcs(*comm, out, success1); \
949 Teuchos::broadcast(*comm, fts, 0, n, x3); \
950 Kokkos::deep_copy(h_x3, x3); \
953 Teuchos::broadcast(*comm, fts, 0, n, h_x3); \
954 bool success2 = checkFadArrays( \
955 h_x, h_x3, std::string(#FAD)+" Broadcast FTS", out); \
956 success2 = checkResultOnAllProcs(*comm, out, success2); \
958 success = success1 && success2; \
960 TEUCHOS_UNIT_TEST( FAD##_Comm_Kokkos_##Device, Fad_SumAll ) { \
961 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
962 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
967 int num_proc = comm->getSize(); \
968 ValueTypeSerializer<int,FadType> fts( \
969 rcp(new ValueTypeSerializer<int,double>), p); \
971 typedef Kokkos::View<FadType*,Device> ViewType; \
972 typedef ViewType::HostMirror HostViewType; \
973 ViewType x("x",n,p+1), sums("sums",n,p+1), \
974 sums2("sums2",n,p+1), sums3("sums3",n,p+1); \
975 HostViewType h_x = Kokkos::create_mirror_view(x); \
976 HostViewType h_sums = Kokkos::create_mirror_view(sums); \
977 HostViewType h_sums2 = Kokkos::create_mirror_view(sums2); \
978 HostViewType h_sums3 = Kokkos::create_mirror_view(sums3); \
979 for (int i=0; i<n; i++) { \
980 h_x[i] = FadType(p, 1.0*(i+1)); \
981 for (int j=0; j<p; j++) \
982 h_x[i].fastAccessDx(j) = 2.0*(i+1); \
984 for (int i=0; i<n; i++) { \
985 h_sums[i] = FadType(p, 1.0*(i+1)*num_proc); \
986 for (int j=0; j<p; j++) \
987 h_sums[i].fastAccessDx(j) = 2.0*(i+1)*num_proc; \
989 for (int i=0; i<n; i++) { \
990 h_sums2[i] = FadType(p, 0.0); \
992 Kokkos::deep_copy(x, h_x); \
993 Kokkos::deep_copy(sums, h_sums); \
994 Kokkos::deep_copy(sums2, h_sums2); \
999 const bool accessible = \
1000 Kokkos::Impl::MemorySpaceAccess< \
1001 Kokkos::HostSpace, \
1002 typename Device::memory_space >::accessible; \
1004 Teuchos::reduceAll(*comm, Teuchos::REDUCE_SUM, n, x, sums2); \
1005 Kokkos::deep_copy(h_sums2, sums2); \
1008 Teuchos::reduceAll(*comm, Teuchos::REDUCE_SUM, n, h_x, h_sums2); \
1009 bool success1 = checkFadArrays( \
1010 h_sums, h_sums2, std::string(#FAD)+" Sum All", out); \
1011 success1 = checkResultOnAllProcs(*comm, out, success1); \
1014 Teuchos::reduceAll(*comm, fts, Teuchos::REDUCE_SUM, n, x, sums3); \
1015 Kokkos::deep_copy(h_sums3, sums3); \
1018 Teuchos::reduceAll(*comm, fts, Teuchos::REDUCE_SUM, n, h_x, h_sums3); \
1019 bool success2 = checkFadArrays( \
1020 h_sums, h_sums3, std::string(#FAD)+" Sum All FTS", out); \
1021 success2 = checkResultOnAllProcs(*comm, out, success2); \
1022 success = success1 && success2; \
1025 TEUCHOS_UNIT_TEST( FAD##_Comm_Kokkos_##Device, Fad_MaxAll ) { \
1026 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
1027 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
1032 int rank = comm->getRank(); \
1033 int num_proc = comm->getSize(); \
1034 ValueTypeSerializer<int,FadType> fts( \
1035 rcp(new ValueTypeSerializer<int,double>), p); \
1037 typedef Kokkos::View<FadType*,Device> ViewType; \
1038 typedef ViewType::HostMirror HostViewType; \
1039 ViewType x("x",n,p+1), maxs("maxs",n,p+1), \
1040 maxs2("maxs2",n,p+1), maxs3("maxs3",n,p+1); \
1041 HostViewType h_x = Kokkos::create_mirror_view(x); \
1042 HostViewType h_maxs = Kokkos::create_mirror_view(maxs); \
1043 HostViewType h_maxs2 = Kokkos::create_mirror_view(maxs2); \
1044 HostViewType h_maxs3 = Kokkos::create_mirror_view(maxs3); \
1045 for (int i=0; i<n; i++) { \
1046 h_x[i] = FadType(p, 1.0*(i+1)*(rank+1)); \
1047 for (int j=0; j<p; j++) \
1048 h_x[i].fastAccessDx(j) = 2.0*(i+1)*(rank+1); \
1050 for (int i=0; i<n; i++) { \
1051 h_maxs[i] = FadType(p, 1.0*(i+1)*num_proc); \
1052 for (int j=0; j<p; j++) \
1053 h_maxs[i].fastAccessDx(j) = 2.0*(i+1)*num_proc; \
1055 for (int i=0; i<n; i++) { \
1056 h_maxs2[i] = FadType(p, 0.0); \
1058 Kokkos::deep_copy(x, h_x); \
1059 Kokkos::deep_copy(maxs, h_maxs); \
1060 Kokkos::deep_copy(maxs2, h_maxs2); \
1065 const bool accessible = \
1066 Kokkos::Impl::MemorySpaceAccess< \
1067 Kokkos::HostSpace, \
1068 typename Device::memory_space >::accessible; \
1070 Teuchos::reduceAll(*comm, Teuchos::REDUCE_MAX, n, x, maxs2); \
1071 Kokkos::deep_copy(h_maxs2, maxs2); \
1074 Teuchos::reduceAll(*comm, Teuchos::REDUCE_MAX, n, h_x, h_maxs2); \
1075 bool success1 = checkFadArrays( \
1076 h_maxs, h_maxs2, std::string(#FAD)+" Max All", out); \
1077 success1 = checkResultOnAllProcs(*comm, out, success1); \
1080 Teuchos::reduceAll(*comm, fts, Teuchos::REDUCE_MAX, n, x, maxs3); \
1081 Kokkos::deep_copy(h_maxs3, maxs3); \
1084 Teuchos::reduceAll(*comm, fts, Teuchos::REDUCE_MAX, n, h_x, h_maxs3); \
1085 bool success2 = checkFadArrays( \
1086 h_maxs, h_maxs3, std::string(#FAD)+" Max All FTS", out); \
1087 success2 = checkResultOnAllProcs(*comm, out, success2); \
1088 success = success1 && success2; \
1091 TEUCHOS_UNIT_TEST( FAD##_Comm_Kokkos_##Device, Fad_MinAll ) { \
1092 Teuchos::RCP<const Teuchos::Comm<Ordinal> > \
1093 comm = Teuchos::DefaultComm<Ordinal>::getComm(); \
1098 int rank = comm->getRank(); \
1099 ValueTypeSerializer<int,FadType> fts( \
1100 rcp(new ValueTypeSerializer<int,double>), p); \
1102 typedef Kokkos::View<FadType*,Device> ViewType; \
1103 typedef ViewType::HostMirror HostViewType; \
1104 ViewType x("x",n,p+1), mins("mins",n,p+1), \
1105 mins2("mins2",n,p+1), mins3("mins3",n,p+1); \
1106 HostViewType h_x = Kokkos::create_mirror_view(x); \
1107 HostViewType h_mins = Kokkos::create_mirror_view(mins); \
1108 HostViewType h_mins2 = Kokkos::create_mirror_view(mins2); \
1109 HostViewType h_mins3 = Kokkos::create_mirror_view(mins3); \
1110 for (int i=0; i<n; i++) { \
1111 h_x[i] = FadType(p, 1.0*(i+1)*(rank+1)); \
1112 for (int j=0; j<p; j++) \
1113 h_x[i].fastAccessDx(j) = 2.0*(i+1)*(rank+1); \
1115 for (int i=0; i<n; i++) { \
1116 h_mins[i] = FadType(p, 1.0*(i+1)); \
1117 for (int j=0; j<p; j++) \
1118 h_mins[i].fastAccessDx(j) = 2.0*(i+1); \
1120 for (int i=0; i<n; i++) { \
1121 h_mins2[i] = FadType(p, 0.0); \
1123 Kokkos::deep_copy(x, h_x); \
1124 Kokkos::deep_copy(mins, h_mins); \
1125 Kokkos::deep_copy(mins2, h_mins2); \
1130 const bool accessible = \
1131 Kokkos::Impl::MemorySpaceAccess< \
1132 Kokkos::HostSpace, \
1133 typename Device::memory_space >::accessible; \
1135 Teuchos::reduceAll(*comm, Teuchos::REDUCE_MIN, n, x, mins2); \
1136 Kokkos::deep_copy(h_mins2, mins2); \
1139 Teuchos::reduceAll(*comm, Teuchos::REDUCE_MIN, n, h_x, h_mins2); \
1140 bool success1 = checkFadArrays( \
1141 h_mins, h_mins2, std::string(#FAD)+" Min All", out); \
1142 success1 = checkResultOnAllProcs(*comm, out, success1); \
1145 Teuchos::reduceAll(*comm, fts, Teuchos::REDUCE_MIN, n, x, mins3); \
1146 Kokkos::deep_copy(h_mins3, mins3); \
1149 Teuchos::reduceAll(*comm, fts, Teuchos::REDUCE_MIN, n, h_x, h_mins3); \
1150 bool success2 = checkFadArrays( \
1151 h_mins, h_mins3, std::string(#FAD)+" Min All FTS", out); \
1152 success2 = checkResultOnAllProcs(*comm, out, success2); \
1153 success = success1 && success2; \
1157 #ifdef KOKKOS_ENABLE_OPENMP
1158 #define FAD_KOKKOS_COMM_TESTS_OPENMP(FadType, FAD) \
1159 using Kokkos::OpenMP; \
1160 FAD_KOKKOS_COMM_TESTS_DEV(FadType, FAD, OpenMP)
1162 #define FAD_KOKKOS_COMM_TESTS_OPENMP(FadType, FAD)
1165 #ifdef KOKKOS_ENABLE_THREADS
1166 #define FAD_KOKKOS_COMM_TESTS_THREADS(FadType, FAD) \
1167 using Kokkos::Threads; \
1168 FAD_KOKKOS_COMM_TESTS_DEV(FadType, FAD, Threads)
1170 #define FAD_KOKKOS_COMM_TESTS_THREADS(FadType, FAD)
1173 #ifdef KOKKOS_ENABLE_CUDA
1174 #define FAD_KOKKOS_COMM_TESTS_CUDA(FadType, FAD) \
1175 using Kokkos::Cuda; \
1176 FAD_KOKKOS_COMM_TESTS_DEV(FadType, FAD, Cuda)
1178 #define FAD_KOKKOS_COMM_TESTS_CUDA(FadType, FAD)
1181 #ifdef KOKKOS_ENABLE_HIP
1182 #define FAD_KOKKOS_COMM_TESTS_HIP(FadType, FAD) \
1183 using Kokkos::HIP; \
1184 FAD_KOKKOS_COMM_TESTS_DEV(FadType, FAD, HIP)
1186 #define FAD_KOKKOS_COMM_TESTS_HIP(FadType, FAD)
1189 #ifdef KOKKOS_ENABLE_SERIAL
1190 #define FAD_KOKKOS_COMM_TESTS_SERIAL(FadType, FAD) \
1191 using Kokkos::Serial; \
1192 FAD_KOKKOS_COMM_TESTS_DEV(FadType, FAD, Serial)
1194 #define FAD_KOKKOS_COMM_TESTS_SERIAL(FadType, FAD)
1197 #define FAD_KOKKOS_COMM_TESTS(FadType, FAD) \
1198 FAD_KOKKOS_COMM_TESTS_OPENMP(FadType, FAD) \
1199 FAD_KOKKOS_COMM_TESTS_THREADS(FadType, FAD) \
1200 FAD_KOKKOS_COMM_TESTS_CUDA(FadType, FAD) \
1201 FAD_KOKKOS_COMM_TESTS_SERIAL(FadType, FAD)
1205 #define FAD_KOKKOS_COMM_TESTS(FadType, FAD)
1209 #define FAD_COMM_TESTS(FadType, FAD) \
1210 FAD_BASE_COMM_TESTS(FadType, FAD)
Sacado::Fad::DFad< double > FadType
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
bool checkFadArrays(const ArrayType &x, const ArrayType &x2, const std::string &tag, Teuchos::FancyOStream &out)
bool checkResultOnAllProcs(const Teuchos::Comm< Ordinal > &comm, Teuchos::FancyOStream &out, const bool result)
static SACADO_INLINE_FUNCTION bool eval(const T &x, const T &y)