Belle II Software development
tools.h
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#pragma once
10
11#include <vector>
12#include <iostream>
13#include <TString.h>
14#include <TRandom.h>
15
16#include <Eigen/Dense>
17
18//If compiled within BASF2
19#ifdef _PACKAGE_
20#include <framework/logging/Logger.h>
21#else
22#ifndef B2FATAL
23#define B2FATAL(arg) { std::cout << arg << std::endl; std::exit(1);}
24#define B2ASSERT(str, cond) { if((cond) == false) {std::cout << __FILE__ <<", line "<<__LINE__ << std::endl << "" << #cond << " failed" << std::endl; std::exit(1);} }
25#endif
26#endif
27
28
29
30
31namespace Belle2 {
38 inline TString rn() {return Form("%d", gRandom->Integer(1000000000)); }
39
41 inline std::vector<std::vector<double>> merge(std::vector<std::vector<std::vector<double>>> toMerge)
42 {
43 std::vector<std::vector<double>> allVecs;
44 for (const auto& v : toMerge)
45 allVecs.insert(allVecs.end(), v.begin(), v.end());
46 return allVecs;
47 }
48
49
51 inline Eigen::VectorXd vec2vec(std::vector<double> vec)
52 {
53 Eigen::VectorXd v(vec.size());
54 for (unsigned i = 0; i < vec.size(); ++i) {
55 v[i] = vec[i];
56 }
57 return v;
58 }
59
61 inline std::vector<double> vec2vec(Eigen::VectorXd v)
62 {
63 std::vector<double> vNew(v.rows());
64 for (int i = 0; i < v.rows(); ++i)
65 vNew[i] = v(i);
66 return vNew;
67 }
68
69
70
72 inline Eigen::MatrixXd vecs2mat(std::vector<std::vector<double>> vecs)
73 {
74 Eigen::MatrixXd m(vecs[0].size(), vecs.size());
75 for (unsigned i = 0; i < vecs[0].size(); ++i)
76 for (unsigned j = 0; j < vecs.size(); ++j) {
77 m(i, j) = vecs[j][i];
78 }
79 return m;
80 }
81
83 inline std::vector<double> getRangeLin(int nVals, double xMin, double xMax)
84 {
85 B2ASSERT("At least one value in the spline required", nVals >= 1);
86 if (nVals == 1) return {};
87 std::vector<double> v(nVals);
88 for (int i = 0; i < nVals; ++i)
89 v[i] = xMin + i * (xMax - xMin) / (nVals - 1);
90 return v;
91 }
92
94 inline std::vector<double> getRangeZero(int nVals, double xMin, double xMax)
95 {
96 B2ASSERT("At least one value in the spline required", nVals >= 1);
97 if (nVals == 1) return {};
98 std::vector<double> v(nVals - 1);
99 for (int i = 1; i < nVals; ++i)
100 v[i - 1] = xMin + i * (xMax - xMin) / (nVals);
101 return v;
102 }
103
104
106 inline std::vector<double> slice(std::vector<double> v, unsigned ind, unsigned n)
107 {
108 std::vector<double> vNew;
109 for (unsigned i = ind; i < ind + n && i < v.size(); ++i)
110 vNew.push_back(v[i]);
111 return vNew;
112 }
113
115 inline double eval(const std::vector<double>& spl, const std::vector<double>& vals, double x)
116 {
117 int order = -1;
118 if (spl.size() == 0)
119 order = 0;
120 else if (spl.size() == vals.size() - 1)
121 order = 0;
122 else if (spl.size() == vals.size())
123 order = 1;
124 else {
125 B2FATAL("Unknown order of spline");
126 }
127 B2ASSERT("Spline order should be zero or one", order == 0 || order == 1);
128
129 if (order == 1) {
130 B2ASSERT("Linear spline only meaningful for two or more nodes", spl.size() >= 2);
131 B2ASSERT("As nodes as values in lin. spline", spl.size() == vals.size());
132
133 if (x <= spl[0]) return vals[0];
134 if (x >= spl.back()) return vals.back();
135
136 // binary search for position
137 int i1 = lower_bound(spl.begin(), spl.end(), x) - spl.begin() - 1;
138
139 if (!(spl[i1] <= x && x <= spl[i1 + 1])) {
140 B2FATAL("Wrong place founded : " << spl[i1] << " " << x << " " << spl[i1 + 1]);
141 }
142
143 // Linear interpolation between neighbouring nodes
144 double v = ((spl[i1 + 1] - x) * vals[i1] + (x - spl[i1]) * vals[i1 + 1]) / (spl[i1 + 1] - spl[i1]);
145 return v;
146 } else if (order == 0) { //zero order polynomial
147 B2ASSERT("#values vs #nodes in zero-order spline", spl.size() + 1 == vals.size());
148 if (vals.size() == 1) {
149 return vals[0];
150 } else {
151 double res = vals[0]; //by default value from lowest node
152 for (unsigned i = 0; i < spl.size(); ++i) {
153 if (spl[i] <= x) res = vals[i + 1];
154 else break;
155 }
156 return res;
157 }
158 }
159 return -99;
160 }
161
162
165 struct Spline {
166 std::vector<double> nodes;
167 std::vector<double> vals;
168 std::vector<double> errs;
169
171 double val(double x) const {return eval(nodes, vals, x);}
172
174 double err(double x) const {return eval(nodes, errs, x);}
175
177 double center() const
178 {
179 if (nodes.size() == 0) //dummy situation
180 return 0;
181 else if (nodes.size() % 2 == 1)
182 return nodes[nodes.size() / 2];
183 else
184 return (nodes[nodes.size() / 2 - 1] + nodes[nodes.size() / 2]) / 2;
185 }
186
188 void print(TString tag = "")
189 {
190 std::cout << tag << " : nodes: ";
191 for (auto n : nodes) {
192 std::cout << n << " ";
193 }
194 std::cout << std::endl;
195 std::cout << tag << " : vals: ";
196 for (auto n : vals) {
197 std::cout << n << " ";
198 }
199 std::cout << std::endl;
200
201 std::cout << tag << " : errs: ";
202 for (auto n : errs) {
203 std::cout << n << " ";
204 }
205 std::cout << std::endl;
206 }
207
208 Spline() : nodes{}, vals{0}, errs{0} {}
209
210 };
211
212
214}
215
Eigen::VectorXd vec2vec(std::vector< double > vec)
std vector -> ROOT vector
Definition: tools.h:51
std::vector< Atom > slice(std::vector< Atom > vec, int s, int e)
Slice the vector to contain only elements with indexes s .. e (included)
Definition: Splitter.h:85
double eval(const std::vector< double > &spl, const std::vector< double > &vals, double x)
Evaluate spline (zero order or first order) in point x.
Definition: tools.h:115
std::vector< std::vector< double > > merge(std::vector< std::vector< std::vector< double > > > toMerge)
merge { vector<double> a, vector<double> b} into {a, b}
Definition: tools.h:41
std::vector< double > getRangeZero(int nVals, double xMin, double xMax)
Equidistant range between xMin and xMax for spline of the zero order.
Definition: tools.h:94
std::vector< double > getRangeLin(int nVals, double xMin, double xMax)
Equidistant range between xMin and xMax for spline of the first order.
Definition: tools.h:83
Eigen::MatrixXd vecs2mat(std::vector< std::vector< double > > vecs)
merge columns (from std::vectors) into ROOT matrix
Definition: tools.h:72
TString rn()
Get random string.
Definition: tools.h:38
Abstract base class for different kinds of events.
Spline structure for zero-order & linear splines.
Definition: tools.h:165
double val(double x) const
get value of spline at point x
Definition: tools.h:171
std::vector< double > nodes
vector of spline nodes
Definition: tools.h:166
std::vector< double > errs
vector of spline errors
Definition: tools.h:168
std::vector< double > vals
vector of spline values
Definition: tools.h:167
void print(TString tag="")
print the spline
Definition: tools.h:188
double center() const
Get center of the spline domain.
Definition: tools.h:177
double err(double x) const
get error of spline at point x
Definition: tools.h:174