Belle II Software development
CreateLogics Class Reference
Inheritance diagram for CreateLogics:

Public Member Functions

def event (self)
 
def terminate (self)
 

Static Public Attributes

ROOT file = ROOT.TFile(argvs[2], 'recreate')
 the output file
 
ROOT tgrl = ROOT.TTree('tgrl', 'tree with GRL_Logic')
 the tree in the output file
 
array ntrk_2dfinder_t = array('i', [-1])
 #2d finder tracks
 
array ntrk_2dfitter_t = array('i', [-1])
 #2d fitter tracks
 
array ntrk_3dfitter_t = array('i', [-1])
 #3d fitter tracks
 
array ntrk_NN_t = array('i', [-1])
 number of NN tracks
 
array ntrk_2Dmatch_t = array('i', [-1])
 #2d matched tracks
 
array ntrk_3Dmatch_t = array('i', [-1])
 #3d matched tracks
 
array max_deltphi_2dfinder_t = array('f', [0.0])
 max phi angle between two 2d finder tracks
 
array cpair_t = array('i', 8 * [0])
 number of cluster pairs with different energy threshold
 
int ncomp_clu = 3
 number of array components
 
array etot_t = array('f', [0.0])
 the total deposited cluster energy in ecl
 
array ncluster_1000b_t = array('i', [-1])
 number of ecl cluster with threshold >1.0GeV, exclude TC ID 1,2, 17
 
array ncluster_2000e_t = array('i', [-1])
 
array ncluster_t = array('i', [-1])
 
array ncluster_neutral_t = array('i', [-1])
 
array max_cluster_neutral_t = array('f', ncomp_clu * [0.0])
 [energy, theta, phi] of the largest energetic ecl neutral cluster
 
array max_cms_cluster_neutral_t = array('f', ncomp_clu * [0.0])
 [energy, theta, phi] of the largest energetic ecl neutral cluster in CMS
 
array max_cluster_t = array('f', ncomp_clu * [0.0])
 [energy, theta, phi] of the largest energetic ecl cluster
 
array max_cms_cluster_t = array('f', ncomp_clu * [0.0])
 [energy, theta, phi] of the largest energetic ecl cluster in CMS
 
array smax_cluster_neutral_t = array('f', ncomp_clu * [0.0])
 [energy, theta, phi] of the secondary energetic ecl neutral cluster
 
array smax_cms_cluster_neutral_t = array('f', ncomp_clu * [0.0])
 [energy, theta, phi] of the secondary energetic ecl neutral cluster in CMS
 
array smax_cluster_t = array('f', ncomp_clu * [0.0])
 [energy, theta, phi] of the secondary energetic ecl cluster
 
array smax_cms_cluster_t = array('f', ncomp_clu * [0.0])
 [energy, theta, phi] of the secondary energetic ecl cluster in CMS
 
array max_deltphi_cluster_t = array('f', [0.0])
 max delt phi angle between two clusters
 
array time_cluster_t = array('f', 100 * [-99999.])
 the timing of clusters
 
array ntrk_klm_t = array('i', [-1])
 number of KLM tracks
 
array nhit_klm_t = array('i', [-1])
 number of KLM hits
 
array npair_tc_t = array('i', [-1])
 
array npair_cc_t = array('i', [-1])
 
int nbha_var = 5
 number of array components
 
array bhabha_t = array('i', [0])
 bhabha veto logic, 1: bhabha, 0: non bhabha
 
array sbhabha_t = array('i', [0])
 bhabha with single track veto logic, 1: bhabha, 0: non bhabha
 
array eclbhabha_t = array('i', [0])
 eclbhabha veto logic, 1: bhabha, 0: non bhabha
 
array bhabha_var_t = array('f', nbha_var * [0.0])
 variables used in bhabha logic, [Delt_theta, Delt_phi, E1, E2, E1+E2]
 
array eclbhabha_var_t = array('f', nbha_var * [0.0])
 variables used in eclbhabha logic, [Delt_theta, Delt_phi, E1, E2, E1+E2]
 

Detailed Description

This module is to calculate some variables which are useful for the trigger development

Definition at line 361 of file ReadTSIMInfo.py.

Member Function Documentation

