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