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 
33 #pragma once
34 
35 #include <gtsam/base/Manifold.h>
37 #include <gtsam/basis/Basis.h>
38 
39 namespace gtsam {
40 
46 class GTSAM_EXPORT Chebyshev2 : public Basis<Chebyshev2> {
47  public:
49 
51  using Parameters = Eigen::Matrix<double, /*Nx1*/ -1, 1>;
52  using DiffMatrix = Eigen::Matrix<double, /*NxN*/ -1, -1>;
53 
55  static double Point(size_t N, int j) {
56  assert(j >= 0 && size_t(j) < N);
57  const double dtheta = M_PI / (N > 1 ? (N - 1) : 1);
58  // We add -PI so that we get values ordered from -1 to +1
59  // sin(- M_PI_2 + dtheta*j); also works
60  return cos(-M_PI + dtheta * j);
61  }
62 
64  static double Point(size_t N, int j, double a, double b) {
65  assert(j >= 0 && size_t(j) < N);
66  const double dtheta = M_PI / (N - 1);
67  // We add -PI so that we get values ordered from -1 to +1
68  return a + (b - a) * (1. + cos(-M_PI + dtheta * j)) / 2;
69  }
70 
72  static Vector Points(size_t N) {
73  Vector points(N);
74  for (size_t j = 0; j < N; j++) points(j) = Point(N, j);
75  return points;
76  }
77 
79  static Vector Points(size_t N, double a, double b) {
80  Vector points = Points(N);
81  const double T1 = (a + b) / 2, T2 = (b - a) / 2;
82  points = T1 + (T2 * points).array();
83  return points;
84  }
85 
94  static Weights CalculateWeights(size_t N, double x, double a = -1,
95  double b = 1);
96 
101  static Weights DerivativeWeights(size_t N, double x, double a = -1,
102  double b = 1);
103 
108  static DiffMatrix DifferentiationMatrix(size_t N, double a = -1,
109  double b = 1);
110 
129  static Weights IntegrationWeights(size_t N, double a = -1, double b = 1);
130 
134  template <size_t M>
135  static Matrix matrix(std::function<Eigen::Matrix<double, M, 1>(double)> f,
136  size_t N, double a = -1, double b = 1) {
137  Matrix Xmat(M, N);
138  for (size_t j = 0; j < N; j++) {
139  Xmat.col(j) = f(Point(N, j, a, b));
140  }
141  return Xmat;
142  }
143 }; // \ Chebyshev2
144 
145 } // namespace gtsam
int array[24]
static double Point(size_t N, int j)
Specific Chebyshev point.
Definition: Chebyshev2.h:55
Matrix< RealScalar, Dynamic, Dynamic > M
Definition: bench_gemm.cpp:51
Jet< T, N > cos(const Jet< T, N > &f)
Definition: jet.h:426
static Vector Points(size_t N)
All Chebyshev points.
Definition: Chebyshev2.h:72
Eigen::MatrixXd Matrix
Definition: base/Matrix.h:39
#define M_PI
Definition: main.h:106
#define EIGEN_MAKE_ALIGNED_OPERATOR_NEW
Definition: Memory.h:841
static const Pose3 T2(Rot3::Rodrigues(0.3, 0.2, 0.1), P2)
#define N
Definition: gksort.c:12
static Vector Points(size_t N, double a, double b)
All Chebyshev points, within [a,b] interval.
Definition: Chebyshev2.h:79
Eigen::VectorXd Vector
Definition: Vector.h:38
Compute an interpolating basis.
Base class and basic functions for Manifold types.
static Matrix matrix(std::function< Eigen::Matrix< double, M, 1 >(double)> f, size_t N, double a=-1, double b=1)
Definition: Chebyshev2.h:135
Point2(* f)(const Point3 &, OptionalJacobian< 2, 3 >)
const G & b
Definition: Group.h:86
static double Point(size_t N, int j, double a, double b)
Specific Chebyshev point, within [a,b] interval.
Definition: Chebyshev2.h:64
traits
Definition: chartTesting.h:28
static const Similarity3 T1(R, Point3(3.5, -8.2, 4.2), 1)
Special class for optional Jacobian arguments.
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
std::ptrdiff_t j


gtsam
Author(s):
autogenerated on Tue Jul 4 2023 02:34:01