Sacado Package Browser (Single Doxygen Collection) Version of the Day
Loading...
Searching...
No Matches
Sacado_Fad_Exp_ExprAssign.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_EXPRASSIGN_HPP
31#define SACADO_FAD_EXP_EXPRASSIGN_HPP
32
33namespace Sacado {
34
35 namespace Fad {
36 namespace Exp {
37
38#ifndef SACADO_FAD_DERIV_LOOP
39#if defined(SACADO_VIEW_CUDA_HIERARCHICAL_DFAD) && !defined(SACADO_DISABLE_CUDA_IN_KOKKOS) && ( defined(__CUDA_ARCH__) || defined(__HIP_DEVICE_COMPILE__) )
40#define SACADO_FAD_DERIV_LOOP(I,SZ) for (int I=threadIdx.x; I<SZ; I+=blockDim.x)
41#else
42#define SACADO_FAD_DERIV_LOOP(I,SZ) for (int I=0; I<SZ; ++I)
43#endif
44#endif
45
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__) || defined(__HIP_DEVICE_COMPILE__) )
48#define SACADO_FAD_THREAD_SINGLE if (threadIdx.x == 0)
49#else
50#define SACADO_FAD_THREAD_SINGLE /* */
51#endif
52#endif
53
55
62 template <typename DstType, typename Enabled = void>
63 class ExprAssign {
64 public:
65
67 typedef typename DstType::value_type value_type;
68
70 template <typename SrcType>
72 static void assign_equal(DstType& dst, const SrcType& x)
73 {
74 const int xsz = x.size();
75
76 if (xsz != dst.size())
77 dst.resizeAndZero(xsz);
78
79 const int sz = dst.size();
80
81 // For ViewStorage, the resize above may not in fact resize the
82 // derivative array, so it is possible that sz != xsz at this point.
83 // The only valid use case here is sz > xsz == 0, so we use sz in the
84 // assignment below
85
86 if (sz) {
87 if (x.hasFastAccess()) {
89 dst.fastAccessDx(i) = x.fastAccessDx(i);
90 }
91 else
93 dst.fastAccessDx(i) = x.dx(i);
94 }
95
96 dst.val() = x.val();
97 }
98
100 template <typename SrcType>
102 static void assign_plus_equal(DstType& dst, const SrcType& x)
103 {
104 const int xsz = x.size(), sz = dst.size();
105
106#if defined(SACADO_DEBUG) && !defined(__CUDA_ARCH__ ) && !defined(__HIP_DEVICE_COMPILE__ )
107 if ((xsz != sz) && (xsz != 0) && (sz != 0))
108 throw "Fad Error: Attempt to assign with incompatible sizes";
109#endif
110
111 if (xsz) {
112 if (sz) {
113 if (x.hasFastAccess())
115 dst.fastAccessDx(i) += x.fastAccessDx(i);
116 else
117 for (int i=0; i<sz; ++i)
118 dst.fastAccessDx(i) += x.dx(i);
119 }
120 else {
121 dst.resizeAndZero(xsz);
122 if (x.hasFastAccess())
124 dst.fastAccessDx(i) = x.fastAccessDx(i);
125 else
127 dst.fastAccessDx(i) = x.dx(i);
128 }
129 }
130
131 dst.val() += x.val();
132 }
133
135 template <typename SrcType>
137 static void assign_minus_equal(DstType& dst, const SrcType& x)
138 {
139 const int xsz = x.size(), sz = dst.size();
140
141#if defined(SACADO_DEBUG) && !defined(__CUDA_ARCH__ ) && !defined(__HIP_DEVICE_COMPILE__ )
142 if ((xsz != sz) && (xsz != 0) && (sz != 0))
143 throw "Fad Error: Attempt to assign with incompatible sizes";
144#endif
145
146 if (xsz) {
147 if (sz) {
148 if (x.hasFastAccess())
150 dst.fastAccessDx(i) -= x.fastAccessDx(i);
151 else
153 dst.fastAccessDx(i) -= x.dx(i);
154 }
155 else {
156 dst.resizeAndZero(xsz);
157 if (x.hasFastAccess())
159 dst.fastAccessDx(i) = -x.fastAccessDx(i);
160 else
162 dst.fastAccessDx(i) = -x.dx(i);
163 }
164 }
165
166 dst.val() -= x.val();
167 }
168
170 template <typename SrcType>
172 static void assign_times_equal(DstType& dst, const SrcType& x)
173 {
174 const int xsz = x.size(), sz = dst.size();
175 const value_type xval = x.val();
176 const value_type v = dst.val();
177
178#if defined(SACADO_DEBUG) && !defined(__CUDA_ARCH__ ) && !defined(__HIP_DEVICE_COMPILE__ )
179 if ((xsz != sz) && (xsz != 0) && (sz != 0))
180 throw "Fad Error: Attempt to assign with incompatible sizes";
181#endif
182
183 if (xsz) {
184 if (sz) {
185 if (x.hasFastAccess())
187 dst.fastAccessDx(i) = v*x.fastAccessDx(i) + dst.fastAccessDx(i)*xval;
188 else
190 dst.fastAccessDx(i) = v*x.dx(i) + dst.fastAccessDx(i)*xval;
191 }
192 else {
193 dst.resizeAndZero(xsz);
194 if (x.hasFastAccess())
196 dst.fastAccessDx(i) = v*x.fastAccessDx(i);
197 else
199 dst.fastAccessDx(i) = v*x.dx(i);
200 }
201 }
202 else {
203 if (sz) {
205 dst.fastAccessDx(i) *= xval;
206 }
207 }
208
209 dst.val() *= xval;
210 }
211
213 template <typename SrcType>
215 static void assign_divide_equal(DstType& dst, const SrcType& x)
216 {
217 const int xsz = x.size(), sz = dst.size();
218 const value_type xval = x.val();
219 const value_type v = dst.val();
220
221#if defined(SACADO_DEBUG) && !defined(__CUDA_ARCH__ ) && !defined(__HIP_DEVICE_COMPILE__ )
222 if ((xsz != sz) && (xsz != 0) && (sz != 0))
223 throw "Fad Error: Attempt to assign with incompatible sizes";
224#endif
225
226 if (xsz) {
227 const value_type xval2 = xval*xval;
228 if (sz) {
229 if (x.hasFastAccess())
231 dst.fastAccessDx(i) =
232 ( dst.fastAccessDx(i)*xval - v*x.fastAccessDx(i) ) / xval2;
233 else
235 dst.fastAccessDx(i) =
236 ( dst.fastAccessDx(i)*xval - v*x.dx(i) ) / xval2;
237 }
238 else {
239 dst.resizeAndZero(xsz);
240 if (x.hasFastAccess())
242 dst.fastAccessDx(i) = - v*x.fastAccessDx(i) / xval2;
243 else
245 dst.fastAccessDx(i) = -v*x.dx(i) / xval2;
246 }
247 }
248 else {
249 if (sz) {
251 dst.fastAccessDx(i) /= xval;
252 }
253 }
254
255 dst.val() /= xval;
256 }
257
258 };
259
261
266 template <typename DstType>
267 class ExprAssign<DstType,
268 typename std::enable_if<Sacado::IsStaticallySized<DstType>::value>::type> {
269 public:
270
272 typedef typename DstType::value_type value_type;
273
275 template <typename SrcType>
277 static void assign_equal(DstType& dst, const SrcType& x)
278 {
279 const int sz = dst.size();
281 dst.fastAccessDx(i) = x.fastAccessDx(i);
282 dst.val() = x.val();
283 }
284
286 template <typename SrcType>
288 static void assign_plus_equal(DstType& dst, const SrcType& x)
289 {
290 const int sz = dst.size();
292 dst.fastAccessDx(i) += x.fastAccessDx(i);
293 dst.val() += x.val();
294 }
295
297 template <typename SrcType>
299 static void assign_minus_equal(DstType& dst, const SrcType& x)
300 {
301 const int sz = dst.size();
303 dst.fastAccessDx(i) -= x.fastAccessDx(i);
304 dst.val() -= x.val();
305 }
306
308 template <typename SrcType>
310 static void assign_times_equal(DstType& dst, const SrcType& x)
311 {
312 const int sz = dst.size();
313 const value_type xval = x.val();
314 const value_type v = dst.val();
316 dst.fastAccessDx(i) = v*x.fastAccessDx(i) + dst.fastAccessDx(i)*xval;
317 dst.val() *= xval;
318 }
319
321 template <typename SrcType>
323 static void assign_divide_equal(DstType& dst, const SrcType& x)
324 {
325 const int sz = dst.size();
326 const value_type xval = x.val();
327 const value_type xval2 = xval*xval;
328 const value_type v = dst.val();
330 dst.fastAccessDx(i) =
331 ( dst.fastAccessDx(i)*xval - v*x.fastAccessDx(i) )/ xval2;
332 dst.val() /= xval;
333 }
334
335 };
336
337 } // namespace Exp
338 } // namespace Fad
339
340} // namespace Sacado
341
342#endif // SACADO_FAD_EXP_EXPRASSIGN_HPP
#define SACADO_INLINE_FUNCTION
#define SACADO_FAD_DERIV_LOOP(I, SZ)
static SACADO_INLINE_FUNCTION void assign_plus_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_divide_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_times_equal(DstType &dst, const SrcType &x)
Implementation of dst *= x.
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_plus_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.
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.