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

Public Member Functions

def set_descr_shifter (self, histogram, description, check)
 
def set_descr_expert (self, histogram, description, check)
 
def set_style (self, histogram, xtitle, ytitle)
 
def get_relative (self, hist1, hist2, title, particle, trgbit)
 
def initialize (self)
 
def beginRun (self)
 
def event (self)
 
def endRun (self)
 
def terminate (self)
 

Public Attributes

 tfile
 output file
 
 Nevent
 number of events
 
 hist_inbit
 validation histogram
 
 hist_outbit
 validation histogram
 
 h_E_ECL
 validation histogram
 
 h_Esum_ECL
 validation histogram
 
 h_theta_ECL
 validation histogram
 
 h_thetaID_ECL
 validation histogram
 
 h_phi_ECL
 validation histogram
 
 h_sector_BKLM
 validation histogram
 
 h_sector_EKLM
 validation histogram
 
 hist_inbit_expert
 validation histogram for experts
 
 hist_outbit_expert
 validation histogram for experts
 
 d_w
 validation histogram
 
 d_w_2
 validation histogram
 
 d_phi
 validation histogram
 
 d_phi_2
 validation histogram
 
 d_phi_3
 validation histogram
 
 d_z0_3d
 validation histogram
 
 d_z0_nn
 validation histogram
 
 d_E_ECL
 validation histogram
 

Detailed Description

Make validation histograms for trg ecl/cdc/klm

Definition at line 358 of file TRGValidation.py.

Member Function Documentation

◆ beginRun()

def beginRun (   self)
reset all histograms at the begin of a new run

Definition at line 621 of file TRGValidation.py.

621 def beginRun(self):
622 '''
623 reset all histograms at the begin of a new run
624 '''
625
626 self.hist_inbit.Reset()
627 self.hist_outbit.Reset()
628 self.hist_inbit_expert.Reset()
629 self.hist_outbit_expert.Reset()
630 self.h_Esum_ECL.Reset()
631 self.h_E_ECL.Reset()
632 self.h_theta_ECL.Reset()
633 self.h_thetaID_ECL.Reset()
634 self.h_phi_ECL.Reset()
635

◆ endRun()

def endRun (   self)
end of run

Definition at line 688 of file TRGValidation.py.

688 def endRun(self):
689 '''
690 end of run
691 '''
692
693 print("end")
694

◆ event()

def event (   self)
event loop

Definition at line 636 of file TRGValidation.py.

636 def event(self):
637 '''
638 event loop
639 '''
640
641 trg_summary = Belle2.PyStoreObj("TRGSummary")
642 clusters = Belle2.PyStoreArray("TRGECLClusters")
643 klmSummary = Belle2.PyStoreObj("KLMTrgSummary")
644
645 for bit in moniInbits:
646 binIndex = inputBits.index(bit["name"])
647 bitIndex = bit["index"]
648 if trg_summary.testInput(bitIndex):
649 self.hist_inbit.Fill(binIndex + 0.5)
650 for bit in moniOutbits:
651 binIndex = outputBits.index(bit["name"])
652 bitIndex = bit["index"]
653 if trg_summary.testFtdl(bitIndex):
654 self.hist_outbit.Fill(binIndex + 0.5)
655
656 for bit in moniInbits_expert:
657 binIndex = inputBits_expert.index(bit["name"])
658 bitIndex = bit["index"]
659 if trg_summary.testInput(bitIndex):
660 self.hist_inbit_expert.Fill(binIndex + 0.5)
661 for bit in moniOutbits_expert:
662 binIndex = outputBits_expert.index(bit["name"])
663 bitIndex = bit["index"]
664 if trg_summary.testFtdl(bitIndex):
665 self.hist_outbit_expert.Fill(binIndex + 0.5)
666
667 etot = 0
668 for cluster in clusters:
669 x = cluster.getPositionX()
670 y = cluster.getPositionY()
671 z = cluster.getPositionZ()
672 e = cluster.getEnergyDep()
673 self.h_E_ECL.Fill(e)
674 etot += e
675 vec = ROOT.Math.XYZVector(x, y, z)
676 self.h_theta_ECL.Fill(vec.Theta() * TMath.RadToDeg())
677 self.h_thetaID_ECL.Fill(cluster.getMaxTCId())
678 self.h_phi_ECL.Fill(vec.Phi() * TMath.RadToDeg())
679
680 if etot > 0:
681 self.h_Esum_ECL.Fill(etot)
682
683 self.h_sector_BKLM.Fill(klmSummary.getBKLM_n_trg_sectors())
684 self.h_sector_EKLM.Fill(klmSummary.getEKLM_n_trg_sectors())
685
686 self.Nevent = self.Nevent + 1
687
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72
a (simplified) python wrapper for StoreObjPtr.
Definition: PyStoreObj.h:67

