10 #ifndef SACADO_FAD_EXP_EXPRASSIGN_HPP
11 #define SACADO_FAD_EXP_EXPRASSIGN_HPP
18 #ifndef SACADO_FAD_DERIV_LOOP
19 #if defined(SACADO_VIEW_CUDA_HIERARCHICAL_DFAD) && !defined(SACADO_DISABLE_CUDA_IN_KOKKOS) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
20 #define SACADO_FAD_DERIV_LOOP(I,SZ) for (int I=threadIdx.x; I<SZ; I+=blockDim.x)
22 #define SACADO_FAD_DERIV_LOOP(I,SZ) for (int I=0; I<SZ; ++I)
26 #ifndef SACADO_FAD_THREAD_SINGLE
27 #if (defined(SACADO_VIEW_CUDA_HIERARCHICAL) || defined(SACADO_VIEW_CUDA_HIERARCHICAL_DFAD)) && !defined(SACADO_DISABLE_CUDA_IN_KOKKOS) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
28 #define SACADO_FAD_THREAD_SINGLE if (threadIdx.x == 0)
30 #define SACADO_FAD_THREAD_SINGLE
42 template <
typename DstType,
typename Enabled =
void>
50 template <
typename SrcType>
54 const int xsz = x.size();
56 if (xsz != dst.size())
57 dst.resizeAndZero(xsz);
59 const int sz = dst.size();
67 if (x.hasFastAccess()) {
69 dst.fastAccessDx(
i) = x.fastAccessDx(
i);
73 dst.fastAccessDx(
i) = x.dx(
i);
80 template <
typename SrcType>
84 const int xsz = x.size(), sz = dst.size();
86 #if defined(SACADO_DEBUG) && !defined(__CUDA_ARCH__ ) && !defined(__HIP_DEVICE_COMPILE__ )
87 if ((xsz != sz) && (xsz != 0) && (sz != 0))
88 throw "Fad Error: Attempt to assign with incompatible sizes";
93 if (x.hasFastAccess())
95 dst.fastAccessDx(
i) += x.fastAccessDx(
i);
97 for (
int i=0;
i<sz; ++
i)
98 dst.fastAccessDx(
i) += x.dx(
i);
101 dst.resizeAndZero(xsz);
102 if (x.hasFastAccess())
104 dst.fastAccessDx(
i) = x.fastAccessDx(
i);
107 dst.fastAccessDx(
i) = x.dx(
i);
111 dst.val() += x.val();
115 template <
typename SrcType>
119 const int xsz = x.size(), sz = dst.size();
121 #if defined(SACADO_DEBUG) && !defined(__CUDA_ARCH__ ) && !defined(__HIP_DEVICE_COMPILE__ )
122 if ((xsz != sz) && (xsz != 0) && (sz != 0))
123 throw "Fad Error: Attempt to assign with incompatible sizes";
128 if (x.hasFastAccess())
130 dst.fastAccessDx(
i) -= x.fastAccessDx(
i);
133 dst.fastAccessDx(
i) -= x.dx(
i);
136 dst.resizeAndZero(xsz);
137 if (x.hasFastAccess())
139 dst.fastAccessDx(
i) = -x.fastAccessDx(
i);
142 dst.fastAccessDx(
i) = -x.dx(
i);
146 dst.val() -= x.val();
150 template <
typename SrcType>
154 const int xsz = x.size(), sz = dst.size();
158 #if defined(SACADO_DEBUG) && !defined(__CUDA_ARCH__ ) && !defined(__HIP_DEVICE_COMPILE__ )
159 if ((xsz != sz) && (xsz != 0) && (sz != 0))
160 throw "Fad Error: Attempt to assign with incompatible sizes";
165 if (x.hasFastAccess())
167 dst.fastAccessDx(
i) = v*x.fastAccessDx(
i) + dst.fastAccessDx(
i)*xval;
170 dst.fastAccessDx(
i) = v*x.dx(
i) + dst.fastAccessDx(
i)*xval;
173 dst.resizeAndZero(xsz);
174 if (x.hasFastAccess())
176 dst.fastAccessDx(
i) = v*x.fastAccessDx(
i);
179 dst.fastAccessDx(
i) = v*x.dx(
i);
185 dst.fastAccessDx(
i) *= xval;
193 template <
typename SrcType>
197 const int xsz = x.size(), sz = dst.size();
201 #if defined(SACADO_DEBUG) && !defined(__CUDA_ARCH__ ) && !defined(__HIP_DEVICE_COMPILE__ )
202 if ((xsz != sz) && (xsz != 0) && (sz != 0))
203 throw "Fad Error: Attempt to assign with incompatible sizes";
209 if (x.hasFastAccess())
211 dst.fastAccessDx(
i) =
212 ( dst.fastAccessDx(
i)*xval - v*x.fastAccessDx(
i) ) / xval2;
215 dst.fastAccessDx(
i) =
216 ( dst.fastAccessDx(
i)*xval - v*x.dx(
i) ) / xval2;
219 dst.resizeAndZero(xsz);
220 if (x.hasFastAccess())
222 dst.fastAccessDx(
i) = - v*x.fastAccessDx(
i) / xval2;
225 dst.fastAccessDx(
i) = -v*x.dx(
i) / xval2;
231 dst.fastAccessDx(
i) /= xval;
246 template <
typename DstType>
248 typename std::enable_if<Sacado::IsStaticallySized<DstType>::value>::type> {
255 template <
typename SrcType>
259 const int sz = dst.size();
261 dst.fastAccessDx(
i) = x.fastAccessDx(
i);
266 template <
typename SrcType>
270 const int sz = dst.size();
272 dst.fastAccessDx(
i) += x.fastAccessDx(
i);
273 dst.val() += x.val();
277 template <
typename SrcType>
281 const int sz = dst.size();
283 dst.fastAccessDx(
i) -= x.fastAccessDx(
i);
284 dst.val() -= x.val();
288 template <
typename SrcType>
292 const int sz = dst.size();
296 dst.fastAccessDx(
i) = v*x.fastAccessDx(
i) + dst.fastAccessDx(
i)*xval;
301 template <
typename SrcType>
305 const int sz = dst.size();
310 dst.fastAccessDx(
i) =
311 ( dst.fastAccessDx(
i)*xval - v*x.fastAccessDx(
i) )/ xval2;
322 #endif // SACADO_FAD_EXP_EXPRASSIGN_HPP
Class that implements various forms of expression assignments.
static SACADO_INLINE_FUNCTION void assign_divide_equal(DstType &dst, const SrcType &x)
Implementation of dst /= x.
DstType::value_type value_type
Typename of values.
static SACADO_INLINE_FUNCTION void assign_times_equal(DstType &dst, const SrcType &x)
Implementation of dst *= x.
static SACADO_INLINE_FUNCTION void assign_equal(DstType &dst, const SrcType &x)
Implementation of dst = x.
static SACADO_INLINE_FUNCTION void assign_equal(DstType &dst, const SrcType &x)
Implementation of dst = x.
static SACADO_INLINE_FUNCTION void assign_divide_equal(DstType &dst, const SrcType &x)
Implementation of dst /= x.
#define SACADO_FAD_DERIV_LOOP(I, SZ)
static SACADO_INLINE_FUNCTION void assign_times_equal(DstType &dst, const SrcType &x)
Implementation of dst *= x.
DstType::value_type value_type
Typename of values.
static SACADO_INLINE_FUNCTION void assign_minus_equal(DstType &dst, const SrcType &x)
Implementation of dst -= x.
static SACADO_INLINE_FUNCTION void assign_minus_equal(DstType &dst, const SrcType &x)
Implementation of dst -= x.
static SACADO_INLINE_FUNCTION void assign_plus_equal(DstType &dst, const SrcType &x)
Implementation of dst += x.
#define SACADO_INLINE_FUNCTION
static SACADO_INLINE_FUNCTION void assign_plus_equal(DstType &dst, const SrcType &x)
Implementation of dst += x.