Belle II Software  light-2205-abys
flavorTagger.py
1 #!/usr/bin/env python3
2 
3 
10 
11 # ************* Flavor Tagging ************
12 # * This script is needed to train *
13 # * and to test the flavor tagger. *
14 # ********************************************
15 
16 from basf2 import B2INFO, B2FATAL, B2WARNING
17 import basf2
18 import basf2_mva
19 import inspect
20 import modularAnalysis as ma
21 import variables
22 from variables import utils
23 import os
24 import glob
25 import b2bii
26 import collections
27 
28 
29 def getBelleOrBelle2():
30  """
31  Gets the global ModeCode.
32  """
33  if b2bii.isB2BII():
34  return 'Belle'
35  else:
36  return 'Belle2'
37 
38 
39 def setInteractionWithDatabase(downloadFromDatabaseIfNotFound=False, uploadToDatabaseAfterTraining=False):
40  """
41  Sets the interaction with the database: download trained weight files or upload weight files after training.
42  """
43 
44  global downloadFlag
45  global uploadFlag
46 
47  downloadFlag = downloadFromDatabaseIfNotFound
48  uploadFlag = uploadToDatabaseAfterTraining
49 
50 
51 # Default list of aliases that should be used to save the flavor tagging information using VariablesToNtuple
52 flavor_tagging = ['FBDT_qrCombined', 'FANN_qrCombined', 'qrMC', 'mcFlavorOfOtherB',
53  'qpElectron', 'hasTrueTargetElectron', 'isRightCategoryElectron',
54  'qpIntermediateElectron', 'hasTrueTargetIntermediateElectron', 'isRightCategoryIntermediateElectron',
55  'qpMuon', 'hasTrueTargetMuon', 'isRightCategoryMuon',
56  'qpIntermediateMuon', 'hasTrueTargetIntermediateMuon', 'isRightCategoryIntermediateMuon',
57  'qpKinLepton', 'hasTrueTargetKinLepton', 'isRightCategoryKinLepton',
58  'qpIntermediateKinLepton', 'hasTrueTargetIntermediateKinLepton', 'isRightCategoryIntermediateKinLepton',
59  'qpKaon', 'hasTrueTargetKaon', 'isRightCategoryKaon',
60  'qpSlowPion', 'hasTrueTargetSlowPion', 'isRightCategorySlowPion',
61  'qpFastHadron', 'hasTrueTargetFastHadron', 'isRightCategoryFastHadron',
62  'qpLambda', 'hasTrueTargetLambda', 'isRightCategoryLambda',
63  'qpFSC', 'hasTrueTargetFSC', 'isRightCategoryFSC',
64  'qpMaximumPstar', 'hasTrueTargetMaximumPstar', 'isRightCategoryMaximumPstar',
65  'qpKaonPion', 'hasTrueTargetKaonPion', 'isRightCategoryKaonPion']
66 
67 
68 def add_default_FlavorTagger_aliases():
69  """
70  This function adds the default aliases for flavor tagging variables
71  and defines the collection of flavor tagging variables.
72  """
73 
74  variables.variables.addAlias('FBDT_qrCombined', 'qrOutput(FBDT)')
75  variables.variables.addAlias('FANN_qrCombined', 'qrOutput(FANN)')
76  variables.variables.addAlias('qrMC', 'isRelatedRestOfEventB0Flavor')
77 
78  for iCategory in AvailableCategories:
79  aliasForQp = 'qp' + iCategory
80  aliasForTrueTarget = 'hasTrueTarget' + iCategory
81  aliasForIsRightCategory = 'isRightCategory' + iCategory
82  variables.variables.addAlias(aliasForQp, 'qpCategory(' + iCategory + ')')
83  variables.variables.addAlias(aliasForTrueTarget, 'hasTrueTargets(' + iCategory + ')')
84  variables.variables.addAlias(aliasForIsRightCategory, 'isTrueFTCategory(' + iCategory + ')')
85 
86  utils.add_collection(flavor_tagging, 'flavor_tagging')
87 
88 
89 def set_FlavorTagger_pid_aliases():
90  """
91  This function adds the pid aliases needed by the flavor tagger.
92  """
93  variables.variables.addAlias('eid_TOP', 'ifNANgiveX(pidPairProbabilityExpert(11, 211, TOP), 0.5)')
94  variables.variables.addAlias('eid_ARICH', 'ifNANgiveX(pidPairProbabilityExpert(11, 211, ARICH), 0.5)')
95  variables.variables.addAlias('eid_ECL', 'ifNANgiveX(pidPairProbabilityExpert(11, 211, ECL), 0.5)')
96 
97  variables.variables.addAlias('muid_TOP', 'ifNANgiveX(pidPairProbabilityExpert(13, 211, TOP), 0.5)')
98  variables.variables.addAlias('muid_ARICH', 'ifNANgiveX(pidPairProbabilityExpert(13, 211, ARICH), 0.5)')
99  variables.variables.addAlias('muid_KLM', 'ifNANgiveX(pidPairProbabilityExpert(13, 211, KLM), 0.5)')
100 
101  variables.variables.addAlias('piid_TOP', 'ifNANgiveX(pidPairProbabilityExpert(211, 321, TOP), 0.5)')
102  variables.variables.addAlias('piid_ARICH', 'ifNANgiveX(pidPairProbabilityExpert(211, 321, ARICH), 0.5)')
103 
104  variables.variables.addAlias('Kid_TOP', 'ifNANgiveX(pidPairProbabilityExpert(321, 211, TOP), 0.5)')
105  variables.variables.addAlias('Kid_ARICH', 'ifNANgiveX(pidPairProbabilityExpert(321, 211, ARICH), 0.5)')
106 
107  if getBelleOrBelle2() == "Belle":
108  variables.variables.addAlias('eid_dEdx', 'ifNANgiveX(pidPairProbabilityExpert(11, 211, CDC, SVD), 0.5)')
109  variables.variables.addAlias('muid_dEdx', 'ifNANgiveX(pidPairProbabilityExpert(13, 211, CDC, SVD), 0.5)')
110  variables.variables.addAlias('piid_dEdx', 'ifNANgiveX(pidPairProbabilityExpert(211, 321, CDC, SVD), 0.5)')
111  variables.variables.addAlias('pi_vs_edEdxid', 'ifNANgiveX(pidPairProbabilityExpert(211, 11, CDC, SVD), 0.5)')
112  variables.variables.addAlias('Kid_dEdx', 'ifNANgiveX(pidPairProbabilityExpert(321, 211, CDC, SVD), 0.5)')
113  else:
114  # Removed SVD PID for Belle II MC and data as it is absent in release 4.
115  variables.variables.addAlias('eid_dEdx', 'ifNANgiveX(pidPairProbabilityExpert(11, 211, CDC), 0.5)')
116  variables.variables.addAlias('muid_dEdx', 'ifNANgiveX(pidPairProbabilityExpert(13, 211, CDC), 0.5)')
117  variables.variables.addAlias('piid_dEdx', 'ifNANgiveX(pidPairProbabilityExpert(211, 321, CDC), 0.5)')
118  variables.variables.addAlias('pi_vs_edEdxid', 'ifNANgiveX(pidPairProbabilityExpert(211, 11, CDC), 0.5)')
119  variables.variables.addAlias('Kid_dEdx', 'ifNANgiveX(pidPairProbabilityExpert(321, 211, CDC), 0.5)')
120 
121 
122 def setInputVariablesWithMask(maskName='all'):
123  """
124  Set aliases for input variables with ROE mask.
125  """
126  variables.variables.addAlias('pMissTag_withMask', 'pMissTag('+maskName+')')
127  variables.variables.addAlias('cosTPTO_withMask', 'cosTPTO('+maskName+')')
128  variables.variables.addAlias('ptTracksRoe_withMask', 'ptTracksRoe('+maskName+')')
129  variables.variables.addAlias('pt2TracksRoe_withMask', 'pt2TracksRoe('+maskName+')')
130  variables.variables.addAlias('ptTracksRoe_withMask', 'ptTracksRoe('+maskName+')')
131 
132 
133 def getFastBDTCategories():
134  '''
135  Helper function for getting the FastBDT categories.
136  It's necessary for removing top-level ROOT imports.
137  '''
138  import ROOT # noqa
139  fastBDTCategories = basf2_mva.FastBDTOptions()
140  fastBDTCategories.m_nTrees = 500
141  fastBDTCategories.m_nCuts = 8
142  fastBDTCategories.m_nLevels = 3
143  fastBDTCategories.m_shrinkage = 0.10
144  fastBDTCategories.m_randRatio = 0.5
145  return fastBDTCategories
146 
147 
148 def getFastBDTCombiner():
149  '''
150  Helper function for getting the FastBDT combiner.
151  It's necessary for removing top-level ROOT imports.
152  '''
153  import ROOT # noqa
154  fastBDTCombiner = basf2_mva.FastBDTOptions()
155  fastBDTCombiner.m_nTrees = 500
156  fastBDTCombiner.m_nCuts = 8
157  fastBDTCombiner.m_nLevels = 3
158  fastBDTCombiner.m_shrinkage = 0.10
159  fastBDTCombiner.m_randRatio = 0.5
160  return fastBDTCombiner
161 
162 
163 def getMlpFANNCombiner():
164  '''
165  Helper function for getting the MLP FANN combiner.
166  It's necessary for removing top-level ROOT imports.
167  '''
168  import ROOT # noqa
169  mlpFANNCombiner = basf2_mva.FANNOptions()
170  mlpFANNCombiner.m_max_epochs = 10000
171  mlpFANNCombiner.m_hidden_layers_architecture = "3*N"
172  mlpFANNCombiner.m_hidden_activiation_function = "FANN_SIGMOID_SYMMETRIC"
173  mlpFANNCombiner.m_output_activiation_function = "FANN_SIGMOID_SYMMETRIC"
174  mlpFANNCombiner.m_error_function = "FANN_ERRORFUNC_LINEAR"
175  mlpFANNCombiner.m_training_method = "FANN_TRAIN_RPROP"
176  mlpFANNCombiner.m_validation_fraction = 0.5
177  mlpFANNCombiner.m_random_seeds = 10
178  mlpFANNCombiner.m_test_rate = 500
179  mlpFANNCombiner.m_number_of_threads = 8
180  mlpFANNCombiner.m_scale_features = True
181  mlpFANNCombiner.m_scale_target = False
182  # mlpFANNCombiner.m_scale_target = True
183  return mlpFANNCombiner
184 
185 
186 # SignalFraction: FBDT feature
187 # For smooth output set to -1, this will break the calibration.
188 # For correct calibration set to -2, leads to peaky combiner output.
189 signalFraction = -2
190 
191 # Maximal number of events to train each method
192 maxEventsNumber = 0 # 0 takes all the sampled events. The number in the past was 500000
193 
194 
195 FTCategoryParameters = collections.namedtuple('FTCategoryParameters',
196  ['particleList', 'trackName', 'eventName', 'variableName', 'code'])
197 # Definition of all available categories, 'standard category name':
198 # ['ParticleList', 'trackLevel category name', 'eventLevel category name',
199 # 'combinerLevel variable name', 'category code']
200 AvailableCategories = {
201  'Electron':
202  FTCategoryParameters('e+:inRoe', 'Electron', 'Electron',
203  'QpOf(e+:inRoe, isRightCategory(Electron), isRightCategory(Electron))',
204  0),
205  'IntermediateElectron':
206  FTCategoryParameters('e+:inRoe', 'IntermediateElectron', 'IntermediateElectron',
207  'QpOf(e+:inRoe, isRightCategory(IntermediateElectron), isRightCategory(IntermediateElectron))',
208  1),
209  'Muon':
210  FTCategoryParameters('mu+:inRoe', 'Muon', 'Muon',
211  'QpOf(mu+:inRoe, isRightCategory(Muon), isRightCategory(Muon))',
212  2),
213  'IntermediateMuon':
214  FTCategoryParameters('mu+:inRoe', 'IntermediateMuon', 'IntermediateMuon',
215  'QpOf(mu+:inRoe, isRightCategory(IntermediateMuon), isRightCategory(IntermediateMuon))',
216  3),
217  'KinLepton':
218  FTCategoryParameters('mu+:inRoe', 'KinLepton', 'KinLepton',
219  'QpOf(mu+:inRoe, isRightCategory(KinLepton), isRightCategory(KinLepton))',
220  4),
221  'IntermediateKinLepton':
222  FTCategoryParameters('mu+:inRoe', 'IntermediateKinLepton', 'IntermediateKinLepton',
223  'QpOf(mu+:inRoe, isRightCategory(IntermediateKinLepton), isRightCategory(IntermediateKinLepton))',
224  5),
225  'Kaon':
226  FTCategoryParameters('K+:inRoe', 'Kaon', 'Kaon',
227  'weightedQpOf(K+:inRoe, isRightCategory(Kaon), isRightCategory(Kaon))',
228  6),
229  'SlowPion':
230  FTCategoryParameters('pi+:inRoe', 'SlowPion', 'SlowPion',
231  'QpOf(pi+:inRoe, isRightCategory(SlowPion), isRightCategory(SlowPion))',
232  7),
233  'FastHadron':
234  FTCategoryParameters('pi+:inRoe', 'FastHadron', 'FastHadron',
235  'QpOf(pi+:inRoe, isRightCategory(FastHadron), isRightCategory(FastHadron))',
236  8),
237  'Lambda':
238  FTCategoryParameters('Lambda0:inRoe', 'Lambda', 'Lambda',
239  'weightedQpOf(Lambda0:inRoe, isRightCategory(Lambda), isRightCategory(Lambda))',
240  9),
241  'FSC':
242  FTCategoryParameters('pi+:inRoe', 'SlowPion', 'FSC',
243  'QpOf(pi+:inRoe, isRightCategory(FSC), isRightCategory(SlowPion))',
244  10),
245  'MaximumPstar':
246  FTCategoryParameters('pi+:inRoe', 'MaximumPstar', 'MaximumPstar',
247  'QpOf(pi+:inRoe, isRightCategory(MaximumPstar), isRightCategory(MaximumPstar))',
248  11),
249  'KaonPion':
250  FTCategoryParameters('K+:inRoe', 'Kaon', 'KaonPion',
251  'QpOf(K+:inRoe, isRightCategory(KaonPion), isRightCategory(Kaon))',
252  12),
253 }
254 
255 
256 def getTrainingVariables(category=None):
257  """
258  Helper function to get training variables.
259 
260  NOTE: This function is not called the Expert mode. It is not necessary to be consistent with variables list of weight files.
261  """
262 
263  KId = {'Belle': 'ifNANgiveX(atcPIDBelle(3,2), 0.5)', 'Belle2': 'kaonID'}
264  muId = {'Belle': 'muIDBelle', 'Belle2': 'muonID'}
265  eId = {'Belle': 'eIDBelle', 'Belle2': 'electronID'}
266 
267  variables = []
268  if category == 'Electron' or category == 'IntermediateElectron':
269  variables = ['useCMSFrame(p)',
270  'useCMSFrame(pt)',
271  'p',
272  'pt',
273  'cosTheta',
274  eId[getBelleOrBelle2()],
275  'eid_TOP',
276  'eid_ARICH',
277  'eid_ECL',
278  'BtagToWBosonVariables(recoilMassSqrd)',
279  'BtagToWBosonVariables(pMissCMS)',
280  'BtagToWBosonVariables(cosThetaMissCMS)',
281  'BtagToWBosonVariables(EW90)',
282  'cosTPTO',
283  'chiProb',
284  ]
285  if getBelleOrBelle2() == "Belle":
286  variables.append('eid_dEdx')
287  variables.append('ImpactXY')
288  variables.append('distance')
289 
290  elif category == 'Muon' or category == 'IntermediateMuon':
291  variables = ['useCMSFrame(p)',
292  'useCMSFrame(pt)',
293  'p',
294  'pt',
295  'cosTheta',
296  muId[getBelleOrBelle2()],
297  'muid_TOP',
298  'muid_ARICH',
299  'muid_KLM',
300  'BtagToWBosonVariables(recoilMassSqrd)',
301  'BtagToWBosonVariables(pMissCMS)',
302  'BtagToWBosonVariables(cosThetaMissCMS)',
303  'BtagToWBosonVariables(EW90)',
304  'cosTPTO',
305  ]
306  if getBelleOrBelle2() == "Belle":
307  variables.append('muid_dEdx')
308  variables.append('ImpactXY')
309  variables.append('distance')
310  variables.append('chiProb')
311 
312  elif category == 'KinLepton' or category == 'IntermediateKinLepton':
313  variables = ['useCMSFrame(p)',
314  'useCMSFrame(pt)',
315  'p',
316  'pt',
317  'cosTheta',
318  muId[getBelleOrBelle2()],
319  'muid_TOP',
320  'muid_ARICH',
321  'muid_KLM',
322  eId[getBelleOrBelle2()],
323  'eid_TOP',
324  'eid_ARICH',
325  'eid_ECL',
326  'BtagToWBosonVariables(recoilMassSqrd)',
327  'BtagToWBosonVariables(pMissCMS)',
328  'BtagToWBosonVariables(cosThetaMissCMS)',
329  'BtagToWBosonVariables(EW90)',
330  'cosTPTO',
331  ]
332  if getBelleOrBelle2() == "Belle":
333  variables.append('eid_dEdx')
334  variables.append('muid_dEdx')
335  variables.append('ImpactXY')
336  variables.append('distance')
337  variables.append('chiProb')
338 
339  elif category == 'Kaon':
340  variables = ['useCMSFrame(p)',
341  'useCMSFrame(pt)',
342  'p',
343  'pt',
344  'cosTheta',
345  KId[getBelleOrBelle2()],
346  'Kid_dEdx',
347  'Kid_TOP',
348  'Kid_ARICH',
349  'NumberOfKShortsInRoe',
350  'ptTracksRoe',
351  'BtagToWBosonVariables(recoilMassSqrd)',
352  'BtagToWBosonVariables(pMissCMS)',
353  'BtagToWBosonVariables(cosThetaMissCMS)',
354  'BtagToWBosonVariables(EW90)',
355  'cosTPTO',
356  'chiProb',
357  ]
358  if getBelleOrBelle2() == "Belle":
359  variables.append('ImpactXY')
360  variables.append('distance')
361 
362  elif category == 'SlowPion':
363  variables = ['useCMSFrame(p)',
364  'useCMSFrame(pt)',
365  'cosTheta',
366  'p',
367  'pt',
368  'pionID',
369  'piid_TOP',
370  'piid_ARICH',
371  'pi_vs_edEdxid',
372  KId[getBelleOrBelle2()],
373  'Kid_dEdx',
374  'Kid_TOP',
375  'Kid_ARICH',
376  'NumberOfKShortsInRoe',
377  'ptTracksRoe',
378  eId[getBelleOrBelle2()],
379  'BtagToWBosonVariables(recoilMassSqrd)',
380  'BtagToWBosonVariables(EW90)',
381  'BtagToWBosonVariables(cosThetaMissCMS)',
382  'BtagToWBosonVariables(pMissCMS)',
383  'cosTPTO',
384  ]
385  if getBelleOrBelle2() == "Belle":
386  variables.append('piid_dEdx')
387  variables.append('ImpactXY')
388  variables.append('distance')
389  variables.append('chiProb')
390 
391  elif category == 'FastHadron':
392  variables = ['useCMSFrame(p)',
393  'useCMSFrame(pt)',
394  'cosTheta',
395  'p',
396  'pt',
397  'pionID',
398  'piid_dEdx',
399  'piid_TOP',
400  'piid_ARICH',
401  'pi_vs_edEdxid',
402  KId[getBelleOrBelle2()],
403  'Kid_dEdx',
404  'Kid_TOP',
405  'Kid_ARICH',
406  'NumberOfKShortsInRoe',
407  'ptTracksRoe',
408  eId[getBelleOrBelle2()],
409  'BtagToWBosonVariables(recoilMassSqrd)',
410  'BtagToWBosonVariables(EW90)',
411  'BtagToWBosonVariables(cosThetaMissCMS)',
412  'cosTPTO',
413  ]
414  if getBelleOrBelle2() == "Belle":
415  variables.append('BtagToWBosonVariables(pMissCMS)')
416  variables.append('ImpactXY')
417  variables.append('distance')
418  variables.append('chiProb')
419 
420  elif category == 'Lambda':
421  variables = ['lambdaFlavor',
422  'NumberOfKShortsInRoe',
423  'M',
424  'cosAngleBetweenMomentumAndVertexVector',
425  'lambdaZError',
426  'daughter(0,p)',
427  'daughter(0,useCMSFrame(p))',
428  'daughter(1,p)',
429  'daughter(1,useCMSFrame(p))',
430  'useCMSFrame(p)',
431  'p',
432  'chiProb',
433  ]
434  if getBelleOrBelle2() == "Belle2":
435  variables.append('daughter(1,protonID)') # protonID always 0 in B2BII check in future
436  variables.append('daughter(0,pionID)') # not very powerful in B2BII
437  else:
438  variables.append('distance')
439 
440  elif category == 'MaximumPstar':
441  variables = ['useCMSFrame(p)',
442  'useCMSFrame(pt)',
443  'p',
444  'pt',
445  'cosTPTO',
446  ]
447  if getBelleOrBelle2() == "Belle2":
448  variables.append('ImpactXY')
449  variables.append('distance')
450 
451  elif category == 'FSC':
452  variables = ['useCMSFrame(p)',
453  'cosTPTO',
454  KId[getBelleOrBelle2()],
455  'FSCVariables(pFastCMS)',
456  'FSCVariables(cosSlowFast)',
457  'FSCVariables(cosTPTOFast)',
458  'FSCVariables(SlowFastHaveOpositeCharges)',
459  ]
460  elif category == 'KaonPion':
461  variables = ['extraInfo(isRightCategory(Kaon))',
462  'HighestProbInCat(pi+:inRoe, isRightCategory(SlowPion))',
463  'KaonPionVariables(cosKaonPion)',
464  'KaonPionVariables(HaveOpositeCharges)',
465  KId[getBelleOrBelle2()]
466  ]
467 
468  return variables
469 
470 
471 def FillParticleLists(maskName='all', categories=None, path=None):
472  """
473  Fills the particle Lists for all categories.
474  """
475 
476  from vertex import kFit
477  readyParticleLists = []
478 
479  if categories is None:
480  categories = []
481 
482  trackCut = 'isInRestOfEvent > 0.5 and passesROEMask(' + maskName + ') > 0.5 and p >= 0'
483 
484  for category in categories:
485  particleList = AvailableCategories[category].particleList
486 
487  if particleList in readyParticleLists:
488  continue
489 
490  # Select particles in ROE for different categories according to mass hypothesis.
491  if particleList == 'Lambda0:inRoe':
492  if 'pi+:inRoe' not in readyParticleLists:
493  ma.fillParticleList('pi+:inRoe', trackCut, path=path)
494  readyParticleLists.append('pi+:inRoe')
495 
496  ma.fillParticleList('p+:inRoe', trackCut, path=path)
497  ma.reconstructDecay(particleList + ' -> pi-:inRoe p+:inRoe', '1.00<=M<=1.23', False, path=path)
498  kFit(particleList, 0.01, path=path)
499  ma.matchMCTruth(particleList, path=path)
500  readyParticleLists.append(particleList)
501 
502  else:
503  # Filling particle list for actual category
504  ma.fillParticleList(particleList, trackCut, path=path)
505  readyParticleLists.append(particleList)
506 
507  # Additional particleLists for K_S0
508  if getBelleOrBelle2() == 'Belle':
509  ma.cutAndCopyList('K_S0:inRoe', 'K_S0:mdst', 'extraInfo(ksnbStandard) == 1 and isInRestOfEvent == 1', path=path)
510  else:
511  if 'pi+:inRoe' not in readyParticleLists:
512  ma.fillParticleList('pi+:inRoe', trackCut, path=path)
513  ma.reconstructDecay('K_S0:inRoe -> pi+:inRoe pi-:inRoe', '0.40<=M<=0.60', False, path=path)
514  kFit('K_S0:inRoe', 0.01, path=path)
515 
516  # Apply BDT-based LID
517  if getBelleOrBelle2() == 'Belle2':
518  default_list_for_lid_BDT = ['e+:inRoe', 'mu+:inRoe']
519  list_for_lid_BDT = []
520 
521  for particleList in default_list_for_lid_BDT:
522  if particleList in readyParticleLists:
523  list_for_lid_BDT.append(particleList)
524 
525  if list_for_lid_BDT: # empty check
526  ma.applyChargedPidMVA(particleLists=list_for_lid_BDT, path=path,
527  trainingMode=0, # binary
528  binaryHypoPDGCodes=(11, 211)) # e vs pi
529  ma.applyChargedPidMVA(particleLists=list_for_lid_BDT, path=path,
530  trainingMode=0, # binary
531  binaryHypoPDGCodes=(13, 211)) # mu vs pi
532  ma.applyChargedPidMVA(particleLists=list_for_lid_BDT, path=path,
533  trainingMode=1) # Multiclass
534 
535 
536 def eventLevel(mode='Expert', weightFiles='B2JpsiKs_mu', categories=None, path=None):
537  """
538  Samples data for training or tests all categories all categories at event level.
539  """
540 
541  from basf2 import create_path
542  from basf2 import register_module
543 
544  B2INFO('EVENT LEVEL')
545 
546  ReadyMethods = 0
547 
548  # Each category has its own Path in order to be skipped if the corresponding particle list is empty
549  identifiersExtraInfosDict = dict()
550  identifiersExtraInfosKaonPion = []
551 
552  if categories is None:
553  categories = []
554 
555  for category in categories:
556  particleList = AvailableCategories[category].particleList
557 
558  methodPrefixEventLevel = "FlavorTagger_" + getBelleOrBelle2() + "_" + weightFiles + 'EventLevel' + category + 'FBDT'
559  identifierEventLevel = methodPrefixEventLevel
560  targetVariable = 'isRightCategory(' + category + ')'
561  extraInfoName = targetVariable
562 
563  if mode == 'Expert':
564 
565  if downloadFlag or useOnlyLocalFlag:
566  identifierEventLevel = filesDirectory + '/' + methodPrefixEventLevel + '_1.root'
567 
568  if downloadFlag:
569  if not os.path.isfile(identifierEventLevel):
570  basf2_mva.download(methodPrefixEventLevel, identifierEventLevel)
571  if not os.path.isfile(identifierEventLevel):
572  B2FATAL('Flavor Tagger: Weight file ' + identifierEventLevel +
573  ' was not downloaded from Database. Please check the buildOrRevision name. Stopped')
574 
575  if useOnlyLocalFlag:
576  if not os.path.isfile(identifierEventLevel):
577  B2FATAL('Flavor Tagger: ' + particleList + ' Eventlevel was not trained. Weight file ' +
578  identifierEventLevel + ' was not found. Stopped')
579 
580  B2INFO('flavorTagger: MVAExpert ' + methodPrefixEventLevel + ' ready.')
581 
582  elif mode == 'Sampler':
583 
584  identifierEventLevel = filesDirectory + '/' + methodPrefixEventLevel + '_1.root'
585  if os.path.isfile(identifierEventLevel):
586  B2INFO('flavorTagger: MVAExpert ' + methodPrefixEventLevel + ' ready.')
587 
588  if 'KaonPion' in categories:
589  methodPrefixEventLevelKaonPion = "FlavorTagger_" + getBelleOrBelle2() + \
590  "_" + weightFiles + 'EventLevelKaonPionFBDT'
591  identifierEventLevelKaonPion = filesDirectory + '/' + methodPrefixEventLevelKaonPion + '_1.root'
592  if not os.path.isfile(identifierEventLevelKaonPion):
593  # Slow Pion and Kaon categories are used if Kaon-Pion is lacking for
594  # sampling. The others are not needed and skipped
595  if category != "SlowPion" and category != "Kaon":
596  continue
597 
598  if mode == 'Expert' or (mode == 'Sampler' and os.path.isfile(identifierEventLevel)):
599 
600  B2INFO('flavorTagger: Applying MVAExpert ' + methodPrefixEventLevel + '.')
601 
602  if category == 'KaonPion':
603  identifiersExtraInfosKaonPion.append((extraInfoName, identifierEventLevel))
604  elif particleList not in identifiersExtraInfosDict:
605  identifiersExtraInfosDict[particleList] = [(extraInfoName, identifierEventLevel)]
606  else:
607  identifiersExtraInfosDict[particleList].append((extraInfoName, identifierEventLevel))
608 
609  ReadyMethods += 1
610 
611  # Each category has its own Path in order to be skipped if the corresponding particle list is empty
612  for particleList in identifiersExtraInfosDict:
613  eventLevelPath = create_path()
614  SkipEmptyParticleList = register_module("SkimFilter")
615  SkipEmptyParticleList.set_name('SkimFilter_EventLevel_' + particleList)
616  SkipEmptyParticleList.param('particleLists', particleList)
617  SkipEmptyParticleList.if_true(eventLevelPath, basf2.AfterConditionPath.CONTINUE)
618  path.add_module(SkipEmptyParticleList)
619 
620  mvaMultipleExperts = register_module('MVAMultipleExperts')
621  mvaMultipleExperts.set_name('MVAMultipleExperts_EventLevel_' + particleList)
622  mvaMultipleExperts.param('listNames', [particleList])
623  mvaMultipleExperts.param('extraInfoNames', [row[0] for row in identifiersExtraInfosDict[particleList]])
624  mvaMultipleExperts.param('signalFraction', signalFraction)
625  mvaMultipleExperts.param('identifiers', [row[1] for row in identifiersExtraInfosDict[particleList]])
626  eventLevelPath.add_module(mvaMultipleExperts)
627 
628  if 'KaonPion' in categories and len(identifiersExtraInfosKaonPion) != 0:
629  eventLevelKaonPionPath = create_path()
630  SkipEmptyParticleList = register_module("SkimFilter")
631  SkipEmptyParticleList.set_name('SkimFilter_' + 'K+:inRoe')
632  SkipEmptyParticleList.param('particleLists', 'K+:inRoe')
633  SkipEmptyParticleList.if_true(eventLevelKaonPionPath, basf2.AfterConditionPath.CONTINUE)
634  path.add_module(SkipEmptyParticleList)
635 
636  mvaExpertKaonPion = register_module("MVAExpert")
637  mvaExpertKaonPion.set_name('MVAExpert_KaonPion_' + 'K+:inRoe')
638  mvaExpertKaonPion.param('listNames', ['K+:inRoe'])
639  mvaExpertKaonPion.param('extraInfoName', identifiersExtraInfosKaonPion[0][0])
640  mvaExpertKaonPion.param('signalFraction', signalFraction)
641  mvaExpertKaonPion.param('identifier', identifiersExtraInfosKaonPion[0][1])
642 
643  eventLevelKaonPionPath.add_module(mvaExpertKaonPion)
644 
645  if mode == 'Sampler':
646 
647  for category in categories:
648  particleList = AvailableCategories[category].particleList
649 
650  methodPrefixEventLevel = "FlavorTagger_" + getBelleOrBelle2() + "_" + weightFiles + 'EventLevel' + category + 'FBDT'
651  identifierEventLevel = filesDirectory + '/' + methodPrefixEventLevel + '_1.root'
652  targetVariable = 'isRightCategory(' + category + ')'
653 
654  if not os.path.isfile(identifierEventLevel):
655 
656  if category == 'KaonPion':
657  methodPrefixEventLevelSlowPion = "FlavorTagger_" + getBelleOrBelle2() + \
658  "_" + weightFiles + 'EventLevelSlowPionFBDT'
659  identifierEventLevelSlowPion = filesDirectory + '/' + methodPrefixEventLevelSlowPion + '_1.root'
660  if not os.path.isfile(identifierEventLevelSlowPion):
661  B2INFO("Flavor Tagger: event level weight file for the Slow Pion category is absent." +
662  "It is required to sample the training information for the KaonPion category." +
663  "An additional sampling step will be needed after the following training step.")
664  continue
665 
666  B2INFO('flavorTagger: file ' + filesDirectory + '/' +
667  methodPrefixEventLevel + "sampled" + fileId + '.root will be saved.')
668 
669  ma.applyCuts(particleList, 'isRightCategory(mcAssociated) > 0', path)
670  eventLevelpath = create_path()
671  SkipEmptyParticleList = register_module("SkimFilter")
672  SkipEmptyParticleList.set_name('SkimFilter_EventLevel' + category)
673  SkipEmptyParticleList.param('particleLists', particleList)
674  SkipEmptyParticleList.if_true(eventLevelpath, basf2.AfterConditionPath.CONTINUE)
675  path.add_module(SkipEmptyParticleList)
676 
677  ntuple = register_module('VariablesToNtuple')
678  ntuple.param('fileName', filesDirectory + '/' + methodPrefixEventLevel + "sampled" + fileId + ".root")
679  ntuple.param('treeName', methodPrefixEventLevel + "_tree")
680  variablesToBeSaved = getTrainingVariables(category) + [targetVariable, 'ancestorHasWhichFlavor',
681  'isSignal', 'mcPDG', 'mcErrors', 'genMotherPDG',
682  'nMCMatches', 'B0mcErrors']
683  if category != 'KaonPion' and category != 'FSC':
684  variablesToBeSaved = variablesToBeSaved + \
685  ['extraInfo(isRightTrack(' + category + '))',
686  'hasHighestProbInCat(' + particleList + ', isRightTrack(' + category + '))']
687  ntuple.param('variables', variablesToBeSaved)
688  ntuple.param('particleList', particleList)
689  eventLevelpath.add_module(ntuple)
690 
691  if ReadyMethods != len(categories):
692  return False
693  else:
694  return True
695 
696 
697 def eventLevelTeacher(weightFiles='B2JpsiKs_mu', categories=None):
698  """
699  Trains all categories at event level.
700  """
701 
702  B2INFO('EVENT LEVEL TEACHER')
703 
704  ReadyMethods = 0
705 
706  if categories is None:
707  categories = []
708 
709  for category in categories:
710  methodPrefixEventLevel = "FlavorTagger_" + getBelleOrBelle2() + "_" + weightFiles + 'EventLevel' + category + 'FBDT'
711  targetVariable = 'isRightCategory(' + category + ')'
712  weightFile = filesDirectory + '/' + methodPrefixEventLevel + "_1.root"
713 
714  if os.path.isfile(weightFile):
715  ReadyMethods += 1
716  continue
717 
718  sampledFilesList = glob.glob(filesDirectory + '/' + methodPrefixEventLevel + 'sampled*.root')
719  if len(sampledFilesList) == 0:
720  B2INFO('flavorTagger: eventLevelTeacher did not find any ' + methodPrefixEventLevel +
721  ".root" + ' file. Please run the flavorTagger in "Sampler" mode afterwards.')
722 
723  else:
724  B2INFO('flavorTagger: MVA Teacher training' + methodPrefixEventLevel + ' .')
725  trainingOptionsEventLevel = basf2_mva.GeneralOptions()
726  trainingOptionsEventLevel.m_datafiles = basf2_mva.vector(*sampledFilesList)
727  trainingOptionsEventLevel.m_treename = methodPrefixEventLevel + "_tree"
728  trainingOptionsEventLevel.m_identifier = weightFile
729  trainingOptionsEventLevel.m_variables = basf2_mva.vector(*getTrainingVariables(category))
730  trainingOptionsEventLevel.m_target_variable = targetVariable
731  trainingOptionsEventLevel.m_max_events = maxEventsNumber
732 
733  basf2_mva.teacher(trainingOptionsEventLevel, getFastBDTCategories())
734 
735  if uploadFlag:
736  basf2_mva.upload(weightFile, methodPrefixEventLevel)
737 
738  if ReadyMethods != len(categories):
739  return False
740  else:
741  return True
742 
743 
744 def combinerLevel(mode='Expert', weightFiles='B2JpsiKs_mu', categories=None,
745  variablesCombinerLevel=None, categoriesCombinationCode=None, path=None):
746  """
747  Samples the input data or tests the combiner according to the selected categories.
748  """
749 
750  B2INFO('COMBINER LEVEL')
751 
752  if categories is None:
753  categories = []
754  if variablesCombinerLevel is None:
755  variablesCombinerLevel = []
756 
757  B2INFO("Flavor Tagger: Required Combiner for Categories:")
758  for category in categories:
759  B2INFO(category)
760 
761  B2INFO("Flavor Tagger: which corresponds to a weight file with categories combination code " + categoriesCombinationCode)
762 
763  methodPrefixCombinerLevel = "FlavorTagger_" + getBelleOrBelle2() + "_" + weightFiles + 'Combiner' \
764  + categoriesCombinationCode
765 
766  if mode == 'Sampler':
767 
768  if os.path.isfile(filesDirectory + '/' + methodPrefixCombinerLevel + 'FBDT' + '_1.root') or \
769  os.path.isfile(filesDirectory + '/' + methodPrefixCombinerLevel + 'FANN' + '_1.root'):
770  B2FATAL('flavorTagger: File' + methodPrefixCombinerLevel + 'FBDT' + "_1.root" + ' or ' + methodPrefixCombinerLevel +
771  'FANN' + '_1.root found. Please run the "Expert" mode or delete the file if a new sampling is desired.')
772 
773  B2INFO('flavorTagger: Sampling Data on Combiner Level. File' +
774  methodPrefixCombinerLevel + ".root" + ' will be saved')
775 
776  ntuple = basf2.register_module('VariablesToNtuple')
777  ntuple.param('fileName', filesDirectory + '/' + methodPrefixCombinerLevel + "sampled" + fileId + ".root")
778  ntuple.param('treeName', methodPrefixCombinerLevel + 'FBDT' + "_tree")
779  ntuple.param('variables', variablesCombinerLevel + ['qrCombined'])
780  ntuple.param('particleList', "")
781  path.add_module(ntuple)
782 
783  if mode == 'Expert':
784 
785  # Check if weight files are ready
786  if TMVAfbdt:
787  identifierFBDT = methodPrefixCombinerLevel + 'FBDT'
788  if downloadFlag or useOnlyLocalFlag:
789  identifierFBDT = filesDirectory + '/' + methodPrefixCombinerLevel + 'FBDT' + '_1.root'
790 
791  if downloadFlag:
792  if not os.path.isfile(identifierFBDT):
793  basf2_mva.download(methodPrefixCombinerLevel + 'FBDT', identifierFBDT)
794  if not os.path.isfile(identifierFBDT):
795  B2FATAL('Flavor Tagger: Weight file ' + identifierFBDT +
796  ' was not downloaded from Database. Please check the buildOrRevision name. Stopped')
797 
798  if useOnlyLocalFlag:
799  if not os.path.isfile(identifierFBDT):
800  B2FATAL('flavorTagger: Combinerlevel FastBDT was not trained with this combination of categories.' +
801  ' Weight file ' + identifierFBDT + ' not found. Stopped')
802 
803  B2INFO('flavorTagger: Ready to be used with weightFile ' + methodPrefixCombinerLevel + 'FBDT' + '_1.root')
804 
805  if FANNmlp:
806  identifierFANN = methodPrefixCombinerLevel + 'FANN'
807  if downloadFlag or useOnlyLocalFlag:
808  identifierFANN = filesDirectory + '/' + methodPrefixCombinerLevel + 'FANN' + '_1.root'
809 
810  if downloadFlag:
811  if not os.path.isfile(identifierFANN):
812  basf2_mva.download(methodPrefixCombinerLevel + 'FANN', identifierFANN)
813  if not os.path.isfile(identifierFANN):
814  B2FATAL('Flavor Tagger: Weight file ' + identifierFANN +
815  ' was not downloaded from Database. Please check the buildOrRevision name. Stopped')
816  if useOnlyLocalFlag:
817  if not os.path.isfile(identifierFANN):
818  B2FATAL('flavorTagger: Combinerlevel FANNMLP was not trained with this combination of categories. ' +
819  ' Weight file ' + identifierFANN + ' not found. Stopped')
820 
821  B2INFO('flavorTagger: Ready to be used with weightFile ' + methodPrefixCombinerLevel + 'FANN' + '_1.root')
822 
823  # At this stage, all necessary weight files should be ready.
824  # Call MVAExpert or MVAMultipleExperts module.
825  if TMVAfbdt and not FANNmlp:
826  B2INFO('flavorTagger: Apply FBDTMethod ' + methodPrefixCombinerLevel + 'FBDT')
827  path.add_module('MVAExpert', listNames=[], extraInfoName='qrCombined' + 'FBDT', signalFraction=signalFraction,
828  identifier=identifierFBDT)
829 
830  if FANNmlp and not TMVAfbdt:
831  B2INFO('flavorTagger: Apply FANNMethod on combiner level')
832  path.add_module('MVAExpert', listNames=[], extraInfoName='qrCombined' + 'FANN', signalFraction=signalFraction,
833  identifier=identifierFANN)
834 
835  if FANNmlp and TMVAfbdt:
836  B2INFO('flavorTagger: Apply FANNMethod and FBDTMethod on combiner level')
837  mvaMultipleExperts = basf2.register_module('MVAMultipleExperts')
838  mvaMultipleExperts.set_name('MVAMultipleExperts_Combiners')
839  mvaMultipleExperts.param('listNames', [])
840  mvaMultipleExperts.param('extraInfoNames', ['qrCombined' + 'FBDT', 'qrCombined' + 'FANN'])
841  mvaMultipleExperts.param('signalFraction', signalFraction)
842  mvaMultipleExperts.param('identifiers', [identifierFBDT, identifierFANN])
843  path.add_module(mvaMultipleExperts)
844 
845 
846 def combinerLevelTeacher(weightFiles='B2JpsiKs_mu', variablesCombinerLevel=None,
847  categoriesCombinationCode=None):
848  """
849  Trains the combiner according to the selected categories.
850  """
851 
852  B2INFO('COMBINER LEVEL TEACHER')
853 
854  if variablesCombinerLevel is None:
855  variablesCombinerLevel = []
856 
857  methodPrefixCombinerLevel = "FlavorTagger_" + getBelleOrBelle2() + "_" + weightFiles + 'Combiner' \
858  + categoriesCombinationCode
859 
860  sampledFilesList = glob.glob(filesDirectory + '/' + methodPrefixCombinerLevel + 'sampled*.root')
861  if len(sampledFilesList) == 0:
862  B2FATAL('FlavorTagger: combinerLevelTeacher did not find any ' +
863  methodPrefixCombinerLevel + 'sampled*.root file. Please run the flavorTagger in "Sampler" mode.')
864 
865  if TMVAfbdt:
866 
867  if not os.path.isfile(filesDirectory + '/' + methodPrefixCombinerLevel + 'FBDT' + '_1.root'):
868 
869  B2INFO('flavorTagger: MVA Teacher training a FastBDT on Combiner Level')
870 
871  trainingOptionsCombinerLevel = basf2_mva.GeneralOptions()
872  trainingOptionsCombinerLevel.m_datafiles = basf2_mva.vector(*sampledFilesList)
873  trainingOptionsCombinerLevel.m_treename = methodPrefixCombinerLevel + 'FBDT' + "_tree"
874  trainingOptionsCombinerLevel.m_identifier = filesDirectory + '/' + methodPrefixCombinerLevel + 'FBDT' + "_1.root"
875  trainingOptionsCombinerLevel.m_variables = basf2_mva.vector(*variablesCombinerLevel)
876  trainingOptionsCombinerLevel.m_target_variable = 'qrCombined'
877  trainingOptionsCombinerLevel.m_max_events = maxEventsNumber
878 
879  basf2_mva.teacher(trainingOptionsCombinerLevel, getFastBDTCombiner())
880 
881  if uploadFlag:
882  basf2_mva.upload(filesDirectory + '/' + methodPrefixCombinerLevel +
883  'FBDT' + "_1.root", methodPrefixCombinerLevel + 'FBDT')
884 
885  elif FANNmlp and not os.path.isfile(filesDirectory + '/' + methodPrefixCombinerLevel + 'FANN' + '_1.root'):
886 
887  B2INFO('flavorTagger: Combinerlevel FBDT was already trained with this combination of categories. Weight file ' +
888  methodPrefixCombinerLevel + 'FBDT' + '_1.root has been found.')
889 
890  else:
891  B2FATAL('flavorTagger: Combinerlevel was already trained with this combination of categories. Weight files ' +
892  methodPrefixCombinerLevel + 'FBDT' + '_1.root and ' +
893  methodPrefixCombinerLevel + 'FANN' + '_1.root has been found. Please use the "Expert" mode')
894 
895  if FANNmlp:
896 
897  if not os.path.isfile(filesDirectory + '/' + methodPrefixCombinerLevel + 'FANN' + '_1.root'):
898 
899  B2INFO('flavorTagger: MVA Teacher training a FANN MLP on Combiner Level')
900 
901  trainingOptionsCombinerLevel = basf2_mva.GeneralOptions()
902  trainingOptionsCombinerLevel.m_datafiles = basf2_mva.vector(*sampledFilesList)
903  trainingOptionsCombinerLevel.m_treename = methodPrefixCombinerLevel + 'FBDT' + "_tree"
904  trainingOptionsCombinerLevel.m_identifier = filesDirectory + '/' + methodPrefixCombinerLevel + 'FANN' + "_1.root"
905  trainingOptionsCombinerLevel.m_variables = basf2_mva.vector(*variablesCombinerLevel)
906  trainingOptionsCombinerLevel.m_target_variable = 'qrCombined'
907  trainingOptionsCombinerLevel.m_max_events = maxEventsNumber
908 
909  basf2_mva.teacher(trainingOptionsCombinerLevel, getMlpFANNCombiner())
910 
911  if uploadFlag:
912  basf2_mva.upload(filesDirectory + '/' + methodPrefixCombinerLevel +
913  'FANN' + "_1.root", methodPrefixCombinerLevel + 'FANN')
914 
915  elif TMVAfbdt and not os.path.isfile(filesDirectory + '/' + methodPrefixCombinerLevel + 'FBDT' + '_1.root'):
916 
917  B2INFO('flavorTagger: Combinerlevel FBDT was already trained with this combination of categories. Weight file ' +
918  methodPrefixCombinerLevel + 'FANN' + '_1.config has been found.')
919 
920  else:
921  B2FATAL('flavorTagger: Combinerlevel was already trained with this combination of categories. Weight files ' +
922  methodPrefixCombinerLevel + 'FBDT' + '_1.root and ' +
923  methodPrefixCombinerLevel + 'FANN' + '_1.root has been found. Please use the "Expert" mode')
924 
925 
926 def getEventLevelParticleLists(categories=None):
927 
928  if categories is None:
929  categories = []
930 
931  eventLevelParticleLists = []
932 
933  for category in categories:
934  ftCategory = AvailableCategories[category]
935  event_tuple = (ftCategory.particleList, ftCategory.eventName, ftCategory.variableName)
936 
937  if event_tuple not in eventLevelParticleLists:
938  eventLevelParticleLists.append(event_tuple)
939  else:
940  B2FATAL('Flavor Tagger: ' + category + ' has been already given')
941 
942  return eventLevelParticleLists
943 
944 
945 def flavorTagger(
946  particleLists=None,
947  mode='Expert',
948  weightFiles='B2nunubarBGx1',
949  workingDirectory='.',
950  combinerMethods=['TMVA-FBDT', 'FANN-MLP'],
951  categories=[
952  'Electron',
953  'IntermediateElectron',
954  'Muon',
955  'IntermediateMuon',
956  'KinLepton',
957  'IntermediateKinLepton',
958  'Kaon',
959  'SlowPion',
960  'FastHadron',
961  'Lambda',
962  'FSC',
963  'MaximumPstar',
964  'KaonPion'],
965  maskName='all',
966  saveCategoriesInfo=True,
967  useOnlyLocalWeightFiles=False,
968  downloadFromDatabaseIfNotFound=False,
969  uploadToDatabaseAfterTraining=False,
970  samplerFileId='',
971  prefix='',
972  path=None,
973 ):
974  """
975  Defines the whole flavor tagging process for each selected Rest of Event (ROE) built in the steering file.
976  The flavor is predicted by Multivariate Methods trained with Variables and MetaVariables which use
977  Tracks, ECL- and KLMClusters from the corresponding RestOfEvent dataobject.
978  This module can be used to sample the training information, to train and/or to test the flavorTagger.
979 
980  @param particleLists The ROEs for flavor tagging are selected from the given particle lists.
981  @param mode The available modes are
982  ``Expert`` (default), ``Sampler``, and ``Teacher``. In the ``Expert`` mode
983  Flavor Tagging is applied to the analysis,. In the ``Sampler`` mode you save
984  save the variables for training. In the ``Teacher`` mode the FlavorTagger is
985  trained, for this step you do not reconstruct any particle or do any analysis,
986  you just run the flavorTagger alone.
987  @param weightFiles Weight files name. Default=
988  ``B2nunubarBGx1`` (official weight files). If the user self
989  wants to train the FlavorTagger, the weightfiles name should correspond to the
990  analysed CP channel in order to avoid confusions. The default name
991  ``B2nunubarBGx1`` corresponds to
992  :math:`B^0_{\\rm sig}\\to \\nu \\overline{\\nu}`.
993  and ``B2JpsiKs_muBGx1`` to
994  :math:`B^0_{\\rm sig}\\to J/\\psi (\\to \\mu^+ \\mu^-) K_s (\\to \\pi^+ \\pi^-)`.
995  BGx1 stays for events simulated with background.
996  @param workingDirectory Path to the directory containing the FlavorTagging/ folder.
997  @param combinerMethods MVAs for the combiner: ``TMVA-FBDT`` or ``FANN-MLP``. Both used by default.
998  @param categories Categories used for flavor tagging. By default all are used.
999  @param maskName Gets ROE particles from a specified ROE mask.
1000  ``all`` (default): all ROE particles are used.
1001  ``_FTDefaultMask``: tentative mask definition that will be created automatically.
1002  Or one can give any mask name defined before calling this function.
1003  @param saveCategoriesInfo Sets to save information of individual categories.
1004  @param useOnlyLocalWeightFiles [Expert] Uses only locally saved weight files.
1005  @param downloadFromDatabaseIfNotFound [Expert] Weight files are downloaded from
1006  the conditions database if not available in workingDirectory.
1007  @param uploadToDatabaseAfterTraining [Expert] For librarians only: uploads weight files to localdb after training.
1008  @param samplerFileId Identifier to paralellize
1009  sampling. Only used in ``Sampler`` mode. If you are training by yourself and
1010  want to parallelize the sampling, you can run several sampling scripts in
1011  parallel. By changing this parameter you will not overwrite an older sample.
1012  @param prefix Prefix of weight files.
1013  @param path Modules are added to this path
1014 
1015  """
1016 
1017  if mode != 'Sampler' and mode != 'Teacher' and mode != 'Expert':
1018  B2FATAL('flavorTagger: Wrong mode given: The available modes are "Sampler", "Teacher" or "Expert"')
1019 
1020  if len(categories) != len(set(categories)):
1021  dup = [cat for cat in set(categories) if categories.count(cat) > 1]
1022  B2WARNING('Flavor Tagger: There are duplicate elements in the given categories list. '
1023  << 'The following duplicate elements are removed; ' << ', '.join(dup))
1024  categories = list(set(categories))
1025 
1026  if len(categories) < 2:
1027  B2FATAL('Flavor Tagger: Invalid amount of categories. At least two are needed.')
1028  B2FATAL(
1029  'Flavor Tagger: Possible categories are "Electron", "IntermediateElectron", "Muon", "IntermediateMuon", '
1030  '"KinLepton", "IntermediateKinLepton", "Kaon", "SlowPion", "FastHadron",'
1031  '"Lambda", "FSC", "MaximumPstar" or "KaonPion" ')
1032 
1033  for category in categories:
1034  if category not in AvailableCategories:
1035  B2FATAL('Flavor Tagger: ' + category + ' is not a valid category name given')
1036  B2FATAL('Flavor Tagger: Available categories are "Electron", "IntermediateElectron", '
1037  '"Muon", "IntermediateMuon", "KinLepton", "IntermediateKinLepton", "Kaon", "SlowPion", "FastHadron", '
1038  '"Lambda", "FSC", "MaximumPstar" or "KaonPion" ')
1039 
1040  # Directory where the weights of the trained Methods are saved
1041  # workingDirectory = os.environ['BELLE2_LOCAL_DIR'] + '/analysis/data'
1042 
1043  basf2.find_file(workingDirectory)
1044 
1045  global filesDirectory
1046  filesDirectory = workingDirectory + '/FlavorTagging/TrainedMethods'
1047 
1048  if mode == 'Sampler' or (mode == 'Expert' and downloadFromDatabaseIfNotFound):
1049  if not basf2.find_file(workingDirectory + '/FlavorTagging', silent=True):
1050  os.mkdir(workingDirectory + '/FlavorTagging')
1051  os.mkdir(workingDirectory + '/FlavorTagging/TrainedMethods')
1052  elif not basf2.find_file(workingDirectory + '/FlavorTagging/TrainedMethods', silent=True):
1053  os.mkdir(workingDirectory + '/FlavorTagging/TrainedMethods')
1054  filesDirectory = workingDirectory + '/FlavorTagging/TrainedMethods'
1055 
1056  if len(combinerMethods) < 1 or len(combinerMethods) > 2:
1057  B2FATAL('flavorTagger: Invalid list of combinerMethods. The available methods are "TMVA-FBDT" and "FANN-MLP"')
1058 
1059  global FANNmlp
1060  global TMVAfbdt
1061 
1062  FANNmlp = False
1063  TMVAfbdt = False
1064 
1065  for method in combinerMethods:
1066  if method == 'TMVA-FBDT':
1067  TMVAfbdt = True
1068  elif method == 'FANN-MLP':
1069  FANNmlp = True
1070  else:
1071  B2FATAL('flavorTagger: Invalid list of combinerMethods. The available methods are "TMVA-FBDT" and "FANN-MLP"')
1072 
1073  global fileId
1074  fileId = samplerFileId
1075 
1076  global useOnlyLocalFlag
1077  useOnlyLocalFlag = useOnlyLocalWeightFiles
1078 
1079  B2INFO('*** FLAVOR TAGGING ***')
1080  B2INFO(' ')
1081  B2INFO(' Working directory is: ' + filesDirectory)
1082  B2INFO(' ')
1083 
1084  setInteractionWithDatabase(downloadFromDatabaseIfNotFound, uploadToDatabaseAfterTraining)
1085  set_FlavorTagger_pid_aliases()
1086  setInputVariablesWithMask()
1087  weightFiles = prefix + weightFiles
1088 
1089  # Create configuration lists and code-name for given category's list
1090  trackLevelParticleLists = []
1091  eventLevelParticleLists = []
1092  variablesCombinerLevel = []
1093  categoriesCombination = []
1094  categoriesCombinationCode = 'CatCode'
1095  for category in categories:
1096  ftCategory = AvailableCategories[category]
1097 
1098  track_tuple = (ftCategory.particleList, ftCategory.trackName)
1099  event_tuple = (ftCategory.particleList, ftCategory.eventName, ftCategory.variableName)
1100 
1101  if track_tuple not in trackLevelParticleLists and category != 'MaximumPstar':
1102  trackLevelParticleLists.append(track_tuple)
1103 
1104  if event_tuple not in eventLevelParticleLists:
1105  eventLevelParticleLists.append(event_tuple)
1106  variablesCombinerLevel.append(ftCategory.variableName)
1107  categoriesCombination.append(ftCategory.code)
1108  else:
1109  B2FATAL('Flavor Tagger: ' + category + ' has been already given')
1110 
1111  for code in sorted(categoriesCombination):
1112  categoriesCombinationCode = categoriesCombinationCode + '%02d' % code
1113 
1114  # Create default ROE-mask
1115  if maskName == '_FTDefaultMask':
1116  _FTDefaultMask = (
1117  '_FTDefaultMask',
1118  'thetaInCDCAcceptance and dr<1 and abs(dz)<3',
1119  'thetaInCDCAcceptance and clusterNHits>1.5 and [[E>0.08 and clusterReg==1] or [E>0.03 and clusterReg==2] or \
1120  [E>0.06 and clusterReg==3]]')
1121  for name in particleLists:
1122  ma.appendROEMasks(list_name=name, mask_tuples=[_FTDefaultMask], path=path)
1123 
1124  # Start ROE-routine
1125  roe_path = basf2.create_path()
1126  deadEndPath = basf2.create_path()
1127 
1128  if mode == 'Sampler':
1129  # Events containing ROE without B-Meson (but not empty) are discarded for training
1130  ma.signalSideParticleListsFilter(
1131  particleLists,
1132  'nROE_Charged(' + maskName + ', 0) > 0 and abs(qrCombined) == 1',
1133  roe_path,
1134  deadEndPath)
1135 
1136  FillParticleLists(maskName, categories, roe_path)
1137 
1138  if eventLevel(mode, weightFiles, categories, roe_path):
1139  combinerLevel(mode, weightFiles, categories, variablesCombinerLevel, categoriesCombinationCode, roe_path)
1140 
1141  path.for_each('RestOfEvent', 'RestOfEvents', roe_path)
1142 
1143  elif mode == 'Expert':
1144  # If trigger returns 1 jump into empty path skipping further modules in roe_path
1145  # run filter with no cut first to get rid of ROEs that are missing the mask of the signal particle
1146  ma.signalSideParticleListsFilter(particleLists, '', roe_path, deadEndPath)
1147  ma.signalSideParticleListsFilter(particleLists, 'nROE_Charged(' + maskName + ', 0) > 0', roe_path, deadEndPath)
1148 
1149  # Initialization of flavorTaggerInfo dataObject needs to be done in the main path
1150  flavorTaggerInfoBuilder = basf2.register_module('FlavorTaggerInfoBuilder')
1151  path.add_module(flavorTaggerInfoBuilder)
1152 
1153  FillParticleLists(maskName, categories, roe_path)
1154 
1155  if eventLevel(mode, weightFiles, categories, roe_path):
1156  combinerLevel(mode, weightFiles, categories, variablesCombinerLevel, categoriesCombinationCode, roe_path)
1157 
1158  flavorTaggerInfoFiller = basf2.register_module('FlavorTaggerInfoFiller')
1159  flavorTaggerInfoFiller.param('trackLevelParticleLists', trackLevelParticleLists)
1160  flavorTaggerInfoFiller.param('eventLevelParticleLists', eventLevelParticleLists)
1161  flavorTaggerInfoFiller.param('TMVAfbdt', TMVAfbdt)
1162  flavorTaggerInfoFiller.param('FANNmlp', FANNmlp)
1163  flavorTaggerInfoFiller.param('qpCategories', saveCategoriesInfo)
1164  flavorTaggerInfoFiller.param('istrueCategories', saveCategoriesInfo)
1165  flavorTaggerInfoFiller.param('targetProb', False)
1166  flavorTaggerInfoFiller.param('trackPointers', False)
1167  roe_path.add_module(flavorTaggerInfoFiller) # Add FlavorTag Info filler to roe_path
1168  add_default_FlavorTagger_aliases()
1169 
1170  path.for_each('RestOfEvent', 'RestOfEvents', roe_path)
1171 
1172  elif mode == 'Teacher':
1173  if eventLevelTeacher(weightFiles, categories):
1174  combinerLevelTeacher(weightFiles, variablesCombinerLevel, categoriesCombinationCode)
1175 
1176 
1177 if __name__ == '__main__':
1178 
1179  desc_list = []
1180 
1181  function = globals()["flavorTagger"]
1182  signature = inspect.formatargspec(*inspect.getfullargspec(function))
1183  desc_list.append((function.__name__, signature + '\n' + function.__doc__))
1184 
1185  from terminal_utils import Pager
1186  from basf2.utils import pretty_print_description_list
1187  with Pager('Flavor Tagger function accepts the following arguments:'):
1188  pretty_print_description_list(desc_list)
def isB2BII()
Definition: b2bii.py:14