◆ get_relative()

def get_relative (   self,
  hist1,
  hist2,
  title,
  particle,
  trgbit 
)
Get relative ratio between two hists.
:param hist1 numerator
:param hist2 denominator
:param title new histogram title
:param particle particle name
:param trgbit trigger bit

Definition at line 412 of file TRGValidation.py.

412 def get_relative(self, hist1, hist2, title, particle, trgbit):
413 '''
414 Get relative ratio between two hists.
415 :param hist1 numerator
416 :param hist2 denominator
417 :param title new histogram title
418 :param particle particle name
419 :param trgbit trigger bit
420 '''
421
422 bin1 = hist1.GetXaxis().GetNbins()
423 bin2 = hist2.GetXaxis().GetNbins()
424 Xmin1 = hist1.GetXaxis().GetXmin()
425 Xmax1 = hist1.GetXaxis().GetXmax()
426 Xmin2 = hist2.GetXaxis().GetXmin()
427 Xmax2 = hist2.GetXaxis().GetXmax()
428 if bin1 != bin2 or Xmin1 != Xmin2 or Xmax1 != Xmax2:
429 print("warning: two histograms cannot be compared!!")
430 return 0
431
432 teff = TEfficiency(hist1, hist2)
433 htitle = f'h_eff_{title}_{particle}_{trgbit}'
434 heff = TH1F(htitle, f"efficiency histogram of {particle} for {trgbit}", bin1, Xmin1, Xmax1)
435 for ibin in range(1, bin1+1):
436 binc0 = teff.GetEfficiency(ibin)
437 bine0 = (teff.GetEfficiencyErrorUp(ibin) + teff.GetEfficiencyErrorLow(ibin)) / 2.0
438 heff.SetBinContent(ibin, binc0)
439 heff.SetBinError(ibin, bine0)
440 heff.GetYaxis().SetRangeUser(0.0, 1.2)
441 return heff
442

◆ initialize()

def initialize (   self)
initialize class members

Definition at line 443 of file TRGValidation.py.

