Belle II Software  release-08-01-10
B2A400-TreeFit.py
1 #!/usr/bin/env python3
2 
3 
10 
11 
34 
35 import basf2 as b2
36 import modularAnalysis as ma
37 import vertex as vx
38 from stdPi0s import stdPi0s
39 
40 # create path
41 my_path = b2.create_path()
42 
43 # load input ROOT file
44 ma.inputMdst(filename=b2.find_file('B02D0pi0_D02pi0pi0.root', 'examples', False),
45  path=my_path)
46 
47 
48 # use standard final state particle lists
49 #
50 # creates "pi0:eff40_May2020Fit" ParticleList
51 # see Standard Particles section at https://software.belle2.org/
52 stdPi0s(listtype='eff40_May2020Fit', path=my_path)
53 
54 # reconstruct D0 -> pi0 pi0 decay
55 # keep only candidates with 1.7 < M(pi0pi0) < 2.0 GeV
56 ma.reconstructDecay(decayString='D0:pi0pi0 -> pi0:eff40_May2020Fit pi0:eff40_May2020Fit',
57  cut='1.7 < M < 2.0',
58  path=my_path)
59 
60 # reconstruct B0 -> D0 pi0 decay
61 # keep only candidates with Mbc > 5.24 GeV
62 # and -1 < Delta E < 1 GeV
63 ma.reconstructDecay(decayString='B0:all -> D0:pi0pi0 pi0:eff40_May2020Fit',
64  cut='5.24 < Mbc < 5.29 and abs(deltaE) < 1.0',
65  path=my_path)
66 
67 # perform MC matching (MC truth association)
68 ma.matchMCTruth(list_name='B0:all', path=my_path)
69 
70 vx.treeFit(list_name='B0:all',
71  conf_level=0, # 0:keep only fit survivors, -1: keep all candidates; optimise this cut for your need
72  ipConstraint=True,
73  # pins the B0 PRODUCTION vertex to the IP (increases SIG and BKG rejection) use for better vertex resolution
74  updateAllDaughters=True, # update momenta of ALL particles
75  massConstraint=['pi0'], # mass constrain ALL pi0
76  path=my_path)
77 
78 # whatever you are interested in
79 #
80 # see analysis/VariableManager/ for implementation of the vars
81 variables = [
82  'isSignal', # MC truth
83  'chiProb', # pValue of the fit use this in your analysis to reject background
84  "p", # momentum (of the B0)
85  "mcP", # generated momentum
86  "pErr", # momentum uncertainty taking the full px, py, pz covariance into account
87 
88  "daughter(0, E)", # momentum of the D0 that is updated if updateAllDaughters=True.
89  "originalDaughter(0, E)", # momentum of the original D0.
90 ]
91 
92 # safe the output
93 output_name = "B2A400-TreeFit.root"
94 my_path.add_module('VariablesToNtuple',
95  treeName='B0',
96  particleList='B0:all',
97  variables=variables,
98  fileName=output_name)
99 
100 # Process the events
101 b2.process(my_path)
102 
103 # print out the summary
104 print(b2.statistics)