Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-02-20 10:36:59

0001 // ----------------------------------------------------------------------------
0002 // 'FillBHCalCalibrationTupleProcessor.cc'
0003 // Derek Anderson
0004 // 03.02.2023
0005 //
0006 // A simple JANA plugin to compare the
0007 // reconstructed hit and cluster energy
0008 // in the HCal to simulated particles.
0009 // ----------------------------------------------------------------------------
0010 
0011 // C includes
0012 #include <vector>
0013 #include <string>
0014 #include <boost/math/special_functions/sign.hpp>
0015 // eicrecon includes 
0016 #include <services/rootfile/RootFile_service.h>
0017 // user includes
0018 #include "FillBHCalCalibrationTupleProcessor.h"
0019 
0020 // The following just makes this a JANA plugin
0021 extern "C" {
0022   void InitPlugin(JApplication *app) {
0023   InitJANAPlugin(app);
0024   app -> Add(new FillBHCalCalibrationTupleProcessor);
0025   }
0026 }
0027 
0028 
0029 
0030 //-------------------------------------------
0031 // InitWithGlobalRootLock
0032 //-------------------------------------------
0033 void FillBHCalCalibrationTupleProcessor::InitWithGlobalRootLock(){
0034 
0035   // create directory in output file
0036   auto rootfile_svc = GetApplication() -> GetService<RootFile_service>();
0037   auto rootfile     = rootfile_svc     -> GetHistFile();
0038   rootfile -> mkdir("FillBHCalCalibrationTuple") -> cd();
0039 
0040   // initialize bhcal histograms
0041   const unsigned long nNumBin(200);
0042   const unsigned long nChrgBin(6);
0043   const unsigned long nMassBin(1000);
0044   const unsigned long nPhiBin(60);
0045   const unsigned long nEtaBin(40);
0046   const unsigned long nEneBin(200);
0047   const unsigned long nMomBin(200);
0048   const unsigned long nPosTrBin(800);
0049   const unsigned long nPosLoBin(30);
0050   const unsigned long nDiffBin(200);
0051   const unsigned long rNumBin[CONST::NRange]   = {0,      200};
0052   const double        rChrgBin[CONST::NRange]  = {-3.,    3.};
0053   const double        rMassBin[CONST::NRange]  = {0.,     5.};
0054   const double        rPhiBin[CONST::NRange]   = {-3.15,  3.15};
0055   const double        rEtaBin[CONST::NRange]   = {-2.,    2.};
0056   const double        rEneBin[CONST::NRange]   = {0.,     100.};
0057   const double        rMomBin[CONST::NRange]   = {-50.,   50.};
0058   const double        rPosTrBin[CONST::NRange] = {-4000., 4000.};
0059   const double        rPosLoBin[CONST::NRange] = {-3000., 3000.};
0060   const double        rDiffBin[CONST::NRange]  = {-5.,    5.};
0061   // particle histograms
0062   hParChrg                   = new TH1D("hParChrg",     "Gen. Particles", nChrgBin, rChrgBin[0], rChrgBin[1]);
0063   hParMass                   = new TH1D("hParMass",     "Gen. Particles", nMassBin, rMassBin[0], rMassBin[1]);
0064   hParPhi                    = new TH1D("hParPhi",      "Gen. Particles", nPhiBin,  rPhiBin[0],  rPhiBin[1]);
0065   hParEta                    = new TH1D("hParEta",      "Gen. Particles", nEtaBin,  rEtaBin[0],  rEtaBin[1]);
0066   hParEne                    = new TH1D("hParEne",      "Gen. Particles", nEneBin,  rEneBin[0],  rEneBin[1]);
0067   hParMom                    = new TH1D("hParMom",      "Gen. Particles", nEneBin,  rEneBin[0],  rEneBin[1]);
0068   hParMomX                   = new TH1D("hParMomX",     "Gen. Particles", nMomBin,  rMomBin[0],  rMomBin[1]);
0069   hParMomY                   = new TH1D("hParMomY",     "Gen. Particles", nMomBin,  rMomBin[0],  rMomBin[1]);
0070   hParMomZ                   = new TH1D("hParMomZ",     "Gen. Particles", nMomBin,  rMomBin[0],  rMomBin[1]);
0071   hParEtaVsPhi               = new TH2D("hParEtaVsPhi", "Gen. Particles", nPhiBin,  rPhiBin[0],  rPhiBin[1], nEtaBin, rEtaBin[0], rEtaBin[1]);
0072   // reco. bhcal hit histograms
0073   hHCalRecHitPhi             = new TH1D("hHCalRecHitPhi",      "Barrel HCal", nPhiBin,   rPhiBin[0],   rPhiBin[1]);
0074   hHCalRecHitEta             = new TH1D("hHCalRecHitEta",      "Barrel HCal", nEtaBin,   rEtaBin[0],   rEtaBin[1]);
0075   hHCalRecHitEne             = new TH1D("hHCalRecHitEne",      "Barrel HCal", nEneBin,   rEneBin[0],   rEneBin[1]);
0076   hHCalRecHitPosZ            = new TH1D("hHCalRecHitPosZ",     "Barrel HCal", nPosLoBin, rPosLoBin[0], rPosLoBin[1]);
0077   hHCalRecHitParDiff         = new TH1D("hHCalRecHitParDiff",  "Barrel HCal", nDiffBin,  rDiffBin[0],  rDiffBin[1]);
0078   hHCalRecHitPosYvsX         = new TH2D("hHCalRecHitPosYvsX",  "Barrel HCal", nPosTrBin, rPosTrBin[0], rPosTrBin[1], nPosTrBin, rPosTrBin[0], rPosTrBin[1]);
0079   hHCalRecHitEtaVsPhi        = new TH2D("hHCalRecHitEtaVsPhi", "Barrel HCal", nPhiBin,   rPhiBin[0],   rPhiBin[1],   nEtaBin,   rEtaBin[0],   rEtaBin[1]);
0080   hHCalRecHitVsParEne        = new TH2D("hHCalRecHitVsParEne", "Barrel HCal", nEneBin,   rEneBin[0],   rEneBin[1],   nEneBin,   rEneBin[0],   rEneBin[1]);
0081   // bhcal cluster hit histograms
0082   hHCalClustHitPhi           = new TH1D("hHCalClustHitPhi",      "Barrel HCal", nPhiBin,   rPhiBin[0],   rPhiBin[1]);
0083   hHCalClustHitEta           = new TH1D("hHCalClustHitEta",      "Barrel HCal", nEtaBin,   rEtaBin[0],   rEtaBin[1]);
0084   hHCalClustHitEne           = new TH1D("hHCalClustHitEne",      "Barrel HCal", nEneBin,   rEneBin[0],   rEneBin[1]);
0085   hHCalClustHitPosZ          = new TH1D("hHCalClustHitPosZ",     "Barrel HCal", nPosLoBin, rPosLoBin[0], rPosLoBin[1]);
0086   hHCalClustHitParDiff       = new TH1D("hHCalClustHitParDiff",  "Barrel HCal", nDiffBin,  rDiffBin[0],  rDiffBin[1]);
0087   hHCalClustHitPosYvsX       = new TH2D("hHCalClustHitPosYvsX",  "Barrel HCal", nPosTrBin, rPosTrBin[0], rPosTrBin[1], nPosTrBin, rPosTrBin[0], rPosTrBin[1]);
0088   hHCalClustHitEtaVsPhi      = new TH2D("hHCalClustHitEtaVsPhi", "Barrel HCal", nPhiBin,   rPhiBin[0],   rPhiBin[1],   nEtaBin,   rEtaBin[0],   rEtaBin[1]);
0089   hHCalClustHitVsParEne      = new TH2D("hHCalClustHitVsParEne", "Barrel HCal", nEneBin,   rEneBin[0],   rEneBin[1],   nEneBin,   rEneBin[0],   rEneBin[1]);
0090   // reco. bhcal cluster histograms
0091   hHCalClustPhi              = new TH1D("hHCalClustPhi",      "Barrel HCal", nPhiBin,   rPhiBin[0],   rPhiBin[1]);
0092   hHCalClustEta              = new TH1D("hHCalClustEta",      "Barrel HCal", nEtaBin,   rEtaBin[0],   rEtaBin[1]);
0093   hHCalClustEne              = new TH1D("hHCalClustEne",      "Barrel HCal", nEneBin,   rEneBin[0],   rEneBin[1]);
0094   hHCalClustPosZ             = new TH1D("hHCalClustPosZ",     "Barrel HCal", nPosLoBin, rPosLoBin[0], rPosLoBin[1]);
0095   hHCalClustNumHit           = new TH1I("hHCalClustNumHit",   "Barrel HCal", nNumBin,   rNumBin[0],   rNumBin[1]);
0096   hHCalClustParDiff          = new TH1D("hHCalClustParDiff",  "Barrel HCal", nDiffBin,  rDiffBin[0],  rDiffBin[1]);
0097   hHCalClustPosYvsX          = new TH2D("hHCalClustPosYvsX",  "Barrel HCal", nPosTrBin, rPosTrBin[0], rPosTrBin[1], nPosTrBin, rPosTrBin[0], rPosTrBin[1]);
0098   hHCalClustEtaVsPhi         = new TH2D("hHCalClustEtaVsPhi", "Barrel HCal", nPhiBin,   rPhiBin[0],   rPhiBin[1],   nEtaBin,   rEtaBin[0],   rEtaBin[1]);
0099   hHCalClustVsParEne         = new TH2D("hHCalClustVsParEne", "Barrel HCal", nEneBin,   rEneBin[0],   rEneBin[1],   nEneBin,   rEneBin[0],   rEneBin[1]);
0100   // bhcal cluster hit histograms
0101   hHCalTruClustHitPhi        = new TH1D("hHCalTruClustHitPhi",      "Barrel HCal", nPhiBin,   rPhiBin[0],   rPhiBin[1]);
0102   hHCalTruClustHitEta        = new TH1D("hHCalTruClustHitEta",      "Barrel HCal", nEtaBin,   rEtaBin[0],   rEtaBin[1]);
0103   hHCalTruClustHitEne        = new TH1D("hHCalTruClustHitEne",      "Barrel HCal", nEneBin,   rEneBin[0],   rEneBin[1]);
0104   hHCalTruClustHitPosZ       = new TH1D("hHCalTruClustHitPosZ",     "Barrel HCal", nPosLoBin, rPosLoBin[0], rPosLoBin[1]);
0105   hHCalTruClustHitParDiff    = new TH1D("hHCalTruClustHitParDiff",  "Barrel HCal", nDiffBin,  rDiffBin[0],  rDiffBin[1]);
0106   hHCalTruClustHitPosYvsX    = new TH2D("hHCalTruClustHitPosYvsX",  "Barrel HCal", nPosTrBin, rPosTrBin[0], rPosTrBin[1], nPosTrBin, rPosTrBin[0], rPosTrBin[1]);
0107   hHCalTruClustHitEtaVsPhi   = new TH2D("hHCalTruClustHitEtaVsPhi", "Barrel HCal", nPhiBin,   rPhiBin[0],   rPhiBin[1],   nEtaBin,   rEtaBin[0],   rEtaBin[1]);
0108   hHCalTruClustHitVsParEne   = new TH2D("hHCalTruClustHitVsParEne", "Barrel HCal", nEneBin,   rEneBin[0],   rEneBin[1],   nEneBin,   rEneBin[0],   rEneBin[1]);
0109   // truth bhcal cluster histograms
0110   hHCalTruClustPhi           = new TH1D("hHCalTruClustPhi",      "Barrel HCal", nPhiBin,   rPhiBin[0],   rPhiBin[1]);
0111   hHCalTruClustEta           = new TH1D("hHCalTruClustEta",      "Barrel HCal", nEtaBin,   rEtaBin[0],   rEtaBin[1]);
0112   hHCalTruClustEne           = new TH1D("hHCalTruClustEne",      "Barrel HCal", nEneBin,   rEneBin[0],   rEneBin[1]);
0113   hHCalTruClustPosZ          = new TH1D("hHCalTruClustPosZ",     "Barrel HCal", nPosLoBin, rPosLoBin[0], rPosLoBin[1]);
0114   hHCalTruClustNumHit        = new TH1I("hHCalTruClustNumHit",   "Barrel HCal", nNumBin,   rNumBin[0],   rNumBin[1]);
0115   hHCalTruClustParDiff       = new TH1D("hHCalTruClustParDiff",  "Barrel HCal", nDiffBin,  rDiffBin[0],  rDiffBin[1]);
0116   hHCalTruClustPosYvsX       = new TH2D("hHCalTruClustPosYvsX",  "Barrel HCal", nPosTrBin, rPosTrBin[0], rPosTrBin[1], nPosTrBin, rPosTrBin[0], rPosTrBin[1]);
0117   hHCalTruClustEtaVsPhi      = new TH2D("hHCalTruClustEtaVsPhi", "Barrel HCal", nPhiBin,   rPhiBin[0],   rPhiBin[1],   nEtaBin,   rEtaBin[0],   rEtaBin[1]);
0118   hHCalTruClustVsParEne      = new TH2D("hHCalTruClustVsParEne", "Barrel HCal", nEneBin,   rEneBin[0],   rEneBin[1],   nEneBin,   rEneBin[0],   rEneBin[1]);
0119   // bhcal general event-wise histograms
0120   hEvtHCalNumPar             = new TH1I("hEvtHCalNumPar", "Barrel HCal", nNumBin, rNumBin[0], rNumBin[1]);
0121   // bhcal hit event-wise histograms
0122   hEvtHCalNumHit             = new TH1I("hEvtHCalNumHit",      "Barrel HCal", nNumBin,  rNumBin[0],  rNumBin[1]);
0123   hEvtHCalSumHitEne          = new TH1D("hEvtHCalSumHitEne",   "Barrel HCal", nEneBin,  rEneBin[0],  rEneBin[1]);
0124   hEvtHCalSumHitDiff         = new TH1D("hEvtHCalSumHitDiff",  "Barrel HCal", nDiffBin, rDiffBin[0], rDiffBin[1]);
0125   hEvtHCalSumHitVsPar        = new TH2D("hEvtHCalSumHitVsPar", "Barrel HCal", nEneBin,  rEneBin[0],  rEneBin[1], nEneBin, rEneBin[0], rEneBin[1]);
0126   // bhcal cluster event-wise histograms
0127   hEvtHCalNumClust           = new TH1I("hEvtHCalNumClust",      "Barrel HCal", nNumBin,  rNumBin[0],  rNumBin[1]);
0128   hEvtHCalSumClustEne        = new TH1D("hEvtHCalSumClustEne",   "Barrel HCal", nEneBin,  rEneBin[0],  rEneBin[1]);
0129   hEvtHCalSumClustDiff       = new TH1D("hEvtHCalSumClustDiff",  "Barrel HCal", nDiffBin, rDiffBin[0], rDiffBin[1]);
0130   hEvtHCalNumClustVsHit      = new TH2I("hEvtHCalNumClustVsHit", "Barrel HCal", nNumBin,  rNumBin[0],  rNumBin[1], nNumBin, rNumBin[0], rNumBin[1]);
0131   hEvtHCalSumClustVsPar      = new TH2D("hEvtHCalSumClustVsPar", "Barrel HCal", nEneBin,  rEneBin[0],  rEneBin[1], nEneBin, rEneBin[0], rEneBin[1]);
0132   // bhcal lead cluster event-wise histograms
0133   hEvtHCalLeadClustNumHit    = new TH1I("hEvtHCalLeadClustNumHit", "Barrel HCal", nNumBin,  rNumBin[0],  rNumBin[1]);
0134   hEvtHCalLeadClustEne       = new TH1D("hEvtHCalLeadClustEne",    "Barrel HCal", nEneBin,  rEneBin[0],  rEneBin[1]);
0135   hEvtHCalLeadClustDiff      = new TH1D("hEvtHCalLeadClustDiff",   "Barrel HCal", nDiffBin, rDiffBin[0], rDiffBin[1]);
0136   hEvtHCalLeadClustVsPar     = new TH2D("hEvtHCalLeadClustVsPar",  "Barrel HCal", nEneBin,  rEneBin[0],  rEneBin[1], nEneBin, rEneBin[0], rEneBin[1]);
0137   // bhcal truth cluster event-wise histograms
0138   hEvtHCalNumTruClust        = new TH1I("hEvtHCalNumTruClust",        "Barrel HCal", nNumBin,  rNumBin[0],  rNumBin[1]);
0139   hEvtHCalSumTruClustEne     = new TH1D("hEvtHCalSumTruClustEne",     "Barrel HCal", nEneBin,  rEneBin[0],  rEneBin[1]);
0140   hEvtHCalSumTruClustDiff    = new TH1D("hEvtHCalSumTruClustDiff",    "Barrel HCal", nDiffBin, rDiffBin[0], rDiffBin[1]);
0141   hEvtHCalNumTruClustVsClust = new TH2I("hEvtHCalNumTruClustVsClust", "Barrel HCal", nNumBin,  rNumBin[0],  rNumBin[1], nNumBin, rNumBin[0], rNumBin[1]);
0142   hEvtHCalSumTruClustVsPar   = new TH2D("hEvtHCalSumTruClustVsPar",   "Barrel HCal", nEneBin,  rEneBin[0],  rEneBin[1], nEneBin, rEneBin[0], rEneBin[1]);
0143   // bhcal truth lead cluster event-wise histograms
0144   hEvtHCalLeadTruClustNumHit = new TH1I("hEvtHCalLeadTruClustNumHit", "Barrel HCal", nNumBin,  rNumBin[0],  rNumBin[1]);
0145   hEvtHCalLeadTruClustEne    = new TH1D("hEvtHCalLeadTruClustEne",    "Barrel HCal", nEneBin,  rEneBin[0],  rEneBin[1]);
0146   hEvtHCalLeadTruClustDiff   = new TH1D("hEvtHCalLeadTruClustDiff",   "Barrel HCal", nDiffBin, rDiffBin[0], rDiffBin[1]);
0147   hEvtHCalLeadTruClustVsPar  = new TH2D("hEvtHCalLeadTruClustVsPar",  "Barrel HCal", nEneBin,  rEneBin[0],  rEneBin[1], nEneBin, rEneBin[0], rEneBin[1]);
0148   // bhcal particle errors
0149   hParChrg                   -> Sumw2();
0150   hParMass                   -> Sumw2();
0151   hParPhi                    -> Sumw2();
0152   hParEta                    -> Sumw2();
0153   hParEne                    -> Sumw2();
0154   hParMom                    -> Sumw2();
0155   hParMomX                   -> Sumw2();
0156   hParMomY                   -> Sumw2();
0157   hParMomZ                   -> Sumw2();
0158   hParEtaVsPhi               -> Sumw2();
0159   // bhcal reco. hit errors
0160   hHCalRecHitPhi             -> Sumw2();
0161   hHCalRecHitEta             -> Sumw2();
0162   hHCalRecHitEne             -> Sumw2();
0163   hHCalRecHitPosZ            -> Sumw2();
0164   hHCalRecHitParDiff         -> Sumw2();
0165   hHCalRecHitPosYvsX         -> Sumw2();
0166   hHCalRecHitEtaVsPhi        -> Sumw2();
0167   hHCalRecHitVsParEne        -> Sumw2();
0168   // bhcal reco. cluster hit errors
0169   hHCalClustHitPhi           -> Sumw2();
0170   hHCalClustHitEta           -> Sumw2();
0171   hHCalClustHitEne           -> Sumw2();
0172   hHCalClustHitPosZ          -> Sumw2();
0173   hHCalClustHitParDiff       -> Sumw2();
0174   hHCalClustHitPosYvsX       -> Sumw2();
0175   hHCalClustHitEtaVsPhi      -> Sumw2();
0176   hHCalClustHitVsParEne      -> Sumw2();
0177   // bhcal reco. cluster errors
0178   hHCalClustPhi              -> Sumw2();
0179   hHCalClustEta              -> Sumw2();
0180   hHCalClustEne              -> Sumw2();
0181   hHCalClustPosZ             -> Sumw2();
0182   hHCalClustNumHit           -> Sumw2();
0183   hHCalClustParDiff          -> Sumw2();
0184   hHCalClustPosYvsX          -> Sumw2();
0185   hHCalClustEtaVsPhi         -> Sumw2();
0186   hHCalClustVsParEne         -> Sumw2();
0187   // bhcal truth cluster hit errors
0188   hHCalTruClustHitPhi        -> Sumw2();
0189   hHCalTruClustHitEta        -> Sumw2();
0190   hHCalTruClustHitEne        -> Sumw2();
0191   hHCalTruClustHitPosZ       -> Sumw2();
0192   hHCalTruClustHitParDiff    -> Sumw2();
0193   hHCalTruClustHitPosYvsX    -> Sumw2();
0194   hHCalTruClustHitEtaVsPhi   -> Sumw2();
0195   hHCalTruClustHitVsParEne   -> Sumw2();
0196   // bhcal truth cluster errors
0197   hHCalTruClustPhi           -> Sumw2();
0198   hHCalTruClustEta           -> Sumw2();
0199   hHCalTruClustEne           -> Sumw2();
0200   hHCalTruClustPosZ          -> Sumw2();
0201   hHCalTruClustNumHit        -> Sumw2();
0202   hHCalTruClustParDiff       -> Sumw2();
0203   hHCalTruClustPosYvsX       -> Sumw2();
0204   hHCalTruClustEtaVsPhi      -> Sumw2();
0205   hHCalTruClustVsParEne      -> Sumw2();
0206   // bhcal general event-wise errors
0207   hEvtHCalNumPar             -> Sumw2();
0208   // bhcal hit event-wise errors
0209   hEvtHCalNumHit             -> Sumw2();
0210   hEvtHCalSumHitEne          -> Sumw2();
0211   hEvtHCalSumHitDiff         -> Sumw2();
0212   hEvtHCalSumHitVsPar        -> Sumw2();
0213   // bhcal cluster event-wise errors
0214   hEvtHCalNumClust           -> Sumw2();
0215   hEvtHCalSumClustEne        -> Sumw2();
0216   hEvtHCalSumClustDiff       -> Sumw2();
0217   hEvtHCalNumClustVsHit      -> Sumw2();
0218   hEvtHCalSumClustVsPar      -> Sumw2();
0219   // bhcal lead cluster event-wise errors
0220   hEvtHCalLeadClustNumHit    -> Sumw2();
0221   hEvtHCalLeadClustEne       -> Sumw2();
0222   hEvtHCalLeadClustDiff      -> Sumw2();
0223   hEvtHCalLeadClustVsPar     -> Sumw2();
0224   // bhcal truth cluster event-wise errors
0225   hEvtHCalNumTruClust        -> Sumw2();
0226   hEvtHCalSumTruClustEne     -> Sumw2();
0227   hEvtHCalSumTruClustDiff    -> Sumw2();
0228   hEvtHCalNumTruClustVsClust -> Sumw2();
0229   // bhcal truth lead cluster event-wise errors
0230   hEvtHCalLeadTruClustNumHit -> Sumw2();
0231   hEvtHCalLeadTruClustEne    -> Sumw2();
0232   hEvtHCalLeadTruClustDiff   -> Sumw2();
0233   hEvtHCalLeadTruClustVsPar  -> Sumw2();
0234 
0235   // initialize becal histograms
0236   const unsigned long nLayerBin(50);
0237   const unsigned long rLayerBin[CONST::NRange] = {0, 50};
0238   // scifi reconstructed hit histograms
0239   hSciFiRecHitNLayer         = new TH1I("hSciFiRecHitNLayer",      "Barrel SciFi", nLayerBin, rLayerBin[0], rLayerBin[1]);
0240   hSciFiRecHitPhi            = new TH1D("hSciFiRecHitPhi",         "Barrel SciFi", nPhiBin,   rPhiBin[0],   rPhiBin[1]);
0241   hSciFiRecHitEta            = new TH1D("hSciFiRecHitEta",         "Barrel SciFi", nEtaBin,   rEtaBin[0],   rEtaBin[1]);
0242   hSciFiRecHitEne            = new TH1D("hSciFiRecHitEne",         "Barrel SciFi", nEneBin,   rEneBin[0],   rEneBin[1]);
0243   hSciFiRecHitPosZ           = new TH1D("hSciFiRecHitPosZ",        "Barrel SciFi", nPosLoBin, rPosLoBin[0], rPosLoBin[1]);
0244   hSciFiRecHitParDiff        = new TH1D("hSciFiRecHitParDiff",     "Barrel SciFi", nDiffBin,  rDiffBin[0],  rDiffBin[1]);
0245   hSciFiRecHitPosYvsX        = new TH2D("hSciFiRecHitPosYvsX",     "Barrel SciFi", nPosTrBin, rPosTrBin[0], rPosTrBin[1], nPosTrBin, rPosTrBin[0], rPosTrBin[1]);
0246   hSciFiRecHitEtaVsPhi       = new TH2D("hSciFiRecHitEtaVsPhi",    "Barrel SciFi", nPhiBin,   rPhiBin[0],   rPhiBin[1],   nEtaBin,   rEtaBin[0],   rEtaBin[1]);
0247   hSciFiRecHitVsParEne       = new TH2D("hSciFiRecHitVsParEne",    "Barrel SciFi", nEneBin,   rEneBin[0],   rEneBin[1],   nEneBin,   rEneBin[0],   rEneBin[1]);
0248   hSciFiRecHitEneVsNLayer    = new TH2D("hSciFiRecHitEneVsNLayer", "Barrel SciFi", nLayerBin, rLayerBin[0], rLayerBin[1], nEneBin,   rEneBin[0],   rEneBin[1]);
0249   // imaging reconstructed hit histograms
0250   hImageRecHitNLayer         = new TH1I("hImageRecHitNLayer",      "Barrel Image", nLayerBin, rLayerBin[0], rLayerBin[1]);
0251   hImageRecHitPhi            = new TH1D("hImageRecHitPhi",         "Barrel Image", nPhiBin,   rPhiBin[0],   rPhiBin[1]);
0252   hImageRecHitEta            = new TH1D("hImageRecHitEta",         "Barrel Image", nEtaBin,   rEtaBin[0],   rEtaBin[1]);
0253   hImageRecHitEne            = new TH1D("hImageRecHitEne",         "Barrel Image", nEneBin,   rEneBin[0],   rEneBin[1]);
0254   hImageRecHitPosZ           = new TH1D("hImageRecHitPosZ",        "Barrel Image", nPosLoBin, rPosLoBin[0], rPosLoBin[1]);
0255   hImageRecHitParDiff        = new TH1D("hImageRecHitParDiff",     "Barrel Image", nDiffBin,  rDiffBin[0],  rDiffBin[1]);
0256   hImageRecHitPosYvsX        = new TH2D("hImageRecHitPosYvsX",     "Barrel Image", nPosTrBin, rPosTrBin[0], rPosTrBin[1], nPosTrBin, rPosTrBin[0], rPosTrBin[1]);
0257   hImageRecHitEtaVsPhi       = new TH2D("hImageRecHitEtaVsPhi",    "Barrel Image", nPhiBin,   rPhiBin[0],   rPhiBin[1],   nEtaBin,   rEtaBin[0],   rEtaBin[1]);
0258   hImageRecHitVsParEne       = new TH2D("hImageRecHitVsParEne",    "Barrel Image", nEneBin,   rEneBin[0],   rEneBin[1],   nEneBin,   rEneBin[0],   rEneBin[1]);
0259   hImageRecHitEneVsNLayer    = new TH2D("hImageRecHitEneVsNLayer", "Barrel Image", nLayerBin, rLayerBin[0], rLayerBin[1], nEneBin,   rEneBin[0],   rEneBin[1]);
0260   // reco. bemc cluster histograms
0261   hECalClustPhi              = new TH1D("hECalClustPhi",      "Barrel ECal", nPhiBin,   rPhiBin[0],   rPhiBin[1]);
0262   hECalClustEta              = new TH1D("hECalClustEta",      "Barrel ECal", nEtaBin,   rEtaBin[0],   rEtaBin[1]);
0263   hECalClustEne              = new TH1D("hECalClustEne",      "Barrel ECal", nEneBin,   rEneBin[0],   rEneBin[1]);
0264   hECalClustPosZ             = new TH1D("hECalClustPosZ",     "Barrel ECal", nPosLoBin, rPosLoBin[0], rPosLoBin[1]);
0265   hECalClustNumHit           = new TH1I("hECalClustNumHit",   "Barrel ECal", nNumBin,   rNumBin[0],   rNumBin[1]);
0266   hECalClustParDiff          = new TH1D("hECalClustParDiff",  "Barrel ECal", nDiffBin,  rDiffBin[0],  rDiffBin[1]);
0267   hECalClustPosYvsX          = new TH2D("hECalClustPosYvsX",  "Barrel ECal", nPosTrBin, rPosTrBin[0], rPosTrBin[1], nPosTrBin, rPosTrBin[0], rPosTrBin[1]);
0268   hECalClustEtaVsPhi         = new TH2D("hECalClustEtaVsPhi", "Barrel ECal", nPhiBin,   rPhiBin[0],   rPhiBin[1],   nEtaBin,   rEtaBin[0],   rEtaBin[1]);
0269   hECalClustVsParEne         = new TH2D("hECalClustVsParEne", "Barrel ECal", nEneBin,   rEneBin[0],   rEneBin[1],   nEneBin,   rEneBin[0],   rEneBin[1]);
0270   // scifi hit event-wise histogram
0271   hEvtSciFiSumEne            = new TH1D("hEvtSciFiSumEne",          "Barrel SciFi", nEneBin,   rEneBin[0],   rEneBin[1]);
0272   hEvtSciFiSumEneVsNLayer    = new TH2D("hEvtSciFiSumEneVsNLayer",  "Barrel SciFi", nLayerBin, rLayerBin[0], rLayerBin[1], nEneBin, rEneBin[0], rEneBin[1]);
0273   hEvtSciFiVsHCalHitSumEne   = new TH2D("hEvtSciFiVsHCalHitSumEne", "Barrel SciFi", nEneBin,   rEneBin[0],   rEneBin[1],   nEneBin, rEneBin[0], rEneBin[1]);
0274   // imaging hit event-wise histogram
0275   hEvtImageSumEne            = new TH1D("hEvtImageSumEne",          "Barrel Image", nEneBin,   rEneBin[0],   rEneBin[1]);
0276   hEvtImageSumEneVsNLayer    = new TH2D("hEvtImageSumEneVsNLayer",  "Barrel Image", nLayerBin, rLayerBin[0], rLayerBin[1], nEneBin, rEneBin[0], rEneBin[1]);
0277   hEvtImageVsHCalHitSumEne   = new TH2D("hEvtImageVsHCalHitSumEne", "Barrel Image", nEneBin,   rEneBin[0],   rEneBin[1],   nEneBin, rEneBin[0], rEneBin[1]);
0278   // bemc cluster event-wise histograms
0279   hEvtECalNumClust           = new TH1I("hEvtECalNumClust",          "Barrel ECal", nNumBin,  rNumBin[0],  rNumBin[1]);
0280   hEvtECalSumClustEne        = new TH1D("hEvtECalSumClustEne",       "Barrel ECal", nEneBin,  rEneBin[0],  rEneBin[1]);
0281   hEvtECalSumClustDiff       = new TH1D("hEvtECalSumClustDiff",      "Barrel ECal", nDiffBin, rDiffBin[0], rDiffBin[1]);
0282   hEvtECalSumClustVsPar      = new TH2D("hEvtECalSumClustVsPar",     "Barrel ECal", nEneBin,  rEneBin[0],  rEneBin[1], nEneBin, rEneBin[0], rEneBin[1]);
0283   hEvtECalVsHCalSumClustEne  = new TH2D("hEvtECalVsHCalSumClustEne", "Barrel ECal", nEneBin,  rEneBin[0],  rEneBin[1], nEneBin, rEneBin[0], rEneBin[1]);
0284   // bemc lead cluster event-wise histograms
0285   hEvtECalLeadClustNumHit    = new TH1I("hEvtECalLeadClustNumHit",    "Barrel ECal", nNumBin,  rNumBin[0],  rNumBin[1]);
0286   hEvtECalLeadClustEne       = new TH1D("hEvtECalLeadClustEne",       "Barrel ECal", nEneBin,  rEneBin[0],  rEneBin[1]);
0287   hEvtECalLeadClustDiff      = new TH1D("hEvtECalLeadClustDiff",      "Barrel ECal", nDiffBin, rDiffBin[0], rDiffBin[1]);
0288   hEvtECalLeadClustVsPar     = new TH2D("hEvtECalLeadClustVsPar",     "Barrel ECal", nEneBin,  rEneBin[0],  rEneBin[1], nEneBin, rEneBin[0], rEneBin[1]);
0289   hEvtECalVsHCalLeadClustEne = new TH2D("hEvtECalVsHCalLeadClustEne", "Barrel ECal", nEneBin,  rEneBin[0],  rEneBin[1], nEneBin, rEneBin[0], rEneBin[1]);
0290   // scifi reconstructed hit errors
0291   hSciFiRecHitNLayer         -> Sumw2();
0292   hSciFiRecHitPhi            -> Sumw2();
0293   hSciFiRecHitEta            -> Sumw2();
0294   hSciFiRecHitEne            -> Sumw2();
0295   hSciFiRecHitPosZ           -> Sumw2();
0296   hSciFiRecHitParDiff        -> Sumw2();
0297   hSciFiRecHitPosYvsX        -> Sumw2();
0298   hSciFiRecHitEtaVsPhi       -> Sumw2();
0299   hSciFiRecHitVsParEne       -> Sumw2();
0300   hSciFiRecHitEneVsNLayer    -> Sumw2();
0301   // imaging reconstructed hit errors
0302   hImageRecHitNLayer         -> Sumw2();
0303   hImageRecHitPhi            -> Sumw2();
0304   hImageRecHitEta            -> Sumw2();
0305   hImageRecHitEne            -> Sumw2();
0306   hImageRecHitPosZ           -> Sumw2();
0307   hImageRecHitParDiff        -> Sumw2();
0308   hImageRecHitPosYvsX        -> Sumw2();
0309   hImageRecHitEtaVsPhi       -> Sumw2();
0310   hImageRecHitVsParEne       -> Sumw2();
0311   hImageRecHitEneVsNLayer    -> Sumw2();
0312   // bemc reconstructed cluster errors
0313   hECalClustEta              -> Sumw2();
0314   hECalClustPhi              -> Sumw2();
0315   hECalClustEne              -> Sumw2();
0316   hECalClustPosZ             -> Sumw2();
0317   hECalClustNumHit           -> Sumw2();
0318   hECalClustParDiff          -> Sumw2();
0319   hECalClustPosYvsX          -> Sumw2();
0320   hECalClustEtaVsPhi         -> Sumw2();
0321   hECalClustVsParEne         -> Sumw2();
0322   // scifi hit event-wise histogram
0323   hEvtSciFiSumEne            -> Sumw2();
0324   hEvtSciFiSumEneVsNLayer    -> Sumw2();
0325   hEvtSciFiVsHCalHitSumEne   -> Sumw2();
0326   // imaging hit event-wise histogram
0327   hEvtImageSumEne            -> Sumw2();
0328   hEvtImageSumEneVsNLayer    -> Sumw2();
0329   hEvtImageVsHCalHitSumEne   -> Sumw2();
0330   // bemc cluster event-wise errors
0331   hEvtECalNumClust           -> Sumw2();
0332   hEvtECalSumClustEne        -> Sumw2();
0333   hEvtECalSumClustDiff       -> Sumw2();
0334   hEvtECalSumClustVsPar      -> Sumw2();
0335   hEvtECalVsHCalLeadClustEne -> Sumw2();
0336   // bemc lead cluster event-wise errors
0337   hEvtECalLeadClustNumHit    -> Sumw2();
0338   hEvtECalLeadClustEne       -> Sumw2();
0339   hEvtECalLeadClustDiff      -> Sumw2();
0340   hEvtECalLeadClustVsPar     -> Sumw2();
0341   hEvtECalVsHCalLeadClustEne -> Sumw2();
0342 
0343   // calibration variables
0344   const std::vector<std::string> vecCalibVars = {
0345     "ePar",
0346     "fracParVsLeadBHCal",
0347     "fracParVsLeadBEMC",
0348     "fracParVsSumBHCal",
0349     "fracParVsSumBEMC",
0350     "fracLeadBHCalVsBEMC",
0351     "fracSumBHCalVsBEMC",
0352     "eLeadBHCal",
0353     "eLeadBEMC",
0354     "eSumBHCal",
0355     "eSumBEMC",
0356     "diffLeadBHCal",
0357     "diffLeadBEMC",
0358     "diffSumBHCal",
0359     "diffSumBEMC",
0360     "nHitsLeadBHCal",
0361     "nHitsLeadBEMC",
0362     "nClustBHCal",
0363     "nClustBEMC",
0364     "hLeadBHCal",
0365     "hLeadBEMC",
0366     "fLeadBHCal",
0367     "fLeadBEMC",
0368     "eLeadImage",
0369     "eSumImage",
0370     "eLeadSciFi",
0371     "eSumSciFi",
0372     "nClustImage",
0373     "nClustSciFi",
0374     "hLeadImage",
0375     "hLeadSciFi",
0376     "fLeadImage",
0377     "fLeadSciFi",
0378     "eSumSciFiLayer1",
0379     "eSumSciFiLayer2",
0380     "eSumSciFiLayer3",
0381     "eSumSciFiLayer4",
0382     "eSumSciFiLayer5",
0383     "eSumSciFiLayer6",
0384     "eSumSciFiLayer7",
0385     "eSumSciFiLayer8",
0386     "eSumSciFiLayer9",
0387     "eSumSciFiLayer10",
0388     "eSumSciFiLayer11",
0389     "eSumSciFiLayer12",
0390     "eSumImageLayer1",
0391     "eSumImageLayer2",
0392     "eSumImageLayer3",
0393     "eSumImageLayer4",
0394     "eSumImageLayer5",
0395     "eSumImageLayer6"
0396   };
0397 
0398   // convert to ntuple argument
0399   std::string argCalibVars("");
0400   for (size_t iCalibVar = 0; iCalibVar < vecCalibVars.size(); iCalibVar++) {
0401     argCalibVars.append(vecCalibVars[iCalibVar]);
0402     if ((iCalibVar + 1) != vecCalibVars.size()) {
0403       argCalibVars.append(":");
0404     }
0405   }
0406 
0407   // ntuple for calibration
0408   ntForCalibration = new TNtuple("ntForCalibration", "For Calibration", argCalibVars.c_str());
0409   return;
0410 
0411 }  // end 'InitWithGlobalRootLock()'
0412 
0413 
0414 
0415 
0416 //-------------------------------------------
0417 // ProcessSequential
0418 //-------------------------------------------
0419 void FillBHCalCalibrationTupleProcessor::ProcessSequential(const std::shared_ptr<const JEvent>& event) {
0420 
0421   // clear array for ntuple
0422   for (size_t iCalibVar = 0; iCalibVar < NCalibVars; iCalibVar++) {
0423     varsForCalibration[iCalibVar] = 0.;
0424   }
0425 
0426   // hit and cluster sums
0427   double eHCalHitSum(0.);
0428   double eHCalClustSum(0.);
0429   double eECalClustSum(0.);
0430   double eImageClustSum(0.);
0431   double eSciFiClustSum(0.);
0432   double eTruHCalClustSum(0.);
0433 
0434   // sum bhcal hit energy
0435   for (auto bhCalHit : bhcalRecHits()) {
0436     eHCalHitSum += bhCalHit -> getEnergy();
0437   }  // end 1st bhcal hit loop
0438 
0439   // if hit sum is 0, skip event
0440   const bool isHCalHitSumNonzero = (eHCalHitSum > 0.);
0441   if (!isHCalHitSumNonzero) {
0442     return;
0443   }
0444 
0445   // MC particle properties
0446   float  cMcPar(0.);
0447   double mMcPar(0.);
0448   double fMcPar(0.);
0449   double hMcPar(0.);
0450   double eMcPar(0.);
0451   double pTotMcPar(0.);
0452   double pMcPar[NComp] = {0., 0., 0.};
0453 
0454   // particle loop
0455   unsigned long nPar(0);
0456   for (auto par : genParticles()) {
0457 
0458     // grab particle properties
0459     const auto typePar = par -> getType();
0460     const auto cPar    = par -> getCharge();
0461     const auto mPar    = par -> getMass();
0462     const auto ePar    = par -> getEnergy();
0463     const auto pParX   = par -> getMomentum().x;
0464     const auto pParY   = par -> getMomentum().y;
0465     const auto pParZ   = par -> getMomentum().z;
0466     const auto pPar    = std::sqrt((pParX * pParX) + (pParY * pParY) + (pParZ * pParZ));
0467     const auto fPar    = std::atan(pParY / pParX);
0468     const auto hPar    = std::atanh(pParZ / pPar);
0469 
0470     // select MC particle
0471     const bool isMcParticle = (typePar == 1);
0472     if (isMcParticle) {
0473       cMcPar    = cPar;
0474       mMcPar    = mPar;
0475       fMcPar    = fPar;
0476       hMcPar    = hPar;
0477       eMcPar    = ePar;
0478       pMcPar[0] = pParX;
0479       pMcPar[1] = pParY;
0480       pMcPar[2] = pParZ;
0481       pTotMcPar = pPar;
0482     }
0483     ++nPar;
0484   }  // end particle loop
0485 
0486   // fill particle histograms
0487   hParChrg     -> Fill(cMcPar);
0488   hParMass     -> Fill(mMcPar);
0489   hParPhi      -> Fill(fMcPar);
0490   hParEta      -> Fill(hMcPar);
0491   hParEne      -> Fill(eMcPar);
0492   hParMom      -> Fill(pTotMcPar);
0493   hParMomX     -> Fill(pMcPar[0]);
0494   hParMomY     -> Fill(pMcPar[1]);
0495   hParMomZ     -> Fill(pMcPar[2]);
0496   hParEtaVsPhi -> Fill(fMcPar, hMcPar);
0497 
0498   // reco. bhcal hit loop
0499   unsigned long nHCalHit(0);
0500   for (auto bhCalHit : bhcalRecHits()) {
0501 
0502     // grab hit properties
0503     const auto rHCalHitX   = bhCalHit -> getPosition().x;
0504     const auto rHCalHitY   = bhCalHit -> getPosition().y;
0505     const auto rHCalHitZ   = bhCalHit -> getPosition().z;
0506     const auto eHCalHit    = bhCalHit -> getEnergy();
0507     const auto rHCalHitS   = std::sqrt((rHCalHitX * rHCalHitX) + (rHCalHitY * rHCalHitY));
0508     const auto rHCalHitR   = std::sqrt((rHCalHitS * rHCalHitS) + (rHCalHitZ * rHCalHitZ));
0509     const auto fHCalHit    = boost::math::sign(rHCalHitY) * acos(rHCalHitX / rHCalHitS);
0510     const auto tHCalHit    = std::acos(rHCalHitZ / rHCalHitR);
0511     const auto hHCalHit    = (-1.) * std::log(std::atan(tHCalHit / 2.));
0512     const auto diffHCalHit = (eHCalHit - eMcPar) / eMcPar;
0513 
0514     // fill hit histograms and increment sums/counters
0515     hHCalRecHitPhi      -> Fill(fHCalHit);
0516     hHCalRecHitEta      -> Fill(hHCalHit);
0517     hHCalRecHitEne      -> Fill(eHCalHit);
0518     hHCalRecHitPosZ     -> Fill(rHCalHitZ);
0519     hHCalRecHitParDiff  -> Fill(diffHCalHit);
0520     hHCalRecHitPosYvsX  -> Fill(rHCalHitX, rHCalHitY);
0521     hHCalRecHitEtaVsPhi -> Fill(fHCalHit, hHCalHit);
0522     hHCalRecHitVsParEne -> Fill(eMcPar, eHCalHit);
0523     ++nHCalHit;
0524   }  // end 2nd bhcal hit loop
0525 
0526   // for highest energy bhcal clusters
0527   int    iLeadHCalClust(-1);
0528   int    iLeadTruHCalClust(-1);
0529   int    nHitLeadHCalClust(-1);
0530   int    nHitLeadTruHCalClust(-1);
0531   double hLeadHCalClust(-999.);
0532   double fLeadHCalClust(-999.);
0533   double eLeadHCalClust(-999.);
0534   double eLeadTruHCalClust(-999.);
0535   double diffLeadHCalClust(-999.);
0536   double diffLeadTruHCalClust(-999.);
0537 
0538   // get protoclusters
0539   auto bhCalProtoClusters = event -> Get<edm4eic::ProtoCluster>("HcalBarrelIslandProtoClusters");
0540 
0541   // reco. bhcal cluster loop
0542   unsigned long iHCalClust(0);
0543   unsigned long nHCalProto(0);
0544   unsigned long nHCalClust(0);
0545   for (auto bhCalClust : bhcalClusters()) {
0546 
0547     // grab cluster properties
0548     const auto rHCalClustX   = bhCalClust -> getPosition().x;
0549     const auto rHCalClustY   = bhCalClust -> getPosition().y;
0550     const auto rHCalClustZ   = bhCalClust -> getPosition().z;
0551     const auto eHCalClust    = bhCalClust -> getEnergy();
0552     const auto nHitHCalClust = bhCalClust -> getNhits();
0553     const auto tHCalClust    = bhCalClust -> getIntrinsicTheta();
0554     const auto diffHCalClust = (eHCalClust - eMcPar) / eMcPar;
0555 
0556     // calculate cluster eta, phi
0557     // FIXME update to XYZvectors
0558     const TVector3 vecPosition(rHCalClustX, rHCalClustY, rHCalClustZ);
0559     const auto     hHCalClust = vecPosition.Eta();
0560     const auto     fHCalClust = vecPosition.Phi();
0561     
0562     // loop over protoclusters
0563     unsigned long iHCalProto(0);
0564     unsigned long nProtoHits(0);
0565     for (auto bhCalProto : bhCalProtoClusters) {
0566 
0567       // check if proto index is same as cluster index
0568       // FIXME: this might not be the correct way to match reco and proto clusters
0569       const bool isSameIndex = (iHCalProto == iHCalClust);
0570       if (!isSameIndex) continue;
0571 
0572       // loop over hits
0573       nProtoHits = bhCalProto -> hits_size();
0574       for (uint32_t iProtoHit = 0; iProtoHit < nProtoHits; iProtoHit++) {
0575 
0576         // get hit
0577         const auto bhCalProtoHit = bhCalProto -> getHits(iProtoHit);
0578 
0579         // grab hit properties
0580         const auto rHCalProtoHitX   = bhCalProtoHit.getPosition().x;
0581         const auto rHCalProtoHitY   = bhCalProtoHit.getPosition().y;
0582         const auto rHCalProtoHitZ   = bhCalProtoHit.getPosition().z;
0583         const auto eHCalProtoHit    = bhCalProtoHit.getEnergy();
0584         const auto rHCalProtoHitS   = std::sqrt((rHCalProtoHitX * rHCalProtoHitX) + (rHCalProtoHitY * rHCalProtoHitY));
0585         const auto rHCalProtoHitR   = std::sqrt((rHCalProtoHitS * rHCalProtoHitS) + (rHCalProtoHitZ * rHCalProtoHitZ));
0586         const auto fHCalProtoHit    = boost::math::sign(rHCalProtoHitY) * acos(rHCalProtoHitX / rHCalProtoHitS);
0587         const auto tHCalProtoHit    = std::acos(rHCalProtoHitZ / rHCalProtoHitR);
0588         const auto hHCalProtoHit    = (-1.) * std::log(std::atan(tHCalProtoHit / 2.));
0589         const auto diffHCalProtoHit = (eHCalProtoHit - eMcPar) / eMcPar;
0590 
0591         // fill hit histograms and increment sums/counters
0592         hHCalClustHitPhi      -> Fill(fHCalProtoHit);
0593         hHCalClustHitEta      -> Fill(hHCalProtoHit);
0594         hHCalClustHitEne      -> Fill(eHCalProtoHit);
0595         hHCalClustHitPosZ     -> Fill(rHCalProtoHitZ);
0596         hHCalClustHitParDiff  -> Fill(diffHCalProtoHit);
0597         hHCalClustHitPosYvsX  -> Fill(rHCalProtoHitX, rHCalProtoHitY);
0598         hHCalClustHitEtaVsPhi -> Fill(fHCalProtoHit, hHCalProtoHit);
0599         hHCalClustHitVsParEne -> Fill(eMcPar, eHCalProtoHit);
0600       }
0601       ++nHCalProto;
0602     }  // end protocluster loop
0603 
0604     // fill cluster histograms and increment counters
0605     hHCalClustPhi      -> Fill(fHCalClust);
0606     hHCalClustEta      -> Fill(hHCalClust);
0607     hHCalClustEne      -> Fill(eHCalClust);
0608     hHCalClustPosZ     -> Fill(rHCalClustZ);
0609     hHCalClustNumHit   -> Fill(nHitHCalClust);
0610     hHCalClustParDiff  -> Fill(diffHCalClust);
0611     hHCalClustPosYvsX  -> Fill(rHCalClustX, rHCalClustY);
0612     hHCalClustEtaVsPhi -> Fill(fHCalClust, hHCalClust);
0613     hHCalClustVsParEne -> Fill(eMcPar, eHCalClust);
0614     eHCalClustSum += eHCalClust;
0615     ++nHCalClust;
0616     ++iHCalClust;
0617 
0618     // select leading cluster
0619     const bool isBiggerEne = (eHCalClust > eLeadHCalClust);
0620     if (isBiggerEne) {
0621       iLeadHCalClust    = iHCalClust;
0622       nHitLeadHCalClust = nHitHCalClust;
0623       hLeadHCalClust    = hHCalClust;
0624       fLeadHCalClust    = fHCalClust;
0625       eLeadHCalClust    = eHCalClust;
0626       diffLeadHCalClust = diffHCalClust;
0627     }
0628   }  // end reco. bhcal cluster loop
0629 
0630   // get truth protoclusters
0631   auto bhCalTruProtoClusters = event -> Get<edm4eic::ProtoCluster>("HcalBarrelTruthProtoClusters");
0632 
0633   // true bhcal cluster loop
0634   unsigned long iTruHCalClust(0);
0635   unsigned long nTruHCalProto(0);
0636   unsigned long nTruHCalClust(0);
0637   for (auto truthHCalClust : bhcalTruthClusters()) {
0638 
0639     // grab cluster properties
0640     const auto rTruHCalClustX   = truthHCalClust -> getPosition().x;
0641     const auto rTruHCalClustY   = truthHCalClust -> getPosition().y;
0642     const auto rTruHCalClustZ   = truthHCalClust -> getPosition().z;
0643     const auto eTruHCalClust    = truthHCalClust -> getEnergy();
0644     const auto nHitTruHCalClust = truthHCalClust -> getNhits();
0645     const auto fTruHCalClust    = truthHCalClust -> getIntrinsicPhi();
0646     const auto tTruHCalClust    = truthHCalClust -> getIntrinsicTheta();
0647     const auto hTruHCalClust    = (-1.) * std::log(std::atan(tTruHCalClust / 2.));
0648     const auto diffTruHCalClust = (eTruHCalClust - eMcPar) / eMcPar;
0649     
0650     // loop over protoclusters
0651     unsigned long iTruHCalProto(0);
0652     unsigned long nTruProtoHits(0);
0653     for (auto bhCalTruProto : bhCalTruProtoClusters) {
0654 
0655       // check if truth proto index is same as truth cluster index
0656       // FIXME: this might not be the correct way to match reco and proto clusters
0657       const bool isSameIndex = (iTruHCalProto == iTruHCalClust);
0658       if (!isSameIndex) continue;
0659 
0660       // loop over hits
0661       nTruProtoHits = bhCalTruProto -> hits_size();
0662       for (uint32_t iTruProtoHit = 0; iTruProtoHit < nTruProtoHits; iTruProtoHit++) {
0663 
0664         // get hit
0665         const auto bhCalTruProtoHit = bhCalTruProto -> getHits(iTruProtoHit);
0666 
0667         // grab hit properties
0668         const auto rTruHCalProtoHitX   = bhCalTruProtoHit.getPosition().x;
0669         const auto rTruHCalProtoHitY   = bhCalTruProtoHit.getPosition().y;
0670         const auto rTruHCalProtoHitZ   = bhCalTruProtoHit.getPosition().z;
0671         const auto eTruHCalProtoHit    = bhCalTruProtoHit.getEnergy();
0672         const auto rTruHCalProtoHitS   = std::sqrt((rTruHCalProtoHitX * rTruHCalProtoHitX) + (rTruHCalProtoHitY * rTruHCalProtoHitY));
0673         const auto rTruHCalProtoHitR   = std::sqrt((rTruHCalProtoHitS * rTruHCalProtoHitS) + (rTruHCalProtoHitZ * rTruHCalProtoHitZ));
0674         const auto fTruHCalProtoHit    = boost::math::sign(rTruHCalProtoHitY) * acos(rTruHCalProtoHitX / rTruHCalProtoHitS);
0675         const auto tTruHCalProtoHit    = std::acos(rTruHCalProtoHitZ / rTruHCalProtoHitR);
0676         const auto hTruHCalProtoHit    = (-1.) * std::log(std::atan(tTruHCalProtoHit / 2.));
0677         const auto diffTruHCalProtoHit = (eTruHCalProtoHit - eMcPar) / eMcPar;
0678 
0679         // fill hit histograms and increment sums/counters
0680         hHCalTruClustHitPhi      -> Fill(fTruHCalProtoHit);
0681         hHCalTruClustHitEta      -> Fill(hTruHCalProtoHit);
0682         hHCalTruClustHitEne      -> Fill(eTruHCalProtoHit);
0683         hHCalTruClustHitPosZ     -> Fill(rTruHCalProtoHitZ);
0684         hHCalTruClustHitParDiff  -> Fill(diffTruHCalProtoHit);
0685         hHCalTruClustHitPosYvsX  -> Fill(rTruHCalProtoHitX, rTruHCalProtoHitY);
0686         hHCalTruClustHitEtaVsPhi -> Fill(fTruHCalProtoHit, hTruHCalProtoHit);
0687         hHCalTruClustHitVsParEne -> Fill(eMcPar, eTruHCalProtoHit);
0688       }
0689       ++nTruHCalProto;
0690     }  // end protocluster loop
0691 
0692     // fill cluster histograms and increment counters
0693     hHCalTruClustPhi      -> Fill(fTruHCalClust);
0694     hHCalTruClustEta      -> Fill(hTruHCalClust);
0695     hHCalTruClustEne      -> Fill(eTruHCalClust);
0696     hHCalTruClustPosZ     -> Fill(rTruHCalClustZ);
0697     hHCalTruClustNumHit   -> Fill(nHitTruHCalClust);
0698     hHCalTruClustParDiff  -> Fill(diffTruHCalClust);
0699     hHCalTruClustPosYvsX  -> Fill(rTruHCalClustX, rTruHCalClustY);
0700     hHCalTruClustEtaVsPhi -> Fill(fTruHCalClust, hTruHCalClust);
0701     hHCalTruClustVsParEne -> Fill(eMcPar, eTruHCalClust);
0702     eTruHCalClustSum += eTruHCalClust;
0703     ++nTruHCalClust;
0704 
0705     // select leading cluster
0706     const bool isBiggerEne = (eTruHCalClust > eLeadTruHCalClust);
0707     if (isBiggerEne) {
0708       iLeadTruHCalClust    = iTruHCalClust;
0709       nHitLeadTruHCalClust = nHitTruHCalClust;
0710       eLeadTruHCalClust    = eTruHCalClust;
0711       diffLeadTruHCalClust = diffTruHCalClust;
0712     }
0713     ++iTruHCalClust;
0714   }  // end true bhcal cluster loop
0715 
0716   // for scifi/image hit sums
0717   double eSciFiHitSum(0.);
0718   double eImageHitSum(0.);
0719 
0720   double eSciFiHitSumVsNLayer[CONST::NSciFiLayer];
0721   double eImageHitSumVsNLayer[CONST::NImageLayer];
0722   for (size_t iSciFi = 0; iSciFi < CONST::NSciFiLayer; iSciFi++) {
0723     eSciFiHitSumVsNLayer[iSciFi] = 0.;
0724   }
0725   for (size_t iImage = 0; iImage < CONST::NImageLayer; iImage++) {
0726     eImageHitSumVsNLayer[iImage] = 0.;
0727   }
0728 
0729   // reco. scifi hit loop
0730   unsigned long nSciFiHit(0);
0731   for (auto scifiHit : scifiRecHits()) {
0732 
0733     // grab hit properties
0734     const auto nLayerSciFi  = scifiHit -> getLayer();
0735     const auto iLayerSciFi  = nLayerSciFi - 1;
0736     const auto rSciFiHitX   = scifiHit -> getPosition().x;
0737     const auto rSciFiHitY   = scifiHit -> getPosition().y;
0738     const auto rSciFiHitZ   = scifiHit -> getPosition().z;
0739     const auto eSciFiHit    = scifiHit -> getEnergy();
0740     const auto rSciFiHitS   = std::sqrt((rSciFiHitX * rSciFiHitX) + (rSciFiHitY * rSciFiHitY));
0741     const auto rSciFiHitR   = std::sqrt((rSciFiHitS * rSciFiHitS) + (rSciFiHitZ * rSciFiHitZ));
0742     const auto fSciFiHit    = boost::math::sign(rSciFiHitY) * acos(rSciFiHitX / rSciFiHitS);
0743     const auto tSciFiHit    = std::acos(rSciFiHitZ / rSciFiHitR);
0744     const auto hSciFiHit    = (-1.) * std::log(std::atan(tSciFiHit / 2.));
0745     const auto diffSciFiHit = (eSciFiHit - eMcPar) / eMcPar;
0746 
0747     // fill hit histograms
0748     hSciFiRecHitNLayer      -> Fill(nLayerSciFi);
0749     hSciFiRecHitPhi         -> Fill(fSciFiHit);
0750     hSciFiRecHitEta         -> Fill(hSciFiHit);
0751     hSciFiRecHitEne         -> Fill(eSciFiHit);
0752     hSciFiRecHitPosZ        -> Fill(rSciFiHitZ);
0753     hSciFiRecHitParDiff     -> Fill(diffSciFiHit);
0754     hSciFiRecHitPosYvsX     -> Fill(rSciFiHitX,  rSciFiHitY);
0755     hSciFiRecHitEtaVsPhi    -> Fill(fSciFiHit,   hSciFiHit);
0756     hSciFiRecHitVsParEne    -> Fill(eMcPar,      eSciFiHit);
0757     hSciFiRecHitEneVsNLayer -> Fill(nLayerSciFi, eSciFiHit);
0758 
0759     // increment sums/counters
0760     eSciFiHitSumVsNLayer[iLayerSciFi] += eSciFiHit;
0761     eSciFiHitSum                      += eSciFiHit;
0762     ++nSciFiHit;
0763   }  // end scifi hit loop
0764 
0765   // reco. image hit loop
0766   unsigned long nImageHit(0);
0767   for (auto imageHit : imageRecHits()) {
0768 
0769     // grab hit properties
0770     const auto nLayerImage  = imageHit -> getLayer();
0771     const auto iLayerImage  = nLayerImage - 1;
0772     const auto rImageHitX   = imageHit -> getPosition().x;
0773     const auto rImageHitY   = imageHit -> getPosition().y;
0774     const auto rImageHitZ   = imageHit -> getPosition().z;
0775     const auto eImageHit    = imageHit -> getEnergy();
0776     const auto rImageHitS   = std::sqrt((rImageHitX * rImageHitX) + (rImageHitY * rImageHitY));
0777     const auto rImageHitR   = std::sqrt((rImageHitS * rImageHitS) + (rImageHitZ * rImageHitZ));
0778     const auto fImageHit    = boost::math::sign(rImageHitY) * acos(rImageHitX / rImageHitS);
0779     const auto tImageHit    = std::acos(rImageHitZ / rImageHitR);
0780     const auto hImageHit    = (-1.) * std::log(std::atan(tImageHit / 2.));
0781     const auto diffImageHit = (eImageHit - eMcPar) / eMcPar;
0782 
0783     // fill hit histograms
0784     hImageRecHitNLayer      -> Fill(nLayerImage);
0785     hImageRecHitPhi         -> Fill(fImageHit);
0786     hImageRecHitEta         -> Fill(hImageHit);
0787     hImageRecHitEne         -> Fill(eImageHit);
0788     hImageRecHitPosZ        -> Fill(rImageHitZ);
0789     hImageRecHitParDiff     -> Fill(diffImageHit);
0790     hImageRecHitPosYvsX     -> Fill(rImageHitX,  rImageHitY);
0791     hImageRecHitEtaVsPhi    -> Fill(fImageHit,   hImageHit);
0792     hImageRecHitVsParEne    -> Fill(eMcPar,      eImageHit);
0793     hImageRecHitEneVsNLayer -> Fill(nLayerImage, eImageHit);
0794 
0795     // increment sums/counters
0796     eImageHitSumVsNLayer[iLayerImage] += eImageHit;
0797     eImageHitSum                      += eImageHit;
0798     ++nImageHit;
0799   }  // end scifi hit loop
0800 
0801   // for highest energy bemc clusters
0802   int    iLeadECalClust(-1);
0803   int    iLeadSciFiClust(-1);
0804   int    iLeadImageClust(-1);
0805   int    nHitLeadECalClust(-1);
0806   int    nHitLeadImageClust(-1);
0807   int    nHitLeadSciFiClust(-1);
0808   double hLeadECalClust(-999.);
0809   double hLeadImageClust(-999.);
0810   double hLeadSciFiClust(-999.);
0811   double fLeadECalClust(-999.);
0812   double fLeadImageClust(-999.);
0813   double fLeadSciFiClust(-999.);
0814   double eLeadECalClust(-999.);
0815   double eLeadImageClust(0.);
0816   double eLeadSciFiClust(0.);
0817   double diffLeadECalClust(-999.);
0818   double diffLeadImageClust(-999.);
0819   double diffLeadSciFiClust(-999.);
0820 
0821   // reco. bemc cluster loop
0822   unsigned long iECalClust(0);
0823   unsigned long iImageClust(0);
0824   unsigned long iSciFiClust(0);
0825   unsigned long nECalClust(0);
0826   unsigned long nImageClust(0);
0827   unsigned long nSciFiClust(0);
0828   for (auto bemcClust : bemcClusters()) {
0829 
0830     // grab cluster properties
0831     const auto rECalClustX   = bemcClust -> getPosition().x;
0832     const auto rECalClustY   = bemcClust -> getPosition().y;
0833     const auto rECalClustZ   = bemcClust -> getPosition().z;
0834     const auto eECalClust    = bemcClust -> getEnergy();
0835     const auto nHitECalClust = bemcClust -> getNhits();
0836     const auto tECalClust    = bemcClust -> getIntrinsicTheta();
0837     const auto diffECalClust = (eECalClust - eMcPar) / eMcPar;
0838 
0839     // calculate cluster eta, phi
0840     // FIXME update to XYZvectors
0841     const TVector3 vecPosition(rECalClustX, rECalClustY, rECalClustZ);
0842     const auto     hECalClust = vecPosition.Eta();
0843     const auto     fECalClust = vecPosition.Phi();
0844 
0845     // fill cluster histograms and increment counters
0846     hECalClustPhi      -> Fill(fECalClust);
0847     hECalClustEta      -> Fill(hECalClust);
0848     hECalClustEne      -> Fill(eECalClust);
0849     hECalClustPosZ     -> Fill(rECalClustZ);
0850     hECalClustNumHit   -> Fill(nHitECalClust);
0851     hECalClustParDiff  -> Fill(diffECalClust);
0852     hECalClustPosYvsX  -> Fill(rECalClustX, rECalClustY);
0853     hECalClustEtaVsPhi -> Fill(fECalClust, hECalClust);
0854     hECalClustVsParEne -> Fill(eMcPar, eECalClust);
0855     eECalClustSum += eECalClust;
0856     ++nECalClust;
0857     ++iECalClust;
0858 
0859     // select leading cluster
0860     const bool isBiggerEne = (eECalClust > eLeadECalClust);
0861     if (isBiggerEne) {
0862       iLeadECalClust    = iECalClust;
0863       nHitLeadECalClust = nHitECalClust;
0864       hLeadECalClust    = hECalClust;
0865       fLeadECalClust    = fECalClust;
0866       eLeadECalClust    = eECalClust;
0867       diffLeadECalClust = diffECalClust;
0868     }
0869   }  // end reco. bemc cluster loop
0870 
0871   // loop over scifi clusters
0872   for (auto scifiClust : scifiClusters()) {
0873 
0874     // grab cluster properties
0875     const auto rSciFiClustX   = scifiClust -> getPosition().x;
0876     const auto rSciFiClustY   = scifiClust -> getPosition().y;
0877     const auto rSciFiClustZ   = scifiClust -> getPosition().z;
0878     const auto eSciFiClust    = scifiClust -> getEnergy();
0879     const auto nHitSciFiClust = scifiClust -> getNhits();
0880     const auto tSciFiClust    = scifiClust -> getIntrinsicTheta();
0881     const auto diffSciFiClust = (eSciFiClust - eMcPar) / eMcPar;
0882 
0883     // calculate cluster eta, phi
0884     // FIXME update to XYZvectors
0885     const TVector3 vecPosition(rSciFiClustX, rSciFiClustY, rSciFiClustZ);
0886     const auto     hSciFiClust = vecPosition.Eta();
0887     const auto     fSciFiClust = vecPosition.Phi();
0888 
0889     // increment counters
0890     eSciFiClustSum += eSciFiClust;
0891     ++nSciFiClust;
0892     ++iSciFiClust;
0893 
0894     // select leading cluster
0895     const bool isBiggerEne = (eSciFiClust > eLeadSciFiClust);
0896     if (isBiggerEne) {
0897       iLeadSciFiClust    = iSciFiClust;
0898       nHitLeadSciFiClust = nHitSciFiClust;
0899       hLeadSciFiClust    = hSciFiClust;
0900       fLeadSciFiClust    = fSciFiClust;
0901       eLeadSciFiClust    = eSciFiClust;
0902       diffLeadSciFiClust = diffSciFiClust;
0903     }
0904   }  // end scifi cluster loop
0905 
0906   // loop over imaging clusters
0907   for (auto imageClust : imageClusters()) {
0908 
0909     // grab cluster properties
0910     const auto rImageClustX   = imageClust -> getPosition().x;
0911     const auto rImageClustY   = imageClust -> getPosition().y;
0912     const auto rImageClustZ   = imageClust -> getPosition().z;
0913     const auto eImageClust    = imageClust -> getEnergy();
0914     const auto nHitImageClust = imageClust -> getNhits();
0915     const auto tImageClust    = imageClust -> getIntrinsicTheta();
0916     const auto diffImageClust = (eImageClust - eMcPar) / eMcPar;
0917 
0918     // calculate cluster eta, phi
0919     // FIXME update to XYZvectors
0920     const TVector3 vecPosition(rImageClustX, rImageClustY, rImageClustZ);
0921     const auto     hImageClust = vecPosition.Eta();
0922     const auto     fImageClust = vecPosition.Phi();
0923 
0924     // increment counters
0925     eImageClustSum += eImageClust;
0926     ++nImageClust;
0927     ++iImageClust;
0928 
0929     // select leading cluster
0930     const bool isBiggerEne = (eImageClust > eLeadImageClust);
0931     if (isBiggerEne) {
0932       iLeadImageClust    = iImageClust;
0933       nHitLeadImageClust = nHitImageClust;
0934       hLeadImageClust    = hImageClust;
0935       fLeadImageClust    = fImageClust;
0936       eLeadImageClust    = eImageClust;
0937       diffLeadImageClust = diffImageClust;
0938     }
0939   }  // end imaging cluster loop
0940 
0941   // do event-wise calculations
0942   const auto fracParVsLeadHCal   = eLeadHCalClust / eMcPar;
0943   const auto fracParVsLeadECal   = eLeadECalClust / eMcPar;
0944   const auto fracParVsSumHCal    = eHCalClustSum / eMcPar;
0945   const auto fracParVsSumECal    = eECalClustSum / eMcPar;
0946   const auto fracLeadHCalVsECal  = eLeadECalClust / (eLeadHCalClust + eLeadECalClust);
0947   const auto fracSumHCalVsECal   = eECalClustSum / (eHCalClustSum + eECalClustSum);
0948   const auto diffHCalHitSum      = (eHCalHitSum - eMcPar) / eMcPar;
0949   const auto diffHCalClustSum    = (eHCalClustSum - eMcPar) / eMcPar;
0950   const auto diffECalClustSum    = (eECalClustSum - eMcPar) / eMcPar;
0951   const auto diffTruHCalClustSum = (eTruHCalClustSum - eMcPar) / eMcPar;
0952 
0953   // fill general event-wise bhcal histograms
0954   hEvtHCalNumPar             -> Fill(nPar);
0955   // fill hit event-wise bhcal histograms
0956   hEvtHCalNumHit             -> Fill(nHCalHit);
0957   hEvtHCalSumHitEne          -> Fill(eHCalHitSum);
0958   hEvtHCalSumHitDiff         -> Fill(diffHCalHitSum);
0959   hEvtHCalSumHitVsPar        -> Fill(eMcPar, eHCalHitSum);
0960   // fill cluster event-wise bhcal histograms
0961   hEvtHCalNumClust           -> Fill(nHCalClust);
0962   hEvtHCalSumClustEne        -> Fill(eHCalClustSum);
0963   hEvtHCalSumClustDiff       -> Fill(diffHCalClustSum);
0964   hEvtHCalNumClustVsHit      -> Fill(nHCalHit, nHCalClust);
0965   hEvtHCalSumClustVsPar      -> Fill(eMcPar,   eHCalClustSum);
0966   // fill lead cluster event-wise bhcal histograms
0967   hEvtHCalLeadClustNumHit    -> Fill(nHitLeadHCalClust);
0968   hEvtHCalLeadClustEne       -> Fill(eLeadHCalClust);
0969   hEvtHCalLeadClustDiff      -> Fill(diffLeadHCalClust);
0970   hEvtHCalLeadClustVsPar     -> Fill(eMcPar, eLeadHCalClust);
0971   // fill truth cluster event-wise bhcal histograms
0972   hEvtHCalNumTruClust        -> Fill(nTruHCalClust);
0973   hEvtHCalSumTruClustEne     -> Fill(eTruHCalClustSum);
0974   hEvtHCalSumTruClustDiff    -> Fill(diffTruHCalClustSum);
0975   hEvtHCalNumTruClustVsClust -> Fill(nHCalClust, nTruHCalClust);
0976   hEvtHCalSumTruClustVsPar   -> Fill(eMcPar,     eTruHCalClustSum);
0977   // fill lead truth cluster event-wise bhcal histograms
0978   hEvtHCalLeadTruClustNumHit -> Fill(nHitLeadTruHCalClust);
0979   hEvtHCalLeadTruClustEne    -> Fill(eLeadTruHCalClust);
0980   hEvtHCalLeadTruClustDiff   -> Fill(diffLeadTruHCalClust);
0981   hEvtHCalLeadTruClustVsPar  -> Fill(eMcPar, eLeadTruHCalClust);
0982 
0983   // fill hit event-wise scifi histograms
0984   hEvtSciFiSumEne          -> Fill(eSciFiHitSum);
0985   hEvtSciFiVsHCalHitSumEne -> Fill(eHCalHitSum, eSciFiHitSum);
0986   for (size_t iSciFi = 0; iSciFi < CONST::NSciFiLayer; iSciFi++) {
0987     hEvtSciFiSumEneVsNLayer -> Fill(iSciFi + 1, eSciFiHitSumVsNLayer[iSciFi]);
0988   }
0989 
0990   // fill hit event-wise image histograms
0991   hEvtImageSumEne          -> Fill(eImageHitSum);
0992   hEvtImageVsHCalHitSumEne -> Fill(eHCalHitSum, eImageHitSum);
0993   for (size_t iImage = 0; iImage < CONST::NImageLayer; iImage++) {
0994     hEvtImageSumEneVsNLayer -> Fill(iImage + 1, eImageHitSumVsNLayer[iImage]);
0995   }
0996 
0997   // fill cluster event-wise bhcal histograms
0998   hEvtECalNumClust          -> Fill(nECalClust);
0999   hEvtECalSumClustEne       -> Fill(eECalClustSum);
1000   hEvtECalSumClustDiff      -> Fill(diffECalClustSum);
1001   hEvtECalSumClustVsPar     -> Fill(eMcPar,        eECalClustSum);
1002   hEvtECalVsHCalSumClustEne -> Fill(eHCalClustSum, eECalClustSum);
1003   // fill lead cluster event-wise bhcal histograms
1004   hEvtECalLeadClustNumHit    -> Fill(nHitLeadECalClust);
1005   hEvtECalLeadClustEne       -> Fill(eLeadECalClust);
1006   hEvtECalLeadClustDiff      -> Fill(diffLeadECalClust);
1007   hEvtECalLeadClustVsPar     -> Fill(eMcPar,         eLeadECalClust);
1008   hEvtECalVsHCalLeadClustEne -> Fill(eLeadHCalClust, eLeadECalClust);
1009 
1010   // set variables for calibration tuple
1011   varsForCalibration[0]  = (Float_t) eMcPar;
1012   varsForCalibration[1]  = (Float_t) fracParVsLeadHCal;
1013   varsForCalibration[2]  = (Float_t) fracParVsLeadECal;
1014   varsForCalibration[3]  = (Float_t) fracParVsSumHCal;
1015   varsForCalibration[4]  = (Float_t) fracParVsSumECal;
1016   varsForCalibration[5]  = (Float_t) fracLeadHCalVsECal;
1017   varsForCalibration[6]  = (Float_t) fracSumHCalVsECal;
1018   varsForCalibration[7]  = (Float_t) eLeadHCalClust;
1019   varsForCalibration[8]  = (Float_t) eLeadECalClust;
1020   varsForCalibration[9]  = (Float_t) eHCalClustSum;
1021   varsForCalibration[10] = (Float_t) eECalClustSum;
1022   varsForCalibration[11] = (Float_t) diffLeadHCalClust;
1023   varsForCalibration[12] = (Float_t) diffLeadECalClust;
1024   varsForCalibration[13] = (Float_t) diffHCalClustSum;
1025   varsForCalibration[14] = (Float_t) diffECalClustSum;
1026   varsForCalibration[15] = (Float_t) nHitLeadHCalClust;
1027   varsForCalibration[16] = (Float_t) nHitLeadECalClust;
1028   varsForCalibration[17] = (Float_t) nHCalClust;
1029   varsForCalibration[18] = (Float_t) nECalClust;
1030   varsForCalibration[19] = (Float_t) hLeadHCalClust;
1031   varsForCalibration[20] = (Float_t) hLeadECalClust;
1032   varsForCalibration[21] = (Float_t) fLeadHCalClust;
1033   varsForCalibration[22] = (Float_t) fLeadECalClust;
1034   varsForCalibration[23] = (Float_t) eLeadImageClust;
1035   varsForCalibration[24] = (Float_t) eImageClustSum;
1036   varsForCalibration[25] = (Float_t) eLeadSciFiClust;
1037   varsForCalibration[26] = (Float_t) eSciFiClustSum;
1038   varsForCalibration[27] = (Float_t) nSciFiClust;
1039   varsForCalibration[28] = (Float_t) nImageClust;
1040   varsForCalibration[29] = (Float_t) hLeadImageClust;
1041   varsForCalibration[30] = (Float_t) hLeadSciFiClust;
1042   varsForCalibration[31] = (Float_t) fLeadImageClust;
1043   varsForCalibration[32] = (Float_t) fLeadSciFiClust;
1044   varsForCalibration[33] = (Float_t) eSciFiHitSumVsNLayer[0];
1045   varsForCalibration[34] = (Float_t) eSciFiHitSumVsNLayer[1];
1046   varsForCalibration[35] = (Float_t) eSciFiHitSumVsNLayer[2];
1047   varsForCalibration[36] = (Float_t) eSciFiHitSumVsNLayer[3];
1048   varsForCalibration[37] = (Float_t) eSciFiHitSumVsNLayer[4];
1049   varsForCalibration[38] = (Float_t) eSciFiHitSumVsNLayer[5];
1050   varsForCalibration[39] = (Float_t) eSciFiHitSumVsNLayer[6];
1051   varsForCalibration[40] = (Float_t) eSciFiHitSumVsNLayer[7];
1052   varsForCalibration[41] = (Float_t) eSciFiHitSumVsNLayer[8];
1053   varsForCalibration[42] = (Float_t) eSciFiHitSumVsNLayer[9];
1054   varsForCalibration[43] = (Float_t) eSciFiHitSumVsNLayer[10];
1055   varsForCalibration[44] = (Float_t) eSciFiHitSumVsNLayer[11];
1056   varsForCalibration[45] = (Float_t) eImageHitSumVsNLayer[0];
1057   varsForCalibration[46] = (Float_t) eImageHitSumVsNLayer[1];
1058   varsForCalibration[47] = (Float_t) eImageHitSumVsNLayer[2];
1059   varsForCalibration[48] = (Float_t) eImageHitSumVsNLayer[3];
1060   varsForCalibration[49] = (Float_t) eImageHitSumVsNLayer[4];
1061   varsForCalibration[50] = (Float_t) eImageHitSumVsNLayer[5];
1062 
1063   // fill tuple
1064   ntForCalibration -> Fill(varsForCalibration);
1065   return;
1066 
1067 }  // end 'ProcessSequential(std::shared_ptr<JEvent>&)'
1068 
1069 
1070 
1071 //-------------------------------------------
1072 // FinishWithGlobalRootLock
1073 //-------------------------------------------
1074 void FillBHCalCalibrationTupleProcessor::FinishWithGlobalRootLock() {
1075 
1076   // generic axis titles
1077   const TString sCount("counts");
1078 
1079   // particle axis titles
1080   const TString sMass("m_{par} [GeV/c^{2}]");
1081   const TString sCharge("charge");
1082   const TString sPhiPar("#varphi_{par}");
1083   const TString sEtaPar("#eta_{Par}");
1084   const TString sEnePar("E_{par} [GeV]");
1085   const TString sMomPar("p_{par} [GeV/c]");
1086   const TString sMomParX("p_{x, par} [GeV/c]");
1087   const TString sMomParY("p_{y, par} [GeV/c]");
1088   const TString sMomParZ("p_{z, par} [GeV/c]");
1089   const TString sNumParEvt("N_{par} per event");
1090 
1091   // hit axis titles
1092   const TString sPosHitX("x_{hit} [mm]");
1093   const TString sPosHitY("y_{hit} [mm]");
1094   const TString sPosHitZ("z_{hit} [mm]");
1095   const TString sPhiHit("#varphi_{hit}");
1096   const TString sEtaHit("#eta_{hit}");
1097   const TString sEneHit("e_{hit} [GeV]");
1098   const TString sEneHitSum("E^{sum}_{hit} = #Sigmae_{hit} [GeV]");
1099   const TString sEneHitDiff("#Deltae_{hit} / e_{hit} = (e_{hit} - E_{par}) / e_{hit} [GeV]");
1100   const TString sEneHitSumDiff("#DeltaE^{sum}_{hit} / E^{sum}_{hit} = (E^{sum}_{hit} - E_{par}) / E^{sum}_{hit} [GeV]");
1101   const TString sNumHitEvt("N_{hit} per event");
1102 
1103   // reco. cluster axis titles
1104   const TString sPosClustX("x_{clust} [mm]");
1105   const TString sPosClustY("y_{clust} [mm]");
1106   const TString sPosClustZ("z_{clust} [mm]");
1107   const TString sEneClust("e_{clust} [GeV]");
1108   const TString sPhiClust("#varphi_{clust}");
1109   const TString sEtaClust("#eta_{clust}");
1110   const TString sEneClustSum("E^{sum}_{clust} = #Sigmae_{clust} [GeV]");
1111   const TString sEneClustDiff("#Deltae_{clust} / e_{clust} = (e_{clust} - E_{par}) / e_{clust} [GeV]");
1112   const TString sEneClustLead("E^{lead}_{clust} [GeV]");
1113   const TString sEneClustSumDiff("#DeltaE^{sum}_{clust} / E^{sum}_{clust} = (E^{sum}_{clust} - E_{par}) / E^{sum}_{clust} [GeV]");
1114   const TString sEneClustLeadDiff("#DeltaE^{lead}_{clust} / E^{lead}_{clust} = (E^{lead}_{clust} - E_{par}) / E^{lead}_{clust} [GeV]");
1115   const TString sNumHitClust("N_{hit} per cluster");
1116   const TString sNumClustEvt("N_{clust} per event");
1117 
1118   // truth cluster axis titles
1119   const TString sPosTruClustX("x_{truth clust} [mm]");
1120   const TString sPosTruClustY("y_{truth clust} [mm]");
1121   const TString sPosTruClustZ("z_{truth clust} [mm]");
1122   const TString sPhiTruClust("#varphi^{truth}_{clust}");
1123   const TString sEtaTruClust("#eta^{truth}_{clust}");
1124   const TString sEneTruClust("e^{truth}_{clust} [GeV]");
1125   const TString sEneTruClustDiff("#Deltae^{truth}_{clust} / e^{truth}_{clust} / (e^{truth}_{clust} - E_{par}) / e^{truth}_{clust} [GeV]");
1126   const TString sEneTruClustSum("E^{sum/truth}_{clust} = #Sigmae^{truth}_{clust} [GeV]");
1127   const TString sEneTruClustLead("E^{lead/truth}_{clust} [GeV]");
1128   const TString sEneTruClustSumDiff("#DeltaE^{sum/truth}_{clust} / E^{sum/truth}_{clust} = (E^{sum/truth}_{clust} - E_{par}) / E^{sum/truth}_{clust} [GeV]");
1129   const TString sEneTruClustLeadDiff("#DeltaE^{lead/truth}_{clust} / E^{lead/truth}_{clust} = (E^{lead/truth} _{clust} - E_{par}) / E^{lead/truth}_{clust} [GeV]");
1130   const TString sNumHitTruClust("N_{hit} per truth cluster");
1131   const TString sNumTruClustEvt("N_{truth clust} per event");
1132 
1133   // set particle axis titles
1134   hParChrg                   -> GetXaxis() -> SetTitle(sCharge.Data());
1135   hParChrg                   -> GetYaxis() -> SetTitle(sCount.Data());
1136   hParMass                   -> GetXaxis() -> SetTitle(sMass.Data());
1137   hParMass                   -> GetYaxis() -> SetTitle(sCount.Data());
1138   hParPhi                    -> GetXaxis() -> SetTitle(sPhiPar.Data());
1139   hParPhi                    -> GetYaxis() -> SetTitle(sCount.Data());
1140   hParEta                    -> GetXaxis() -> SetTitle(sEtaPar.Data());
1141   hParEta                    -> GetYaxis() -> SetTitle(sCount.Data());
1142   hParEne                    -> GetXaxis() -> SetTitle(sEnePar.Data());
1143   hParEne                    -> GetYaxis() -> SetTitle(sCount.Data());
1144   hParMom                    -> GetXaxis() -> SetTitle(sMomPar.Data());
1145   hParMom                    -> GetYaxis() -> SetTitle(sCount.Data());
1146   hParMomX                   -> GetXaxis() -> SetTitle(sMomParX.Data());
1147   hParMomX                   -> GetYaxis() -> SetTitle(sCount.Data());
1148   hParMomY                   -> GetXaxis() -> SetTitle(sMomParY.Data());
1149   hParMomY                   -> GetYaxis() -> SetTitle(sCount.Data());
1150   hParMomZ                   -> GetXaxis() -> SetTitle(sMomParZ.Data());
1151   hParMomZ                   -> GetYaxis() -> SetTitle(sCount.Data());
1152   hParEtaVsPhi               -> GetXaxis() -> SetTitle(sPhiPar.Data());
1153   hParEtaVsPhi               -> GetYaxis() -> SetTitle(sEtaPar.Data());
1154   hParEtaVsPhi               -> GetZaxis() -> SetTitle(sCount.Data());
1155   // set reco. hit bhcal axis titles
1156   hHCalRecHitPhi             -> GetXaxis() -> SetTitle(sPhiHit.Data());
1157   hHCalRecHitPhi             -> GetYaxis() -> SetTitle(sCount.Data());
1158   hHCalRecHitEta             -> GetXaxis() -> SetTitle(sEtaHit.Data());
1159   hHCalRecHitEta             -> GetYaxis() -> SetTitle(sCount.Data());
1160   hHCalRecHitEne             -> GetXaxis() -> SetTitle(sEneHit.Data());
1161   hHCalRecHitEne             -> GetYaxis() -> SetTitle(sCount.Data());
1162   hHCalRecHitPosZ            -> GetXaxis() -> SetTitle(sPosHitZ.Data());
1163   hHCalRecHitPosZ            -> GetYaxis() -> SetTitle(sCount.Data());
1164   hHCalRecHitParDiff         -> GetXaxis() -> SetTitle(sEneHitDiff.Data());
1165   hHCalRecHitParDiff         -> GetYaxis() -> SetTitle(sCount.Data());
1166   hHCalRecHitPosYvsX         -> GetXaxis() -> SetTitle(sPosHitX.Data());
1167   hHCalRecHitPosYvsX         -> GetYaxis() -> SetTitle(sPosHitY.Data());
1168   hHCalRecHitPosYvsX         -> GetZaxis() -> SetTitle(sCount.Data());
1169   hHCalRecHitEtaVsPhi        -> GetXaxis() -> SetTitle(sPhiHit.Data());
1170   hHCalRecHitEtaVsPhi        -> GetYaxis() -> SetTitle(sEtaHit.Data());
1171   hHCalRecHitEtaVsPhi        -> GetZaxis() -> SetTitle(sCount.Data());
1172   hHCalRecHitVsParEne        -> GetXaxis() -> SetTitle(sEnePar.Data());
1173   hHCalRecHitVsParEne        -> GetYaxis() -> SetTitle(sEneHit.Data());
1174   hHCalRecHitVsParEne        -> GetZaxis() -> SetTitle(sCount.Data());
1175   // set cluster hit bhcal axis titles
1176   hHCalClustHitPhi           -> GetXaxis() -> SetTitle(sPhiHit.Data());
1177   hHCalClustHitPhi           -> GetYaxis() -> SetTitle(sCount.Data());
1178   hHCalClustHitEta           -> GetXaxis() -> SetTitle(sEtaHit.Data());
1179   hHCalClustHitEta           -> GetYaxis() -> SetTitle(sCount.Data());
1180   hHCalClustHitEne           -> GetXaxis() -> SetTitle(sEneHit.Data());
1181   hHCalClustHitEne           -> GetYaxis() -> SetTitle(sCount.Data());
1182   hHCalClustHitPosZ          -> GetXaxis() -> SetTitle(sPosHitZ.Data());
1183   hHCalClustHitPosZ          -> GetYaxis() -> SetTitle(sCount.Data());
1184   hHCalClustHitParDiff       -> GetXaxis() -> SetTitle(sEneHitDiff.Data());
1185   hHCalClustHitParDiff       -> GetYaxis() -> SetTitle(sCount.Data());
1186   hHCalClustHitPosYvsX       -> GetXaxis() -> SetTitle(sPosHitX.Data());
1187   hHCalClustHitPosYvsX       -> GetYaxis() -> SetTitle(sPosHitY.Data());
1188   hHCalClustHitPosYvsX       -> GetZaxis() -> SetTitle(sCount.Data());
1189   hHCalClustHitEtaVsPhi      -> GetXaxis() -> SetTitle(sPhiHit.Data());
1190   hHCalClustHitEtaVsPhi      -> GetYaxis() -> SetTitle(sEtaHit.Data());
1191   hHCalClustHitEtaVsPhi      -> GetZaxis() -> SetTitle(sCount.Data());
1192   hHCalClustHitVsParEne      -> GetXaxis() -> SetTitle(sEnePar.Data());
1193   hHCalClustHitVsParEne      -> GetYaxis() -> SetTitle(sEneHit.Data());
1194   hHCalClustHitVsParEne      -> GetZaxis() -> SetTitle(sCount.Data());
1195   // set reco. cluster bhcal axis titles
1196   hHCalClustPhi              -> GetXaxis() -> SetTitle(sPhiClust.Data());
1197   hHCalClustPhi              -> GetYaxis() -> SetTitle(sCount.Data());
1198   hHCalClustEta              -> GetXaxis() -> SetTitle(sEtaClust.Data());
1199   hHCalClustEta              -> GetYaxis() -> SetTitle(sCount.Data());
1200   hHCalClustEne              -> GetXaxis() -> SetTitle(sEneClust.Data());
1201   hHCalClustEne              -> GetYaxis() -> SetTitle(sCount.Data());
1202   hHCalClustPosZ             -> GetXaxis() -> SetTitle(sPosClustZ.Data());
1203   hHCalClustPosZ             -> GetYaxis() -> SetTitle(sCount.Data());
1204   hHCalClustNumHit           -> GetXaxis() -> SetTitle(sNumHitClust.Data());
1205   hHCalClustNumHit           -> GetYaxis() -> SetTitle(sCount.Data());
1206   hHCalClustParDiff          -> GetXaxis() -> SetTitle(sEneClustDiff.Data());
1207   hHCalClustParDiff          -> GetYaxis() -> SetTitle(sCount.Data());
1208   hHCalClustPosYvsX          -> GetXaxis() -> SetTitle(sPosClustX.Data());
1209   hHCalClustPosYvsX          -> GetYaxis() -> SetTitle(sPosClustY.Data());
1210   hHCalClustPosYvsX          -> GetZaxis() -> SetTitle(sCount.Data());
1211   hHCalClustEtaVsPhi         -> GetXaxis() -> SetTitle(sPhiClust.Data());
1212   hHCalClustEtaVsPhi         -> GetYaxis() -> SetTitle(sEtaClust.Data());
1213   hHCalClustEtaVsPhi         -> GetZaxis() -> SetTitle(sCount.Data());
1214   hHCalClustVsParEne         -> GetXaxis() -> SetTitle(sEnePar.Data());
1215   hHCalClustVsParEne         -> GetYaxis() -> SetTitle(sEneClust.Data());
1216   hHCalClustVsParEne         -> GetZaxis() -> SetTitle(sCount.Data());
1217   // set truth cluster bhcal axis titles
1218   hHCalTruClustPhi           -> GetXaxis() -> SetTitle(sPhiTruClust.Data());
1219   hHCalTruClustPhi           -> GetYaxis() -> SetTitle(sCount.Data());
1220   hHCalTruClustEta           -> GetXaxis() -> SetTitle(sEtaTruClust.Data());
1221   hHCalTruClustEta           -> GetYaxis() -> SetTitle(sCount.Data());
1222   hHCalTruClustEne           -> GetXaxis() -> SetTitle(sEneTruClust.Data());
1223   hHCalTruClustEne           -> GetYaxis() -> SetTitle(sCount.Data());
1224   hHCalTruClustPosZ          -> GetXaxis() -> SetTitle(sPosTruClustZ.Data());
1225   hHCalTruClustPosZ          -> GetYaxis() -> SetTitle(sCount.Data());
1226   hHCalTruClustNumHit        -> GetXaxis() -> SetTitle(sNumHitTruClust.Data());
1227   hHCalTruClustNumHit        -> GetYaxis() -> SetTitle(sCount.Data());
1228   hHCalTruClustParDiff       -> GetXaxis() -> SetTitle(sEneTruClustDiff.Data());
1229   hHCalTruClustParDiff       -> GetYaxis() -> SetTitle(sCount.Data());
1230   hHCalTruClustPosYvsX       -> GetXaxis() -> SetTitle(sPosTruClustX.Data());
1231   hHCalTruClustPosYvsX       -> GetYaxis() -> SetTitle(sPosTruClustY.Data());
1232   hHCalTruClustPosYvsX       -> GetZaxis() -> SetTitle(sCount.Data());
1233   hHCalTruClustEtaVsPhi      -> GetXaxis() -> SetTitle(sPhiTruClust.Data());
1234   hHCalTruClustEtaVsPhi      -> GetYaxis() -> SetTitle(sEtaTruClust.Data());
1235   hHCalTruClustEtaVsPhi      -> GetZaxis() -> SetTitle(sCount.Data());
1236   hHCalTruClustVsParEne      -> GetXaxis() -> SetTitle(sEnePar.Data());
1237   hHCalTruClustVsParEne      -> GetYaxis() -> SetTitle(sEneTruClust.Data());
1238   hHCalTruClustVsParEne      -> GetZaxis() -> SetTitle(sCount.Data());
1239   // set general event-wise bhcal axis titles
1240   hEvtHCalNumPar             -> GetXaxis() -> SetTitle(sNumParEvt.Data());
1241   hEvtHCalNumPar             -> GetYaxis() -> SetTitle(sCount.Data());
1242   // set hit event-wise bhcal axis titles
1243   hEvtHCalNumHit             -> GetXaxis() -> SetTitle(sNumHitEvt.Data());
1244   hEvtHCalNumHit             -> GetYaxis() -> SetTitle(sCount.Data());
1245   hEvtHCalSumHitEne          -> GetXaxis() -> SetTitle(sEneHitSum.Data());
1246   hEvtHCalSumHitEne          -> GetYaxis() -> SetTitle(sCount.Data());
1247   hEvtHCalSumHitDiff         -> GetXaxis() -> SetTitle(sEneHitSumDiff.Data());
1248   hEvtHCalSumHitDiff         -> GetYaxis() -> SetTitle(sCount.Data());
1249   hEvtHCalSumHitVsPar        -> GetXaxis() -> SetTitle(sEnePar.Data());
1250   hEvtHCalSumHitVsPar        -> GetYaxis() -> SetTitle(sEneHitSum.Data());
1251   hEvtHCalSumHitVsPar        -> GetZaxis() -> SetTitle(sCount.Data());
1252   // set cluster event-wise bhcal axis titles
1253   hEvtHCalNumClust           -> GetXaxis() -> SetTitle(sNumClustEvt.Data());
1254   hEvtHCalNumClust           -> GetYaxis() -> SetTitle(sCount.Data());
1255   hEvtHCalSumClustEne        -> GetXaxis() -> SetTitle(sEneClustSum.Data());
1256   hEvtHCalSumClustEne        -> GetYaxis() -> SetTitle(sCount.Data());
1257   hEvtHCalSumClustDiff       -> GetXaxis() -> SetTitle(sEneClustSumDiff.Data());
1258   hEvtHCalSumClustDiff       -> GetYaxis() -> SetTitle(sCount.Data());
1259   hEvtHCalNumClustVsHit      -> GetXaxis() -> SetTitle(sNumHitEvt.Data());
1260   hEvtHCalNumClustVsHit      -> GetYaxis() -> SetTitle(sNumClustEvt.Data());
1261   hEvtHCalNumClustVsHit      -> GetZaxis() -> SetTitle(sCount.Data());
1262   hEvtHCalSumClustVsPar      -> GetXaxis() -> SetTitle(sEnePar.Data());
1263   hEvtHCalSumClustVsPar      -> GetYaxis() -> SetTitle(sEneClustSum.Data());
1264   hEvtHCalSumClustVsPar      -> GetZaxis() -> SetTitle(sCount.Data());
1265   // set lead cluster event-wise bhcal axis titles
1266   hEvtHCalLeadClustNumHit    -> GetXaxis() -> SetTitle(sNumHitClust.Data());
1267   hEvtHCalLeadClustNumHit    -> GetYaxis() -> SetTitle(sCount.Data());
1268   hEvtHCalLeadClustEne       -> GetXaxis() -> SetTitle(sEneClustLead.Data());
1269   hEvtHCalLeadClustEne       -> GetYaxis() -> SetTitle(sCount.Data());
1270   hEvtHCalLeadClustDiff      -> GetXaxis() -> SetTitle(sEneClustLeadDiff.Data());
1271   hEvtHCalLeadClustDiff      -> GetYaxis() -> SetTitle(sCount.Data());
1272   hEvtHCalLeadClustVsPar     -> GetXaxis() -> SetTitle(sEnePar.Data());
1273   hEvtHCalLeadClustVsPar     -> GetYaxis() -> SetTitle(sEneClustLead.Data());
1274   hEvtHCalLeadClustVsPar     -> GetZaxis() -> SetTitle(sCount.Data());
1275   // set truth cluster event-wise bhcal axis titles
1276   hEvtHCalNumTruClust        -> GetXaxis() -> SetTitle(sNumTruClustEvt.Data());
1277   hEvtHCalNumTruClust        -> GetYaxis() -> SetTitle(sCount.Data());
1278   hEvtHCalSumTruClustEne     -> GetXaxis() -> SetTitle(sEneTruClustSum.Data());
1279   hEvtHCalSumTruClustEne     -> GetYaxis() -> SetTitle(sCount.Data());
1280   hEvtHCalSumTruClustDiff    -> GetXaxis() -> SetTitle(sEneTruClustSumDiff.Data());
1281   hEvtHCalSumTruClustDiff    -> GetYaxis() -> SetTitle(sCount.Data());
1282   hEvtHCalSumTruClustVsPar   -> GetXaxis() -> SetTitle(sEnePar.Data());
1283   hEvtHCalSumTruClustVsPar   -> GetYaxis() -> SetTitle(sEneTruClustSum.Data());
1284   hEvtHCalSumTruClustVsPar   -> GetZaxis() -> SetTitle(sCount.Data());
1285   hEvtHCalNumTruClustVsClust -> GetXaxis() -> SetTitle(sNumClustEvt.Data());
1286   hEvtHCalNumTruClustVsClust -> GetYaxis() -> SetTitle(sNumTruClustEvt.Data());
1287   hEvtHCalNumTruClustVsClust -> GetZaxis() -> SetTitle(sCount.Data());
1288   // set truth lead cluster event-wise bhcal axis titles
1289   hEvtHCalLeadTruClustNumHit -> GetXaxis() -> SetTitle(sNumHitClust.Data());
1290   hEvtHCalLeadTruClustNumHit -> GetYaxis() -> SetTitle(sCount.Data());
1291   hEvtHCalLeadTruClustEne    -> GetXaxis() -> SetTitle(sEneTruClustLead.Data());
1292   hEvtHCalLeadTruClustEne    -> GetYaxis() -> SetTitle(sCount.Data());
1293   hEvtHCalLeadTruClustDiff   -> GetXaxis() -> SetTitle(sEneTruClustLeadDiff.Data());
1294   hEvtHCalLeadTruClustDiff   -> GetYaxis() -> SetTitle(sCount.Data());
1295   hEvtHCalLeadTruClustVsPar  -> GetXaxis() -> SetTitle(sEnePar.Data());
1296   hEvtHCalLeadTruClustVsPar  -> GetYaxis() -> SetTitle(sEneTruClustLead.Data());
1297   hEvtHCalLeadTruClustVsPar  -> GetZaxis() -> SetTitle(sCount.Data());
1298 
1299   // TODO add hit histogram axis titles
1300   // set reco. cluster bemc axis titles
1301   hECalClustPhi           -> GetXaxis() -> SetTitle(sPhiClust.Data());
1302   hECalClustPhi           -> GetYaxis() -> SetTitle(sCount.Data());
1303   hECalClustEta           -> GetXaxis() -> SetTitle(sEtaClust.Data());
1304   hECalClustEta           -> GetYaxis() -> SetTitle(sCount.Data());
1305   hECalClustEne           -> GetXaxis() -> SetTitle(sEneClust.Data());
1306   hECalClustEne           -> GetYaxis() -> SetTitle(sCount.Data());
1307   hECalClustPosZ          -> GetXaxis() -> SetTitle(sPosClustZ.Data());
1308   hECalClustPosZ          -> GetYaxis() -> SetTitle(sCount.Data());
1309   hECalClustNumHit        -> GetXaxis() -> SetTitle(sNumHitClust.Data());
1310   hECalClustNumHit        -> GetYaxis() -> SetTitle(sCount.Data());
1311   hECalClustParDiff       -> GetXaxis() -> SetTitle(sEneClustDiff.Data());
1312   hECalClustParDiff       -> GetYaxis() -> SetTitle(sCount.Data());
1313   hECalClustPosYvsX       -> GetXaxis() -> SetTitle(sPosClustX.Data());
1314   hECalClustPosYvsX       -> GetYaxis() -> SetTitle(sPosClustY.Data());
1315   hECalClustPosYvsX       -> GetZaxis() -> SetTitle(sCount.Data());
1316   hECalClustEtaVsPhi      -> GetXaxis() -> SetTitle(sPhiClust.Data());
1317   hECalClustEtaVsPhi      -> GetYaxis() -> SetTitle(sEtaClust.Data());
1318   hECalClustEtaVsPhi      -> GetZaxis() -> SetTitle(sCount.Data());
1319   hECalClustVsParEne      -> GetXaxis() -> SetTitle(sEnePar.Data());
1320   hECalClustVsParEne      -> GetYaxis() -> SetTitle(sEneClust.Data());
1321   hECalClustVsParEne      -> GetZaxis() -> SetTitle(sCount.Data());
1322   // set cluster event-wise bemc axis titles
1323   hEvtECalNumClust        -> GetXaxis() -> SetTitle(sNumClustEvt.Data());
1324   hEvtECalNumClust        -> GetYaxis() -> SetTitle(sCount.Data());
1325   hEvtECalSumClustEne     -> GetXaxis() -> SetTitle(sEneClustSum.Data());
1326   hEvtECalSumClustEne     -> GetYaxis() -> SetTitle(sCount.Data());
1327   hEvtECalSumClustDiff    -> GetXaxis() -> SetTitle(sEneClustSumDiff.Data());
1328   hEvtECalSumClustDiff    -> GetYaxis() -> SetTitle(sCount.Data());
1329   hEvtECalSumClustVsPar   -> GetXaxis() -> SetTitle(sEnePar.Data());
1330   hEvtECalSumClustVsPar   -> GetYaxis() -> SetTitle(sEneClustSum.Data());
1331   hEvtECalSumClustVsPar   -> GetZaxis() -> SetTitle(sCount.Data());
1332   // set lead cluster event-wise bemc axis titles
1333   hEvtECalLeadClustNumHit -> GetXaxis() -> SetTitle(sNumHitClust.Data());
1334   hEvtECalLeadClustNumHit -> GetYaxis() -> SetTitle(sCount.Data());
1335   hEvtECalLeadClustEne    -> GetXaxis() -> SetTitle(sEneClustLead.Data());
1336   hEvtECalLeadClustEne    -> GetYaxis() -> SetTitle(sCount.Data());
1337   hEvtECalLeadClustDiff   -> GetXaxis() -> SetTitle(sEneClustLeadDiff.Data());
1338   hEvtECalLeadClustDiff   -> GetYaxis() -> SetTitle(sCount.Data());
1339   hEvtECalLeadClustVsPar  -> GetXaxis() -> SetTitle(sEnePar.Data());
1340   hEvtECalLeadClustVsPar  -> GetYaxis() -> SetTitle(sEneClustLead.Data());
1341   hEvtECalLeadClustVsPar  -> GetZaxis() -> SetTitle(sCount.Data());
1342   return;
1343 
1344 }  // end 'FinishWithGlobalRootLock()'
1345 
1346 // end ------------------------------------------------------------------------