443 def initialize(self):
444 '''
445 initialize class members
446 '''
447
448 print('start read')
449
450 self.tfile = TFile("./TRGValidation.root", "update")
451
452 self.Nevent = 0
453
454 m_dbinput = Belle2.PyDBObj("TRGGDLDBInputBits")
455 m_dbftdl = Belle2.PyDBObj("TRGGDLDBFTDLBits")
456 for i in range(320):
457 inbitname = m_dbinput.getinbitname(i)
458 outbitname = m_dbftdl.getoutbitname(i)
459 if inbitname in inputBits:
460 moniInbits.append({"name": inbitname, "index": i})
461 if outbitname in outputBits:
462 moniOutbits.append({"name": outbitname, "index": i})
463 if inbitname in inputBits_expert:
464 moniInbits_expert.append({"name": inbitname, "index": i})
465 if outbitname in outputBits_expert:
466 moniOutbits_expert.append({"name": outbitname, "index": i})
467
468 mc_px = "MCParticles.m_momentum_x"
469 mc_py = "MCParticles.m_momentum_y"
470 trk2d_omega = "TRGCDC2DFinderTracks.m_omega"
471 trk2d_phi = "TRGCDC2DFinderTracks.m_phi0"
472
473 ROOT.gROOT.SetBatch(True)
474 ROOT.gStyle.SetOptStat(ROOT.kFALSE)
475
476
477 n_inbit_test = len(inputBits)
478 self.hist_inbit = TH1F("hin", "trigger input bits", n_inbit_test, 0, n_inbit_test)
479 self.hist_inbit.GetXaxis().SetRangeUser(0, n_inbit_test)
480 self.hist_inbit.GetXaxis().SetLabelSize(0.04)
481 self.set_descr_shifter(self.hist_inbit, "trigger input bits for shifter",
482 "Efficiency should not drop significantly for any trigger bit")
483 self.set_style(self.hist_inbit, "", "Efficiency")
484
485
486 n_oubit_test = len(outputBits)
487 self.hist_outbit = TH1F("hout", "trigger output bits", n_oubit_test, 0, n_oubit_test)
488 self.hist_outbit.GetXaxis().SetRangeUser(0, n_oubit_test)
489 self.hist_outbit.GetXaxis().SetLabelSize(0.04)
490 self.set_descr_shifter(
491 self.hist_outbit,
492 "trigger ftdl bits for shifter",
493 "Efficiency should not drop significantly. For some output bits the efficiency is very low, close to zero.")
494 self.set_style(self.hist_outbit, "", "Efficiency")
495
496 for i in range(n_inbit_test):
497 self.hist_inbit.GetXaxis().SetBinLabel(
498 self.hist_inbit.GetXaxis().FindBin(i + 0.5), inputBits[i])
499 for i in range(n_oubit_test):
500 self.hist_outbit.GetXaxis().SetBinLabel(
501 self.hist_outbit.GetXaxis().FindBin(i + 0.5), outputBits[i])
502
503
504 self.h_E_ECL = TH1F("h_E_ECL", "ECL cluster energy [50 MeV]", 100, 0, 5)
505 self.set_descr_shifter(self.h_E_ECL, "TRG ECL cluster energy", "Exponentially falling distribution")
506 self.set_style(self.h_E_ECL, "ECL cluster energy (GeV)", "Events/(50 MeV)")
507
508
509 self.h_Esum_ECL = TH1F("h_Esum_ECL", "sum of ECL cluster energy [50 MeV]", 100, 0, 5)
510 self.set_descr_shifter(self.h_Esum_ECL, "sum of TRG ECL cluster energy", "Peak at 200 MeV with tail")
511 self.set_style(self.h_Esum_ECL, "sum of ECL cluster energy (GeV)", "Events/(50 MeV)")
512
513
514 self.h_theta_ECL = TH1F("h_theta_ECL", "TRG ECL cluster #theta [1.4 degrees]", 128, 0, 180)
515 self.set_descr_shifter(self.h_theta_ECL, "TRG ECL cluster polar angle", "uniform in barrel")
516 self.set_style(self.h_theta_ECL, "TRG ECL cluster polar angle (degree)", "Events/(1.4 degree)")
517
518
519 self.h_thetaID_ECL = TH1F("h_thetaID_ECL", "ECL cluster TC ID", 610, 0, 610)
520 self.set_descr_shifter(self.h_thetaID_ECL, "TRG ECL cluster theta ID", "uniform in barrel")
521 self.set_style(self.h_thetaID_ECL, "ECL cluster TC ID", "Events/(1)")
522
523
524 self.h_phi_ECL = TH1F("h_phi_ECL", "TRG ECL cluster phi [2.0 degrees]", 180, -180, 180)
525 self.set_descr_shifter(self.h_phi_ECL, "TRG ECL cluster phi distribution", "distribute uniformly")
526 self.set_style(self.h_phi_ECL, "ECL cluster #phi (degree) ", "Events/(2.0 degrees)")
527
528
529 self.h_sector_BKLM = TH1F("h_sector_BKLM", "BKLM TRG hit sector", 10, 0, 10)
530 self.set_descr_shifter(self.h_sector_BKLM, "# of BKLM TRG sector", "peak at 0")
531 self.set_style(self.h_sector_BKLM, "# of BKLM TRG sector", "Events/(1)")
532
533
534 self.h_sector_EKLM = TH1F("h_sector_EKLM", "EKLM TRG hit sector", 10, 0, 10)
535 self.set_descr_shifter(self.h_sector_EKLM, "# of EKLM TRG sector", "peak at 0")
536 self.set_style(self.h_sector_EKLM, "# of EKLM TRG sector", "Events/(1)")
537
538
539 n_inbit_expert = len(inputBits_expert)
540 self.hist_inbit_expert = TH1F("hin_expert", "trigger input bits", n_inbit_expert, 0, n_inbit_expert)
541 self.hist_inbit_expert.GetXaxis().SetRangeUser(0, n_inbit_expert)
542 self.set_descr_expert(self.hist_inbit_expert, "# of EKLM TRG sector", "peak at 0")
543 self.set_style(self.hist_inbit_expert, "", "Efficiency")
544
545
546 n_oubit_expert = len(outputBits_expert)
547 self.hist_outbit_expert = TH1F("hout_expert", "trigger output bits", n_oubit_expert, 0, n_oubit_expert)
548 self.hist_outbit_expert.GetXaxis().SetRangeUser(0, n_oubit_expert)
549 self.set_descr_expert(self.hist_outbit_expert, "# of EKLM TRG sector", "peak at 0")
550 self.set_style(self.hist_outbit_expert, "", "Efficiency")
551
552 for i in range(n_inbit_expert):
553 self.hist_inbit_expert.GetXaxis().SetBinLabel(
554 self.hist_inbit_expert.GetXaxis().FindBin(i + 0.5), inputBits_expert[i])
555 for i in range(n_oubit_expert):
556 self.hist_outbit_expert.GetXaxis().SetBinLabel(
557 self.hist_outbit_expert.GetXaxis().FindBin(i + 0.5), outputBits_expert[i])
558
559 mc = "abs(MCParticles.m_pdg)==11&&MCParticles.m_status==11"
560 tree = ROOT.TChain("tree")
561 tree.Add("../TRGValidationGen.root")
562
563
564 self.d_w = TH1F("d_w", "#Deltaw of CDC 2D finder, w = 0.00449/p_{t}", 50, -0.02, 0.02)
565 tree.Draw(f"({trk2d_omega} - 0.00449)/sqrt({mc_px}*{mc_px} + {mc_py}*{mc_py})>> d_w",
566 "MCParticles.m_pdg<0&&" + mc)
567
568 self.d_w_2 = TH1F("d_w_2", "d_w_2", 50, -0.02, 0.02)
569 tree.Draw(f"({trk2d_omega} + 0.00449)/sqrt({mc_px}*{mc_px} + {mc_py}*{mc_py}) >> d_w_2",
570 "MCParticles.m_pdg>0 && " + mc)
571 self.d_w.Add(self.d_w_2)
572 self.set_descr_shifter(self.d_w, "Comparison on w (=0.00449/pt) of a track between CDC 2D finder output and MC.",
573 "A clear peak at 0 with tail.")
574 self.set_style(self.d_w, "#Deltaw", "Events/(0.08)")
575
576
577 self.d_phi = TH1F("d_phi", "#Delta#phi of CDC 2D finder", 50, -0.5, 0.5)
578 tree.Draw(f"{trk2d_phi}-atan({mc_py}/{mc_px})>>d_phi",
579 "MCParticles.m_status==11&&abs(MCParticles.m_pdg)==11",
580 f"fabs({trk2d_phi}-atan({mc_py}/{mc_px}))<{PI}&&" + mc)
581
582
583 self.d_phi_2 = TH1F("d_phi_2", "d_phi_2", 50, -0.5, 0.5)
584 tree.Draw(f"{trk2d_phi}-atan({mc_py}/{mc_px})-{PI}>>d_phi_2",
585 "MCParticles.m_status==11&&abs(MCParticles.m_pdg)==11",
586 f"{trk2d_phi}-atan({mc_py}/{mc_px})>={PI} &&" + mc)
587
588
589 self.d_phi_3 = TH1F("d_phi_3", "d_phi_3", 50, -0.5, 0.5)
590 tree.Draw(f"{trk2d_phi}-atan({mc_py}/{mc_px})+{PI}>>d_phi_2",
591 "MCParticles.m_status==11&&abs(MCParticles.m_pdg)==11",
592 f"{trk2d_phi}-atan({mc_py}/{mc_px})<=-{PI}&&" + mc)
593
594 self.d_phi.Add(self.d_phi_2)
595 self.d_phi.Add(self.d_phi_3)
596 self.set_descr_shifter(self.d_phi, "Comparison on phi_i of a track between CDC 2D finder output and MC.",
597 "A Gaussian peak at 0.")
598 self.set_style(self.d_phi, "#Delta#phi [rad]", "Events/(0.02 rad)")
599
600
601 self.d_z0_3d = TH1F("d_z0_3d", "#Deltaz0 of CDC 3D fitter", 60, -30, 30)
602 tree.Draw("TRGCDC3DFitterTracks.m_z0-MCParticles.m_productionVertex_z>>d_z0_3d", mc)
603 self.set_descr_shifter(self.d_z0_3d, "Comparison on z0 of a track between CDC 2D fitter output and MC.",
604 "A Gaussian peak at 0 with small tail.")
605 self.set_style(self.d_z0_3d, "#Deltaz0 [cm]", "Events/(1 cm)")
606
607
608 self.d_z0_nn = TH1F("d_z0_nn", "#Deltaz0 of CDC Neuro", 60, -30, 30)
609 tree.Draw("TRGCDCNeuroTracks.m_z0-MCParticles.m_productionVertex_z>>d_z0_nn", mc)
610 self.set_descr_shifter(self.d_z0_nn, "Comparison on z0 of a track between CDC Neuro output and MC.",
611 "A Gaussian peak at 0 with small tail.")
612 self.set_style(self.d_z0_nn, "#Deltaz0 [cm]", "Events/(1 cm)")
613
614
615 self.d_E_ECL = TH1F("d_E_ECL", "#DeltaE of ECL clustering", 100, -6, 6)
616 tree.Draw("TRGECLClusters.m_edep-MCParticles.m_energy>>d_E_ECL", mc)
617 self.set_descr_shifter(self.d_E_ECL, "Comparison on deposit energy between ECL cluster output and MC.",
618 "A peak around -0.5 ~ 0 with a tail toward -6.")
619 self.set_style(self.d_E_ECL, "#DeltaE [GeV]", "Events/(0.12 GeV)")
620
Class to access a DBObjPtr from Python.
Definition: PyDBObj.h:50

