Belle II Software  release-06-00-14
059_vertex_fitting.py
1 #!/usr/bin/env python3
2 
3 import sys
4 import basf2 as b2
5 import modularAnalysis as ma
6 import stdV0s
7 import flavorTagger as ft
8 from variables import variables as vm # shorthand for VariableManager
9 import variables.collections as vc
10 import variables.utils as vu
11 import vertex
12 
13 # get input file number from the command line
14 filenumber = sys.argv[1]
15 
16 # set analysis global tag (needed for flavor tagging)
17 b2.conditions.prepend_globaltag("analysis_tools_release-04-02")
18 
19 # create path
20 main = b2.Path()
21 
22 # load input data from mdst/udst file
23 ma.inputMdstList(
24  environmentType="default",
25  filelist=[b2.find_file(f"starterkit/2021/1111540100_eph3_BGx0_{filenumber}.root", "examples")],
26  path=main,
27 )
28 
29 # fill final state particle lists
30 ma.fillParticleList(
31  "e+:uncorrected",
32  "electronID > 0.1 and dr < 0.5 and abs(dz) < 2 and thetaInCDCAcceptance",
33  path=main,
34 )
35 stdV0s.stdKshorts(path=main)
36 
37 # apply Bremsstrahlung correction to electrons
38 vm.addAlias(
39  "goodFWDGamma", "passesCut(clusterReg == 1 and clusterE > 0.075)"
40 )
41 vm.addAlias(
42  "goodBRLGamma", "passesCut(clusterReg == 2 and clusterE > 0.05)"
43 )
44 vm.addAlias(
45  "goodBWDGamma", "passesCut(clusterReg == 3 and clusterE > 0.1)"
46 )
47 vm.addAlias(
48  "goodGamma", "passesCut(goodFWDGamma or goodBRLGamma or goodBWDGamma)"
49 )
50 ma.fillParticleList("gamma:brems", "goodGamma", path=main)
51 ma.correctBrems("e+:corrected", "e+:uncorrected", "gamma:brems", path=main)
52 vm.addAlias("isBremsCorrected", "extraInfo(bremsCorrected)")
53 
54 # combine final state particles to form composite particles
55 ma.reconstructDecay(
56  "J/psi:ee -> e+:corrected e-:corrected ?addbrems",
57  cut="dM < 0.11",
58  path=main,
59 )
60 
61 # perform vertex fit of J/psi candidates
62 vertex.kFit("J/psi:ee", conf_level=0.0, path=main)
63 
64 # combine J/psi and KS candidates to form B0 candidates
65 ma.reconstructDecay(
66  "B0 -> J/psi:ee K_S0:merged",
67  cut="Mbc > 5.2 and abs(deltaE) < 0.15",
68  path=main,
69 )
70 
71 # match reconstructed with MC particles
72 ma.matchMCTruth("B0", path=main)
73 
74 # build the rest of the event
75 ma.buildRestOfEvent("B0", fillWithMostLikely=True, path=main)
76 track_based_cuts = "thetaInCDCAcceptance and pt > 0.075 and dr < 5 and abs(dz) < 10"
77 ecl_based_cuts = "thetaInCDCAcceptance and E > 0.05"
78 roe_mask = ("my_mask", track_based_cuts, ecl_based_cuts)
79 ma.appendROEMasks("B0", [roe_mask], path=main)
80 
81 # call flavor tagging
82 ft.flavorTagger(["B0"], path=main)
83 
84 # fit B vertex on the tag-side
85 vertex.TagV("B0", fitAlgorithm="Rave", path=main)
86 
87 # perform best candidate selection
88 b2.set_random_seed("Belle II StarterKit")
89 ma.rankByHighest("B0", variable="random", numBest=1, path=main)
90 
91 # Create list of variables to save into the output file
92 b_vars = []
93 
94 standard_vars = vc.kinematics + vc.mc_kinematics + vc.mc_truth
95 b_vars += vc.deltae_mbc
96 b_vars += standard_vars
97 
98 # ROE variables
99 roe_kinematics = ["roeE()", "roeM()", "roeP()", "roeMbc()", "roeDeltae()"]
100 roe_multiplicities = [
101  "nROE_Charged()",
102  "nROE_Photons()",
103  "nROE_NeutralHadrons()",
104 ]
105 b_vars += roe_kinematics + roe_multiplicities
106 # Let's also add a version of the ROE variables that includes the mask:
107 for roe_variable in roe_kinematics + roe_multiplicities:
108  # e.g. instead of 'roeE()' (no mask) we want 'roeE(my_mask)'
109  roe_variable_with_mask = roe_variable.replace("()", "(my_mask)")
110  b_vars.append(roe_variable_with_mask)
111 
112 b_vars += ft.flavor_tagging
113 b_vars += vc.tag_vertex + vc.mc_tag_vertex
114 
115 # Variables for final states (electrons, positrons, pions)
116 fs_vars = vc.pid + vc.track + vc.track_hits + standard_vars
117 b_vars += vu.create_aliases_for_selected(
118  fs_vars + ["isBremsCorrected"],
119  "B0 -> [J/psi -> ^e+ ^e-] K_S0",
120  prefix=["ep", "em"],
121 )
122 b_vars += vu.create_aliases_for_selected(
123  fs_vars, "B0 -> J/psi [K_S0 -> ^pi+ ^pi-]", prefix=["pip", "pim"]
124 )
125 # Variables for J/Psi, KS
126 jpsi_ks_vars = vc.inv_mass + standard_vars
127 jpsi_ks_vars += vc.vertex + vc.mc_vertex
128 b_vars += vu.create_aliases_for_selected(jpsi_ks_vars, "B0 -> ^J/psi ^K_S0")
129 # Add the J/Psi mass calculated with uncorrected electrons:
130 vm.addAlias(
131  "Jpsi_M_uncorrected", "daughter(0, daughterCombination(M,0:0,1:0))"
132 )
133 b_vars += ["Jpsi_M_uncorrected"]
134 # Also add kinematic variables boosted to the center of mass frame (CMS)
135 # for all particles
136 cmskinematics = vu.create_aliases(
137  vc.kinematics, "useCMSFrame({variable})", "CMS"
138 )
139 b_vars += vu.create_aliases_for_selected(
140  cmskinematics, "^B0 -> [^J/psi -> ^e+ ^e-] [^K_S0 -> ^pi+ ^pi-]"
141 )
142 
143 vm.addAlias(
144  "withBremsCorrection",
145  "passesCut(passesCut(ep_isBremsCorrected == 1) or passesCut(em_isBremsCorrected == 1))",
146 )
147 b_vars += ["withBremsCorrection"]
148 
149 # Save variables to an output file (ntuple)
150 ma.variablesToNtuple(
151  "B0",
152  variables=b_vars,
153  filename="Bd2JpsiKS.root",
154  treename="tree",
155  path=main,
156 )
157 
158 # Start the event loop (actually start processing things)
159 b2.process(main)
160 
161 # print out the summary
162 print(b2.statistics)
def stdKshorts(prioritiseV0=True, fitter='TreeFit', path=None)
Definition: stdV0s.py:17
def TagV(list_name, MCassociation='', confidenceLevel=0., trackFindingType="standard_PXD", constraintType="IP", askMCInfo=False, reqPXDHits=0, maskName='', fitAlgorithm='Rave', useTruthInFit=False, useRollBack=False, path=None)
Definition: vertex.py:286
def kFit(list_name, conf_level, fit_type='vertex', constraint='', daughtersUpdate=False, decay_string='', massConstraint=[], smearing=0, path=None)
Definition: vertex.py:121