Belle II Software  release-06-02-00
nodes.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 
10 
11 #include <cassert>
12 #include <vector>
13 #include <iostream>
14 #include <cmath>
15 
16 #include <Eigen/Dense>
17 
18 namespace Belle2 {
24  using namespace std;
25  using Eigen::VectorXd;
26  using Eigen::MatrixXd;
27 
28 
29  // to integrate function between 0 and 1 (more in nodes.h)
30  VectorXd getWeights(int Size)
31  {
32  const int N = Size - 1;
33  assert(N % 2 == 0);
34 
35  vector<vector<double>> coef(Size);
36  for (auto& el : coef) el.resize(Size);
37 
38 
39  for (int k = 0; k <= N / 2; ++k) {
40  coef[2 * k][N] = 1. / N;
41  coef[2 * k][0] = 1. / N ;
42 
43  coef[2 * k][N / 2] = 2. / N * (2 * ((k + 1) % 2) - 1);
44 
45  for (int n = 1; n <= N / 2 - 1; ++n)
46  coef[2 * k][n] = coef[2 * k][N - n] = 2. / N * cos(n * k * M_PI * 2 / N);
47  }
48 
49  VectorXd wgt = VectorXd::Zero(Size);
50 
51 
52  for (int i = 0; i < Size; ++i) {
53  wgt[i] += coef[0][i];
54  wgt[i] += coef[N][i] / (1. - N * N);
55  for (int k = 1; k <= N / 2 - 1; ++k) {
56  double w = 2. / (1 - 4 * k * k);
57  wgt[i] += w * coef[2 * k][i];
58  }
59 
60  wgt[i] *= 0.5; //for interval (0,1)
61  }
62  return wgt;
63  }
64 
65  // vector of nodes to integrate function between 0 and 1 (more in nodes.h)
66  VectorXd getNodes(int Size)
67  {
68  assert((Size - 1) % 2 == 0);
69  VectorXd xi = VectorXd::Zero(Size);
70  for (int i = 0; i < Size; ++i) {
71  double Cos = cos(i / (Size - 1.) * M_PI);
72  xi[i] = (1 - Cos) / 2;
73  }
74  return xi;
75  }
76 
77 // Evaluate Cheb pols up to Size at point x, x is between 0 and 1 (more in nodes.h)
78  VectorXd getPols(int Size, double x)
79  {
80  VectorXd pol(Size);
81  double C = 2 * (2 * x - 1);
82 
83  if (Size >= 1) pol[0] = 1;
84  if (Size >= 2) pol[1] = C / 2;
85 
86  for (int i = 2; i < Size; ++i)
87  pol[i] = C * pol[i - 1] - pol[i - 2];
88  return pol;
89  }
90 
91 
97  VectorXd getPolsSum(int Size, VectorXd x)
98  {
99  assert(Size > 2);
100 
101  VectorXd polSum(Size);
102 
103  VectorXd pol0 = 0 * x.array() + 1;
104  VectorXd pol1 = 2 * x.array() - 1;
105  VectorXd C = 2 * pol1;
106 
107  VectorXd pol2(x.size());
108  for (int i = 2; i < Size; ++i) {
109  polSum(i - 2) = pol0.sum();
110 
111  pol2 = C.array() * pol1.array() - pol0.array();
112 
113  pol0 = pol1;
114  pol1 = pol2;
115  }
116 
117  polSum(Size - 2) = pol0.sum();
118  polSum(Size - 1) = pol1.sum();
119 
120  return polSum;
121  }
122 
123 
124 
125 
126 
127 // Transformation matrix between Cheb nodes and Cheb. Coefficients (more in nodes.h)
128  MatrixXd getCoefs(int Size, bool isInverse = false)
129  {
130  const int N = Size - 1;
131  assert(N % 2 == 0);
132 
133  MatrixXd coef(Size, Size);
134 
135  double mul = 1;
136  double C = 1. / N;
137  if (isInverse == true) {C = 1. / 2; }
138 
139  for (int k = 0; k <= N; ++k) {
140  if (!isInverse) {
141  coef(k, N) = C;
142  coef(k, 0) = C * (k % 2 == 1 ? -1 : 1);
143  } else {
144  mul = k % 2 == 1 ? -1 : 1;
145  coef(N - k, N) = C * mul;
146  coef(N - k, 0) = C ;
147  }
148 
149  for (int n = 1; n <= N - 1; ++n) {
150  double el = cos(n * k * M_PI / N) * 2.*C * mul;
151  if (!isInverse) coef(k, N - n) = el;
152  else coef(N - k, N - n) = el;
153  }
154  }
155 
156  return coef;
157  }
158 
159 
160 
161 
162 // getCoef with better normalization of the borders (more in nodes.h)
163  MatrixXd getCoefsCheb(int oldSize)
164  {
165  auto coef = getCoefs(oldSize);
166 
167  coef.row(0) *= 0.5;
168  coef.row(coef.rows() - 1) *= 0.5;
169 
170  return coef;
171  }
172 
173 
174 
175 
177  double evalPol(const VectorXd& polCoef, double x)
178  {
179  VectorXd pols = getPols(polCoef.size(), x);
180 
181  double s = pols.dot(polCoef);
182 
183  return s;
184  }
185 
186 
190  VectorXd interpol(const VectorXd& xi, double x)
191  {
192  double Norm = (xi[xi.size() - 1] - xi[0]) / 2;
193  VectorXd coefs(xi.size());
194  for (int i = 0; i < xi.size(); ++i) {
195  double num = 1, den = 1;
196  for (int j = 0; j < xi.size(); ++j)
197  if (j != i) {
198  num *= (x - xi(j)) / Norm;
199  den *= (xi(i) - xi(j)) / Norm;
200  }
201  coefs(i) = num / den;
202  }
203  return coefs;
204  }
205 
206 
211  double interpol(VectorXd xi, VectorXd vals, double x)
212  {
213  VectorXd coefs = interpol(xi, x);
214  return coefs.dot(vals);
215  }
216 
218 }
Eigen::VectorXd interpol(const Eigen::VectorXd &xi, double x)
Get Interpolation vector k_i for point x from the function values at points xi (polynomial interpolat...
Eigen::VectorXd getWeights(int Size)
Get the vector of weights to calculate the integral over the Chebyshev nodes The nodes are by definit...
Definition: nodes.cc:30
Eigen::VectorXd getPols(int Size, double x)
Evaluate Chebyshev polynomials up to Size at point x It returns a vector of the P_i(x) for i=0....
Definition: nodes.cc:78
Eigen::MatrixXd getCoefsCheb(int oldSize)
Transformation matrix between Cheb nodes and Cheb coefficients with better normalization of the borde...
Definition: nodes.cc:163
Eigen::VectorXd getPolsSum(int Size, Eigen::VectorXd x)
Calculate the Chebyshev polynomials of order i=0..Size-1 at points given in vector x_j and sum it ove...
Eigen::MatrixXd getCoefs(int Size, bool isInverse=false)
Transformation matrix between Cheb nodes and coefficients of the Cheb polynomials Notice,...
Definition: nodes.cc:128
Eigen::VectorXd getNodes(int Size)
Get the vector of positions of the Chebyshev nodes The nodes are by definition between 0 and 1,...
Definition: nodes.cc:66
double evalPol(const Eigen::VectorXd &polCoef, double x)
Evaluate Cheb.
Abstract base class for different kinds of events.