◆ set_descr_expert()

def set_descr_expert (   self,
  histogram,
  description,
  check 
)
Sets description, check and contact to validation histogram.
:param histogram validation histogram
:param description description text
:param check information on what to check in comparison with reference

Definition at line 380 of file TRGValidation.py.

380 def set_descr_expert(self, histogram, description, check):
381 '''
382 Sets description, check and contact to validation histogram.
383 :param histogram validation histogram
384 :param description description text
385 :param check information on what to check in comparison with reference
386 '''
387
388 descr = TNamed("Description", description)
389 histogram.GetListOfFunctions().Add(descr)
390 Check = TNamed("Check", check)
391 histogram.GetListOfFunctions().Add(Check)
392 contact = TNamed("Contact", "yinjh@nankai.edu.cn")
393 histogram.GetListOfFunctions().Add(contact)
394 Meta = TNamed("MetaOptions", "expert,pvalue-warn=0.02,pvalue-error=0.01")
395 histogram.GetListOfFunctions().Add(Meta)
396

◆ set_descr_shifter()

def set_descr_shifter (   self,
  histogram,
  description,
  check 
)
Sets description, check and contact to validation histogram.
:param histogram validation histogram
:param description description text
:param check information on what to check in comparison with reference

Definition at line 363 of file TRGValidation.py.

