30 #ifndef SACADO_FAD_EXP_EXPRASSIGN_HPP
31 #define SACADO_FAD_EXP_EXPRASSIGN_HPP
38 #ifndef SACADO_FAD_DERIV_LOOP
39 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL_DFAD) && !defined(SACADO_DISABLE_CUDA_IN_KOKKOS) && defined(__CUDA_ARCH__)
40 #define SACADO_FAD_DERIV_LOOP(I,SZ) for (int I=threadIdx.x; I<SZ; I+=blockDim.x)
42 #define SACADO_FAD_DERIV_LOOP(I,SZ) for (int I=0; I<SZ; ++I)
46 #ifndef SACADO_FAD_THREAD_SINGLE
47 #if (defined(SACADO_VIEW_CUDA_HIERARCHICAL) || defined(SACADO_VIEW_CUDA_HIERARCHICAL_DFAD)) && !defined(SACADO_DISABLE_CUDA_IN_KOKKOS) && defined(__CUDA_ARCH__)
48 #define SACADO_FAD_THREAD_SINGLE if (threadIdx.x == 0)
50 #define SACADO_FAD_THREAD_SINGLE
62 template <
typename DstType,
typename Enabled =
void>
70 template <
typename SrcType>
74 const int xsz = x.size();
76 if (xsz != dst.size())
77 dst.resizeAndZero(xsz);
79 const int sz = dst.size();
87 if (x.hasFastAccess()) {
89 dst.fastAccessDx(i) = x.fastAccessDx(i);
93 dst.fastAccessDx(i) = x.dx(i);
100 template <
typename SrcType>
104 const int xsz = x.size(), sz = dst.size();
106 #if defined(SACADO_DEBUG) && !defined(__CUDA_ARCH__ )
107 if ((xsz != sz) && (xsz != 0) && (sz != 0))
108 throw "Fad Error: Attempt to assign with incompatible sizes";
113 if (x.hasFastAccess())
115 dst.fastAccessDx(i) += x.fastAccessDx(i);
117 for (
int i=0; i<sz; ++i)
118 dst.fastAccessDx(i) += x.dx(i);
121 dst.resizeAndZero(xsz);
122 if (x.hasFastAccess())
124 dst.fastAccessDx(i) = x.fastAccessDx(i);
127 dst.fastAccessDx(i) = x.dx(i);
131 dst.val() += x.val();
135 template <
typename SrcType>
139 const int xsz = x.size(), sz = dst.size();
141 #if defined(SACADO_DEBUG) && !defined(__CUDA_ARCH__ )
142 if ((xsz != sz) && (xsz != 0) && (sz != 0))
143 throw "Fad Error: Attempt to assign with incompatible sizes";
148 if (x.hasFastAccess())
150 dst.fastAccessDx(i) -= x.fastAccessDx(i);
153 dst.fastAccessDx(i) -= x.dx(i);
156 dst.resizeAndZero(xsz);
157 if (x.hasFastAccess())
159 dst.fastAccessDx(i) = -x.fastAccessDx(i);
162 dst.fastAccessDx(i) = -x.dx(i);
166 dst.val() -= x.val();
170 template <
typename SrcType>
174 const int xsz = x.size(), sz = dst.size();
178 #if defined(SACADO_DEBUG) && !defined(__CUDA_ARCH__ )
179 if ((xsz != sz) && (xsz != 0) && (sz != 0))
180 throw "Fad Error: Attempt to assign with incompatible sizes";
185 if (x.hasFastAccess())
187 dst.fastAccessDx(i) = v*x.fastAccessDx(i) + dst.fastAccessDx(i)*xval;
190 dst.fastAccessDx(i) = v*x.dx(i) + dst.fastAccessDx(i)*xval;
193 dst.resizeAndZero(xsz);
194 if (x.hasFastAccess())
196 dst.fastAccessDx(i) = v*x.fastAccessDx(i);
199 dst.fastAccessDx(i) = v*x.dx(i);
205 dst.fastAccessDx(i) *= xval;
213 template <
typename SrcType>
217 const int xsz = x.size(), sz = dst.size();
221 #if defined(SACADO_DEBUG) && !defined(__CUDA_ARCH__ )
222 if ((xsz != sz) && (xsz != 0) && (sz != 0))
223 throw "Fad Error: Attempt to assign with incompatible sizes";
229 if (x.hasFastAccess())
231 dst.fastAccessDx(i) =
232 ( dst.fastAccessDx(i)*xval - v*x.fastAccessDx(i) ) / xval2;
235 dst.fastAccessDx(i) =
236 ( dst.fastAccessDx(i)*xval - v*x.dx(i) ) / xval2;
239 dst.resizeAndZero(xsz);
240 if (x.hasFastAccess())
242 dst.fastAccessDx(i) = - v*x.fastAccessDx(i) / xval2;
245 dst.fastAccessDx(i) = -v*x.dx(i) / xval2;
251 dst.fastAccessDx(i) /= xval;
266 template <
typename DstType>
268 typename std::enable_if<Sacado::IsStaticallySized<DstType>::value>::type> {
275 template <
typename SrcType>
279 const int sz = dst.size();
281 dst.fastAccessDx(i) = x.fastAccessDx(i);
286 template <
typename SrcType>
290 const int sz = dst.size();
292 dst.fastAccessDx(i) += x.fastAccessDx(i);
293 dst.val() += x.val();
297 template <
typename SrcType>
301 const int sz = dst.size();
303 dst.fastAccessDx(i) -= x.fastAccessDx(i);
304 dst.val() -= x.val();
308 template <
typename SrcType>
312 const int sz = dst.size();
316 dst.fastAccessDx(i) = v*x.fastAccessDx(i) + dst.fastAccessDx(i)*xval;
321 template <
typename SrcType>
325 const int sz = dst.size();
330 dst.fastAccessDx(i) =
331 ( dst.fastAccessDx(i)*xval - v*x.fastAccessDx(i) )/ xval2;
342 #endif // SACADO_FAD_EXP_EXPRASSIGN_HPP
Class that implements various forms of expression assignments.
DstType::value_type value_type
Typename of values.
static KOKKOS_INLINE_FUNCTION void assign_times_equal(DstType &dst, const SrcType &x)
Implementation of dst *= x.
static KOKKOS_INLINE_FUNCTION void assign_minus_equal(DstType &dst, const SrcType &x)
Implementation of dst -= x.
static KOKKOS_INLINE_FUNCTION void assign_equal(DstType &dst, const SrcType &x)
Implementation of dst = x.
static KOKKOS_INLINE_FUNCTION void assign_equal(DstType &dst, const SrcType &x)
Implementation of dst = x.
#define KOKKOS_INLINE_FUNCTION
static KOKKOS_INLINE_FUNCTION void assign_divide_equal(DstType &dst, const SrcType &x)
Implementation of dst /= x.
static KOKKOS_INLINE_FUNCTION void assign_plus_equal(DstType &dst, const SrcType &x)
Implementation of dst += x.
static KOKKOS_INLINE_FUNCTION void assign_times_equal(DstType &dst, const SrcType &x)
Implementation of dst *= x.
#define SACADO_FAD_DERIV_LOOP(I, SZ)
DstType::value_type value_type
Typename of values.
static KOKKOS_INLINE_FUNCTION void assign_minus_equal(DstType &dst, const SrcType &x)
Implementation of dst -= x.
static KOKKOS_INLINE_FUNCTION void assign_plus_equal(DstType &dst, const SrcType &x)
Implementation of dst += x.
static KOKKOS_INLINE_FUNCTION void assign_divide_equal(DstType &dst, const SrcType &x)
Implementation of dst /= x.