Belle II Software development
treeFit_beam_momentum_constrain.py
1#!/usr/bin/env python3
2
3
10import unittest
11import tempfile
12import basf2
13import b2test_utils
14import modularAnalysis as ma
15from variables import variables as vm
16from ROOT import TFile
17
18
19basf2.set_random_seed('aSeed')
20
21
22class TestTreeFits(unittest.TestCase):
23 """The unit test case for TreeFitter"""
24
25 def testFit(self):
26 """Run the test fit"""
27
28 testFile = tempfile.NamedTemporaryFile()
29
30 main = basf2.create_path()
31
32 ma.inputMdst(b2test_utils.require_file('analysis/tests/150_noBKG_DtoPiNuNu.root'), path=main)
33 ma.fillParticleList('pi+:a', 'pionID > 0.5', path=main)
34 ma.fillParticleList('K+:a', 'kaonID > 0.5', path=main)
35
36 ma.reconstructDecay(
37 "D-:tag -> K+:a pi-:a pi-:a",
38 cut='1.859 < M < 1.879 and useCMSFrame(p)>2.0',
39 path=main
40 )
41
42 ma.reconstructDecay(
43 "D+:sig -> pi+:a",
44 '',
45 path=main
46 )
47
48 ma.reconstructDecay('Z0:rec -> D+:sig D-:tag pi-:a pi+:a', '', path=main)
49 ma.matchMCTruth('Z0:rec', path=main)
50
51 ma.setBeamConstrainedMomentum(
52 'Z0:rec',
53 '^Z0:rec -> D+:sig D-:tag pi-:a pi+:a',
54 'Z0:rec -> D+:sig D-:tag pi-:a pi+:a',
55 path=main
56 )
57
58 ma.setBeamConstrainedMomentum(
59 'Z0:rec',
60 'Z0:rec -> ^D+:sig D-:tag pi-:a pi+:a',
61 'Z0:rec -> D+:sig ^D-:tag ^pi-:a ^pi+:a',
62 path=main
63 )
64
65 conf = 0
66 main.add_module('TreeFitter',
67 particleList='Z0:rec',
68 confidenceLevel=conf,
69 massConstraintList=[],
70 massConstraintListParticlename=[],
71 expertBeamConstraintPDG=23,
72 treatAsInvisible='Z0:rec -> ^D+:sig D-:tag pi-:a pi+:a',
73 ipConstraint=True,
74 updateAllDaughters=True)
75
76 ma.applyCuts('Z0:rec',
77 "1.6<daughter(0,M)<2.2",
78 path=main
79 )
80
81 ma.rankByHighest("Z0:rec",
82 "chiProb",
83 numBest=1,
84 path=main)
85
86 vm.addAlias("signalD", "daughter(0,isSignalAcceptMissingNeutrino)")
87
88 ntupler = basf2.register_module('VariablesToNtuple')
89 ntupler.param('fileName', testFile.name)
90 ntupler.param('variables', ['chiProb', 'M', 'signalD'])
91 ntupler.param('particleList', 'Z0:rec')
92 main.add_module(ntupler)
93
94 basf2.process(main)
95
96 ntuplefile = TFile(testFile.name)
97 ntuple = ntuplefile.Get('ntuple')
98
99 self.assertFalse(ntuple.GetEntries() == 0, "Ntuple is empty.")
100
101 allBkg = ntuple.GetEntries("signalD == 0")
102 allSig = ntuple.GetEntries("signalD > 0")
103
104 truePositives = ntuple.GetEntries("(chiProb > 0) && (signalD > 0)")
105 falsePositives = ntuple.GetEntries("(chiProb > 0) && (signalD == 0)")
106
107 mustBeZero = ntuple.GetEntries(f"(chiProb < {conf})")
108
109 print(f"True fit survivors: {truePositives} out of {allSig} true candidates")
110 print(f"False fit survivors: {falsePositives} out of {allBkg} false candidates")
111
112 self.assertFalse(truePositives == 0, "No signal survived the fit.")
113
114 self.assertTrue(falsePositives < 1, "Background rejection too small.")
115
116 self.assertTrue(truePositives > 8, "Signal rejection too high")
117 self.assertFalse(mustBeZero, f"We should have dropped all candidates with confidence level less than {conf}.")
118
119 print("Test passed, cleaning up.")
120
121
122if __name__ == '__main__':
124 unittest.main()
def require_file(filename, data_type="", py_case=None)
Definition: __init__.py:54
def clean_working_directory()
Definition: __init__.py:189