gtsam 4.1.1
gtsam
Chebyshev2.h
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
34#pragma once
35
36#include <gtsam/base/Manifold.h>
38#include <gtsam/basis/Basis.h>
39
40#include <boost/function.hpp>
41
42namespace gtsam {
43
49class GTSAM_EXPORT Chebyshev2 : public Basis<Chebyshev2> {
50 public:
51 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
52
53 using Base = Basis<Chebyshev2>;
54 using Parameters = Eigen::Matrix<double, /*Nx1*/ -1, 1>;
55 using DiffMatrix = Eigen::Matrix<double, /*NxN*/ -1, -1>;
56
58 static double Point(size_t N, int j) {
59 assert(j >= 0 && size_t(j) < N);
60 const double dtheta = M_PI / (N > 1 ? (N - 1) : 1);
61 // We add -PI so that we get values ordered from -1 to +1
62 // sin(- M_PI_2 + dtheta*j); also works
63 return cos(-M_PI + dtheta * j);
64 }
65
67 static double Point(size_t N, int j, double a, double b) {
68 assert(j >= 0 && size_t(j) < N);
69 const double dtheta = M_PI / (N - 1);
70 // We add -PI so that we get values ordered from -1 to +1
71 return a + (b - a) * (1. + cos(-M_PI + dtheta * j)) / 2;
72 }
73
75 static Vector Points(size_t N) {
76 Vector points(N);
77 for (size_t j = 0; j < N; j++) points(j) = Point(N, j);
78 return points;
79 }
80
82 static Vector Points(size_t N, double a, double b) {
83 Vector points = Points(N);
84 const double T1 = (a + b) / 2, T2 = (b - a) / 2;
85 points = T1 + (T2 * points).array();
86 return points;
87 }
88
97 static Weights CalculateWeights(size_t N, double x, double a = -1,
98 double b = 1);
99
104 static Weights DerivativeWeights(size_t N, double x, double a = -1,
105 double b = 1);
106
111 static DiffMatrix DifferentiationMatrix(size_t N, double a = -1,
112 double b = 1);
113
132 static Weights IntegrationWeights(size_t N, double a = -1, double b = 1);
133
137 template <size_t M>
138 static Matrix matrix(boost::function<Eigen::Matrix<double, M, 1>(double)> f,
139 size_t N, double a = -1, double b = 1) {
140 Matrix Xmat(M, N);
141 for (size_t j = 0; j < N; j++) {
142 Xmat.col(j) = f(Point(N, j, a, b));
143 }
144 return Xmat;
145 }
146}; // \ Chebyshev2
147
148} // namespace gtsam
Special class for optional Jacobian arguments.
Base class and basic functions for Manifold types.
Compute an interpolating basis.
Global functions in a separate testing namespace.
Definition: chartTesting.h:28
CRTP Base class for function bases.
Definition: Basis.h:95
Chebyshev Interpolation on Chebyshev points of the second kind Note that N here, the number of points...
Definition: Chebyshev2.h:49
static Matrix matrix(boost::function< Eigen::Matrix< double, M, 1 >(double)> f, size_t N, double a=-1, double b=1)
Create matrix of values at Chebyshev points given vector-valued function.
Definition: Chebyshev2.h:138
static double Point(size_t N, int j)
Specific Chebyshev point.
Definition: Chebyshev2.h:58
static double Point(size_t N, int j, double a, double b)
Specific Chebyshev point, within [a,b] interval.
Definition: Chebyshev2.h:67
static Vector Points(size_t N)
All Chebyshev points.
Definition: Chebyshev2.h:75
static Vector Points(size_t N, double a, double b)
All Chebyshev points, within [a,b] interval.
Definition: Chebyshev2.h:82