363 def set_descr_shifter(self, histogram, description, check):
364 '''
365 Sets description, check and contact to validation histogram.
366 :param histogram validation histogram
367 :param description description text
368 :param check information on what to check in comparison with reference
369 '''
370
371 descr = TNamed("Description", description)
372 histogram.GetListOfFunctions().Add(descr)
373 Check = TNamed("Check", check)
374 histogram.GetListOfFunctions().Add(Check)
375 contact = TNamed("Contact", "yinjh@nankai.edu.cn")
376 histogram.GetListOfFunctions().Add(contact)
377 Meta = TNamed("MetaOptions", "shifter,pvalue-warn=0.02,pvalue-error=0.01")
378 histogram.GetListOfFunctions().Add(Meta)
379

◆ set_style()

def set_style (   self,
  histogram,
  xtitle,
  ytitle 
)
Sets x-y titles, and sets histogram style.
:param histogram validation histogram
:param xtitle X-axis title
:param xtitle Y-axis title

Definition at line 397 of file TRGValidation.py.

397 def set_style(self, histogram, xtitle, ytitle):
398 '''
399 Sets x-y titles, and sets histogram style.
400 :param histogram validation histogram
401 :param xtitle X-axis title
402 :param xtitle Y-axis title
403 '''
404
405 histogram.GetXaxis().SetTitle(xtitle)
406 histogram.GetYaxis().SetTitle(ytitle)
407 histogram.GetYaxis().SetTitleOffset(1.4)
408 histogram.GetXaxis().SetTitleSize(0.045)
409 histogram.GetYaxis().SetLabelSize(0.020)
410 histogram.SetLineColor(ROOT.kBlack)
411

