Chebyshev2.cpp
Go to the documentation of this file.
1 /* ----------------------------------------------------------------------------
2 
3  * GTSAM Copyright 2010, Georgia Tech Research Corporation,
4  * Atlanta, Georgia 30332-0415
5  * All Rights Reserved
6  * Authors: Frank Dellaert, et al. (see THANKS for the full author list)
7 
8  * See LICENSE for the license information
9 
10  * -------------------------------------------------------------------------- */
11 
19 #include <gtsam/basis/Chebyshev2.h>
20 
21 namespace gtsam {
22 
23 double Chebyshev2::Point(size_t N, int j, double a, double b) {
24  assert(j >= 0 && size_t(j) < N);
25  const double dtheta = M_PI / (N > 1 ? (N - 1) : 1);
26  // We add -PI so that we get values ordered from -1 to +1
27  // sin(-M_PI_2 + dtheta*j); also works
28  return a + (b - a) * (1. + cos(-M_PI + dtheta * j)) / 2;
29 }
30 
31 Weights Chebyshev2::CalculateWeights(size_t N, double x, double a, double b) {
32  // Allocate space for weights
33  Weights weights(N);
34 
35  // We start by getting distances from x to all Chebyshev points
36  // as well as getting smallest distance
37  Weights distances(N);
38 
39  for (size_t j = 0; j < N; j++) {
40  const double dj =
41  x - Point(N, j, a, b); // only thing that depends on [a,b]
42 
43  if (std::abs(dj) < 1e-12) {
44  // exceptional case: x coincides with a Chebyshev point
45  weights.setZero();
46  weights(j) = 1;
47  return weights;
48  }
49  distances(j) = dj;
50  }
51 
52  // Beginning of interval, j = 0, x(0) = a
53  weights(0) = 0.5 / distances(0);
54 
55  // All intermediate points j=1:N-2
56  double d = weights(0), s = -1; // changes sign s at every iteration
57  for (size_t j = 1; j < N - 1; j++, s = -s) {
58  weights(j) = s / distances(j);
59  d += weights(j);
60  }
61 
62  // End of interval, j = N-1, x(N-1) = b
63  weights(N - 1) = 0.5 * s / distances(N - 1);
64  d += weights(N - 1);
65 
66  // normalize
67  return weights / d;
68 }
69 
70 Weights Chebyshev2::DerivativeWeights(size_t N, double x, double a, double b) {
71  // Allocate space for weights
72  Weights weightDerivatives(N);
73 
74  // toggle variable so we don't need to use `pow` for -1
75  double t = -1;
76 
77  // We start by getting distances from x to all Chebyshev points
78  // as well as getting smallest distance
79  Weights distances(N);
80 
81  for (size_t j = 0; j < N; j++) {
82  const double dj =
83  x - Point(N, j, a, b); // only thing that depends on [a,b]
84  if (std::abs(dj) < 1e-12) {
85  // exceptional case: x coincides with a Chebyshev point
86  weightDerivatives.setZero();
87  // compute the jth row of the differentiation matrix for this point
88  double cj = (j == 0 || j == N - 1) ? 2. : 1.;
89  for (size_t k = 0; k < N; k++) {
90  if (j == 0 && k == 0) {
91  // we reverse the sign since we order the cheb points from -1 to 1
92  weightDerivatives(k) = -(cj * (N - 1) * (N - 1) + 1) / 6.0;
93  } else if (j == N - 1 && k == N - 1) {
94  // we reverse the sign since we order the cheb points from -1 to 1
95  weightDerivatives(k) = (cj * (N - 1) * (N - 1) + 1) / 6.0;
96  } else if (k == j) {
97  double xj = Point(N, j);
98  double xj2 = xj * xj;
99  weightDerivatives(k) = -0.5 * xj / (1 - xj2);
100  } else {
101  double xj = Point(N, j);
102  double xk = Point(N, k);
103  double ck = (k == 0 || k == N - 1) ? 2. : 1.;
104  t = ((j + k) % 2) == 0 ? 1 : -1;
105  weightDerivatives(k) = (cj / ck) * t / (xj - xk);
106  }
107  }
108  return 2 * weightDerivatives / (b - a);
109  }
110  distances(j) = dj;
111  }
112 
113  // This section of code computes the derivative of
114  // the Barycentric Interpolation weights formula by applying
115  // the chain rule on the original formula.
116 
117  // g and k are multiplier terms which represent the derivatives of
118  // the numerator and denominator
119  double g = 0, k = 0;
120  double w = 1;
121 
122  for (size_t j = 0; j < N; j++) {
123  if (j == 0 || j == N - 1) {
124  w = 0.5;
125  } else {
126  w = 1.0;
127  }
128 
129  t = (j % 2 == 0) ? 1 : -1;
130 
131  double c = t / distances(j);
132  g += w * c;
133  k += (w * c / distances(j));
134  }
135 
136  double s = 1; // changes sign s at every iteration
137  double g2 = g * g;
138 
139  for (size_t j = 0; j < N; j++, s = -s) {
140  // Beginning of interval, j = 0, x0 = -1.0 and end of interval, j = N-1,
141  // x0 = 1.0
142  if (j == 0 || j == N - 1) {
143  w = 0.5;
144  } else {
145  // All intermediate points j=1:N-2
146  w = 1.0;
147  }
148  weightDerivatives(j) = (w * -s / (g * distances(j) * distances(j))) -
149  (w * -s * k / (g2 * distances(j)));
150  }
151 
152  return weightDerivatives;
153 }
154 
156  double b) {
157  DiffMatrix D(N, N);
158  if (N == 1) {
159  D(0, 0) = 1;
160  return D;
161  }
162 
163  // toggle variable so we don't need to use `pow` for -1
164  double t = -1;
165 
166  for (size_t i = 0; i < N; i++) {
167  double xi = Point(N, i);
168  double ci = (i == 0 || i == N - 1) ? 2. : 1.;
169  for (size_t j = 0; j < N; j++) {
170  if (i == 0 && j == 0) {
171  // we reverse the sign since we order the cheb points from -1 to 1
172  D(i, j) = -(ci * (N - 1) * (N - 1) + 1) / 6.0;
173  } else if (i == N - 1 && j == N - 1) {
174  // we reverse the sign since we order the cheb points from -1 to 1
175  D(i, j) = (ci * (N - 1) * (N - 1) + 1) / 6.0;
176  } else if (i == j) {
177  double xi2 = xi * xi;
178  D(i, j) = -xi / (2 * (1 - xi2));
179  } else {
180  double xj = Point(N, j);
181  double cj = (j == 0 || j == N - 1) ? 2. : 1.;
182  t = ((i + j) % 2) == 0 ? 1 : -1;
183  D(i, j) = (ci / cj) * t / (xi - xj);
184  }
185  }
186  }
187  // scale the matrix to the range
188  return D / ((b - a) / 2.0);
189 }
190 
191 Weights Chebyshev2::IntegrationWeights(size_t N, double a, double b) {
192  // Allocate space for weights
193  Weights weights(N);
194  size_t K = N - 1, // number of intervals between N points
195  K2 = K * K;
196  weights(0) = 0.5 * (b - a) / (K2 + K % 2 - 1);
197  weights(N - 1) = weights(0);
198 
199  size_t last_k = K / 2 + K % 2 - 1;
200 
201  for (size_t i = 1; i <= N - 2; ++i) {
202  double theta = i * M_PI / K;
203  weights(i) = (K % 2 == 0) ? 1 - cos(K * theta) / (K2 - 1) : 1;
204 
205  for (size_t k = 1; k <= last_k; ++k)
206  weights(i) -= 2 * cos(2 * k * theta) / (4 * k * k - 1);
207  weights(i) *= (b - a) / K;
208  }
209 
210  return weights;
211 }
212 
213 } // namespace gtsam
w
RowVector3d w
Definition: Matrix_resize_int.cpp:3
Chebyshev2.h
Pseudo-spectral parameterization for Chebyshev polynomials of the second kind.
D
MatrixXcd D
Definition: EigenSolver_EigenSolver_MatrixType.cpp:14
s
RealScalar s
Definition: level1_cplx_impl.h:126
e
Array< double, 1, 3 > e(1./3., 0.5, 2.)
d
static const double d[K][N]
Definition: igam.h:11
screwPose2::xi
Vector xi
Definition: testPose2.cpp:148
c
Scalar Scalar * c
Definition: benchVecAdd.cpp:17
x
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
Definition: gnuplot_common_settings.hh:12
K2
static const Cal3Bundler K2(625, 1e-3, 1e-3)
ceres::cos
Jet< T, N > cos(const Jet< T, N > &f)
Definition: jet.h:426
j
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2
gtsam::Chebyshev2::DifferentiationMatrix
static DiffMatrix DifferentiationMatrix(size_t N, double a=-1, double b=1)
Definition: Chebyshev2.cpp:155
g2
Pose3 g2(g1.expmap(h *V1_g1))
gtsam::Chebyshev2::Point
static double Point(size_t N, int j, double a=-1, double b=1)
Specific Chebyshev point, within [a,b] interval. Default interval is [-1, 1].
Definition: Chebyshev2.cpp:23
g
void g(const string &key, int i)
Definition: testBTree.cpp:41
gtsam::b
const G & b
Definition: Group.h:79
gtsam::Chebyshev2::CalculateWeights
static Weights CalculateWeights(size_t N, double x, double a=-1, double b=1)
Definition: Chebyshev2.cpp:31
a
ArrayXXi a
Definition: Array_initializer_list_23_cxx11.cpp:1
gtsam
traits
Definition: SFMdata.h:40
K
#define K
Definition: igam.h:8
gtsam::Chebyshev2::IntegrationWeights
static Weights IntegrationWeights(size_t N, double a=-1, double b=1)
Definition: Chebyshev2.cpp:191
Eigen::Matrix< double, 1, -1 >
abs
#define abs(x)
Definition: datatypes.h:17
N
#define N
Definition: igam.h:9
Eigen::PlainObjectBase::setZero
EIGEN_DEVICE_FUNC Derived & setZero(Index size)
Definition: CwiseNullaryOp.h:562
M_PI
#define M_PI
Definition: mconf.h:117
align_3::t
Point2 t(10, 10)
i
int i
Definition: BiCGSTAB_step_by_step.cpp:9
gtsam::Chebyshev2::DerivativeWeights
static Weights DerivativeWeights(size_t N, double x, double a=-1, double b=1)
Definition: Chebyshev2.cpp:70


gtsam
Author(s):
autogenerated on Sun Dec 22 2024 04:11:14