Belle II Software development
BoostVectorStandAlone.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#include <iostream>
10#include <vector>
11#include <tuple>
12#include <TTree.h>
13#include <Eigen/Dense>
14
15//if compiled within BASF2
16#ifdef _PACKAGE_
17#include <reconstruction/calibration/BeamSpotBoostInvMass/BoostVectorStandAlone.h>
18#include <reconstruction/calibration/BeamSpotBoostInvMass/Splitter.h>
19#include <reconstruction/calibration/BeamSpotBoostInvMass/tools.h>
20#include <reconstruction/calibration/BeamSpotBoostInvMass/minimizer.h>
21#else
22#include <BoostVectorStandAlone.h>
23#include <Splitter.h>
24#include <tools.h>
25#include <minimizer.h>
26#endif
27
28using Eigen::VectorXd;
29using Eigen::Vector3d;
30using Eigen::MatrixXd;
31using Eigen::Matrix3d;
32
33
34namespace Belle2::BoostVectorCalib {
35
37 double median(double* v, int n)
38 {
39 B2ASSERT("At least 3 points to get median", n > 2);
40
41 if (n % 2 == 1) { // odd number of elements
42 std::nth_element(v, v + n / 2, v + n);
43 return v[n / 2];
44 } else { //even
45 std::nth_element(v, v + n / 2, v + n);
46 double v1 = v[n / 2];
47
48 std::nth_element(v, v + n / 2 - 1, v + n);
49 double v2 = v[n / 2 - 1];
50
51 return (v1 + v2) / 2;
52 }
53 }
54
55
56
57
58
59 // Read events from TTree to std::vector
60 std::vector<Event> getEvents(TTree* tr)
61 {
62
63 std::vector<Event> events;
64 events.reserve(tr->GetEntries());
65
66 Event evt;
67
68 tr->SetBranchAddress("run", &evt.run);
69 tr->SetBranchAddress("exp", &evt.exp);
70 tr->SetBranchAddress("event", &evt.evtNo);
71
72 B2Vector3D* p0 = nullptr;
73 B2Vector3D* p1 = nullptr;
74
75 tr->SetBranchAddress("mu0_p", &p0);
76 tr->SetBranchAddress("mu1_p", &p1);
77
78 tr->SetBranchAddress("mu0_pid", &evt.mu0.pid);
79 tr->SetBranchAddress("mu1_pid", &evt.mu1.pid);
80
81
82 tr->SetBranchAddress("time", &evt.t); //time in hours
83
84
85 for (int i = 0; i < tr->GetEntries(); ++i) {
86 tr->GetEntry(i);
87
88 evt.mu0.p = *p0;
89 evt.mu1.p = *p1;
90
91 evt.nBootStrap = 1;
92 evt.isSig = true;
93 events.push_back(evt);
94 }
95
96 //sort by time
97 sort(events.begin(), events.end(), [](Event e1, Event e2) {return e1.t < e2.t;});
98
99
100 return events;
101 }
102
104 struct FunBoost {
105 MatrixXd mat;
106 VectorXd res;
107
109 double operator()(double c, double s)
110 {
111 // c and s are slopes of 2D linear function
112 Vector3d pars(-c, -s, -1);
113 res = mat * pars; // calculate residuals
114 res = res.array().square();
115
116 return median(res.data(), res.rows());
117 }
118 };
119
120
122 MatrixXd toMat(const std::vector<std::vector<double>>& vecs)
123 {
124 MatrixXd mat(vecs[0].size(), vecs.size());
125
126 for (unsigned k = 0; k < vecs.size(); ++k)
127 for (unsigned i = 0; i < vecs[k].size(); ++i)
128 mat(i, k) = vecs[k][i];
129 return mat;
130 }
131
132
134 std::vector<double> getRapidities(std::vector<B2Vector3D> vecs, std::vector<double> boostDir)
135 {
136 B2Vector3D boost(boostDir[0], boostDir[1], 1);
137 boost = boost.Unit();
138
139 double th0 = vecs[0].Angle(boost);
140 double th1 = vecs[1].Angle(boost);
141
142 const double mL = 105.6583745e-3; //muon mass [GeV]
143
144 double pMu0 = vecs[0].Mag();
145 double pMu1 = vecs[1].Mag();
146
147 double C0 = 1. / sqrt(pow(mL / pMu0, 2) + 1);
148 double C1 = 1. / sqrt(pow(mL / pMu1, 2) + 1);
149
150 double y0 = atanh(cos(th0) * C0);
151 double y1 = atanh(cos(th1) * C1);
152
153 return {y0, y1};
154 }
155
157 std::vector<double> fitBoostFast(const std::vector<Event>& evts)
158 {
159 std::vector<double> vCos, vSin, vData;
160
161 for (const auto& e : evts) {
162 if (e.nBootStrap * e.isSig == 0) continue;
163 std::vector<B2Vector3D> vecs = {e.mu0.p, e.mu1.p};
164
165
166 B2Vector3D n = vecs[0].Cross(vecs[1]);
167 double phi = n.Phi();
168 double angle = M_PI / 2 - n.Theta();
169 auto res = std::make_pair(phi, tan(angle));
170
171 for (int i = 0; i < e.nBootStrap * e.isSig; ++i) {
172 vCos.push_back(cos(res.first));
173 vSin.push_back(sin(res.first));
174 vData.push_back(res.second);
175 }
176 }
177
178
179 FunBoost fun;
180
181 fun.mat = toMat({vCos, vSin, vData});
182 fun.res.resize(vCos.size());
183
184 auto res = getMinimum(fun, 130e-3, 170e-3, -20e-3, 20e-3);
185
186
187 return res;
188 }
189
190
198 std::vector<Event> filter(const std::vector<Event>& evts, const std::vector<double>& boostDir, double pidCut, double rapCut)
199 {
200 std::vector<Event> evtsF;
201
202 for (const auto& e : evts) {
203 for (int i = 0; i < e.nBootStrap * e.isSig; ++i) {
204
205 //Keep only muons
206 if (e.mu0.pid < pidCut || e.mu1.pid < pidCut) {
207 continue;
208 }
209
210 std::vector<B2Vector3D> vecs = {e.mu0.p, e.mu1.p};
211 std::vector<double> raps = getRapidities(vecs, boostDir);
212
213 double yDiff = abs(raps[0] - raps[1]) / 2;
214
215 // Exclude events with fwd/bwd muons where det acc is limited
216 if (yDiff > rapCut) continue;
217
218 evtsF.push_back(e);
219 }
220 }
221
222 return evtsF;
223 }
224
225
226
228 double fitBoostMagnitude(const std::vector<Event>& evts, const std::vector<double>& boostDir)
229 {
230
231 std::vector<double> yAvgVec;
232
233 for (auto e : evts) {
234 if (e.nBootStrap * e.isSig == 0) continue;
235 std::vector<B2Vector3D> vecs = {e.mu0.p, e.mu1.p};
236
237 std::vector<double> raps = getRapidities(vecs, boostDir);
238
239 double yAvg = (raps[0] + raps[1]) / 2.;
240 for (int i = 0; i < e.nBootStrap * e.isSig; ++i) {
241 yAvgVec.push_back(yAvg);
242 }
243 }
244
245 return median(yAvgVec.data(), yAvgVec.size());
246 }
247
248
250 struct VectorVar {
251 std::vector<Vector3d> vecs;
252
254 void add(Vector3d v) { vecs.push_back(v); }
255
257 Vector3d getNom() const { return vecs[0]; }
258
260 MatrixXd getCov() const
261 {
262 Matrix3d res = Matrix3d::Zero();
263 for (unsigned i = 1; i < vecs.size(); ++i) {
264 Vector3d diff = vecs[i] - vecs[0];
265 res += diff * diff.transpose();
266 }
267 if (vecs.size() > 1)
268 res *= 1. / (vecs.size() - 1);
269 return res;
270
271 }
272
273
274 };
275
276
278 B2Vector3D getBoostVector(const std::vector<Event>& evts)
279 {
280 // Get the direction of the boost vector (tanXZ and tanYZ)
281 std::vector<double> boostDir = fitBoostFast(evts);
282
283 // Get the rapidity of the boost vector
284 double yMag = fitBoostMagnitude(evts, boostDir);
285
286 //from rapidity to velocity
287 double beta = tanh(yMag);
288
289 // Vector pointing in the boost vector direction
290 B2Vector3D bVec = B2Vector3D(boostDir[0], boostDir[1], 1);
291
292 // Adding proper magnitude of the vector (in speed of light units)
293 bVec = beta * bVec.Unit();
294
295 return bVec;
296 }
297
299 std::pair<Vector3d, Matrix3d> getBoostAndError(std::vector<Event> evts)
300 {
301 evts = filter(evts, {151.986e-3 /*TanNomAngle*/, 0}, 0.9/*muon pid*/, 1.0 /*rap cut*/);
302
303 VectorVar var;
304
305 const int nBoost = 10;
306 for (int i = 0; i < nBoost; ++i) {
307 if (i != 0)
308 for (auto& e : evts)
309 e.nBootStrap = gRandom->Poisson(1);
310
311 B2Vector3D b = getBoostVector(evts);
312 Vector3d boost(b.X(), b.Y(), b.Z());
313 var.add(boost);
314 }
315 return std::make_pair(var.getNom(), var.getCov());
316
317 }
318
320 std::vector<std::vector<Event>> separate(const std::vector<Event>& evts, const std::vector<double>& splitPoints)
321 {
322 std::vector<std::vector<Event>> evtsOut(splitPoints.size() + 1);
323
324 for (const auto& ev : evts) {
325 for (int i = 0; i < int(splitPoints.size()) - 1; ++i) {
326 if (splitPoints[i] < ev.t && ev.t < splitPoints[i + 1]) {
327 evtsOut[i + 1].push_back(ev);
328 break;
329 }
330 }
331
332 if (splitPoints.size() >= 1) {
333 if (ev.t < splitPoints.front()) evtsOut.front().push_back(ev);
334 if (ev.t > splitPoints.back()) evtsOut.back().push_back(ev);
335 } else {
336 evtsOut[0].push_back(ev);
337 }
338 }
339 return evtsOut;
340 }
341
342
343
344
345
346 // Returns tuple with the boost vector parameters
347 // cppcheck-suppress passedByValue
348 std::tuple<std::vector<VectorXd>, std::vector<MatrixXd>, MatrixXd> runBoostVectorAnalysis(std::vector<Event> evts,
349 const std::vector<double>& splitPoints)
350 {
351 int n = splitPoints.size() + 1;
352 std::vector<VectorXd> boostVec(n);
353 std::vector<MatrixXd> boostVecUnc(n);
354 MatrixXd boostVecSpred;
355
356 std::vector<std::vector<Event>> evtsSep = separate(evts, splitPoints);
357
358
359 for (int i = 0; i < n; ++i) {
360
361 std::tie(boostVec[i], boostVecUnc[i]) = getBoostAndError(evtsSep[i]);
362
363
364 // TanThetaXZ of boostVector [mrad]
365 double aX = 1e3 * boostVec[i][0] / boostVec[i][2];
366 // TanThetaYZ of boostVector [mrad]
367 double aY = 1e3 * boostVec[i][1] / boostVec[i][2];
368
369 // TanThetaXZ unc of boostVector [mrad]
370 double eX = 1e3 * sqrt(boostVecUnc[i](0, 0)) / boostVec[i](2);
371 // TanThetaYZ unc of boostVector [mrad]
372 double eY = 1e3 * sqrt(boostVecUnc[i](1, 1)) / boostVec[i](2);
373
374 B2INFO(evtsSep[i][0].run << " " << i << " " << evtsSep[i].size() << " : " << aX << "+-" << eX << " " << aY << "+-" << eY <<
375 " : " << 1e3 * boostVec[i].norm());
376 }
377
378 // Spread is currently not calculated
379 boostVecSpred = MatrixXd::Zero(3, 3);
380
381 return std::make_tuple(boostVec, boostVecUnc, boostVecSpred);
382 }
383
384}
B2Vector3< DataType > Unit() const
Unit vector parallel to this.
Definition: B2Vector3.h:269
B2Vector3< double > B2Vector3D
typedef for common usage with double
Definition: B2Vector3.h:516
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
std::vector< double > getMinimum(std::function< double(double, double)> fun, double xMin, double xMax, double yMin, double yMax)
Get minimum of 2D function in the rectangular domain defined by xMin,xMax & yMin,yMax.
Definition: minimizer.h:44
const std::vector< double > v2
MATLAB generated random vector.
const std::vector< double > v1
MATLAB generated random vector.
Functor to minimize median of residuals^2.
double operator()(double c, double s)
get median of residuals^2
VectorXd res
vector with residuals^2
MatrixXd mat
Matrix of the linear system.
Structure to store all bootstrap replicas.
MatrixXd getCov() const
Get cov matrix with unc.
Vector3d getNom() const
Get the nominal result.
std::vector< Vector3d > vecs
Vector of replicas.