◆ terminate()

def terminate (   self)
write histograms

Definition at line 695 of file TRGValidation.py.

695 def terminate(self):
696 '''
697 write histograms
698 '''
699
700 bits = ['ty_0', 'ehigh', 'clst_0']
701 part = 'e'
702 h_gen_p = self.tfile.Get(f"{part}_all/p")
703 h_gen_E = self.tfile.Get(f"{part}_all/clusterE")
704 h_gen_theta = self.tfile.Get(f"{part}_all/theta")
705 h_gen_phi = self.tfile.Get(f"{part}_all/phi")
706
707 for bit in bits:
708 h_p = self.tfile.Get(f"{part}_{bit}/p")
709 h_E = self.tfile.Get(f"{part}_{bit}/clusterE")
710 h_theta = self.tfile.Get(f"{part}_{bit}/theta")
711 h_phi = self.tfile.Get(f"{part}_{bit}/phi")
712
713 h_eff_p = self.get_relative(h_p, h_gen_p, 'p', part, bit)
714 self.set_style(h_eff_p, "Momentum", "Efficiency")
715 self.set_descr_expert(h_eff_p, "Momentum dependent efficiency for experts", "")
716
717 h_eff_E = self.get_relative(h_E, h_gen_E, 'E', part, bit)
718 self.set_style(h_eff_E, "Energy", "Efficiency")
719 if bit == 'ehigh':
720 self.set_descr_shifter(h_eff_E, "Energy dependent efficiency for shifter",
721 "Efficiency around 40% below 1 GeV, then jump up to 100%")
722 elif bit == 'clst_0':
723 self.set_descr_shifter(h_eff_E, "Energy dependent efficiency for shifter",
724 "Turn-on curve from 40% efficiency at 0.5 GeV to nearly 100% above 1.5 GeV")
725 else:
726 self.set_descr_expert(h_eff_E, "Energy dependent efficiency for experts", "")
727
728 h_eff_theta = self.get_relative(h_theta, h_gen_theta, 'theta', part, bit)
729 self.set_style(h_eff_theta, "Polar angle", "Efficiency")
730 self.set_descr_expert(h_eff_theta, "polar angle dependent efficiency for experts", "")
731
732 h_eff_phi = self.get_relative(h_phi, h_gen_phi, 'phi', part, bit)
733 self.set_style(h_eff_phi, "phi angle", "Efficiency")
734 self.set_descr_expert(h_eff_phi, "azimuth angle dependent efficiency for experts", "")
735
736 h_eff_p.Write()
737 h_eff_E.Write()
738 h_eff_theta.Write()
739 h_eff_phi.Write()
740
741 bits = ['ty_0', 'klm_0', "klm_hit"]
742 part = "mu"
743 h_gen_p = self.tfile.Get(f"{part}_all/p")
744 h_gen_E = self.tfile.Get(f"{part}_all/clusterE")
745 h_gen_theta = self.tfile.Get(f"{part}_all/theta")
746 h_gen_phi = self.tfile.Get(f"{part}_all/phi")
747
748 for bit in bits:
749 h_p = self.tfile.Get(f"{part}_{bit}/p")
750 h_theta = self.tfile.Get(f"{part}_{bit}/theta")
751 h_phi = self.tfile.Get(f"{part}_{bit}/phi")
752
753 h_eff_p = self.get_relative(h_p, h_gen_p, 'p', part, bit)
754 self.set_style(h_eff_p, "Momentum", "Efficiency")
755 if bit == 'ty_0':
756 self.set_descr_shifter(h_eff_p, "Momentum dependent efficiency for shifter",
757 "Efficiency should be above 90% for the whole momentum range")
758 elif bit == 'klm_0':
759 self.set_descr_shifter(h_eff_p, "Momentum dependent efficiency for shifter",
760 "Efficiency should rise to about 90% for momenta greater than 1 GeV")
761 else:
762 self.set_descr_expert(h_eff_p, "Momentum dependent efficiency for experts", "")
763
764 h_eff_theta = self.get_relative(h_theta, h_gen_theta, 'theta', part, bit)
765 self.set_style(h_eff_theta, "Polar angle", "Efficiency")
766 self.set_descr_expert(h_eff_theta, "polar angle dependent efficiency for experts", "")
767
768 h_eff_phi = self.get_relative(h_phi, h_gen_phi, 'phi', part, bit)
769 self.set_style(h_eff_phi, "phi angle", "Efficiency")
770 self.set_descr_expert(h_eff_phi, "azimuth angle dependent efficiency for experts", "")
771
772 h_eff_p.Write()
773 h_eff_theta.Write()
774 h_eff_phi.Write()
775
776 self.hist_inbit.Scale(1./self.Nevent)
777 self.hist_outbit.Scale(1.0/self.Nevent)
778 self.hist_inbit_expert.Scale(1./self.Nevent)
779 self.hist_outbit_expert.Scale(1.0/self.Nevent)
780
781 self.hist_inbit.Write()
782 self.hist_outbit.Write()
783 self.hist_inbit_expert.Write()
784 self.hist_outbit_expert.Write()
785 self.h_Esum_ECL.Write()
786 self.h_E_ECL.Write()
787 self.h_theta_ECL.Write()
788 self.h_thetaID_ECL.Write()
789 self.h_phi_ECL.Write()
790 self.h_sector_BKLM.Write()
791 self.h_sector_EKLM.Write()
792 self.d_w.Write()
793 self.d_phi.Write()
794 self.d_z0_3d.Write()
795 self.d_z0_nn.Write()
796 self.d_E_ECL.Write()
797 self.tfile.Close()
798

