Belle II Software light-2406-ragdoll
V0DaughterTrackVariables.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// Own header.
10#include <analysis/variables/V0DaughterTrackVariables.h>
11
12// include VariableManager
13#include <analysis/VariableManager/Manager.h>
14
15#include <analysis/dataobjects/Particle.h>
16#include <analysis/variables/TrackVariables.h>
17
18// framework - DataStore
19#include <framework/dataobjects/Helix.h>
20
21// dataobjects from the MDST
22#include <mdst/dataobjects/Track.h>
23#include <mdst/dataobjects/MCParticle.h>
24#include <mdst/dataobjects/TrackFitResult.h>
25#include <mdst/dataobjects/HitPatternCDC.h>
26#include <mdst/dataobjects/HitPatternVXD.h>
27
28// framework aux
29#include <framework/logging/Logger.h>
30
31#include <algorithm>
32#include <TVector2.h>
33#include <cmath>
34
35using namespace std;
36
37namespace Belle2 {
42 namespace Variable {
43
44 // helper function to get the helix parameters of the V0 daughter tracks
45 // Not registered in variable manager
46 double getV0DaughterTrackDetNHits(const Particle* particle, const double daughterID, const Const::EDetector& det)
47 {
48 auto daughter = particle->getDaughter(daughterID);
49 return trackNHits(daughter, det);
50 }
51
52 double v0DaughterTrackNCDCHits(const Particle* part, const std::vector<double>& daughterID)
53 {
54 return getV0DaughterTrackDetNHits(part, daughterID[0], Const::EDetector::CDC);
55 }
56
57 double v0DaughterTrackNSVDHits(const Particle* part, const std::vector<double>& daughterID)
58 {
59 return getV0DaughterTrackDetNHits(part, daughterID[0], Const::EDetector::SVD);
60 }
61
62 double v0DaughterTrackNPXDHits(const Particle* part, const std::vector<double>& daughterID)
63 {
64 return getV0DaughterTrackDetNHits(part, daughterID[0], Const::EDetector::PXD);
65 }
66
67 double v0DaughterTrackNVXDHits(const Particle* part, const std::vector<double>& daughterID)
68 {
69 return v0DaughterTrackNPXDHits(part, daughterID) + v0DaughterTrackNSVDHits(part, daughterID);
70 }
71
72 double v0DaughterTrackNRemovedHits(const Particle* part, const std::vector<double>& daughterID)
73 {
74 // returns 0 if called for non-V0 particles
75 if (part->getParticleSource() != Particle::EParticleSourceObject::c_V0)
76 return 0;
77 auto daughter = part->getDaughter(daughterID[0]);
78 const Track* track = daughter->getTrack();
79 if (!track) return Const::doubleNaN;
80 const TrackFitResult* trackFit = track->getTrackFitResultWithClosestMass(Const::ChargedStable(abs(
81 daughter->getPDGCode())));
82 if (!trackFit) return Const::doubleNaN;
83 double nHitsBeforeRemoval = trackFit->getHitPatternCDC().getNHits()
84 + trackFit->getHitPatternVXD().getNSVDHits()
85 + trackFit->getHitPatternVXD().getNPXDHits();
86 double nHitsAfterRemoval = trackNVXDHits(daughter) + trackNCDCHits(daughter);
87 return nHitsBeforeRemoval - nHitsAfterRemoval;
88 }
89
90 double v0DaughterTrackFirstSVDLayer(const Particle* part, const std::vector<double>& daughterID)
91 {
92 auto daughter = part->getDaughter(daughterID[0]);
93 return trackFirstSVDLayer(daughter);
94 }
95
96 double v0DaughterTrackFirstPXDLayer(const Particle* part, const std::vector<double>& daughterID)
97 {
98 auto daughter = part->getDaughter(daughterID[0]);
99 return trackFirstPXDLayer(daughter);
100 }
101
102 double v0DaughterTrackFirstCDCLayer(const Particle* part, const std::vector<double>& daughterID)
103 {
104 auto daughter = part->getDaughter(daughterID[0]);
105 return trackFirstCDCLayer(daughter);
106 }
107
108 double v0DaughterTrackLastCDCLayer(const Particle* part, const std::vector<double>& daughterID)
109 {
110 auto daughter = part->getDaughter(daughterID[0]);
111 return trackLastCDCLayer(daughter);
112 }
113
114 double v0DaughterTrackPValue(const Particle* part, const std::vector<double>& daughterID)
115 {
116 auto daughter = part->getDaughter(daughterID[0]);
117 return trackPValue(daughter);
118 }
119
120 double v0DaughterTrackD0(const Particle* part, const std::vector<double>& daughterID)
121 {
122 auto daughter = part->getDaughter(daughterID[0]);
123 return trackD0(daughter);
124 }
125
126 double v0DaughterTrackPhi0(const Particle* part, const std::vector<double>& daughterID)
127 {
128 auto daughter = part->getDaughter(daughterID[0]);
129 return trackPhi0(daughter);
130 }
131
132 double v0DaughterTrackOmega(const Particle* part, const std::vector<double>& daughterID)
133 {
134 auto daughter = part->getDaughter(daughterID[0]);
135 return trackOmega(daughter);
136 }
137
138 double v0DaughterTrackZ0(const Particle* part, const std::vector<double>& daughterID)
139 {
140 auto daughter = part->getDaughter(daughterID[0]);
141 return trackZ0(daughter);
142 }
143
144 double v0DaughterTrackTanLambda(const Particle* part, const std::vector<double>& daughterID)
145 {
146 auto daughter = part->getDaughter(daughterID[0]);
147 return trackTanLambda(daughter);
148 }
149
150 double v0DaughterTrackD0Error(const Particle* part, const std::vector<double>& daughterID)
151 {
152 auto daughter = part->getDaughter(daughterID[0]);
153 return trackD0Error(daughter);
154 }
155
156 double v0DaughterTrackPhi0Error(const Particle* part, const std::vector<double>& daughterID)
157 {
158 auto daughter = part->getDaughter(daughterID[0]);
159 return trackPhi0Error(daughter);
160 }
161
162 double v0DaughterTrackOmegaError(const Particle* part, const std::vector<double>& daughterID)
163 {
164 auto daughter = part->getDaughter(daughterID[0]);
165 return trackOmegaError(daughter);
166 }
167
168 double v0DaughterTrackZ0Error(const Particle* part, const std::vector<double>& daughterID)
169 {
170 auto daughter = part->getDaughter(daughterID[0]);
171 return trackZ0Error(daughter);
172 }
173
174 double v0DaughterTrackTanLambdaError(const Particle* part, const std::vector<double>& daughterID)
175 {
176 auto daughter = part->getDaughter(daughterID[0]);
177 return trackTanLambdaError(daughter);
178 }
179
180 double v0DaughterD0(const Particle* particle, const std::vector<double>& daughterID)
181 {
182 if (!particle)
183 return Const::doubleNaN;
184
185 ROOT::Math::XYZVector v0Vertex = particle->getVertex();
186
187 const Particle* daug = particle->getDaughter(daughterID[0]);
188
189 const TrackFitResult* trackFit = daug->getTrackFitResult();
190 if (!trackFit) return Const::doubleNaN;
191
192 UncertainHelix helix = trackFit->getUncertainHelix();
193 helix.passiveMoveBy(v0Vertex);
194
195 return helix.getD0();
196 }
197
198 double v0DaughterD0Diff(const Particle* particle)
199 {
200 return v0DaughterD0(particle, {0}) - v0DaughterD0(particle, {1});
201 }
202
203 double v0DaughterZ0(const Particle* particle, const std::vector<double>& daughterID)
204 {
205 if (!particle)
206 return Const::doubleNaN;
207
208 ROOT::Math::XYZVector v0Vertex = particle->getVertex();
209
210 const Particle* daug = particle->getDaughter(daughterID[0]);
211
212 const TrackFitResult* trackFit = daug->getTrackFitResult();
213 if (!trackFit) return Const::doubleNaN;
214
215 UncertainHelix helix = trackFit->getUncertainHelix();
216 helix.passiveMoveBy(v0Vertex);
217
218 return helix.getZ0();
219 }
220
221 double v0DaughterZ0Diff(const Particle* particle)
222 {
223 return v0DaughterZ0(particle, {0}) - v0DaughterZ0(particle, {1});
224 }
225
226 // helper function to get pull of the helix parameters of the V0 daughter tracks with the true vertex as the pivot.
227 // Not registered in variable manager
228 double getHelixParameterPullOfV0DaughterWithTrueVertexAsPivotAtIndex(const Particle* particle, const double daughterID,
229 const int tauIndex)
230 {
231 if (!particle) { return Const::doubleNaN; }
232
233 const int dID = int(std::lround(daughterID));
234 if (not(dID == 0 || dID == 1)) { return Const::doubleNaN; }
235
236 const MCParticle* mcparticle_v0 = particle->getMCParticle();
237 if (!mcparticle_v0) { return Const::doubleNaN; }
238
239 if (!(particle->getDaughter(dID))) { return Const::doubleNaN; }
240
241 const MCParticle* mcparticle = particle->getDaughter(dID)->getMCParticle();
242 if (!mcparticle) { return Const::doubleNaN; }
243
244 const TrackFitResult* trackFit = particle->getDaughter(dID)->getTrackFitResult();
245 if (!trackFit) { return Const::doubleNaN; }
246
247 // MC information
248 const ROOT::Math::XYZVector mcProdVertex = mcparticle->getVertex();
249 const ROOT::Math::XYZVector mcMomentum = mcparticle->getMomentum();
250 const double mcParticleCharge = mcparticle->getCharge();
251 const double BzAtProdVertex = BFieldManager::getFieldInTesla(mcProdVertex).Z();
252 Helix mcHelix = Helix(mcProdVertex, mcMomentum, mcParticleCharge, BzAtProdVertex);
253 mcHelix.passiveMoveBy(mcProdVertex);
254 const std::vector<double> mcHelixPars = { mcHelix.getD0(), mcHelix.getPhi0(), mcHelix.getOmega(),
255 mcHelix.getZ0(), mcHelix.getTanLambda()
256 };
257
258 // measured information (from the reconstructed particle)
259 UncertainHelix measHelix = trackFit->getUncertainHelix();
260 measHelix.passiveMoveBy(mcProdVertex);
261 const TMatrixDSym measCovariance = measHelix.getCovariance();
262 const std::vector<double> measHelixPars = {measHelix.getD0(), measHelix.getPhi0(), measHelix.getOmega(),
263 measHelix.getZ0(), measHelix.getTanLambda()
264 };
265 const std::vector<double> measErrSquare = {measCovariance[0][0], measCovariance[1][1], measCovariance[2][2],
266 measCovariance[3][3], measCovariance[4][4]
267 };
268
269 if (measErrSquare.at(tauIndex) > 0)
270 return (mcHelixPars.at(tauIndex) - measHelixPars.at(tauIndex)) / std::sqrt(measErrSquare.at(tauIndex));
271 else
272 return Const::doubleNaN;
273 }
274
275 double v0DaughterHelixWithTrueVertexAsPivotD0Pull(const Particle* part, const std::vector<double>& daughterID)
276 {
277 return getHelixParameterPullOfV0DaughterWithTrueVertexAsPivotAtIndex(part, daughterID[0], 0);
278 }
279
280 double v0DaughterHelixWithTrueVertexAsPivotPhi0Pull(const Particle* part, const std::vector<double>& daughterID)
281 {
282 return getHelixParameterPullOfV0DaughterWithTrueVertexAsPivotAtIndex(part, daughterID[0], 1);
283 }
284
285 double v0DaughterHelixWithTrueVertexAsPivotOmegaPull(const Particle* part, const std::vector<double>& daughterID)
286 {
287 return getHelixParameterPullOfV0DaughterWithTrueVertexAsPivotAtIndex(part, daughterID[0], 2);
288 }
289
290 double v0DaughterHelixWithTrueVertexAsPivotZ0Pull(const Particle* part, const std::vector<double>& daughterID)
291 {
292 return getHelixParameterPullOfV0DaughterWithTrueVertexAsPivotAtIndex(part, daughterID[0], 3);
293 }
294
295 double v0DaughterHelixWithTrueVertexAsPivotTanLambdaPull(const Particle* part, const std::vector<double>& daughterID)
296 {
297 return getHelixParameterPullOfV0DaughterWithTrueVertexAsPivotAtIndex(part, daughterID[0], 4);
298 }
299
300 double v0DaughterHelixWithOriginAsPivotD0Pull(const Particle* part, const std::vector<double>& daughterID)
301 {
302 auto daughter = part->getDaughter(daughterID[0]);
303 return getHelixD0Pull(daughter);
304 }
305
306 double v0DaughterHelixWithOriginAsPivotPhi0Pull(const Particle* part, const std::vector<double>& daughterID)
307 {
308 auto daughter = part->getDaughter(daughterID[0]);
309 return getHelixPhi0Pull(daughter);
310 }
311
312 double v0DaughterHelixWithOriginAsPivotOmegaPull(const Particle* part, const std::vector<double>& daughterID)
313 {
314 auto daughter = part->getDaughter(daughterID[0]);
315 return getHelixOmegaPull(daughter);
316 }
317
318 double v0DaughterHelixWithOriginAsPivotZ0Pull(const Particle* part, const std::vector<double>& daughterID)
319 {
320 auto daughter = part->getDaughter(daughterID[0]);
321 return getHelixZ0Pull(daughter);
322 }
323
324 double v0DaughterHelixWithOriginAsPivotTanLambdaPull(const Particle* part, const std::vector<double>& daughterID)
325 {
326 auto daughter = part->getDaughter(daughterID[0]);
327 return getHelixTanLambdaPull(daughter);
328 }
329
330 double v0DaughterTrackParam5AtIPPerigee(const Particle* part, const std::vector<double>& params)
331 {
332 auto daughter = part->getDaughter(params[0]);
333 if (!daughter) {
334 return Const::doubleNaN;
335 }
336 auto trackFit = daughter->getTrackFitResult();
337 if (!trackFit) {
338 return Const::doubleNaN;
339 }
340
341 const int paramID = int(std::lround(params[1]));
342 if (not(0 <= paramID && paramID < 5))
343 return Const::doubleNaN;
344
345 std::vector<float> tau = trackFit->getTau();
346 return tau[paramID];
347 }
348
349 double v0DaughterTrackParamCov5x5AtIPPerigee(const Particle* part, const std::vector<double>& params)
350 {
351 auto daughter = part->getDaughter(params[0]);
352 if (!daughter) {
353 return Const::doubleNaN;
354 }
355 auto trackFit = daughter->getTrackFitResult();
356 if (!trackFit) {
357 return Const::doubleNaN;
358 }
359
360 const int paramID = int(std::lround(params[1]));
361 if (not(0 <= paramID && paramID < 15))
362 return Const::doubleNaN;
363
364 std::vector<float> cov = trackFit->getCov();
365 return cov[paramID];
366 }
367
368 int convertedPhotonErrorChecks(const Particle* gamma, const std::vector<double>& daughterIndices)
369 {
370 //Check that exactly two daughter indices are provided
371 if (daughterIndices.size() != 2) {
372 B2ERROR("Invalid number of daughter indices. Please specify exactly two valid daughter indices.");
373 return -1;
374 }
375
376 //Check that there are at least (r+1) daughters where r is the bigger of the two indices provided
377 int daughterIndex1 = int(daughterIndices[0]);
378 int daughterIndex2 = int(daughterIndices[1]);
379 if (int(gamma->getNDaughters()) <= std::max(daughterIndex1, daughterIndex2)) {
380 B2ERROR("Invalid daughter indices provided. Particle does not have that many daughters.");
381 return -1;
382 }
383
384 //Check that there exists tracks associated with the daughter indices provided
385 if (!gamma->getDaughter(daughterIndex1)->getTrack()) {
386 B2ERROR("There is no track associated with daughter index " << daughterIndex1);
387 return -1;
388 }
389 if (!gamma->getDaughter(daughterIndex2)->getTrack()) {
390 B2ERROR("There is no track associated with daughter index " << daughterIndex2);
391 return -1;
392 }
393
394 //Check whether tracks used to calculate variable has been reconstructed as electrons/positrons or not (INCONSEQUENTIAL)
395 if (fabs(gamma->getDaughter(daughterIndex1)->getPDGCode()) != 11) {
396 B2INFO("The first track provided has not been reconstructed as an electron/positron. It has PDG code " << gamma->getDaughter(
397 daughterIndex1)->getPDGCode() << ". However, this is still fully admissible.");
398 }
399 if (fabs(gamma->getDaughter(daughterIndex2)->getPDGCode()) != 11) {
400 B2INFO("The second track provided has not been reconstructed as an electron/positron. It has PDG code " << gamma->getDaughter(
401 daughterIndex1)->getPDGCode() << ".However, this is still fully admissible.");
402 }
403
404 return 0;
405 }
406
407
408 int convertedPhotonLoadHelixParams(const Particle* gamma, int daughterIndex1, int daughterIndex2, double& Phi01, double& D01,
409 double& Omega1, double& Z01, double& TanLambda1, double& Phi02, double& D02, double& Omega2, double& Z02,
410 double& TanLambda2)
411 {
412 //Get helix parameters
413 //Electron/track 1
414 Helix e1Helix = gamma->getDaughter(daughterIndex1)->getTrackFitResult()->getHelix();
415
416 Phi01 = e1Helix.getPhi0();
417 D01 = e1Helix.getD0() ;
418 Omega1 = e1Helix.getOmega();
419 Z01 = e1Helix.getZ0();
420 TanLambda1 = e1Helix.getTanLambda();
421
422 //Electron/track 2
423 Helix e2Helix = gamma->getDaughter(daughterIndex2)->getTrackFitResult()->getHelix();
424
425 Phi02 = e2Helix.getPhi0();
426 D02 = e2Helix.getD0() ;
427 Omega2 = e2Helix.getOmega();
428 Z02 = e2Helix.getZ0();
429 TanLambda2 = e2Helix.getTanLambda();
430
431 //Check if either track has zero curvature
432 if (Omega1 == 0) {return -1;}
433 else if (Omega2 == 0) {return -2;}
434 else {return 0;}
435
436 }
437
438 double convertedPhotonInvariantMass(const Particle* gamma, const std::vector<double>& daughterIndices)
439 {
440 //Do basic checks
441 int errFlag = convertedPhotonErrorChecks(gamma, daughterIndices);
442 if (errFlag == -1) {return Const::doubleNaN;}
443
444 //Load helix parameters
445 double Phi01, D01, Omega1, Z01, TanLambda1, Phi02, D02, Omega2, Z02, TanLambda2;
446 int daughterIndex1 = int(daughterIndices[0]);
447 int daughterIndex2 = int(daughterIndices[1]);
448 convertedPhotonLoadHelixParams(gamma, daughterIndex1, daughterIndex2, Phi01, D01, Omega1, Z01, TanLambda1, Phi02, D02,
449 Omega2, Z02, TanLambda2);
450
451
452 //Calculating invariant mass
453 //Sine and cosine Lambda
454 double sinlam1 = TanLambda1 / sqrt(1 + (TanLambda1 * TanLambda1));
455 double coslam1 = 1 / sqrt(1 + (TanLambda1 * TanLambda1));
456 double sinlam2 = TanLambda2 / sqrt(1 + (TanLambda2 * TanLambda2));
457 double coslam2 = 1 / sqrt(1 + (TanLambda2 * TanLambda2));
458
459 //Transverse and longitudinal momentum components; energy with electron mass hypothesis
460 //electron 1
461 double p1 = gamma->getDaughter(daughterIndex1)->getMomentumMagnitude();
462 double pt1 = p1 * coslam1, pz1 = p1 * sinlam1;
463 double e1 = sqrt((p1 * p1) + (Const::electronMass * Const::electronMass));
464 //electron 2
465 double p2 = gamma->getDaughter(daughterIndex2)->getMomentumMagnitude();
466 double pt2 = p2 * coslam2, pz2 = p2 * sinlam2;
467 double e2 = sqrt((p2 * p2) + (Const::electronMass * Const::electronMass));
468
469 //Invariant mass of the two track system
470 double vtxMass = sqrt(pow(e1 + e2, 2.0) - pow(pt1 + pt2, 2.0) - pow(pz1 + pz2, 2.0));
471 return vtxMass;
472 }
473
474 double convertedPhotonDelTanLambda(const Particle* gamma, const std::vector<double>& daughterIndices)
475 {
476 //Do basic checks
477 int errFlag = convertedPhotonErrorChecks(gamma, daughterIndices);
478 if (errFlag == -1) {return Const::doubleNaN;}
479
480 //Load helix parameters
481 double Phi01, D01, Omega1, Z01, TanLambda1, Phi02, D02, Omega2, Z02, TanLambda2;
482 int daughterIndex1 = int(daughterIndices[0]);
483 int daughterIndex2 = int(daughterIndices[1]);
484 convertedPhotonLoadHelixParams(gamma, daughterIndex1, daughterIndex2, Phi01, D01, Omega1, Z01, TanLambda1, Phi02, D02,
485 Omega2, Z02, TanLambda2);
486
487
488
489 //Delta-TanLambda
490 return (TanLambda2 - TanLambda1);
491 }
492
493 double convertedPhotonDelR(const Particle* gamma, const std::vector<double>& daughterIndices)
494 {
495 //Do basic checks
496 int errFlag = convertedPhotonErrorChecks(gamma, daughterIndices);
497 if (errFlag == -1) {return Const::doubleNaN;}
498
499 //Load helix parameters
500 double Phi01, D01, Omega1, Z01, TanLambda1, Phi02, D02, Omega2, Z02, TanLambda2;
501 int daughterIndex1 = int(daughterIndices[0]);
502 int daughterIndex2 = int(daughterIndices[1]);
503 errFlag = convertedPhotonLoadHelixParams(gamma, daughterIndex1, daughterIndex2, Phi01, D01, Omega1, Z01, TanLambda1, Phi02, D02,
504 Omega2, Z02, TanLambda2);
505 if (errFlag == -1) {
506 B2ERROR("First track provided has curvature zero. Calculation of convertedPhotonDelR failed.");
507 return Const::doubleNaN;
508 }
509 if (errFlag == -2) {
510 B2ERROR("Second track provided has curvature zero. Calculation of convertedPhotonDelR failed.");
511 return Const::doubleNaN;
512 }
513
514 //Delta-R
515 double radius1 = 1 / Omega1;
516 double radius2 = 1 / Omega2;
517
518 TVector2 center1((radius1 + D01) * sin(Phi01), -1 * (radius1 + D01) * cos(Phi01));
519 TVector2 center2((radius2 + D02) * sin(Phi02), -1 * (radius2 + D02) * cos(Phi02));
520 TVector2 cenDiff = center1 - center2;
521
522 double delR = fabs(radius1) + fabs(radius2) - cenDiff.Mod();
523 return delR;
524 }
525
526 std::pair<double, double> convertedPhotonZ1Z2(const Particle* gamma, const std::vector<double>& daughterIndices)
527 {
528 //Do basic checks
529 int errFlag = convertedPhotonErrorChecks(gamma, daughterIndices);
530 if (errFlag == -1) {return std::pair<double, double>(Const::doubleNaN, Const::doubleNaN);}
531
532 //Load helix parameters
533 double Phi01, D01, Omega1, Z01, TanLambda1, Phi02, D02, Omega2, Z02, TanLambda2;
534 int daughterIndex1 = int(daughterIndices[0]);
535 int daughterIndex2 = int(daughterIndices[1]);
536 errFlag = convertedPhotonLoadHelixParams(gamma, daughterIndex1, daughterIndex2, Phi01, D01, Omega1, Z01, TanLambda1, Phi02, D02,
537 Omega2, Z02, TanLambda2);
538 if (errFlag == -1) {
539 B2ERROR("First track provided has curvature zero. Calculation of convertedPhotonZ1Z2 failed.");
540 return std::pair<double, double>(Const::doubleNaN, Const::doubleNaN);
541 }
542 if (errFlag == -2) {
543 B2ERROR("Second track provided has curvature zero. Calculation of convertedPhotonZ1Z2 failed.");
544 return std::pair<double, double>(Const::doubleNaN, Const::doubleNaN);
545 }
546
547 //Delta-Z
548 //Radial unit vectors
549 double radius1 = 1 / Omega1;
550 double radius2 = 1 / Omega2;
551
552 TVector2 center1((radius1 + D01) * sin(Phi01), -1 * (radius1 + D01) * cos(Phi01));
553 TVector2 center2((radius2 + D02) * sin(Phi02), -1 * (radius2 + D02) * cos(Phi02));
554
555 TVector2 n1 = center1 - center2; n1 = n1.Unit();
556 TVector2 n2 = -1 * n1;
557 n1 = copysign(1.0, Omega1) * n1;
558 n2 = copysign(1.0, Omega2) * n2;
559
560 //Getting running parameter phi at nominal vertex
561 double phiN1 = atan2(n1.X(), -n1.Y());
562 double phiN2 = atan2(n2.X(), -n2.Y());
563 double Phi01Intersect = phiN1 - Phi01;
564 double Phi02Intersect = phiN2 - Phi02;
565
566 double z1 = Z01 - (radius1 * TanLambda1 * Phi01Intersect);
567 double z2 = Z02 - (radius2 * TanLambda2 * Phi02Intersect);
568 std::pair<double, double> z1z2(z1, z2);
569 return z1z2;
570 }
571
572 double convertedPhotonDelZ(const Particle* gamma, const std::vector<double>& daughterIndices)
573 {
574 std::pair<double, double> z1z2 = convertedPhotonZ1Z2(gamma, daughterIndices);
575 double z1 = z1z2.first; double z2 = z1z2.second;
576 return (z1 - z2);
577 }
578
579 double convertedPhotonZ(const Particle* gamma, const std::vector<double>& daughterIndices)
580 {
581 std::pair<double, double> z1z2 = convertedPhotonZ1Z2(gamma, daughterIndices);
582 double z1 = z1z2.first; double z2 = z1z2.second;
583 return (z1 + z2) * 0.5;
584 }
585
586 TVector2 convertedPhotonXY(const Particle* gamma, const std::vector<double>& daughterIndices)
587 {
588 //Do basic checks
589 int errFlag = convertedPhotonErrorChecks(gamma, daughterIndices);
590 if (errFlag == -1) {return TVector2(Const::doubleNaN, Const::doubleNaN);}
591
592 //Load helix parameters
593 double Phi01, D01, Omega1, Z01, TanLambda1, Phi02, D02, Omega2, Z02, TanLambda2;
594 int daughterIndex1 = int(daughterIndices[0]);
595 int daughterIndex2 = int(daughterIndices[1]);
596 errFlag = convertedPhotonLoadHelixParams(gamma, daughterIndex1, daughterIndex2, Phi01, D01, Omega1, Z01, TanLambda1, Phi02, D02,
597 Omega2, Z02, TanLambda2);
598 if (errFlag == -1) {
599 B2ERROR("First track provided has curvature zero. Calculation of convertedPhotonXY failed.");
600 return TVector2(Const::doubleNaN, Const::doubleNaN);
601 }
602 if (errFlag == -2) {
603 B2ERROR("Second track provided has curvature zero. Calculation of convertedPhotonXY failed.");
604 return TVector2(Const::doubleNaN, Const::doubleNaN);
605 }
606
607 //Radial unit vectors
608 double radius1 = 1 / Omega1;
609 double radius2 = 1 / Omega2;
610
611 TVector2 center1((radius1 + D01) * sin(Phi01), -1 * (radius1 + D01) * cos(Phi01));
612 TVector2 center2((radius2 + D02) * sin(Phi02), -1 * (radius2 + D02) * cos(Phi02));
613 TVector2 cenDiff = center2 - center1;
614 double delR = fabs(radius1) + fabs(radius2) - cenDiff.Mod();
615
616 //Calculate transverse vertex
617 TVector2 n1 = cenDiff.Unit();
618 TVector2 vtxXY = center1 + ((fabs(radius1) - (delR / 2)) * n1);
619 return vtxXY;
620 }
621
622 double convertedPhotonX(const Particle* gamma, const std::vector<double>& daughterIndices)
623 {
624 auto vtxXY = convertedPhotonXY(gamma, daughterIndices);
625 return vtxXY.X();
626 }
627
628 double convertedPhotonY(const Particle* gamma, const std::vector<double>& daughterIndices)
629 {
630 auto vtxXY = convertedPhotonXY(gamma, daughterIndices);
631 return vtxXY.Y();
632 }
633
634 double convertedPhotonRho(const Particle* gamma, const std::vector<double>& daughterIndices)
635 {
636 auto vtxXY = convertedPhotonXY(gamma, daughterIndices);
637 return vtxXY.Mod();
638 }
639
640 B2Vector3D convertedPhoton3Momentum(const Particle* gamma, const std::vector<double>& daughterIndices)
641 {
642 //Do basic checks
643 int errFlag = convertedPhotonErrorChecks(gamma, daughterIndices);
644 if (errFlag == -1) {return B2Vector3D(Const::doubleNaN, Const::doubleNaN, Const::doubleNaN);}
645
646 //Load helix parameters
647 double Phi01, D01, Omega1, Z01, TanLambda1, Phi02, D02, Omega2, Z02, TanLambda2;
648 int daughterIndex1 = int(daughterIndices[0]);
649 int daughterIndex2 = int(daughterIndices[1]);
650 errFlag = convertedPhotonLoadHelixParams(gamma, daughterIndex1, daughterIndex2, Phi01, D01, Omega1, Z01, TanLambda1, Phi02, D02,
651 Omega2, Z02,
652 TanLambda2);
653 if (errFlag == -1) {
654 B2ERROR("First track provided has curvature zero. Calculation of convertedPhoton3Momentum failed.");
656 }
657 if (errFlag == -2) {
658 B2ERROR("Second track provided has curvature zero. Calculation of convertedPhoton3Momentum failed.");
660 }
661
662 //Delta-Z
663 //Radial unit vectors
664 double radius1 = 1 / Omega1;
665 double radius2 = 1 / Omega2;
666
667 TVector2 center1((radius1 + D01) * sin(Phi01), -1 * (radius1 + D01) * cos(Phi01));
668 TVector2 center2((radius2 + D02) * sin(Phi02), -1 * (radius2 + D02) * cos(Phi02));
669 TVector2 n1 = center1 - center2; n1 = n1.Unit();
670 TVector2 n2 = -1 * n1;
671 n1 = copysign(1.0, Omega1) * n1;
672 n2 = copysign(1.0, Omega2) * n2;
673
674 //Getting running parameter phi at nominal vertex
675 double phiN1 = atan2(n1.X(), -n1.Y());
676 double phiN2 = atan2(n2.X(), -n2.Y());
677
678 //Sine and cosine Lambda
679 double sinlam1 = TanLambda1 / sqrt(1 + (TanLambda1 * TanLambda1));
680 double coslam1 = 1 / sqrt(1 + (TanLambda1 * TanLambda1));
681 double sinlam2 = TanLambda2 / sqrt(1 + (TanLambda2 * TanLambda2));
682 double coslam2 = 1 / sqrt(1 + (TanLambda2 * TanLambda2));
683
684 //Photon 3-momentum
685 double p1 = gamma->getDaughter(daughterIndex1)->getMomentumMagnitude();
686 B2Vector3D e1Momentum(coslam1 * cos(phiN1), coslam1 * sin(phiN1), sinlam1);
687 double p2 = gamma->getDaughter(daughterIndex2)->getMomentumMagnitude();
688 B2Vector3D e2Momentum(coslam2 * cos(phiN2), coslam2 * sin(phiN2), sinlam2);
689 B2Vector3D gammaMomentum = (e1Momentum * p1) + (e2Momentum * p2);
690
691 return gammaMomentum;
692 }
693
694 double convertedPhotonPx(const Particle* gamma, const std::vector<double>& daughterIndices)
695 {
696 auto gammaMomentum = convertedPhoton3Momentum(gamma, daughterIndices);
697 return gammaMomentum.Px();
698 }
699
700 double convertedPhotonPy(const Particle* gamma, const std::vector<double>& daughterIndices)
701 {
702 auto gammaMomentum = convertedPhoton3Momentum(gamma, daughterIndices);
703 return gammaMomentum.Py();
704 }
705
706 double convertedPhotonPz(const Particle* gamma, const std::vector<double>& daughterIndices)
707 {
708 auto gammaMomentum = convertedPhoton3Momentum(gamma, daughterIndices);
709 return gammaMomentum.Pz();
710 }
711
712 int v0DaughtersShareInnermostHit(const Particle* part)
713 {
714 if (!part)
715 return -1;
716 auto daughterPlus = part->getDaughter(0);
717 auto daughterMinus = part->getDaughter(1);
718 if (!daughterPlus || !daughterMinus)
719 return -1;
720 auto trackFitPlus = daughterPlus->getTrackFitResult();
721 auto trackFitMinus = daughterMinus->getTrackFitResult();
722 if (!trackFitPlus || !trackFitMinus)
723 return -1;
724 int flagPlus = trackFitPlus->getHitPatternVXD().getInnermostHitShareStatus();
725 int flagMinus = trackFitMinus->getHitPatternVXD().getInnermostHitShareStatus();
726 if (flagPlus != flagMinus)
727 return -1;
728 return flagPlus;
729 }
730
731 bool v0DaughtersShareInnermostUHit(const Particle* part)
732 {
733 return ((v0DaughtersShareInnermostHit(part) / 2) == 1);
734 }
735
736 bool v0DaughtersShareInnermostVHit(const Particle* part)
737 {
738 return ((v0DaughtersShareInnermostHit(part) % 2) == 1);
739 }
740
741 VARIABLE_GROUP("V0Daughter");
742
743 REGISTER_VARIABLE("v0DaughterNCDCHits(i)", v0DaughterTrackNCDCHits, "Number of CDC hits associated to the i-th daughter track");
744 MAKE_DEPRECATED("v0DaughterNCDCHits", false, "light-2104-poseidon", R"DOC(
745 The same value can be calculated with the more generic variable `nCDCHits`,
746 so replace the current call with ``daughter(i, nCDCHits)``.)DOC");
747 REGISTER_VARIABLE("v0DaughterNSVDHits(i)", v0DaughterTrackNSVDHits, "Number of SVD hits associated to the i-th daughter track");
748 MAKE_DEPRECATED("v0DaughterNSVDHits", false, "light-2104-poseidon", R"DOC(
749 The same value can be calculated with the more generic variable `nSVDHits`,
750 so replace the current call with ``daughter(i, nSVDHits)``.)DOC");
751 REGISTER_VARIABLE("v0DaughterNPXDHits(i)", v0DaughterTrackNPXDHits, "Number of PXD hits associated to the i-th daughter track");
752 MAKE_DEPRECATED("v0DaughterNPXDHits", false, "light-2104-poseidon", R"DOC(
753 The same value can be calculated with the more generic variable `nPXDHits`,
754 so replace the current call with ``daughter(i, nPXDHits)``.)DOC");
755 REGISTER_VARIABLE("v0DaughterNVXDHits(i)", v0DaughterTrackNVXDHits, "Number of PXD+SVD hits associated to the i-th daughter track");
756 MAKE_DEPRECATED("v0DaughterNVXDHits", false, "light-2104-poseidon", R"DOC(
757 The same value can be calculated with the more generic variable `nVXDHits`,
758 so replace the current call with ``daughter(i, nVXDHits)``.)DOC");
759 REGISTER_VARIABLE("v0DaughterNRemovedHits(i)", v0DaughterTrackNRemovedHits,
760 "The number of the i-th daughter track hits removed in V0Finder. Returns 0 if called for something other than V0 daughters.");
761 REGISTER_VARIABLE("v0DaughterFirstSVDLayer(i)", v0DaughterTrackFirstSVDLayer,
762 "First activated SVD layer associated to the i-th daughter track");
763 MAKE_DEPRECATED("v0DaughterFirstSVDLayer", false, "light-2104-poseidon", R"DOC(
764 The same value can be calculated with the more generic variable `firstSVDLayer`,
765 so replace the current call with ``daughter(i, firstSVDLayer)``.)DOC");
766 REGISTER_VARIABLE("v0DaughterFirstPXDLayer(i)", v0DaughterTrackFirstPXDLayer,
767 "First activated PXD layer associated to the i-th daughter track");
768 MAKE_DEPRECATED("v0DaughterFirstPXDLayer", false, "light-2104-poseidon", R"DOC(
769 The same value can be calculated with the more generic variable `firstPXDLayer`,
770 so replace the current call with ``daughter(i, firstPXDLayer)``.)DOC");
771 REGISTER_VARIABLE("v0DaughterFirstCDCLayer(i)", v0DaughterTrackFirstCDCLayer,
772 "First activated CDC layer associated to the i-th daughter track");
773 MAKE_DEPRECATED("v0DaughterFirstCDCLayer", false, "light-2104-poseidon", R"DOC(
774 The same value can be calculated with the more generic variable `firstCDCLayer`,
775 so replace the current call with ``daughter(i, firstCDCLayer)``.)DOC");
776 REGISTER_VARIABLE("v0DaughterLastCDCLayer(i)", v0DaughterTrackLastCDCLayer,
777 "Last CDC layer associated to the i-th daughter track");
778 MAKE_DEPRECATED("v0DaughterLastCDCLayer", false, "light-2104-poseidon", R"DOC(
779 The same value can be calculated with the more generic variable `lastCDCLayer`,
780 so replace the current call with ``daughter(i, lastCDCLayer)``.)DOC");
781 REGISTER_VARIABLE("v0DaughterPValue(i)", v0DaughterTrackPValue,
782 "chi2 probalility of the i-th daughter track fit");
783 MAKE_DEPRECATED("v0DaughterPValue", false, "light-2104-poseidon", R"DOC(
784 The same value can be calculated with the more generic variable `pValue`,
785 so replace the current call with ``daughter(i, pValue)``.)DOC");
787 REGISTER_VARIABLE("v0DaughterD0(i)", v0DaughterTrackD0, "d0 of the i-th daughter track fit\n\n", "cm");
788 MAKE_DEPRECATED("v0DaughterD0", false, "light-2104-poseidon", R"DOC(
789 The same value can be calculated with the more generic variable `d0`,
790 so replace the current call with ``daughter(i, d0)``.)DOC");
791 REGISTER_VARIABLE("v0DaughterPhi0(i)", v0DaughterTrackPhi0, "phi0 of the i-th daughter track fit\n\n", "rad");
792 MAKE_DEPRECATED("v0DaughterPhi0", false, "light-2104-poseidon", R"DOC(
793 The same value can be calculated with the more generic variable `phi0`,
794 so replace the current call with ``daughter(i, phi0)``.)DOC");
795 REGISTER_VARIABLE("v0DaughterOmega(i)", v0DaughterTrackOmega, "omega of the i-th daughter track fit\n\n",
796 ":math:`\\text{cm}^{-1}`");
797 MAKE_DEPRECATED("v0DaughterOmega", false, "light-2104-poseidon", R"DOC(
798 The same value can be calculated with the more generic variable `omega`,
799 so replace the current call with ``daughter(i, omega)``.)DOC");
800 REGISTER_VARIABLE("v0DaughterZ0(i)", v0DaughterTrackZ0, "z0 of the i-th daughter track fit\n\n", "cm");
801 MAKE_DEPRECATED("v0DaughterZ0", false, "light-2104-poseidon", R"DOC(
802 The same value can be calculated with the more generic variable `z0`,
803 so replace the current call with ``daughter(i, z0)``.)DOC");
804 REGISTER_VARIABLE("v0DaughterTanLambda(i)", v0DaughterTrackTanLambda, "tan(lambda) of the i-th daughter track fit");
805 MAKE_DEPRECATED("v0DaughterTanLambda", false, "light-2104-poseidon", R"DOC(
806 The same value can be calculated with the more generic variable `tanLambda`,
807 so replace the current call with ``daughter(i, tanLambda)``.)DOC");
809 REGISTER_VARIABLE("v0DaughterD0Error(i)", v0DaughterTrackD0Error, "d0 error of the i-th daughter track fit\n\n", "cm");
810 MAKE_DEPRECATED("v0DaughterD0Error", false, "light-2104-poseidon", R"DOC(
811 The same value can be calculated with the more generic variable `d0Err`,
812 so replace the current call with ``daughter(i, d0Err)``.)DOC");
813 REGISTER_VARIABLE("v0DaughterPhi0Error(i)", v0DaughterTrackPhi0Error, "phi0 error of the i-th daughter track fit\n\n", "rad");
814 MAKE_DEPRECATED("v0DaughterPhi0Error", false, "light-2104-poseidon", R"DOC(
815 The same value can be calculated with the more generic variable `phi0Err`,
816 so replace the current call with ``daughter(i, phi0Err)``.)DOC");
817 REGISTER_VARIABLE("v0DaughterOmegaError(i)", v0DaughterTrackOmegaError, "omega error of the i-th daughter track fit\n\n",
818 ":math:`\\text{cm}^{-1}`");
819 MAKE_DEPRECATED("v0DaughterOmegaError", false, "light-2104-poseidon", R"DOC(
820 The same value can be calculated with the more generic variable `omegaErr`,
821 so replace the current call with ``daughter(i, omegaErr)``.)DOC");
822 REGISTER_VARIABLE("v0DaughterZ0Error(i)", v0DaughterTrackZ0Error, "z0 error of the i-th daughter track fit\n\n", "cm");
823 MAKE_DEPRECATED("v0DaughterZ0Error", false, "light-2104-poseidon", R"DOC(
824 The same value can be calculated with the more generic variable `z0Err`,
825 so replace the current call with ``daughter(i, z0Err)``.)DOC");
826 REGISTER_VARIABLE("v0DaughterTanLambdaError(i)", v0DaughterTrackTanLambdaError, "tan(lambda) error of the i-th daughter track fit");
827 MAKE_DEPRECATED("v0DaughterTanLambdaError", false, "light-2104-poseidon", R"DOC(
828 The same value can be calculated with the more generic variable `tanLambdaErr`,
829 so replace the current call with ``daughter(i, tanLambdaErr)``.)DOC");
830
832 REGISTER_VARIABLE("V0d0(id)", v0DaughterD0,
833 "Return the d0 impact parameter of a V0's daughter with daughterID index with the V0 vertex point as a pivot for the track.\n\n",
834 "cm");
835 REGISTER_VARIABLE("V0Deltad0", v0DaughterD0Diff,
836 "Return the difference between d0 impact parameters of V0's daughters with the V0 vertex point as a pivot for the track.\n\n",
837 "cm");
838 REGISTER_VARIABLE("V0z0(id)", v0DaughterZ0,
839 "Return the z0 impact parameter of a V0's daughter with daughterID index with the V0 vertex point as a pivot for the track.\n\n",
840 "cm");
841 REGISTER_VARIABLE("V0Deltaz0", v0DaughterZ0Diff,
842 "Return the difference between z0 impact parameters of V0's daughters with the V0 vertex point as a pivot for the track.\n\n",
843 "cm");
844
846 REGISTER_VARIABLE("v0DaughterD0PullWithTrueVertexAsPivot(i)", v0DaughterHelixWithTrueVertexAsPivotD0Pull,
847 "d0 pull of the i-th daughter track with the true V0 vertex as the track pivot");
848 REGISTER_VARIABLE("v0DaughterPhi0PullWithTrueVertexAsPivot(i)", v0DaughterHelixWithTrueVertexAsPivotPhi0Pull,
849 "phi0 pull of the i-th daughter track with the true V0 vertex as the track pivot");
850 REGISTER_VARIABLE("v0DaughterOmegaPullWithTrueVertexAsPivot(i)", v0DaughterHelixWithTrueVertexAsPivotOmegaPull,
851 "omega pull of the i-th daughter track with the true V0 vertex as the track pivot");
852 REGISTER_VARIABLE("v0DaughterZ0PullWithTrueVertexAsPivot(i)", v0DaughterHelixWithTrueVertexAsPivotZ0Pull,
853 "z0 pull of the i-th daughter track with the true V0 vertex as the track pivot");
854 REGISTER_VARIABLE("v0DaughterTanLambdaPullWithTrueVertexAsPivot(i)", v0DaughterHelixWithTrueVertexAsPivotTanLambdaPull,
855 "tan(lambda) pull of the i-th daughter track with the true V0 vertex as the track pivot");
857 REGISTER_VARIABLE("v0DaughterD0PullWithOriginAsPivot(i)", v0DaughterHelixWithOriginAsPivotD0Pull,
858 "d0 pull of the i-th daughter track with the origin as the track pivot");
859 MAKE_DEPRECATED("v0DaughterD0PullWithOriginAsPivot", false, "light-2104-poseidon", R"DOC(
860 The same value can be calculated with the more generic variable `d0Pull`,
861 so replace the current call with ``daughter(i, d0Pull)``.)DOC");
862 REGISTER_VARIABLE("v0DaughterPhi0PullWithOriginAsPivot(i)", v0DaughterHelixWithOriginAsPivotPhi0Pull,
863 "phi0 pull of the i-th daughter track with the origin as the track pivot");
864 MAKE_DEPRECATED("v0DaughterPhi0PullWithOriginAsPivot", false, "light-2104-poseidon", R"DOC(
865 The same value can be calculated with the more generic variable `phi0Pull`,
866 so replace the current call with ``daughter(i, phi0Pull)``.)DOC");
867 REGISTER_VARIABLE("v0DaughterOmegaPullWithOriginAsPivot(i)", v0DaughterHelixWithOriginAsPivotOmegaPull,
868 "omega pull of the i-th daughter track with the origin as the track pivot");
869 MAKE_DEPRECATED("v0DaughterOmegaPullWithOriginAsPivot", false, "light-2104-poseidon", R"DOC(
870 The same value can be calculated with the more generic variable `omegaPull`,
871 so replace the current call with ``daughter(i, omegaPull)``.)DOC");
872 REGISTER_VARIABLE("v0DaughterZ0PullWithOriginAsPivot(i)", v0DaughterHelixWithOriginAsPivotZ0Pull,
873 "z0 pull of the i-th daughter track with the origin as the track pivot");
874 MAKE_DEPRECATED("v0DaughterZ0PullWithOriginAsPivot", false, "light-2104-poseidon", R"DOC(
875 The same value can be calculated with the more generic variable `z0Pull`,
876 so replace the current call with ``daughter(i, z0Pull)``.)DOC");
877 REGISTER_VARIABLE("v0DaughterTanLambdaPullWithOriginAsPivot(i)", v0DaughterHelixWithOriginAsPivotTanLambdaPull,
878 "tan(lambda) pull of the i-th daughter track with the origin as the track pivot");
879 MAKE_DEPRECATED("v0DaughterTanLambdaPullWithOriginAsPivot", false, "light-2104-poseidon", R"DOC(
880 The same value can be calculated with the more generic variable `tanLambdaPull`,
881 so replace the current call with ``daughter(i, tanLambdaPull)``.)DOC");
883 REGISTER_VARIABLE("v0DaughterTau(i,j)", v0DaughterTrackParam5AtIPPerigee,
884 "j-th track parameter (at IP perigee) of the i-th daughter track. "
885 "j: 0:d0, 1:phi0, 2:omega, 3:z0, 4:tanLambda\n\n", "cm, rad, :math:`\\text{cm}^{-1}`, cm, unitless");
886 REGISTER_VARIABLE("v0DaughterCov(i,j)", v0DaughterTrackParamCov5x5AtIPPerigee,
887 "j-th element of the 15 covariance matrix elements (at IP perigee) of the i-th daughter track. "
888 "(0,0), (0,1) ... (1,1), (1,2) ... (2,2) ..."
889 "index order is: 0:d0, 1:phi0, 2:omega, 3:z0, 4:tanLambda\n\n", "cm, rad, :math:`\\text{cm}^{-1}`, cm, unitless");
891 REGISTER_VARIABLE("convertedPhotonInvariantMass(i,j)", convertedPhotonInvariantMass,
892 "Invariant mass of the i-j daughter system as defined in https://indico.belle2.org/event/3644/contributions/18622/attachments/9401/14443/Photon_vertexin_B2GM.pdf, assuming it's a converted photon\n\n",
893 "GeV/:math:`\\text{c}^2`");
894 REGISTER_VARIABLE("convertedPhotonDelTanLambda(i,j)", convertedPhotonDelTanLambda,
895 "Discriminating variable Delta-TanLambda calculated for daughters (i,j) as defined in https://indico.belle2.org/event/3644/contributions/18622/attachments/9401/14443/Photon_vertexin_B2GM.pdf, assuming it's a converted photon");
896 REGISTER_VARIABLE("convertedPhotonDelR(i,j)", convertedPhotonDelR,
897 "Discriminating variable Delta-R calculated for daughters (i,j) as defined in https://indico.belle2.org/event/3644/contributions/18622/attachments/9401/14443/Photon_vertexin_B2GM.pdf, assuming it's a converted photon\n\n",
898 "cm");
899 REGISTER_VARIABLE("convertedPhotonDelZ(i,j)", convertedPhotonDelZ,
900 "Discriminating variable Delta-Z calculated for daughters (i,j) as defined in https://indico.belle2.org/event/3644/contributions/18622/attachments/9401/14443/Photon_vertexin_B2GM.pdf, assuming it's a converted photon\n\n",
901 "cm");
902 REGISTER_VARIABLE("convertedPhotonX(i,j)", convertedPhotonX,
903 "Estimate of vertex X coordinate calculated for daughters (i,j) as defined in https://indico.belle2.org/event/3644/contributions/18622/attachments/9401/14443/Photon_vertexin_B2GM.pdf, assuming it's a converted photon\n\n",
904 "cm");
905 REGISTER_VARIABLE("convertedPhotonY(i,j)", convertedPhotonY,
906 "Estimate of vertex Y coordinate calculated for daughters (i,j) as defined in https://indico.belle2.org/event/3644/contributions/18622/attachments/9401/14443/Photon_vertexin_B2GM.pdf, assuming it's a converted photon\n\n",
907 "cm");
908 REGISTER_VARIABLE("convertedPhotonZ(i,j)", convertedPhotonZ,
909 "Estimate of vertex Z coordinate calculated for daughters (i,j) as defined in https://indico.belle2.org/event/3644/contributions/18622/attachments/9401/14443/Photon_vertexin_B2GM.pdf, assuming it's a converted photon\n\n",
910 "cm");
911 REGISTER_VARIABLE("convertedPhotonRho(i,j)", convertedPhotonRho,
912 "Estimate of vertex Rho calculated for daughters (i,j) as defined in https://indico.belle2.org/event/3644/contributions/18622/attachments/9401/14443/Photon_vertexin_B2GM.pdf, assuming it's a converted photon\n\n",
913 "cm");
914 REGISTER_VARIABLE("convertedPhotonPx(i,j)", convertedPhotonPx,
915 "Estimate of x-component of photon momentum calculated for daughters (i,j) as defined in https://indico.belle2.org/event/3644/contributions/18622/attachments/9401/14443/Photon_vertexin_B2GM.pdf, assuming it's a converted photon\n\n",
916 "GeV/c");
917 REGISTER_VARIABLE("convertedPhotonPy(i,j)", convertedPhotonPy,
918 "Estimate of y-component of photon momentum calculated for daughters (i,j) as defined in https://indico.belle2.org/event/3644/contributions/18622/attachments/9401/14443/Photon_vertexin_B2GM.pdf, assuming it's a converted photon\n\n",
919 "GeV/c");
920 REGISTER_VARIABLE("convertedPhotonPz(i,j)", convertedPhotonPz,
921 "Estimate of z-component of photon momentum calculated for daughters (i,j) as defined in https://indico.belle2.org/event/3644/contributions/18622/attachments/9401/14443/Photon_vertexin_B2GM.pdf, assuming it's a converted photon\n\n",
922 "GeV/c");
924 REGISTER_VARIABLE("v0DaughtersShare1stHit", v0DaughtersShareInnermostHit,
925 "flag for V0 daughters sharing the first(innermost) VXD hit. 0x1(0x2) bit represents V/z(U/r-phi)-hit share.");
926 REGISTER_VARIABLE("v0DaughtersShare1stUHit", v0DaughtersShareInnermostUHit,
927 "flag for V0 daughters sharing the first(innermost) VXD U-side hit.");
928 REGISTER_VARIABLE("v0DaughtersShare1stVHit", v0DaughtersShareInnermostVHit,
929 "flag for V0 daughters sharing the first(innermost) VXD V-side hit.");
930 }
932}
static ROOT::Math::XYZVector getFieldInTesla(const ROOT::Math::XYZVector &pos)
return the magnetic field at a given position in Tesla.
Definition: BFieldManager.h:61
EDetector
Enum for identifying the detector components (detector and subdetector).
Definition: Const.h:42
static const double electronMass
electron mass
Definition: Const.h:685
static const double doubleNaN
quiet_NaN
Definition: Const.h:703
const MCParticle * getMCParticle() const
Returns the pointer to the MCParticle object that was used to create this Particle (ParticleType == c...
Definition: Particle.cc:942
const Particle * getDaughter(unsigned i) const
Returns a pointer to the i-th daughter particle.
Definition: Particle.cc:631
#define MAKE_DEPRECATED(name, make_fatal, version, description)
Registers a variable as deprecated.
Definition: Manager.h:443
B2Vector3< double > B2Vector3D
typedef for common usage with double
Definition: B2Vector3.h:516
Abstract base class for different kinds of events.
Definition: ClusterUtils.h:24
STL namespace.