◆ event()

def event (   self)
MCPart = Belle2.PyStoreArray('MCParticles')
if len(MCPart)>=4:
    self.eplus_t[0]=Theta_MCTruth(MCPart[2])
    self.eminus_t[0]=Theta_MCTruth(MCPart[3])
    if self.eplus_t[0] > self.eminus_t[0]:
         self.efrd_t[0]=self.eminus_t[0]
         self.ebkd_t[0]=self.eplus_t[0]
    else:
         self.efrd_t[0]=self.eplus_t[0]
         self.ebkd_t[0]=self.eminus_t[0]
count_part = Count_part_inDetector(MCPart)
self.n_par_t[0]=count_part[0]
self.n_par_t[1]=count_part[1]
count_mcpart = Count_mcpart(MCPart)
self.n_mcpar_t[0]=count_mcpart[0]
self.n_mcpar_t[1]=count_mcpart[1]

Definition at line 477 of file ReadTSIMInfo.py.

477 def event(self):
478 """
479 MCPart = Belle2.PyStoreArray('MCParticles')
480 if len(MCPart)>=4:
481 self.eplus_t[0]=Theta_MCTruth(MCPart[2])
482 self.eminus_t[0]=Theta_MCTruth(MCPart[3])
483 if self.eplus_t[0] > self.eminus_t[0]:
484 self.efrd_t[0]=self.eminus_t[0]
485 self.ebkd_t[0]=self.eplus_t[0]
486 else:
487 self.efrd_t[0]=self.eplus_t[0]
488 self.ebkd_t[0]=self.eminus_t[0]
489 count_part = Count_part_inDetector(MCPart)
490 self.n_par_t[0]=count_part[0]
491 self.n_par_t[1]=count_part[1]
492 count_mcpart = Count_mcpart(MCPart)
493 self.n_mcpar_t[0]=count_mcpart[0]
494 self.n_mcpar_t[1]=count_mcpart[1]
495 """
496 trk_2d_finder = Belle2.PyStoreArray('TRG2DFinderTracks')
497 self.ntrk_2dfinder_t[0] = len(trk_2d_finder)
498 trk_2d_fitter = Belle2.PyStoreArray('TRG2DFitterTracks')
499 self.ntrk_2dfitter_t[0] = len(trk_2d_fitter)
500 trk_3d_fitter = Belle2.PyStoreArray('TRG3DFitterTracks')
501 self.ntrk_3dfitter_t[0] = len(trk_3d_fitter)
502 tracks = Belle2.PyStoreArray('TRGNNTracks')
503 self.ntrk_NN_t[0] = len(tracks)
504
505 # clusters
506 clusters_original = Belle2.PyStoreArray('TRGECLClusters')
507 eventtime = Belle2.PyStoreArray('TRGECLTrgs')
508 # get clusters list in time window [-100,100]
509 clusters = Time_Window(clusters_original, eventtime)
510 self.ncluster_t[0] = len(clusters)
511 # get the cluster timing before time window requirement
512 time_clu_list = Time_Cluster(clusters_original, eventtime)
513 for i in range(len(time_clu_list)):
514 if i < 100:
515 self.time_cluster_t[i] = time_clu_list[i]
516
517 clusters_300_cms = Cluster_Threshold(clusters, 0.3, False)
518 clusters_400_cms = Cluster_Threshold(clusters, 0.4, False)
519 clusters_500_cms = Cluster_Threshold(clusters, 0.5, False)
520 clusters_700_cms = Cluster_Threshold(clusters, 0.7, False)
521 clusters_1000_cms = Cluster_Threshold(clusters, 1.0, False)
522 clusters_2000_cms = Cluster_Threshold(clusters, 2.0, False)
523 clusters_2500_cms = Cluster_Threshold(clusters, 2.5, False)
524 self.cpair_t[0] = Back_to_Back(clusters, clusters)
525 self.cpair_t[1] = Back_to_Back(clusters_300_cms, clusters)
526 self.cpair_t[2] = Back_to_Back(clusters_400_cms, clusters)
527 self.cpair_t[3] = Back_to_Back(clusters_500_cms, clusters)
528 self.cpair_t[4] = Back_to_Back(clusters_700_cms, clusters)
529 self.cpair_t[5] = Back_to_Back(clusters_1000_cms, clusters)
530 self.cpair_t[6] = Back_to_Back(clusters_2000_cms, clusters)
531 self.cpair_t[7] = Back_to_Back(clusters_2500_cms, clusters)
532
533 klmtrkcol = Belle2.PyStoreArray('TRGKLMTracks')
534 self.ntrk_klm_t[0] = len(klmtrkcol)
535 klmhitcol = Belle2.PyStoreArray('TRGKLMHits')
536 self.nhit_klm_t[0] = len(klmhitcol)
537
538 matchlist = Belle2.PyStoreArray('TRG3DMatchTracks')
539 self.ntrk_3Dmatch_t[0] = len(matchlist)
540
541 trginfo = Belle2.PyStoreObj('TRGGRLObjects')
542 self.npair_tc_t[0] = trginfo.getNbbTrkCluster()
543 self.npair_cc_t[0] = trginfo.getNbbCluster()
544 self.ncluster_1000b_t[0] = trginfo.getNhighcluster2()
545 self.ncluster_2000e_t[0] = trginfo.getNhighcluster4()
546 self.max_deltphi_2dfinder_t[0] = Max_DeltPhi_trk(trk_2d_finder)
547
548 neutral_clusters = NeutralCluster(clusters)
549 self.ncluster_neutral_t[0] = len(neutral_clusters)
550 max_cluster_neu = Max_Cluster(neutral_clusters, 0)
551 max_cluster = Max_Cluster(clusters, 0)
552 # max_cms_cluster_neu = Max_Cluster(neutral_clusters, 0)
553 # max_cms_cluster = Max_Cluster(clusters, 0)
554 smax_cluster_neu = Max_Cluster(neutral_clusters, 1)
555 smax_cluster = Max_Cluster(clusters, 1)
556 # smax_cms_cluster_neu = Max_Cluster(neutral_clusters, 1)
557 # smax_cms_cluster = Max_Cluster(clusters, 1)
558 for i in range(self.ncomp_clu):
559 self.max_cluster_neutral_t[i] = max_cluster_neu[i]
560 # self.max_cms_cluster_neutral_t[i] = max_cms_cluster_neu[i]
561 self.max_cluster_t[i] = max_cluster[i]
562 # self.max_cms_cluster_t[i] = max_cms_cluster[i]
563 self.smax_cluster_neutral_t[i] = smax_cluster_neu[i]
564 # self.smax_cms_cluster_neutral_t[i] = smax_cms_cluster_neu[i]
565 self.smax_cluster_t[i] = smax_cluster[i]
566 # self.smax_cms_cluster_t[i] = smax_cms_cluster[i]
567
568 self.etot_t[0] = Etot_Cluster(clusters)
569
570 # bhabha
571 bhabhaveto_1 = BhabhaVeto1(matchlist)
572 bha_logic = 0
573 for bha in bhabhaveto_1:
574 if math.fabs(bha[0]) < 50 and math.fabs(bha[0]) > 10 and math.fabs(bha[1]) < 20:
575 if bha[2] > 2.0 and bha[3] > 2.0 and bha[4] > 6.0:
576 if bha[2] > 3.0 or bha[3] > 3.0:
577 if len(trk_2d_finder) == 2:
578 bha_logic = 1
579 self.bhabha_t[0] = bha_logic
580
581 if len(bhabhaveto_1) >= 1:
582 for i in range(self.nbha_var):
583 self.bhabha_var_t[i] = bhabhaveto_1[0][i]
584
585 # eclbhabha
586 eclbhabhaveto = eclBhabhaVeto(clusters)
587 eclbha_logic = 0
588 for eclbha in eclbhabhaveto:
589 if math.fabs(eclbha[0]) < 50 and math.fabs(eclbha[1]) < 50:
590 if eclbha[2] > 2.0 and eclbha[3] > 2.0 and eclbha[4] > 6.0:
591 if eclbha[2] > 3.0 or eclbha[3] > 3.0:
592 eclbha_logic = 1
593 self.eclbhabha_t[0] = eclbha_logic
594
595 # sbhahba
596 sbha_logic = 0
597 if len(trk_2d_finder) == 1:
598 if eclbha_logic == 1:
599 ecol = SBhabhaVeto(matchlist)
600 for i, etr in ecol:
601 if etr[i] > 1.0:
602 sbha_logic = 1
603 self.sbhabha_t[0] = sbha_logic
604
605 self.max_deltphi_cluster_t[0] = Max_DeltPhi_cluster(clusters)
606 if len(eclbhabhaveto) >= 1:
607 for i in range(self.nbha_var):
608 self.eclbhabha_var_t[i] = eclbhabhaveto[0][i]
609 self.tgrl.Fill()
610
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72
a (simplified) python wrapper for StoreObjPtr.
Definition: PyStoreObj.h:67

