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_ExprAssign.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Sacado Package
4 //
5 // Copyright 2006 NTESS and the Sacado contributors.
6 // SPDX-License-Identifier: LGPL-2.1-or-later
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef SACADO_FAD_EXP_EXPRASSIGN_HPP
11 #define SACADO_FAD_EXP_EXPRASSIGN_HPP
12 
13 namespace Sacado {
14 
15  namespace Fad {
16  namespace Exp {
17 
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)
21 #else
22 #define SACADO_FAD_DERIV_LOOP(I,SZ) for (int I=0; I<SZ; ++I)
23 #endif
24 #endif
25 
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)
29 #else
30 #define SACADO_FAD_THREAD_SINGLE /* */
31 #endif
32 #endif
33 
35 
42  template <typename DstType, typename Enabled = void>
43  class ExprAssign {
44  public:
45 
47  typedef typename DstType::value_type value_type;
48 
50  template <typename SrcType>
52  static void assign_equal(DstType& dst, const SrcType& x)
53  {
54  const int xsz = x.size();
55 
56  if (xsz != dst.size())
57  dst.resizeAndZero(xsz);
58 
59  const int sz = dst.size();
60 
61  // For ViewStorage, the resize above may not in fact resize the
62  // derivative array, so it is possible that sz != xsz at this point.
63  // The only valid use case here is sz > xsz == 0, so we use sz in the
64  // assignment below
65 
66  if (sz) {
67  if (x.hasFastAccess()) {
69  dst.fastAccessDx(i) = x.fastAccessDx(i);
70  }
71  else
73  dst.fastAccessDx(i) = x.dx(i);
74  }
75 
76  dst.val() = x.val();
77  }
78 
80  template <typename SrcType>
82  static void assign_plus_equal(DstType& dst, const SrcType& x)
83  {
84  const int xsz = x.size(), sz = dst.size();
85 
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";
89 #endif
90 
91  if (xsz) {
92  if (sz) {
93  if (x.hasFastAccess())
95  dst.fastAccessDx(i) += x.fastAccessDx(i);
96  else
97  for (int i=0; i<sz; ++i)
98  dst.fastAccessDx(i) += x.dx(i);
99  }
100  else {
101  dst.resizeAndZero(xsz);
102  if (x.hasFastAccess())
104  dst.fastAccessDx(i) = x.fastAccessDx(i);
105  else
107  dst.fastAccessDx(i) = x.dx(i);
108  }
109  }
110 
111  dst.val() += x.val();
112  }
113 
115  template <typename SrcType>
117  static void assign_minus_equal(DstType& dst, const SrcType& x)
118  {
119  const int xsz = x.size(), sz = dst.size();
120 
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";
124 #endif
125 
126  if (xsz) {
127  if (sz) {
128  if (x.hasFastAccess())
130  dst.fastAccessDx(i) -= x.fastAccessDx(i);
131  else
133  dst.fastAccessDx(i) -= x.dx(i);
134  }
135  else {
136  dst.resizeAndZero(xsz);
137  if (x.hasFastAccess())
139  dst.fastAccessDx(i) = -x.fastAccessDx(i);
140  else
142  dst.fastAccessDx(i) = -x.dx(i);
143  }
144  }
145 
146  dst.val() -= x.val();
147  }
148 
150  template <typename SrcType>
152  static void assign_times_equal(DstType& dst, const SrcType& x)
153  {
154  const int xsz = x.size(), sz = dst.size();
155  const value_type xval = x.val();
156  const value_type v = dst.val();
157 
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";
161 #endif
162 
163  if (xsz) {
164  if (sz) {
165  if (x.hasFastAccess())
167  dst.fastAccessDx(i) = v*x.fastAccessDx(i) + dst.fastAccessDx(i)*xval;
168  else
170  dst.fastAccessDx(i) = v*x.dx(i) + dst.fastAccessDx(i)*xval;
171  }
172  else {
173  dst.resizeAndZero(xsz);
174  if (x.hasFastAccess())
176  dst.fastAccessDx(i) = v*x.fastAccessDx(i);
177  else
179  dst.fastAccessDx(i) = v*x.dx(i);
180  }
181  }
182  else {
183  if (sz) {
185  dst.fastAccessDx(i) *= xval;
186  }
187  }
188 
189  dst.val() *= xval;
190  }
191 
193  template <typename SrcType>
195  static void assign_divide_equal(DstType& dst, const SrcType& x)
196  {
197  const int xsz = x.size(), sz = dst.size();
198  const value_type xval = x.val();
199  const value_type v = dst.val();
200 
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";
204 #endif
205 
206  if (xsz) {
207  const value_type xval2 = xval*xval;
208  if (sz) {
209  if (x.hasFastAccess())
211  dst.fastAccessDx(i) =
212  ( dst.fastAccessDx(i)*xval - v*x.fastAccessDx(i) ) / xval2;
213  else
215  dst.fastAccessDx(i) =
216  ( dst.fastAccessDx(i)*xval - v*x.dx(i) ) / xval2;
217  }
218  else {
219  dst.resizeAndZero(xsz);
220  if (x.hasFastAccess())
222  dst.fastAccessDx(i) = - v*x.fastAccessDx(i) / xval2;
223  else
225  dst.fastAccessDx(i) = -v*x.dx(i) / xval2;
226  }
227  }
228  else {
229  if (sz) {
231  dst.fastAccessDx(i) /= xval;
232  }
233  }
234 
235  dst.val() /= xval;
236  }
237 
238  };
239 
241 
246  template <typename DstType>
247  class ExprAssign<DstType,
248  typename std::enable_if<Sacado::IsStaticallySized<DstType>::value>::type> {
249  public:
250 
252  typedef typename DstType::value_type value_type;
253 
255  template <typename SrcType>
257  static void assign_equal(DstType& dst, const SrcType& x)
258  {
259  const int sz = dst.size();
261  dst.fastAccessDx(i) = x.fastAccessDx(i);
262  dst.val() = x.val();
263  }
264 
266  template <typename SrcType>
268  static void assign_plus_equal(DstType& dst, const SrcType& x)
269  {
270  const int sz = dst.size();
272  dst.fastAccessDx(i) += x.fastAccessDx(i);
273  dst.val() += x.val();
274  }
275 
277  template <typename SrcType>
279  static void assign_minus_equal(DstType& dst, const SrcType& x)
280  {
281  const int sz = dst.size();
283  dst.fastAccessDx(i) -= x.fastAccessDx(i);
284  dst.val() -= x.val();
285  }
286 
288  template <typename SrcType>
290  static void assign_times_equal(DstType& dst, const SrcType& x)
291  {
292  const int sz = dst.size();
293  const value_type xval = x.val();
294  const value_type v = dst.val();
296  dst.fastAccessDx(i) = v*x.fastAccessDx(i) + dst.fastAccessDx(i)*xval;
297  dst.val() *= xval;
298  }
299 
301  template <typename SrcType>
303  static void assign_divide_equal(DstType& dst, const SrcType& x)
304  {
305  const int sz = dst.size();
306  const value_type xval = x.val();
307  const value_type xval2 = xval*xval;
308  const value_type v = dst.val();
310  dst.fastAccessDx(i) =
311  ( dst.fastAccessDx(i)*xval - v*x.fastAccessDx(i) )/ xval2;
312  dst.val() /= xval;
313  }
314 
315  };
316 
317  } // namespace Exp
318  } // namespace Fad
319 
320 } // namespace Sacado
321 
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.
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.