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>
22#include <BoostVectorStandAlone.h>
34namespace Belle2::BoostVectorCalib {
37 double median(
double* v,
int n)
39 B2ASSERT(
"At least 3 points to get median", n > 2);
42 std::nth_element(v, v + n / 2, v + n);
45 std::nth_element(v, v + n / 2, v + n);
48 std::nth_element(v, v + n / 2 - 1, v + n);
49 double v2 = v[n / 2 - 1];
60 std::vector<Event> getEvents(TTree* tr)
63 std::vector<Event> events;
64 events.reserve(tr->GetEntries());
68 tr->SetBranchAddress(
"run", &evt.run);
69 tr->SetBranchAddress(
"exp", &evt.exp);
70 tr->SetBranchAddress(
"event", &evt.evtNo);
75 tr->SetBranchAddress(
"mu0_p", &p0);
76 tr->SetBranchAddress(
"mu1_p", &p1);
78 tr->SetBranchAddress(
"mu0_pid", &evt.mu0.pid);
79 tr->SetBranchAddress(
"mu1_pid", &evt.mu1.pid);
82 tr->SetBranchAddress(
"time", &evt.t);
85 for (
int i = 0; i < tr->GetEntries(); ++i) {
93 events.push_back(evt);
97 sort(events.begin(), events.end(), [](Event e1, Event e2) {return e1.t < e2.t;});
112 Vector3d pars(-c, -s, -1);
114 res =
res.array().square();
116 return median(
res.data(),
res.rows());
122 MatrixXd toMat(
const std::vector<std::vector<double>>& vecs)
124 MatrixXd mat(vecs[0].size(), vecs.size());
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];
134 std::vector<double> getRapidities(std::vector<B2Vector3D> vecs, std::vector<double> boostDir)
136 B2Vector3D boost(boostDir[0], boostDir[1], 1);
137 boost = boost.Unit();
139 double th0 = vecs[0].Angle(boost);
140 double th1 = vecs[1].Angle(boost);
142 const double mL = 105.6583745e-3;
144 double pMu0 = vecs[0].Mag();
145 double pMu1 = vecs[1].Mag();
147 double C0 = 1. / sqrt(pow(mL / pMu0, 2) + 1);
148 double C1 = 1. / sqrt(pow(mL / pMu1, 2) + 1);
150 double y0 = atanh(cos(th0) * C0);
151 double y1 = atanh(cos(th1) * C1);
157 std::vector<double> fitBoostFast(
const std::vector<Event>& evts)
159 std::vector<double> vCos, vSin, vData;
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};
167 double phi = n.Phi();
168 double angle = M_PI / 2 - n.Theta();
169 auto res = std::make_pair(phi, tan(angle));
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);
181 fun.mat = toMat({vCos, vSin, vData});
182 fun.res.resize(vCos.size());
184 auto res =
getMinimum(fun, 130e-3, 170e-3, -20e-3, 20e-3);
198 std::vector<Event> filter(
const std::vector<Event>& evts,
const std::vector<double>& boostDir,
double pidCut,
double rapCut)
200 std::vector<Event> evtsF;
202 for (
const auto& e : evts) {
203 for (
int i = 0; i < e.nBootStrap * e.isSig; ++i) {
206 if (e.mu0.pid < pidCut || e.mu1.pid < pidCut) {
210 std::vector<B2Vector3D> vecs = {e.mu0.p, e.mu1.p};
211 std::vector<double> raps = getRapidities(vecs, boostDir);
213 double yDiff = abs(raps[0] - raps[1]) / 2;
216 if (yDiff > rapCut)
continue;
228 double fitBoostMagnitude(
const std::vector<Event>& evts,
const std::vector<double>& boostDir)
231 std::vector<double> yAvgVec;
233 for (
auto e : evts) {
234 if (e.nBootStrap * e.isSig == 0)
continue;
235 std::vector<B2Vector3D> vecs = {e.mu0.p, e.mu1.p};
237 std::vector<double> raps = getRapidities(vecs, boostDir);
239 double yAvg = (raps[0] + raps[1]) / 2.;
240 for (
int i = 0; i < e.nBootStrap * e.isSig; ++i) {
241 yAvgVec.push_back(yAvg);
245 return median(yAvgVec.data(), yAvgVec.size());
262 Matrix3d res = Matrix3d::Zero();
263 for (
unsigned i = 1; i <
vecs.size(); ++i) {
265 res += diff * diff.transpose();
268 res *= 1. / (
vecs.size() - 1);
278 B2Vector3D getBoostVector(
const std::vector<Event>& evts)
281 std::vector<double> boostDir = fitBoostFast(evts);
284 double yMag = fitBoostMagnitude(evts, boostDir);
287 double beta = tanh(yMag);
293 bVec = beta * bVec.
Unit();
299 std::pair<Vector3d, Matrix3d> getBoostAndError(std::vector<Event> evts)
301 evts = filter(evts, {151.986e-3 , 0}, 0.9, 1.0 );
305 const int nBoost = 10;
306 for (
int i = 0; i < nBoost; ++i) {
309 e.nBootStrap = gRandom->Poisson(1);
312 Vector3d boost(b.X(), b.Y(), b.Z());
315 return std::make_pair(var.getNom(), var.getCov());
320 std::vector<std::vector<Event>> separate(
const std::vector<Event>& evts,
const std::vector<double>& splitPoints)
322 std::vector<std::vector<Event>> evtsOut(splitPoints.size() + 1);
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);
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);
336 evtsOut[0].push_back(ev);
348 std::tuple<std::vector<VectorXd>, std::vector<MatrixXd>, MatrixXd> runBoostVectorAnalysis(std::vector<Event> evts,
349 const std::vector<double>& splitPoints)
351 int n = splitPoints.size() + 1;
352 std::vector<VectorXd> boostVec(n);
353 std::vector<MatrixXd> boostVecUnc(n);
354 MatrixXd boostVecSpred;
356 std::vector<std::vector<Event>> evtsSep = separate(evts, splitPoints);
359 for (
int i = 0; i < n; ++i) {
361 std::tie(boostVec[i], boostVecUnc[i]) = getBoostAndError(evtsSep[i]);
365 double aX = 1e3 * boostVec[i][0] / boostVec[i][2];
367 double aY = 1e3 * boostVec[i][1] / boostVec[i][2];
370 double eX = 1e3 *
sqrt(boostVecUnc[i](0, 0)) / boostVec[i](2);
372 double eY = 1e3 *
sqrt(boostVecUnc[i](1, 1)) / boostVec[i](2);
374 B2INFO(evtsSep[i][0].run <<
" " << i <<
" " << evtsSep[i].size() <<
" : " << aX <<
"+-" << eX <<
" " << aY <<
"+-" << eY <<
375 " : " << 1e3 * boostVec[i].norm());
379 boostVecSpred = MatrixXd::Zero(3, 3);
381 return std::make_tuple(boostVec, boostVecUnc, boostVecSpred);
B2Vector3< DataType > Unit() const
Unit vector parallel to this.
B2Vector3< double > B2Vector3D
typedef for common usage with double
double sqrt(double a)
sqrt for double
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.
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.
void add(Vector3d v)
Add replica.
std::vector< Vector3d > vecs
Vector of replicas.