◆ terminate()

def terminate (   self)
Write and close the file

Definition at line 611 of file ReadTSIMInfo.py.

611 def terminate(self):
612 """Write and close the file"""
613 self.file.cd()
614 self.file.Write()
615 self.file.Close()
616 PrintBranchDef()
617
618

Member Data Documentation

◆ bhabha_t

array bhabha_t = array('i', [0])
static

bhabha veto logic, 1: bhabha, 0: non bhabha

Definition at line 426 of file ReadTSIMInfo.py.

◆ bhabha_var_t

array bhabha_var_t = array('f', nbha_var * [0.0])
static

variables used in bhabha logic, [Delt_theta, Delt_phi, E1, E2, E1+E2]

Definition at line 432 of file ReadTSIMInfo.py.

◆ cpair_t

array cpair_t = array('i', 8 * [0])
static

number of cluster pairs with different energy threshold

Definition at line 382 of file ReadTSIMInfo.py.

◆ eclbhabha_t

array eclbhabha_t = array('i', [0])
static

eclbhabha veto logic, 1: bhabha, 0: non bhabha

Definition at line 430 of file ReadTSIMInfo.py.

◆ eclbhabha_var_t

array eclbhabha_var_t = array('f', nbha_var * [0.0])
static

variables used in eclbhabha logic, [Delt_theta, Delt_phi, E1, E2, E1+E2]

