forked from RobotLocomotion/drake
-
Notifications
You must be signed in to change notification settings - Fork 0
/
autodiffxd.h
423 lines (359 loc) · 15.7 KB
/
autodiffxd.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
#pragma once
// This file is a modification of Eigen-3.3.3's AutoDiffScalar.h file which is
// available at
// https://bitbucket.org/eigen/eigen/raw/67e894c6cd8f5f1f604b27d37ed47fdf012674ff/unsupported/Eigen/src/AutoDiff/AutoDiffScalar.h
//
// Copyright (C) 2009 Gael Guennebaud <[email protected]>
// Copyright (C) 2017 Drake Authors
//
// This Source Code Form is subject to the terms of the Mozilla
// Public License v. 2.0. If a copy of the MPL was not distributed
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
#ifndef DRAKE_COMMON_AUTODIFF_HEADER
// TODO(soonho-tri): Change to #error.
#warning Do not directly include this file. Include "drake/common/autodiff.h".
#endif
#include <cmath>
#include <ostream>
#include <Eigen/Dense>
namespace Eigen {
#if !defined(DRAKE_DOXYGEN_CXX)
// Explicit template specializations of Eigen::AutoDiffScalar for VectorXd.
//
// AutoDiffScalar tries to call internal::make_coherent to promote empty
// derivatives. However, it fails to do the promotion when an operand is an
// expression tree (i.e. CwiseBinaryOp). Our solution is to provide special
// overloading for VectorXd and change the return types of its operators. With
// this change, the operators evaluate terms immediately and return an
// AutoDiffScalar<VectorXd> instead of expression trees (such as CwiseBinaryOp).
// Eigen's implementation of internal::make_coherent makes use of const_cast in
// order to promote zero sized derivatives. This however interferes badly with
// our caching system and produces unexpected behaviors. See #10971 for details.
// Therefore our implementation stops using internal::make_coherent and treats
// scalars with zero sized derivatives as constants, as it should.
//
// We also provide overloading of math functions for AutoDiffScalar<VectorXd>
// which return AutoDiffScalar<VectorXd> instead of an expression tree.
//
// See https://github.com/RobotLocomotion/drake/issues/6944 for more
// information. See also drake/common/autodiff_overloads.h.
//
// TODO(soonho-tri): Next time when we upgrade Eigen, please check if we still
// need these specializations.
template <>
class AutoDiffScalar<VectorXd>
: public internal::auto_diff_special_op<VectorXd, false> {
public:
typedef internal::auto_diff_special_op<VectorXd, false> Base;
typedef typename internal::remove_all<VectorXd>::type DerType;
typedef typename internal::traits<DerType>::Scalar Scalar;
typedef typename NumTraits<Scalar>::Real Real;
using Base::operator+;
using Base::operator*;
AutoDiffScalar() {}
AutoDiffScalar(const Scalar& value, int nbDer, int derNumber)
: m_value(value), m_derivatives(DerType::Zero(nbDer)) {
m_derivatives.coeffRef(derNumber) = Scalar(1);
}
// NOLINTNEXTLINE(runtime/explicit): Code from Eigen.
AutoDiffScalar(const Real& value) : m_value(value) {
if (m_derivatives.size() > 0) m_derivatives.setZero();
}
AutoDiffScalar(const Scalar& value, const DerType& der)
: m_value(value), m_derivatives(der) {}
template <typename OtherDerType>
AutoDiffScalar(
const AutoDiffScalar<OtherDerType>& other
#ifndef EIGEN_PARSED_BY_DOXYGEN
,
typename internal::enable_if<
internal::is_same<
Scalar, typename internal::traits<typename internal::remove_all<
OtherDerType>::type>::Scalar>::value,
void*>::type = 0
#endif
)
: m_value(other.value()), m_derivatives(other.derivatives()) {
}
friend std::ostream& operator<<(std::ostream& s, const AutoDiffScalar& a) {
return s << a.value();
}
AutoDiffScalar(const AutoDiffScalar& other)
: m_value(other.value()), m_derivatives(other.derivatives()) {}
template <typename OtherDerType>
inline AutoDiffScalar& operator=(const AutoDiffScalar<OtherDerType>& other) {
m_value = other.value();
m_derivatives = other.derivatives();
return *this;
}
inline AutoDiffScalar& operator=(const AutoDiffScalar& other) {
m_value = other.value();
m_derivatives = other.derivatives();
return *this;
}
inline AutoDiffScalar& operator=(const Scalar& other) {
m_value = other;
if (m_derivatives.size() > 0) m_derivatives.setZero();
return *this;
}
inline const Scalar& value() const { return m_value; }
inline Scalar& value() { return m_value; }
inline const DerType& derivatives() const { return m_derivatives; }
inline DerType& derivatives() { return m_derivatives; }
inline bool operator<(const Scalar& other) const { return m_value < other; }
inline bool operator<=(const Scalar& other) const { return m_value <= other; }
inline bool operator>(const Scalar& other) const { return m_value > other; }
inline bool operator>=(const Scalar& other) const { return m_value >= other; }
inline bool operator==(const Scalar& other) const { return m_value == other; }
inline bool operator!=(const Scalar& other) const { return m_value != other; }
friend inline bool operator<(const Scalar& a, const AutoDiffScalar& b) {
return a < b.value();
}
friend inline bool operator<=(const Scalar& a, const AutoDiffScalar& b) {
return a <= b.value();
}
friend inline bool operator>(const Scalar& a, const AutoDiffScalar& b) {
return a > b.value();
}
friend inline bool operator>=(const Scalar& a, const AutoDiffScalar& b) {
return a >= b.value();
}
friend inline bool operator==(const Scalar& a, const AutoDiffScalar& b) {
return a == b.value();
}
friend inline bool operator!=(const Scalar& a, const AutoDiffScalar& b) {
return a != b.value();
}
template <typename OtherDerType>
inline bool operator<(const AutoDiffScalar<OtherDerType>& b) const {
return m_value < b.value();
}
template <typename OtherDerType>
inline bool operator<=(const AutoDiffScalar<OtherDerType>& b) const {
return m_value <= b.value();
}
template <typename OtherDerType>
inline bool operator>(const AutoDiffScalar<OtherDerType>& b) const {
return m_value > b.value();
}
template <typename OtherDerType>
inline bool operator>=(const AutoDiffScalar<OtherDerType>& b) const {
return m_value >= b.value();
}
template <typename OtherDerType>
inline bool operator==(const AutoDiffScalar<OtherDerType>& b) const {
return m_value == b.value();
}
template <typename OtherDerType>
inline bool operator!=(const AutoDiffScalar<OtherDerType>& b) const {
return m_value != b.value();
}
inline const AutoDiffScalar<DerType> operator+(const Scalar& other) const {
return AutoDiffScalar<DerType>(m_value + other, m_derivatives);
}
friend inline const AutoDiffScalar<DerType> operator+(
const Scalar& a, const AutoDiffScalar& b) {
return AutoDiffScalar<DerType>(a + b.value(), b.derivatives());
}
inline AutoDiffScalar& operator+=(const Scalar& other) {
value() += other;
return *this;
}
template <typename OtherDerType>
inline const AutoDiffScalar<DerType> operator+(
const AutoDiffScalar<OtherDerType>& other) const {
const bool has_this_der = m_derivatives.size() > 0;
const bool has_both_der = has_this_der && (other.derivatives().size() > 0);
return MakeAutoDiffScalar(
m_value + other.value(),
has_both_der
? VectorXd(m_derivatives + other.derivatives())
: has_this_der ? m_derivatives : VectorXd(other.derivatives()));
}
template <typename OtherDerType>
inline AutoDiffScalar& operator+=(const AutoDiffScalar<OtherDerType>& other) {
(*this) = (*this) + other;
return *this;
}
inline const AutoDiffScalar<DerType> operator-(const Scalar& b) const {
return AutoDiffScalar<DerType>(m_value - b, m_derivatives);
}
friend inline const AutoDiffScalar<DerType> operator-(
const Scalar& a, const AutoDiffScalar& b) {
return AutoDiffScalar<DerType>(a - b.value(), -b.derivatives());
}
inline AutoDiffScalar& operator-=(const Scalar& other) {
value() -= other;
return *this;
}
template <typename OtherDerType>
inline const AutoDiffScalar<DerType> operator-(
const AutoDiffScalar<OtherDerType>& other) const {
const bool has_this_der = m_derivatives.size() > 0;
const bool has_both_der = has_this_der && (other.derivatives().size() > 0);
return MakeAutoDiffScalar(
m_value - other.value(),
has_both_der
? VectorXd(m_derivatives - other.derivatives())
: has_this_der ? m_derivatives : VectorXd(-other.derivatives()));
}
template <typename OtherDerType>
inline AutoDiffScalar& operator-=(const AutoDiffScalar<OtherDerType>& other) {
*this = *this - other;
return *this;
}
inline const AutoDiffScalar<DerType> operator-() const {
return AutoDiffScalar<DerType>(-m_value, -m_derivatives);
}
inline const AutoDiffScalar<DerType> operator*(const Scalar& other) const {
return MakeAutoDiffScalar(m_value * other, m_derivatives * other);
}
friend inline const AutoDiffScalar<DerType> operator*(
const Scalar& other, const AutoDiffScalar& a) {
return MakeAutoDiffScalar(a.value() * other, a.derivatives() * other);
}
inline const AutoDiffScalar<DerType> operator/(const Scalar& other) const {
return MakeAutoDiffScalar(m_value / other,
(m_derivatives * (Scalar(1) / other)));
}
friend inline const AutoDiffScalar<DerType> operator/(
const Scalar& other, const AutoDiffScalar& a) {
return MakeAutoDiffScalar(
other / a.value(),
a.derivatives() * (Scalar(-other) / (a.value() * a.value())));
}
template <typename OtherDerType>
inline const AutoDiffScalar<DerType> operator/(
const AutoDiffScalar<OtherDerType>& other) const {
const auto& this_der = m_derivatives;
const auto& other_der = other.derivatives();
const bool has_this_der = m_derivatives.size() > 0;
const bool has_both_der = has_this_der && (other.derivatives().size() > 0);
const double scale = 1. / (other.value() * other.value());
return MakeAutoDiffScalar(
m_value / other.value(),
has_both_der ?
VectorXd(this_der * other.value() - other_der * m_value) * scale :
has_this_der ?
VectorXd(this_der * other.value()) * scale :
// has_other_der || has_neither
VectorXd(other_der * -m_value) * scale);
}
template <typename OtherDerType>
inline const AutoDiffScalar<DerType> operator*(
const AutoDiffScalar<OtherDerType>& other) const {
const bool has_this_der = m_derivatives.size() > 0;
const bool has_both_der = has_this_der && (other.derivatives().size() > 0);
return MakeAutoDiffScalar(
m_value * other.value(),
has_both_der ? VectorXd(m_derivatives * other.value() +
other.derivatives() * m_value)
: has_this_der ? VectorXd(m_derivatives * other.value())
: VectorXd(other.derivatives() * m_value));
}
inline AutoDiffScalar& operator*=(const Scalar& other) {
*this = *this * other;
return *this;
}
template <typename OtherDerType>
inline AutoDiffScalar& operator*=(const AutoDiffScalar<OtherDerType>& other) {
*this = *this * other;
return *this;
}
inline AutoDiffScalar& operator/=(const Scalar& other) {
*this = *this / other;
return *this;
}
template <typename OtherDerType>
inline AutoDiffScalar& operator/=(const AutoDiffScalar<OtherDerType>& other) {
*this = *this / other;
return *this;
}
protected:
Scalar m_value;
DerType m_derivatives;
};
#define DRAKE_EIGEN_AUTODIFFXD_DECLARE_GLOBAL_UNARY(FUNC, CODE) \
inline const AutoDiffScalar<VectorXd> FUNC( \
const AutoDiffScalar<VectorXd>& x) { \
EIGEN_UNUSED typedef double Scalar; \
CODE; \
}
DRAKE_EIGEN_AUTODIFFXD_DECLARE_GLOBAL_UNARY(
abs, using std::abs; return Eigen::MakeAutoDiffScalar(
abs(x.value()), x.derivatives() * (x.value() < 0 ? -1 : 1));)
DRAKE_EIGEN_AUTODIFFXD_DECLARE_GLOBAL_UNARY(
abs2, using numext::abs2; return Eigen::MakeAutoDiffScalar(
abs2(x.value()), x.derivatives() * (Scalar(2) * x.value()));)
DRAKE_EIGEN_AUTODIFFXD_DECLARE_GLOBAL_UNARY(
sqrt, using std::sqrt; Scalar sqrtx = sqrt(x.value());
return Eigen::MakeAutoDiffScalar(sqrtx,
x.derivatives() * (Scalar(0.5) / sqrtx));)
DRAKE_EIGEN_AUTODIFFXD_DECLARE_GLOBAL_UNARY(
cos, using std::cos; using std::sin;
return Eigen::MakeAutoDiffScalar(cos(x.value()),
x.derivatives() * (-sin(x.value())));)
DRAKE_EIGEN_AUTODIFFXD_DECLARE_GLOBAL_UNARY(
sin, using std::sin; using std::cos;
return Eigen::MakeAutoDiffScalar(sin(x.value()),
x.derivatives() * cos(x.value()));)
DRAKE_EIGEN_AUTODIFFXD_DECLARE_GLOBAL_UNARY(
exp, using std::exp; Scalar expx = exp(x.value());
return Eigen::MakeAutoDiffScalar(expx, x.derivatives() * expx);)
DRAKE_EIGEN_AUTODIFFXD_DECLARE_GLOBAL_UNARY(
log, using std::log; return Eigen::MakeAutoDiffScalar(
log(x.value()), x.derivatives() * (Scalar(1) / x.value()));)
DRAKE_EIGEN_AUTODIFFXD_DECLARE_GLOBAL_UNARY(
tan, using std::tan; using std::cos; return Eigen::MakeAutoDiffScalar(
tan(x.value()),
x.derivatives() * (Scalar(1) / numext::abs2(cos(x.value()))));)
DRAKE_EIGEN_AUTODIFFXD_DECLARE_GLOBAL_UNARY(
asin, using std::sqrt; using std::asin; return Eigen::MakeAutoDiffScalar(
asin(x.value()),
x.derivatives() * (Scalar(1) / sqrt(1 - numext::abs2(x.value()))));)
DRAKE_EIGEN_AUTODIFFXD_DECLARE_GLOBAL_UNARY(
acos, using std::sqrt; using std::acos; return Eigen::MakeAutoDiffScalar(
acos(x.value()),
x.derivatives() * (Scalar(-1) / sqrt(1 - numext::abs2(x.value()))));)
DRAKE_EIGEN_AUTODIFFXD_DECLARE_GLOBAL_UNARY(
atan, using std::atan; return Eigen::MakeAutoDiffScalar(
atan(x.value()),
x.derivatives() * (Scalar(1) / (1 + x.value() * x.value())));)
DRAKE_EIGEN_AUTODIFFXD_DECLARE_GLOBAL_UNARY(
tanh, using std::cosh; using std::tanh; return Eigen::MakeAutoDiffScalar(
tanh(x.value()),
x.derivatives() * (Scalar(1) / numext::abs2(cosh(x.value()))));)
DRAKE_EIGEN_AUTODIFFXD_DECLARE_GLOBAL_UNARY(
sinh, using std::sinh; using std::cosh;
return Eigen::MakeAutoDiffScalar(sinh(x.value()),
x.derivatives() * cosh(x.value()));)
DRAKE_EIGEN_AUTODIFFXD_DECLARE_GLOBAL_UNARY(
cosh, using std::sinh; using std::cosh;
return Eigen::MakeAutoDiffScalar(cosh(x.value()),
x.derivatives() * sinh(x.value()));)
#undef DRAKE_EIGEN_AUTODIFFXD_DECLARE_GLOBAL_UNARY
// We have this specialization here because the Eigen-3.3.3's atan2
// implementation for AutoDiffScalar does not make a return with properly sized
// derivatives.
inline const AutoDiffScalar<VectorXd> atan2(const AutoDiffScalar<VectorXd>& a,
const AutoDiffScalar<VectorXd>& b) {
const bool has_a_der = a.derivatives().size() > 0;
const bool has_both_der = has_a_der && (b.derivatives().size() > 0);
const double squared_hypot = a.value() * a.value() + b.value() * b.value();
return MakeAutoDiffScalar(
std::atan2(a.value(), b.value()),
VectorXd((has_both_der
? VectorXd(a.derivatives() * b.value() -
a.value() * b.derivatives())
: has_a_der ? VectorXd(a.derivatives() * b.value())
: VectorXd(-a.value() * b.derivatives())) /
squared_hypot));
}
inline const AutoDiffScalar<VectorXd> pow(const AutoDiffScalar<VectorXd>& a,
double b) {
using std::pow;
return MakeAutoDiffScalar(pow(a.value(), b),
a.derivatives() * (b * pow(a.value(), b - 1)));
}
#endif
} // namespace Eigen