Sacado Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Sacado_Fad_Exp_Atomic.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Sacado Package
5 // Copyright (2006) Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // This library is free software; you can redistribute it and/or modify
11 // it under the terms of the GNU Lesser General Public License as
12 // published by the Free Software Foundation; either version 2.1 of the
13 // License, or (at your option) any later version.
14 //
15 // This library is distributed in the hope that it will be useful, but
16 // WITHOUT ANY WARRANTY; without even the implied warranty of
17 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 // Lesser General Public License for more details.
19 //
20 // You should have received a copy of the GNU Lesser General Public
21 // License along with this library; if not, write to the Free Software
22 // Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301
23 // USA
24 // Questions? Contact David M. Gay (dmgay@sandia.gov) or Eric T. Phipps
25 // (etphipp@sandia.gov).
26 //
27 // ***********************************************************************
28 // @HEADER
29 
30 #ifndef SACADO_FAD_EXP_ATOMIC_HPP
31 #define SACADO_FAD_EXP_ATOMIC_HPP
32 
33 #include "Sacado_ConfigDefs.h"
34 #if defined(HAVE_SACADO_KOKKOS)
35 
37 #include "Kokkos_Atomic.hpp"
38 #include "impl/Kokkos_Error.hpp"
39 
40 namespace Sacado {
41 
42  namespace Fad {
43  namespace Exp {
44 
45  // Overload of Kokkos::atomic_add for ViewFad types.
46  template <typename ValT, unsigned sl, unsigned ss, typename U, typename T>
48  void atomic_add(ViewFadPtr<ValT,sl,ss,U> dst, const Expr<T>& xx) {
49  using Kokkos::atomic_add;
50 
51  const typename Expr<T>::derived_type& x = xx.derived();
52 
53  const int xsz = x.size();
54  const int sz = dst->size();
55 
56  // We currently cannot handle resizing since that would need to be
57  // done atomically.
58  if (xsz > sz)
59  Kokkos::abort(
60  "Sacado error: Fad resize within atomic_add() not supported!");
61 
62  if (xsz != sz && sz > 0 && xsz > 0)
63  Kokkos::abort(
64  "Sacado error: Fad assignment of incompatiable sizes!");
65 
66 
67  if (sz > 0 && xsz > 0) {
69  atomic_add(&(dst->fastAccessDx(i)), x.fastAccessDx(i));
70  }
72  atomic_add(&(dst->val()), x.val());
73  }
74 
75  namespace Impl {
76  // Our implementation of Kokkos::atomic_oper_fetch() and
77  // Kokkos::atomic_fetch_oper() for Sacado types on host
78  template <typename Oper, typename DestPtrT, typename ValT, typename T>
79  typename Sacado::BaseExprType< Expr<T> >::type
80  atomic_oper_fetch_host(const Oper& op, DestPtrT dest, ValT* dest_val,
81  const Expr<T>& x)
82  {
83  typedef typename Sacado::BaseExprType< Expr<T> >::type return_type;
84  const typename Expr<T>::derived_type& val = x.derived();
85 
86 #ifdef KOKKOS_INTERNAL_NOT_PARALLEL
87  auto scope = desul::MemoryScopeCaller();
88 #else
89  auto scope = desul::MemoryScopeDevice();
90 #endif
91 
92  while (!desul::Impl::lock_address((void*)dest_val, scope))
93  ;
94  desul::atomic_thread_fence(desul::MemoryOrderAcquire(), scope);
95  return_type return_val = op.apply(*dest, val);
96  *dest = return_val;
97  desul::atomic_thread_fence(desul::MemoryOrderRelease(), scope);
98  desul::Impl::unlock_address((void*)dest_val, scope);
99  return return_val;
100  }
101 
102  template <typename Oper, typename DestPtrT, typename ValT, typename T>
103  typename Sacado::BaseExprType< Expr<T> >::type
104  atomic_fetch_oper_host(const Oper& op, DestPtrT dest, ValT* dest_val,
105  const Expr<T>& x)
106  {
107  typedef typename Sacado::BaseExprType< Expr<T> >::type return_type;
108  const typename Expr<T>::derived_type& val = x.derived();
109 
110 #ifdef KOKKOS_INTERNAL_NOT_PARALLEL
111  auto scope = desul::MemoryScopeCaller();
112 #else
113  auto scope = desul::MemoryScopeDevice();
114 #endif
115 
116  while (!desul::Impl::lock_address((void*)dest_val, scope))
117  ;
118  desul::atomic_thread_fence(desul::MemoryOrderAcquire(), scope);
119  return_type return_val = *dest;
120  *dest = op.apply(return_val, val);
121  desul::atomic_thread_fence(desul::MemoryOrderRelease(), scope);
122  desul::Impl::unlock_address((void*)dest_val, scope);
123  return return_val;
124  }
125 
126  // Helper function to decide if we are using team-based parallelism
127 #if defined(KOKKOS_ENABLE_CUDA) || defined(KOKKOS_ENABLE_HIP)
128  __device__
129  inline bool atomics_use_team() {
130 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL) || defined(SACADO_VIEW_CUDA_HIERARCHICAL_DFAD)
131  // It is not allowed to define SACADO_VIEW_CUDA_HIERARCHICAL or
132  // SACADO_VIEW_CUDA_HIERARCHICAL_DFAD and use Sacado inside a team-based
133  // kernel without Sacado hierarchical parallelism. So use the
134  // team-based version only if blockDim.x > 1 (i.e., a team policy)
135  return (blockDim.x > 1);
136 #else
137  return false;
138 #endif
139  }
140 #endif
141 
142 #if defined(KOKKOS_ENABLE_CUDA)
143 
144  // Our implementation of Kokkos::atomic_oper_fetch() and
145  // Kokkos::atomic_fetch_oper() for Sacado types on device
146  template <typename Oper, typename DestPtrT, typename ValT, typename T>
147  __device__
148  typename Sacado::BaseExprType< Expr<T> >::type
149  atomic_oper_fetch_device(const Oper& op, DestPtrT dest, ValT* dest_val,
150  const Expr<T>& x)
151  {
152  typedef typename Sacado::BaseExprType< Expr<T> >::type return_type;
153  const typename Expr<T>::derived_type& val = x.derived();
154 
155  auto scope = desul::MemoryScopeDevice();
156 
157  if (atomics_use_team()) {
158  int go = 1;
159  while (go) {
160  if (threadIdx.x == 0)
161  go = !desul::Impl::lock_address_cuda((void*)dest_val, scope);
162  go = Kokkos::shfl(go, 0, blockDim.x);
163  }
164  desul::atomic_thread_fence(desul::MemoryOrderAcquire(), scope);
165  return_type return_val = op.apply(*dest, val);
166  *dest = return_val;
167  desul::atomic_thread_fence(desul::MemoryOrderRelease(), scope);
168  if (threadIdx.x == 0)
169  desul::Impl::unlock_address_cuda((void*)dest_val, scope);
170  return return_val;
171  }
172  else {
173  return_type return_val;
174  // This is a way to avoid dead lock in a warp
175  int done = 0;
176  unsigned int mask = __activemask() ;
177  unsigned int active = __ballot_sync(mask, 1);
178  unsigned int done_active = 0;
179  while (active != done_active) {
180  if (!done) {
181  if (desul::Impl::lock_address_cuda((void*)dest_val, scope)) {
182  desul::atomic_thread_fence(desul::MemoryOrderAcquire(), scope);
183  return_val = op.apply(*dest, val);
184  *dest = return_val;
185  desul::atomic_thread_fence(desul::MemoryOrderRelease(), scope);
186  desul::Impl::unlock_address_cuda((void*)dest_val, scope);
187  done = 1;
188  }
189  }
190  done_active = __ballot_sync(mask, done);
191  }
192  return return_val;
193  }
194  }
195 
196  template <typename Oper, typename DestPtrT, typename ValT, typename T>
197  __device__
198  typename Sacado::BaseExprType< Expr<T> >::type
199  atomic_fetch_oper_device(const Oper& op, DestPtrT dest, ValT* dest_val,
200  const Expr<T>& x)
201  {
202  typedef typename Sacado::BaseExprType< Expr<T> >::type return_type;
203  const typename Expr<T>::derived_type& val = x.derived();
204 
205  auto scope = desul::MemoryScopeDevice();
206 
207  if (atomics_use_team()) {
208  int go = 1;
209  while (go) {
210  if (threadIdx.x == 0)
211  go = !desul::Impl::lock_address_cuda((void*)dest_val, scope);
212  go = Kokkos::shfl(go, 0, blockDim.x);
213  }
214  desul::atomic_thread_fence(desul::MemoryOrderAcquire(), scope);
215  return_type return_val = *dest;
216  *dest = op.apply(return_val, val);
217  desul::atomic_thread_fence(desul::MemoryOrderRelease(), scope);
218  if (threadIdx.x == 0)
219  desul::Impl::unlock_address_cuda((void*)dest_val, scope);
220  return return_val;
221  }
222  else {
223  return_type return_val;
224  // This is a way to (hopefully) avoid dead lock in a warp
225  int done = 0;
226  unsigned int mask = __activemask() ;
227  unsigned int active = __ballot_sync(mask, 1);
228  unsigned int done_active = 0;
229  while (active != done_active) {
230  if (!done) {
231  if (desul::Impl::lock_address_cuda((void*)dest_val, scope)) {
232  desul::atomic_thread_fence(desul::MemoryOrderAcquire(), scope);
233  return_val = *dest;
234  *dest = op.apply(return_val, val);
235  desul::atomic_thread_fence(desul::MemoryOrderRelease(), scope);
236  desul::Impl::unlock_address_cuda((void*)dest_val, scope);
237  done = 1;
238  }
239  }
240  done_active = __ballot_sync(mask, done);
241  }
242  return return_val;
243  }
244  }
245 
246 #elif defined(KOKKOS_ENABLE_HIP)
247 
248  // Our implementation of Kokkos::atomic_oper_fetch() and
249  // Kokkos::atomic_fetch_oper() for Sacado types on device
250  template <typename Oper, typename DestPtrT, typename ValT, typename T>
251  __device__
252  typename Sacado::BaseExprType< Expr<T> >::type
253  atomic_oper_fetch_device(const Oper& op, DestPtrT dest, ValT* dest_val,
254  const Expr<T>& x)
255  {
256  typedef typename Sacado::BaseExprType< Expr<T> >::type return_type;
257  const typename Expr<T>::derived_type& val = x.derived();
258 
259  auto scope = desul::MemoryScopeDevice();
260 
261  if (atomics_use_team()) {
262  int go = 1;
263  while (go) {
264  if (threadIdx.x == 0)
265  go = !desul::Impl::lock_address_hip((void*)dest_val, scope);
266  go = Kokkos::shfl(go, 0, blockDim.x);
267  }
268  desul::atomic_thread_fence(desul::MemoryOrderAcquire(), scope);
269  return_type return_val = op.apply(*dest, val);
270  *dest = return_val;
271  desul::atomic_thread_fence(desul::MemoryOrderRelease(), scope);
272  if (threadIdx.x == 0)
273  desul::Impl::unlock_address_hip((void*)dest_val, scope);
274  return return_val;
275  }
276  else {
277  return_type return_val;
278  int done = 0;
279  unsigned int active = __ballot(1);
280  unsigned int done_active = 0;
281  while (active != done_active) {
282  if (!done) {
283  if (desul::Impl::lock_address_hip((void*)dest_val, scope)) {
284  return_val = op.apply(*dest, val);
285  *dest = return_val;
286  desul::Impl::unlock_address_hip((void*)dest_val, scope);
287  done = 1;
288  }
289  }
290  done_active = __ballot(done);
291  }
292  return return_val;
293  }
294  }
295 
296  template <typename Oper, typename DestPtrT, typename ValT, typename T>
297  __device__
298  typename Sacado::BaseExprType< Expr<T> >::type
299  atomic_fetch_oper_device(const Oper& op, DestPtrT dest, ValT* dest_val,
300  const Expr<T>& x)
301  {
302  typedef typename Sacado::BaseExprType< Expr<T> >::type return_type;
303  const typename Expr<T>::derived_type& val = x.derived();
304 
305  auto scope = desul::MemoryScopeDevice();
306 
307  if (atomics_use_team()) {
308  int go = 1;
309  while (go) {
310  if (threadIdx.x == 0)
311  go = !desul::Impl::lock_address_hip((void*)dest_val, scope);
312  go = Kokkos::shfl(go, 0, blockDim.x);
313  }
314  desul::atomic_thread_fence(desul::MemoryOrderAcquire(), scope);
315  return_type return_val = *dest;
316  *dest = op.apply(return_val, val);
317  desul::atomic_thread_fence(desul::MemoryOrderRelease(), scope);
318  if (threadIdx.x == 0)
319  desul::Impl::unlock_address_hip((void*)dest_val, scope);
320  return return_val;
321  }
322  else {
323  return_type return_val;
324  int done = 0;
325  unsigned int active = __ballot(1);
326  unsigned int done_active = 0;
327  while (active != done_active) {
328  if (!done) {
329  if (desul::Impl::lock_address_hip((void*)dest_val, scope)) {
330  return_val = *dest;
331  *dest = op.apply(return_val, val);
332  desul::Impl::unlock_address_hip((void*)dest_val, scope);
333  done = 1;
334  }
335  }
336  done_active = __ballot(done);
337  }
338  return return_val;
339  }
340  }
341 
342 #elif defined(KOKKOS_ENABLE_SYCL)
343 
344  // Our implementation of Kokkos::atomic_oper_fetch() and
345  // Kokkos::atomic_fetch_oper() for Sacado types on device
346  template <typename Oper, typename DestPtrT, typename ValT, typename T>
347  typename Sacado::BaseExprType< Expr<T> >::type
348  atomic_oper_fetch_device(const Oper& op, DestPtrT dest, ValT* dest_val,
349  const Expr<T>& x)
350  {
351  Kokkos::abort("Not implemented!");
352  return {};
353  }
354 
355  template <typename Oper, typename DestPtrT, typename ValT, typename T>
356  typename Sacado::BaseExprType< Expr<T> >::type
357  atomic_fetch_oper_device(const Oper& op, DestPtrT dest, ValT* dest_val,
358  const Expr<T>& x)
359  {
360  Kokkos::abort("Not implemented!");
361  return {};
362  }
363 #endif
364 
365  // Overloads of Kokkos::atomic_oper_fetch/Kokkos::atomic_fetch_oper
366  // for Sacado types
367  template <typename Oper, typename S>
368  SACADO_INLINE_FUNCTION GeneralFad<S>
369  atomic_oper_fetch(const Oper& op, GeneralFad<S>* dest,
370  const GeneralFad<S>& val)
371  {
372  KOKKOS_IF_ON_HOST(return Impl::atomic_oper_fetch_host(op, dest, &(dest->val()), val);)
373  KOKKOS_IF_ON_DEVICE(return Impl::atomic_oper_fetch_device(op, dest, &(dest->val()), val);)
374  }
375  template <typename Oper, typename ValT, unsigned sl, unsigned ss,
376  typename U, typename T>
378  atomic_oper_fetch(const Oper& op, ViewFadPtr<ValT,sl,ss,U> dest,
379  const Expr<T>& val)
380  {
381  KOKKOS_IF_ON_HOST(return Impl::atomic_oper_fetch_host(op, dest, &dest.val(), val);)
382  KOKKOS_IF_ON_DEVICE(return Impl::atomic_oper_fetch_device(op, dest, &dest.val(), val);)
383  }
384 
385  template <typename Oper, typename S>
386  SACADO_INLINE_FUNCTION GeneralFad<S>
387  atomic_fetch_oper(const Oper& op, GeneralFad<S>* dest,
388  const GeneralFad<S>& val)
389  {
390  KOKKOS_IF_ON_HOST(return Impl::atomic_fetch_oper_host(op, dest, &(dest->val()), val);)
391  KOKKOS_IF_ON_DEVICE(return Impl::atomic_fetch_oper_device(op, dest, &(dest->val()), val);)
392  }
393  template <typename Oper, typename ValT, unsigned sl, unsigned ss,
394  typename U, typename T>
396  atomic_fetch_oper(const Oper& op, ViewFadPtr<ValT,sl,ss,U> dest,
397  const Expr<T>& val)
398  {
399  KOKKOS_IF_ON_HOST(return Impl::atomic_fetch_oper_host(op, dest, &dest.val(), val);)
400  KOKKOS_IF_ON_DEVICE(return Impl::atomic_fetch_oper_device(op, dest, &dest.val(), val);)
401  }
402 
403  // Our definition of the various Oper classes to be more type-flexible
404  struct MaxOper {
405  template <class Scalar1, class Scalar2>
406  KOKKOS_FORCEINLINE_FUNCTION
407  static auto apply(const Scalar1& val1, const Scalar2& val2)
408  -> decltype(max(val1,val2))
409  {
410  return max(val1,val2);
411  }
412  };
413  struct MinOper {
414  template <class Scalar1, class Scalar2>
415  KOKKOS_FORCEINLINE_FUNCTION
416  static auto apply(const Scalar1& val1, const Scalar2& val2)
417  -> decltype(min(val1,val2))
418  {
419  return min(val1,val2);
420  }
421  };
422  struct AddOper {
423  template <class Scalar1, class Scalar2>
424  KOKKOS_FORCEINLINE_FUNCTION
425  static auto apply(const Scalar1& val1, const Scalar2& val2)
426  -> decltype(val1+val2)
427  {
428  return val1 + val2;
429  }
430  };
431  struct SubOper {
432  template <class Scalar1, class Scalar2>
433  KOKKOS_FORCEINLINE_FUNCTION
434  static auto apply(const Scalar1& val1, const Scalar2& val2)
435  -> decltype(val1-val2)
436  {
437  return val1 - val2;
438  }
439  };
440  struct MulOper {
441  template <class Scalar1, class Scalar2>
442  KOKKOS_FORCEINLINE_FUNCTION
443  static auto apply(const Scalar1& val1, const Scalar2& val2)
444  -> decltype(val1*val2)
445  {
446  return val1 * val2;
447  }
448  };
449  struct DivOper {
450  template <class Scalar1, class Scalar2>
451  KOKKOS_FORCEINLINE_FUNCTION
452  static auto apply(const Scalar1& val1, const Scalar2& val2)
453  -> decltype(val1/val2)
454  {
455  return val1 / val2;
456  }
457  };
458 
459  } // Impl
460 
461  // Overload of Kokkos::atomic_*_fetch() and Kokkos::atomic_fetch_*()
462  // for Sacado types
463  template <typename S>
464  SACADO_INLINE_FUNCTION GeneralFad<S>
465  atomic_max_fetch(GeneralFad<S>* dest, const GeneralFad<S>& val) {
466  return Impl::atomic_oper_fetch(Impl::MaxOper(), dest, val);
467  }
468  template <typename ValT, unsigned sl, unsigned ss, typename U, typename T>
470  atomic_max_fetch(ViewFadPtr<ValT,sl,ss,U> dest, const Expr<T>& val) {
471  return Impl::atomic_oper_fetch(Impl::MaxOper(), dest, val);
472  }
473  template <typename S>
474  SACADO_INLINE_FUNCTION GeneralFad<S>
475  atomic_min_fetch(GeneralFad<S>* dest, const GeneralFad<S>& val) {
476  return Impl::atomic_oper_fetch(Impl::MinOper(), dest, val);
477  }
478  template <typename ValT, unsigned sl, unsigned ss, typename U, typename T>
480  atomic_min_fetch(ViewFadPtr<ValT,sl,ss,U> dest, const Expr<T>& val) {
481  return Impl::atomic_oper_fetch(Impl::MinOper(), dest, val);
482  }
483  template <typename S>
484  SACADO_INLINE_FUNCTION GeneralFad<S>
485  atomic_add_fetch(GeneralFad<S>* dest, const GeneralFad<S>& val) {
486  return Impl::atomic_oper_fetch(Impl::AddOper(), dest, val);
487  }
488  template <typename ValT, unsigned sl, unsigned ss, typename U, typename T>
490  atomic_add_fetch(ViewFadPtr<ValT,sl,ss,U> dest, const Expr<T>& val) {
491  return Impl::atomic_oper_fetch(Impl::AddOper(), dest, val);
492  }
493  template <typename S>
494  SACADO_INLINE_FUNCTION GeneralFad<S>
495  atomic_sub_fetch(GeneralFad<S>* dest, const GeneralFad<S>& val) {
496  return Impl::atomic_oper_fetch(Impl::SubOper(), dest, val);
497  }
498  template <typename ValT, unsigned sl, unsigned ss, typename U, typename T>
500  atomic_sub_fetch(ViewFadPtr<ValT,sl,ss,U> dest, const Expr<T>& val) {
501  return Impl::atomic_oper_fetch(Impl::SubOper(), dest, val);
502  }
503  template <typename S>
504  SACADO_INLINE_FUNCTION GeneralFad<S>
505  atomic_mul_fetch(GeneralFad<S>* dest, const GeneralFad<S>& val) {
506  return atomic_oper_fetch(Impl::MulOper(), dest, val);
507  }
508  template <typename ValT, unsigned sl, unsigned ss, typename U, typename T>
510  atomic_mul_fetch(ViewFadPtr<ValT,sl,ss,U> dest, const Expr<T>& val) {
511  return Impl::atomic_oper_fetch(Impl::MulOper(), dest, val);
512  }
513  template <typename S>
514  SACADO_INLINE_FUNCTION GeneralFad<S>
515  atomic_div_fetch(GeneralFad<S>* dest, const GeneralFad<S>& val) {
516  return Impl::atomic_oper_fetch(Impl::DivOper(), dest, val);
517  }
518  template <typename ValT, unsigned sl, unsigned ss, typename U, typename T>
520  atomic_div_fetch(ViewFadPtr<ValT,sl,ss,U> dest, const Expr<T>& val) {
521  return Impl::atomic_oper_fetch(Impl::DivOper(), dest, val);
522  }
523 
524  template <typename S>
525  SACADO_INLINE_FUNCTION GeneralFad<S>
526  atomic_fetch_max(GeneralFad<S>* dest, const GeneralFad<S>& val) {
527  return Impl::atomic_fetch_oper(Impl::MaxOper(), dest, val);
528  }
529  template <typename ValT, unsigned sl, unsigned ss, typename U, typename T>
531  atomic_fetch_max(ViewFadPtr<ValT,sl,ss,U> dest, const Expr<T>& val) {
532  return Impl::atomic_fetch_oper(Impl::MaxOper(), dest, val);
533  }
534  template <typename S>
535  SACADO_INLINE_FUNCTION GeneralFad<S>
536  atomic_fetch_min(GeneralFad<S>* dest, const GeneralFad<S>& val) {
537  return Impl::atomic_fetch_oper(Impl::MinOper(), dest, val);
538  }
539  template <typename ValT, unsigned sl, unsigned ss, typename U, typename T>
541  atomic_fetch_min(ViewFadPtr<ValT,sl,ss,U> dest, const Expr<T>& val) {
542  return Impl::atomic_fetch_oper(Impl::MinOper(), dest, val);
543  }
544  template <typename S>
545  SACADO_INLINE_FUNCTION GeneralFad<S>
546  atomic_fetch_add(GeneralFad<S>* dest, const GeneralFad<S>& val) {
547  return Impl::atomic_fetch_oper(Impl::AddOper(), dest, val);
548  }
549  template <typename ValT, unsigned sl, unsigned ss, typename U, typename T>
551  atomic_fetch_add(ViewFadPtr<ValT,sl,ss,U> dest, const Expr<T>& val) {
552  return Impl::atomic_fetch_oper(Impl::AddOper(), dest, val);
553  }
554  template <typename S>
555  SACADO_INLINE_FUNCTION GeneralFad<S>
556  atomic_fetch_sub(GeneralFad<S>* dest, const GeneralFad<S>& val) {
557  return Impl::atomic_fetch_oper(Impl::SubOper(), dest, val);
558  }
559  template <typename ValT, unsigned sl, unsigned ss, typename U, typename T>
561  atomic_fetch_sub(ViewFadPtr<ValT,sl,ss,U> dest, const Expr<T>& val) {
562  return Impl::atomic_fetch_oper(Impl::SubOper(), dest, val);
563  }
564  template <typename S>
565  SACADO_INLINE_FUNCTION GeneralFad<S>
566  atomic_fetch_mul(GeneralFad<S>* dest, const GeneralFad<S>& val) {
567  return Impl::atomic_fetch_oper(Impl::MulOper(), dest, val);
568  }
569  template <typename ValT, unsigned sl, unsigned ss, typename U, typename T>
571  atomic_fetch_mul(ViewFadPtr<ValT,sl,ss,U> dest, const Expr<T>& val) {
572  return Impl::atomic_fetch_oper(Impl::MulOper(), dest, val);
573  }
574  template <typename S>
575  SACADO_INLINE_FUNCTION GeneralFad<S>
576  atomic_fetch_div(GeneralFad<S>* dest, const GeneralFad<S>& val) {
577  return Impl::atomic_fetch_oper(Impl::DivOper(), dest, val);
578  }
579  template <typename ValT, unsigned sl, unsigned ss, typename U, typename T>
581  atomic_fetch_div(ViewFadPtr<ValT,sl,ss,U> dest, const Expr<T>& val) {
582  return Impl::atomic_fetch_oper(Impl::DivOper(), dest, val);
583  }
584 
585  } // namespace Exp
586  } // namespace Fad
587 
588 } // namespace Sacado
589 
590 #endif // HAVE_SACADO_KOKKOS
591 #endif // SACADO_FAD_EXP_VIEWFAD_HPP
#define SACADO_FAD_THREAD_SINGLE
expr val()
#define T
Definition: Sacado_rad.hpp:573
SimpleFad< ValueT > min(const SimpleFad< ValueT > &a, const SimpleFad< ValueT > &b)
#define SACADO_FAD_DERIV_LOOP(I, SZ)
Get the base Fad type from a view/expression.
expr expr expr fastAccessDx(i)) FAD_UNARYOP_MACRO(exp
SimpleFad< ValueT > max(const SimpleFad< ValueT > &a, const SimpleFad< ValueT > &b)
#define SACADO_INLINE_FUNCTION