Definition at line 434 of file ReadTSIMInfo.py.

◆ etot_t

array etot_t = array('f', [0.0])
static

the total deposited cluster energy in ecl

Definition at line 386 of file ReadTSIMInfo.py.

◆ file

ROOT file = ROOT.TFile(argvs[2], 'recreate')
static

the output file

Definition at line 364 of file ReadTSIMInfo.py.

◆ max_cluster_neutral_t

array max_cluster_neutral_t = array('f', ncomp_clu * [0.0])
static

[energy, theta, phi] of the largest energetic ecl neutral cluster

Definition at line 396 of file ReadTSIMInfo.py.

◆ max_cluster_t

array max_cluster_t = array('f', ncomp_clu * [0.0])
static

[energy, theta, phi] of the largest energetic ecl cluster

Definition at line 400 of file ReadTSIMInfo.py.

◆ max_cms_cluster_neutral_t

array max_cms_cluster_neutral_t = array('f', ncomp_clu * [0.0])
static

[energy, theta, phi] of the largest energetic ecl neutral cluster in CMS

Definition at line 398 of file ReadTSIMInfo.py.

◆ max_cms_cluster_t

array max_cms_cluster_t = array('f', ncomp_clu * [0.0])
static

[energy, theta, phi] of the largest energetic ecl cluster in CMS

Definition at line 402 of file ReadTSIMInfo.py.

◆ max_deltphi_2dfinder_t

array max_deltphi_2dfinder_t = array('f', [0.0])
static

max phi angle between two 2d finder tracks

Definition at line 380 of file ReadTSIMInfo.py.

◆ max_deltphi_cluster_t

array max_deltphi_cluster_t = array('f', [0.0])
static

max delt phi angle between two clusters

Definition at line 412 of file ReadTSIMInfo.py.

◆ nbha_var

int nbha_var = 5
static

number of array components

Definition at line 424 of file ReadTSIMInfo.py.

◆ ncluster_1000b_t

array ncluster_1000b_t = array('i', [-1])
static

number of ecl cluster with threshold >1.0GeV, exclude TC ID 1,2, 17