Member Data Documentation

◆ d_E_ECL

d_E_ECL

validation histogram

Definition at line 615 of file TRGValidation.py.

◆ d_phi

d_phi

validation histogram

Definition at line 577 of file TRGValidation.py.

◆ d_phi_2

d_phi_2

validation histogram

Definition at line 583 of file TRGValidation.py.

◆ d_phi_3

d_phi_3

validation histogram

Definition at line 589 of file TRGValidation.py.

◆ d_w

d_w

validation histogram

Definition at line 564 of file TRGValidation.py.

◆ d_w_2

d_w_2

validation histogram

Definition at line 568 of file TRGValidation.py.

◆ d_z0_3d

d_z0_3d

validation histogram

Definition at line 601 of file TRGValidation.py.

◆ d_z0_nn

d_z0_nn

validation histogram

Definition at line 608 of file TRGValidation.py.

◆ h_E_ECL

h_E_ECL

validation histogram

Definition at line 504 of file TRGValidation.py.

◆ h_Esum_ECL

h_Esum_ECL

validation histogram

Definition at line 509 of file TRGValidation.py.

◆ h_phi_ECL

h_phi_ECL

validation histogram

Definition at line 524 of file TRGValidation.py.

◆ h_sector_BKLM

h_sector_BKLM

validation histogram

Definition at line 529 of file TRGValidation.py.

◆ h_sector_EKLM

h_sector_EKLM

validation histogram

Definition at line 534 of file TRGValidation.py.

◆ h_theta_ECL

h_theta_ECL

validation histogram

Definition at line 514 of file TRGValidation.py.

◆ h_thetaID_ECL

h_thetaID_ECL

validation histogram

Definition at line 519 of file TRGValidation.py.

◆ hist_inbit

hist_inbit

validation histogram

Definition at line 478 of file TRGValidation.py.

◆ hist_inbit_expert

hist_inbit_expert

validation histogram for experts

Definition at line 540 of file TRGValidation.py.

◆ hist_outbit

hist_outbit

validation histogram

Definition at line 487 of file TRGValidation.py.

◆ hist_outbit_expert

hist_outbit_expert

validation histogram for experts

Definition at line 547 of file TRGValidation.py.

◆ Nevent

Nevent

number of events

Definition at line 452 of file TRGValidation.py.

◆ tfile

tfile

output file

Definition at line 450 of file TRGValidation.py.


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