Definition at line 388 of file ReadTSIMInfo.py.

◆ ncluster_2000e_t

array ncluster_2000e_t = array('i', [-1])
static

ecl cluster with threshold >2.0GeV in TC ID 1, 17

Definition at line 390 of file ReadTSIMInfo.py.

◆ ncluster_neutral_t

array ncluster_neutral_t = array('i', [-1])
static

ecl cluster w/o associated cdc track

Definition at line 394 of file ReadTSIMInfo.py.

◆ ncluster_t

array ncluster_t = array('i', [-1])
static

ecl clusters

Definition at line 392 of file ReadTSIMInfo.py.

◆ ncomp_clu

int ncomp_clu = 3
static

number of array components

Definition at line 384 of file ReadTSIMInfo.py.

◆ nhit_klm_t

array nhit_klm_t = array('i', [-1])
static

number of KLM hits

Definition at line 418 of file ReadTSIMInfo.py.

◆ npair_cc_t

array npair_cc_t = array('i', [-1])
static

back to back cluster pairs

Definition at line 422 of file ReadTSIMInfo.py.

◆ npair_tc_t

array npair_tc_t = array('i', [-1])
static

back to back track and cluster pairs

Definition at line 420 of file ReadTSIMInfo.py.

◆ ntrk_2dfinder_t

array ntrk_2dfinder_t = array('i', [-1])
static

#2d finder tracks

Definition at line 368 of file ReadTSIMInfo.py.

◆ ntrk_2dfitter_t

array ntrk_2dfitter_t = array('i', [-1])
static

#2d fitter tracks

Definition at line 370 of file ReadTSIMInfo.py.

◆ ntrk_2Dmatch_t

array ntrk_2Dmatch_t = array('i', [-1])
static

#2d matched tracks

Definition at line 376 of file ReadTSIMInfo.py.

◆ ntrk_3dfitter_t

array ntrk_3dfitter_t = array('i', [-1])
static

#3d fitter tracks

Definition at line 372 of file ReadTSIMInfo.py.

◆ ntrk_3Dmatch_t

array ntrk_3Dmatch_t = array('i', [-1])
static

#3d matched tracks

Definition at line 378 of file ReadTSIMInfo.py.

◆ ntrk_klm_t

array ntrk_klm_t = array('i', [-1])
static

number of KLM tracks

Definition at line 416 of file ReadTSIMInfo.py.

◆ ntrk_NN_t

array ntrk_NN_t = array('i', [-1])
static

number of NN tracks

Definition at line 374 of file ReadTSIMInfo.py.

◆ sbhabha_t

array sbhabha_t = array('i', [0])
static

bhabha with single track veto logic, 1: bhabha, 0: non bhabha

Definition at line 428 of file ReadTSIMInfo.py.

◆ smax_cluster_neutral_t

array smax_cluster_neutral_t = array('f', ncomp_clu * [0.0])
static

[energy, theta, phi] of the secondary energetic ecl neutral cluster

Definition at line 404 of file ReadTSIMInfo.py.

◆ smax_cluster_t

array smax_cluster_t = array('f', ncomp_clu * [0.0])
static

[energy, theta, phi] of the secondary energetic ecl cluster

Definition at line 408 of file ReadTSIMInfo.py.

◆ smax_cms_cluster_neutral_t

array smax_cms_cluster_neutral_t = array('f', ncomp_clu * [0.0])
static

[energy, theta, phi] of the secondary energetic ecl neutral cluster in CMS

Definition at line 406 of file ReadTSIMInfo.py.

◆ smax_cms_cluster_t

array smax_cms_cluster_t = array('f', ncomp_clu * [0.0])
static

[energy, theta, phi] of the secondary energetic ecl cluster in CMS

Definition at line 410 of file ReadTSIMInfo.py.

◆ tgrl

ROOT tgrl = ROOT.TTree('tgrl', 'tree with GRL_Logic')
static

the tree in the output file

Definition at line 366 of file ReadTSIMInfo.py.

◆ time_cluster_t

array time_cluster_t = array('f', 100 * [-99999.])
static

the timing of clusters

Definition at line 414 of file ReadTSIMInfo.py.


The documentation for this class was generated from the following file: