Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:59:24

0001 #include "TreeAnalysisRecoXtalsTh.h"
0002 #include <TStyle.h>
0003 #include <TCanvas.h>
0004 
0005 Bool_t FillChain(TChain *chain, const TString &inputFileList);
0006 int main(Int_t argc, Char_t *argv[]) {
0007   if( argc<12 ){
0008     std::cerr << "Please give 10 arguments \n"
0009               << "runList \n"      << "outputFileName \n" 
0010           << "dataType \n"     << "ecalCharIso \n"
0011           << "hcalCharIso \n"  << "ebNeutIso \n"
0012           << "eeNeutIso \n"    << "hhNeutIso \n"
0013           << "nGoodPV \n"      << "L1Seed \n"
0014           << "dRL1Jet \n"      << std::endl;
0015     return -1;
0016   }
0017   
0018   const char *inputFileList = argv[1];
0019   const char *outFileName   = argv[2];
0020   const char *dataType      = argv[3];
0021   const char *ecalCharIso   = argv[4];
0022   const char *hcalCharIso   = argv[5];
0023   const char *ebNeutIso     = argv[6];
0024   const char *eeNeutIso     = argv[7];
0025   const char *hhNeutIso     = argv[8];
0026   const char *nGoodPV       = argv[9];
0027   const char *L1Seed        = argv[10];
0028   const char *dRL1Jet       = argv[11];
0029 
0030   int cut = 0;
0031 
0032   std::cout << "dataType "    << dataType    << std::endl;
0033   std::cout << "ecalCharIso " << ecalCharIso << std::endl;
0034   std::cout << "hcalCharIso " << hcalCharIso << std::endl;
0035   std::cout << "ebNeutIso "   << ebNeutIso   << std::endl;
0036   std::cout << "eeNeutIso "   << eeNeutIso   << std::endl;
0037   std::cout << "hhNeutIso "   << hhNeutIso   << std::endl;
0038   std::cout << "nGoodPV "     << nGoodPV     << std::endl;
0039   std::cout << "L1Seed "      << L1Seed      << std::endl;
0040   std::cout << "dRL1Jet "     << dRL1Jet     << std::endl;
0041 
0042   // Reading Tree                                                        
0043   std::cout << "---------------------" << std::endl;
0044   std::cout << "Reading List of input trees from " << inputFileList << std::endl;
0045   
0046   TChain *chain = new TChain("/isolatedTracksNxN/tree");
0047   if( ! FillChain(chain, inputFileList) ) {
0048     std::cerr << "Cannot get the tree " << std::endl;
0049     return(0);
0050   }
0051   TreeAnalysisRecoXtalsTh tree(chain, outFileName);
0052   
0053   tree.ecalCharIso = ecalCharIso;
0054   tree.hcalCharIso = hcalCharIso;
0055   tree.dataType    = dataType;
0056   tree.ebNeutIso   = atof(ebNeutIso);
0057   tree.eeNeutIso   = atof(eeNeutIso);
0058   tree.hhNeutIso   = atof(hhNeutIso);
0059   tree.GoodPVCut   = atoi(nGoodPV);
0060   tree.L1Seed      = L1Seed;
0061   tree.dRL1Jet     = atof(dRL1Jet);
0062   tree.Loop(cut);
0063 
0064   return 0;
0065 }
0066 
0067 Bool_t FillChain(TChain *chain, const TString &inputFileList) {
0068 
0069   ifstream infile(inputFileList);
0070   std::string buffer;
0071 
0072   if(!infile.is_open()) {
0073     std::cerr << "** ERROR: Can't open '" << inputFileList << "' for input" << std::endl;
0074     return kFALSE;
0075   }
0076   
0077   std::cout << "TreeUtilities : FillChain " << std::endl;
0078   while(1) {
0079     infile >> buffer;
0080     if(!infile.good()) break;
0081     //std::cout << "Adding tree from " << buffer.c_str() << std::endl;
0082     chain->Add(buffer.c_str());
0083   }
0084   std::cout << "No. of Entries in this tree : " << chain->GetEntries() << std::endl;
0085   return kTRUE;
0086 }
0087 
0088 TreeAnalysisRecoXtalsTh::TreeAnalysisRecoXtalsTh(TChain *tree, const char *outFileName) {
0089 
0090   //  double tempgen_TH[22] = { 0.0,  1.0,  2.0,  3.0,  4.0,  
0091   //                5.0,  6.0,  7.0,  8.0,  9.0, 
0092   //                10.0, 12.0, 15.0, 20.0, 25.0, 
0093   //                30.0, 40.0, 50.0, 60.0, 70.0, 80.0, 100};
0094   //{5.0, 6.0, 7.0, 9.0, 11.0, 15.0, 20.0, 30.0}; // from anton
0095 
0096   double tempgen_TH[NPBins+1] = { 0.0,  1.0,  2.0,  3.0,  4.0,  
0097                   5.0,  6.0,  7.0,  9.0, 11.0, 
0098                  15.0, 20.0, 30.0, 50.0, 75.0, 100.0};
0099 
0100 
0101   for (int i=0; i<NPBins+1; i++)  
0102     genPartPBins[i]  = tempgen_TH[i];
0103 
0104   //  double tempgen_Eta[NEtaBins+1] = {0.0, 1.131, 1.653, 2.172};
0105   double tempgen_Eta[NEtaBins+1] = {0.0, 0.2, 0.4, 0.6, 0.8, 
0106                     1.0, 1.2, 1.4, 1.6, 1.8, 
0107                     2.0, 2.2, 2.4};
0108   
0109   for(int i=0; i<NEtaBins+1; i++) genPartEtaBins[i] = tempgen_Eta[i];
0110 
0111   Init(tree);
0112    
0113   BookHistograms(outFileName);
0114 }
0115 
0116 TreeAnalysisRecoXtalsTh::~TreeAnalysisRecoXtalsTh() {
0117   if (!fChain) return;
0118   delete fChain->GetCurrentFile();
0119   
0120   fout->cd();
0121   fout->Write();
0122   fout->Close();   
0123 }
0124 
0125 Int_t TreeAnalysisRecoXtalsTh::Cut(Long64_t entry) {
0126   // This function may be called from Loop.
0127   // returns  1 if entry is accepted.
0128   // returns -1 otherwise.
0129   return 1;
0130 }
0131 
0132 Int_t TreeAnalysisRecoXtalsTh::GetEntry(Long64_t entry) {
0133   
0134   // Read contents of entry.
0135   if (!fChain) return 0;
0136   return fChain->GetEntry(entry);
0137 }
0138 
0139 Long64_t TreeAnalysisRecoXtalsTh::LoadTree(Long64_t entry) {
0140 
0141   // Set the environment to read one entry
0142   if (!fChain) return -5;
0143   Long64_t centry = fChain->LoadTree(entry);
0144   if (centry < 0) return centry;
0145   if (!fChain->InheritsFrom(TChain::Class()))  return centry;
0146   TChain *chain = (TChain*)fChain;
0147   if (chain->GetTreeNumber() != fCurrent) {
0148     fCurrent = chain->GetTreeNumber();
0149     Notify();
0150   }
0151   return centry;
0152 }
0153 
0154 void TreeAnalysisRecoXtalsTh::Init(TChain *tree) {
0155 
0156   // The Init() function is called when the selector needs to initialize
0157   // a new tree or chain. Typically here the branch addresses and branch
0158   // pointers of the tree will be set.
0159   // It is normally not necessary to make changes to the generated
0160   // code, but the routine can be extended by the user if needed.
0161   // Init() will be called many times when running on PROOF
0162   // (once per file to be processed).
0163   
0164 
0165   // Set object pointer
0166   PVx = 0;
0167   PVy = 0;
0168   PVz = 0;
0169   PVisValid = 0;
0170   PVndof = 0;
0171   PVNTracks = 0;
0172   PVNTracksWt = 0;
0173   t_PVTracksSumPt = 0;
0174   t_PVTracksSumPtWt = 0;
0175   t_L1CenJetPt = 0;
0176   t_L1CenJetEta = 0;
0177   t_L1CenJetPhi = 0;
0178   t_L1FwdJetPt = 0;
0179   t_L1FwdJetEta = 0;
0180   t_L1FwdJetPhi = 0;
0181   t_L1TauJetPt = 0;
0182   t_L1TauJetEta = 0;
0183   t_L1TauJetPhi = 0;
0184   t_L1MuonPt = 0;
0185   t_L1MuonEta = 0;
0186   t_L1MuonPhi = 0;
0187   t_L1IsoEMPt = 0;
0188   t_L1IsoEMEta = 0;
0189   t_L1IsoEMPhi = 0;
0190   t_L1NonIsoEMPt = 0;
0191   t_L1NonIsoEMEta = 0;
0192   t_L1NonIsoEMPhi = 0;
0193   t_L1METPt = 0;
0194   t_L1METEta = 0;
0195   t_L1METPhi = 0;
0196   t_jetPt = 0;
0197   t_jetEta = 0;
0198   t_jetPhi = 0;
0199   t_nTrksJetCalo = 0;
0200   t_nTrksJetVtx = 0;
0201   t_trackPAll = 0;
0202   t_trackPhiAll = 0;
0203   t_trackEtaAll = 0;
0204   t_trackPtAll = 0;
0205   t_trackDxyAll = 0;
0206   t_trackDzAll = 0;
0207   t_trackDxyPVAll = 0;
0208   t_trackDzPVAll = 0;
0209   t_trackChiSqAll = 0;
0210   t_trackP = 0;
0211   t_trackPt = 0;
0212   t_trackEta = 0;
0213   t_trackPhi = 0;
0214   t_trackEcalEta = 0;
0215   t_trackEcalPhi = 0;
0216   t_trackHcalEta = 0;
0217   t_trackHcalPhi = 0;
0218   t_trackNOuterHits = 0;
0219   t_NLayersCrossed = 0;
0220   t_trackHitsTOB = 0;
0221   t_trackHitsTEC = 0;
0222   t_trackHitInMissTOB = 0;
0223   t_trackHitInMissTEC = 0;
0224   t_trackHitInMissTIB = 0;
0225   t_trackHitInMissTID = 0;
0226   t_trackHitOutMissTOB = 0;
0227   t_trackHitOutMissTEC = 0;
0228   t_trackHitOutMissTIB = 0;
0229   t_trackHitOutMissTID = 0;
0230   t_trackHitInMeasTOB = 0;
0231   t_trackHitInMeasTEC = 0;
0232   t_trackHitInMeasTIB = 0;
0233   t_trackHitInMeasTID = 0;
0234   t_trackHitOutMeasTOB = 0;
0235   t_trackHitOutMeasTEC = 0;
0236   t_trackHitOutMeasTIB = 0;
0237   t_trackHitOutMeasTID = 0;
0238   t_trackDxy = 0;
0239   t_trackDz = 0;
0240   t_trackDxyPV = 0;
0241   t_trackDzPV = 0;
0242   t_trackChiSq = 0;
0243   t_trackPVIdx = 0;
0244   t_maxNearP31x31 = 0;
0245   t_maxNearP21x21 = 0;
0246   t_ecalSpike11x11 = 0;
0247   t_e7x7 = 0;
0248   t_e9x9 = 0;
0249   t_e11x11 = 0;
0250   t_e15x15 = 0;
0251   t_e7x7_20Sig = 0;
0252   t_e9x9_20Sig = 0;
0253   t_e11x11_20Sig = 0;
0254   t_e15x15_20Sig = 0;
0255   t_maxNearHcalP3x3 = 0;
0256   t_maxNearHcalP5x5 = 0;
0257   t_maxNearHcalP7x7 = 0;
0258   t_h3x3 = 0;
0259   t_h5x5 = 0;
0260   t_h7x7 = 0;
0261   t_infoHcal = 0;
0262 
0263   // Set branch addresses and branch pointers
0264   if (!tree) return;
0265   fChain = tree;
0266   fCurrent = -1;
0267   fChain->SetMakeClass(1);
0268 
0269 
0270   fChain->SetBranchAddress("t_EvtNo", &t_EvtNo, &b_t_EvtNo);
0271   fChain->SetBranchAddress("t_RunNo", &t_RunNo, &b_t_RunNo);
0272   fChain->SetBranchAddress("t_Lumi", &t_Lumi, &b_t_Lumi);
0273   fChain->SetBranchAddress("t_Bunch", &t_Bunch, &b_t_Bunch);
0274   fChain->SetBranchAddress("PVx", &PVx, &b_PVx);
0275   fChain->SetBranchAddress("PVy", &PVy, &b_PVy);
0276   fChain->SetBranchAddress("PVz", &PVz, &b_PVz);
0277   fChain->SetBranchAddress("PVisValid", &PVisValid, &b_PVisValid);
0278   fChain->SetBranchAddress("PVndof", &PVndof, &b_PVndof);
0279   fChain->SetBranchAddress("PVNTracks", &PVNTracks, &b_PVNTracks);
0280   fChain->SetBranchAddress("PVNTracksWt", &PVNTracksWt, &b_PVNTracksWt);
0281   fChain->SetBranchAddress("t_PVTracksSumPt", &t_PVTracksSumPt, &b_t_PVTracksSumPt);
0282   fChain->SetBranchAddress("t_PVTracksSumPtWt", &t_PVTracksSumPtWt, &b_t_PVTracksSumPtWt);
0283   fChain->SetBranchAddress("t_L1Decision", t_L1Decision, &b_t_L1Decision);
0284   fChain->SetBranchAddress("t_L1CenJetPt", &t_L1CenJetPt, &b_t_L1CenJetPt);
0285   fChain->SetBranchAddress("t_L1CenJetEta", &t_L1CenJetEta, &b_t_L1CenJetEta);
0286   fChain->SetBranchAddress("t_L1CenJetPhi", &t_L1CenJetPhi, &b_t_L1CenJetPhi);
0287   fChain->SetBranchAddress("t_L1FwdJetPt", &t_L1FwdJetPt, &b_t_L1FwdJetPt);
0288   fChain->SetBranchAddress("t_L1FwdJetEta", &t_L1FwdJetEta, &b_t_L1FwdJetEta);
0289   fChain->SetBranchAddress("t_L1FwdJetPhi", &t_L1FwdJetPhi, &b_t_L1FwdJetPhi);
0290   fChain->SetBranchAddress("t_L1TauJetPt", &t_L1TauJetPt, &b_t_L1TauJetPt);
0291   fChain->SetBranchAddress("t_L1TauJetEta", &t_L1TauJetEta, &b_t_L1TauJetEta);
0292   fChain->SetBranchAddress("t_L1TauJetPhi", &t_L1TauJetPhi, &b_t_L1TauJetPhi);
0293   fChain->SetBranchAddress("t_L1MuonPt", &t_L1MuonPt, &b_t_L1MuonPt);
0294   fChain->SetBranchAddress("t_L1MuonEta", &t_L1MuonEta, &b_t_L1MuonEta);
0295   fChain->SetBranchAddress("t_L1MuonPhi", &t_L1MuonPhi, &b_t_L1MuonPhi);
0296   fChain->SetBranchAddress("t_L1IsoEMPt", &t_L1IsoEMPt, &b_t_L1IsoEMPt);
0297   fChain->SetBranchAddress("t_L1IsoEMEta", &t_L1IsoEMEta, &b_t_L1IsoEMEta);
0298   fChain->SetBranchAddress("t_L1IsoEMPhi", &t_L1IsoEMPhi, &b_t_L1IsoEMPhi);
0299   fChain->SetBranchAddress("t_L1NonIsoEMPt", &t_L1NonIsoEMPt, &b_t_L1NonIsoEMPt);
0300   fChain->SetBranchAddress("t_L1NonIsoEMEta", &t_L1NonIsoEMEta, &b_t_L1NonIsoEMEta);
0301   fChain->SetBranchAddress("t_L1NonIsoEMPhi", &t_L1NonIsoEMPhi, &b_t_L1NonIsoEMPhi);
0302   fChain->SetBranchAddress("t_L1METPt", &t_L1METPt, &b_t_L1METPt);
0303   fChain->SetBranchAddress("t_L1METEta", &t_L1METEta, &b_t_L1METEta);
0304   fChain->SetBranchAddress("t_L1METPhi", &t_L1METPhi, &b_t_L1METPhi);
0305   fChain->SetBranchAddress("t_jetPt", &t_jetPt, &b_t_jetPt);
0306   fChain->SetBranchAddress("t_jetEta", &t_jetEta, &b_t_jetEta);
0307   fChain->SetBranchAddress("t_jetPhi", &t_jetPhi, &b_t_jetPhi);
0308   fChain->SetBranchAddress("t_nTrksJetCalo", &t_nTrksJetCalo, &b_t_nTrksJetCalo);
0309   fChain->SetBranchAddress("t_nTrksJetVtx", &t_nTrksJetVtx, &b_t_nTrksJetVtx);
0310   fChain->SetBranchAddress("t_trackPAll", &t_trackPAll, &b_t_trackPAll);
0311   fChain->SetBranchAddress("t_trackPhiAll", &t_trackPhiAll, &b_t_trackPhiAll);
0312   fChain->SetBranchAddress("t_trackEtaAll", &t_trackEtaAll, &b_t_trackEtaAll);
0313   fChain->SetBranchAddress("t_trackPtAll", &t_trackPtAll, &b_t_trackPtAll);
0314   fChain->SetBranchAddress("t_trackDxyAll", &t_trackDxyAll, &b_t_trackDxyAll);
0315   fChain->SetBranchAddress("t_trackDzAll", &t_trackDzAll, &b_t_trackDzAll);
0316   fChain->SetBranchAddress("t_trackDxyPVAll", &t_trackDxyPVAll, &b_t_trackDxyPVAll);
0317   fChain->SetBranchAddress("t_trackDzPVAll", &t_trackDzPVAll, &b_t_trackDzPVAll);
0318   fChain->SetBranchAddress("t_trackChiSqAll", &t_trackChiSqAll, &b_t_trackChiSqAll);
0319   fChain->SetBranchAddress("t_trackP", &t_trackP, &b_t_trackP);
0320   fChain->SetBranchAddress("t_trackPt", &t_trackPt, &b_t_trackPt);
0321   fChain->SetBranchAddress("t_trackEta", &t_trackEta, &b_t_trackEta);
0322   fChain->SetBranchAddress("t_trackPhi", &t_trackPhi, &b_t_trackPhi);
0323   fChain->SetBranchAddress("t_trackEcalEta", &t_trackEcalEta, &b_t_trackEcalEta);
0324   fChain->SetBranchAddress("t_trackEcalPhi", &t_trackEcalPhi, &b_t_trackEcalPhi);
0325   fChain->SetBranchAddress("t_trackHcalEta", &t_trackHcalEta, &b_t_trackHcalEta);
0326   fChain->SetBranchAddress("t_trackHcalPhi", &t_trackHcalPhi, &b_t_trackHcalPhi);
0327   fChain->SetBranchAddress("t_trackNOuterHits", &t_trackNOuterHits, &b_t_trackNOuterHits);
0328   fChain->SetBranchAddress("t_NLayersCrossed", &t_NLayersCrossed, &b_t_NLayersCrossed);
0329   fChain->SetBranchAddress("t_trackHitsTOB", &t_trackHitsTOB, &b_t_trackHitsTOB);
0330   fChain->SetBranchAddress("t_trackHitsTEC", &t_trackHitsTEC, &b_t_trackHitsTEC);
0331   fChain->SetBranchAddress("t_trackHitInMissTOB", &t_trackHitInMissTOB, &b_t_trackHitInMissTOB);
0332   fChain->SetBranchAddress("t_trackHitInMissTEC", &t_trackHitInMissTEC, &b_t_trackHitInMissTEC);
0333   fChain->SetBranchAddress("t_trackHitInMissTIB", &t_trackHitInMissTIB, &b_t_trackHitInMissTIB);
0334   fChain->SetBranchAddress("t_trackHitInMissTID", &t_trackHitInMissTID, &b_t_trackHitInMissTID);
0335   fChain->SetBranchAddress("t_trackHitOutMissTOB", &t_trackHitOutMissTOB, &b_t_trackHitOutMissTOB);
0336   fChain->SetBranchAddress("t_trackHitOutMissTEC", &t_trackHitOutMissTEC, &b_t_trackHitOutMissTEC);
0337   fChain->SetBranchAddress("t_trackHitOutMissTIB", &t_trackHitOutMissTIB, &b_t_trackHitOutMissTIB);
0338   fChain->SetBranchAddress("t_trackHitOutMissTID", &t_trackHitOutMissTID, &b_t_trackHitOutMissTID);
0339   fChain->SetBranchAddress("t_trackHitInMeasTOB", &t_trackHitInMeasTOB, &b_t_trackHitInMeasTOB);
0340   fChain->SetBranchAddress("t_trackHitInMeasTEC", &t_trackHitInMeasTEC, &b_t_trackHitInMeasTEC);
0341   fChain->SetBranchAddress("t_trackHitInMeasTIB", &t_trackHitInMeasTIB, &b_t_trackHitInMeasTIB);
0342   fChain->SetBranchAddress("t_trackHitInMeasTID", &t_trackHitInMeasTID, &b_t_trackHitInMeasTID);
0343   fChain->SetBranchAddress("t_trackHitOutMeasTOB", &t_trackHitOutMeasTOB, &b_t_trackHitOutMeasTOB);
0344   fChain->SetBranchAddress("t_trackHitOutMeasTEC", &t_trackHitOutMeasTEC, &b_t_trackHitOutMeasTEC);
0345   fChain->SetBranchAddress("t_trackHitOutMeasTIB", &t_trackHitOutMeasTIB, &b_t_trackHitOutMeasTIB);
0346   fChain->SetBranchAddress("t_trackHitOutMeasTID", &t_trackHitOutMeasTID, &b_t_trackHitOutMeasTID);
0347   fChain->SetBranchAddress("t_trackDxy", &t_trackDxy, &b_t_trackDxy);
0348   fChain->SetBranchAddress("t_trackDz", &t_trackDz, &b_t_trackDz);
0349   fChain->SetBranchAddress("t_trackDxyPV", &t_trackDxyPV, &b_t_trackDxyPV);
0350   fChain->SetBranchAddress("t_trackDzPV", &t_trackDzPV, &b_t_trackDzPV);
0351   fChain->SetBranchAddress("t_trackChiSq", &t_trackChiSq, &b_t_trackChiSq);
0352   fChain->SetBranchAddress("t_trackPVIdx", &t_trackPVIdx, &b_t_trackPVIdx);
0353   fChain->SetBranchAddress("t_maxNearP31x31", &t_maxNearP31x31, &b_t_maxNearP31x31);
0354   fChain->SetBranchAddress("t_maxNearP21x21", &t_maxNearP21x21, &b_t_maxNearP21x21);
0355   fChain->SetBranchAddress("t_ecalSpike11x11", &t_ecalSpike11x11, &b_t_ecalSpike11x11);
0356   fChain->SetBranchAddress("t_e7x7", &t_e7x7, &b_t_e7x7);
0357   fChain->SetBranchAddress("t_e9x9", &t_e9x9, &b_t_e9x9);
0358   fChain->SetBranchAddress("t_e11x11", &t_e11x11, &b_t_e11x11);
0359   fChain->SetBranchAddress("t_e15x15", &t_e15x15, &b_t_e15x15);
0360   fChain->SetBranchAddress("t_e7x7_20Sig", &t_e7x7_20Sig, &b_t_e7x7_20Sig);
0361   fChain->SetBranchAddress("t_e9x9_20Sig", &t_e9x9_20Sig, &b_t_e9x9_20Sig);
0362   fChain->SetBranchAddress("t_e11x11_20Sig", &t_e11x11_20Sig, &b_t_e11x11_20Sig);
0363   fChain->SetBranchAddress("t_e15x15_20Sig", &t_e15x15_20Sig, &b_t_e15x15_20Sig);
0364   fChain->SetBranchAddress("t_maxNearHcalP3x3", &t_maxNearHcalP3x3, &b_t_maxNearHcalP3x3);
0365   fChain->SetBranchAddress("t_maxNearHcalP5x5", &t_maxNearHcalP5x5, &b_t_maxNearHcalP5x5);
0366   fChain->SetBranchAddress("t_maxNearHcalP7x7", &t_maxNearHcalP7x7, &b_t_maxNearHcalP7x7);
0367   fChain->SetBranchAddress("t_h3x3", &t_h3x3, &b_t_h3x3);
0368   fChain->SetBranchAddress("t_h5x5", &t_h5x5, &b_t_h5x5);
0369   fChain->SetBranchAddress("t_h7x7", &t_h7x7, &b_t_h7x7);
0370   fChain->SetBranchAddress("t_infoHcal", &t_infoHcal, &b_t_infoHcal);
0371   fChain->SetBranchAddress("t_nTracks", &t_nTracks, &b_t_nTracks);
0372   Notify();
0373 
0374 }
0375 
0376 void TreeAnalysisRecoXtalsTh::Loop(int cut) {
0377 
0378   if (fChain == 0) return;
0379 
0380   Long64_t nentries = fChain->GetEntriesFast();
0381   std::cout << "No. of Entries in tree " << nentries << std::endl;
0382 
0383   Long64_t nEventsGoodRuns=0, nEventsValidPV=0, nEventsPVTracks=0;
0384 
0385   Long64_t nbytes = 0, nb = 0;
0386   std::map<unsigned int, unsigned int> runEvtList;
0387   std::map<unsigned int, unsigned int> runNTrkList;
0388   std::map<unsigned int, unsigned int> runNIsoTrkList;
0389 
0390   //************Number of goodPV required
0391   int nTrk_trksel_1=0, nTrk_trksel_2=0, nTrk_trksel_3=0, nTrk_trksel_4=0, nTrk_trksel_5=0, nTrk_ecalcharIso=0, nTrk_hcalcharIso=0, nTrk_ecalNeutIso=0, nTrk_hcalNeutIso=0;
0392 
0393   for (Long64_t jentry=0; jentry<nentries;jentry++) {
0394 
0395     // load tree and get current entry
0396     Long64_t ientry = LoadTree(jentry);
0397     if (ientry < 0) break;
0398     nb = fChain->GetEntry(jentry);   nbytes += nb;
0399 
0400     if( !(jentry%1000000) ) {
0401       std::cout << "processing event " << jentry+1 << std::endl;
0402     }
0403 
0404     bool goodRun=false;
0405     goodRun = true;
0406     bool evtSel = true;
0407     if(dataType=="Data" && !goodRun )
0408       evtSel = false;
0409     if( !evtSel ) continue;  //<=================    
0410 
0411     if( runEvtList.find(t_RunNo) != runEvtList.end() ) {
0412       runEvtList[t_RunNo] += 1; 
0413     } else {
0414       runEvtList.insert( std::pair<unsigned int, unsigned int>(t_RunNo,1) );
0415       //      std::cout << "runNo " << t_RunNo <<" "<<runEvtList[t_RunNo]<<std::endl;
0416     }
0417     nEventsGoodRuns++;
0418     
0419     h_NPV_1           ->Fill(PVz->size());
0420 
0421     // I guess the vertex filter works on any PV of the collection
0422     bool anyGoodPV  =false;
0423     bool firstPVGood=false;
0424     int nGoodPV = 0;
0425     int nQltyVtx = 0;
0426 
0427     bool VtxQlty[PVz->size()];
0428 
0429     for(int ipv=0; ipv<PVz->size(); ipv++) {
0430       VtxQlty[ipv] = false;
0431       if (std::abs((*PVz)[ipv])<=25.0 && (*PVndof)[ipv]>4 && 
0432       sqrt((*PVx)[ipv]*(*PVx)[ipv] + (*PVy)[ipv]*(*PVy)[ipv])<=2.0) {
0433     VtxQlty[ipv] = true;
0434     nQltyVtx++;
0435       }
0436       
0437       if ((*PVndof)[ipv]>4) {
0438     anyGoodPV=true;
0439     if(ipv==0) firstPVGood=true;
0440     nGoodPV++;
0441       }
0442     }
0443     if(anyGoodPV)   h_NPV_AnyGoodPV   ->Fill(PVz->size());
0444     if(firstPVGood) h_NPV_FirstGoodPV ->Fill(PVz->size());
0445     
0446     h_nGoodPV->Fill(nGoodPV);
0447     h_nQltyVtx->Fill(nQltyVtx);
0448 
0449     if( firstPVGood ) { 
0450       nEventsValidPV++;
0451       nEventsPVTracks++;
0452       h_NPV_2             ->Fill(PVz->size());
0453       h_PVx_2             ->Fill( (*PVx)[0]             );
0454       h_PVy_2             ->Fill( (*PVy)[0]             );
0455       h_PVr_2             ->Fill( sqrt((*PVx)[0]*(*PVx)[0] + (*PVy)[0]*(*PVy)[0]) );
0456       h_PVz_2             ->Fill( (*PVz)[0]             );
0457     }
0458     
0459     bool pvSel = true;
0460     if(GoodPVCut==4){
0461       if(nGoodPV<4) pvSel = false;
0462     } else if(nGoodPV!=GoodPVCut) pvSel = false;
0463 
0464     if (dataType=="DiPion") pvSel = true;  // no PV selection for DiPion sample
0465     if (!pvSel )  continue;
0466     //================= avoid trigger bias ===================
0467     double maxL1Pt=-1.0, maxL1Eta=999.0, maxL1Phi=999.0;
0468     bool checkL1=false, l1Seed=false;
0469 
0470     if (L1Seed=="L1Jet" || L1Seed=="L1JetL1Tau" || L1Seed=="L1JetL1TauL1EM") {
0471       l1Seed = true;
0472       if (t_L1Decision[15]>0 || t_L1Decision[16]>0 ||  t_L1Decision[17]>0 ||
0473       t_L1Decision[18]>0 || t_L1Decision[19]>0 ||  t_L1Decision[20]>0 ||
0474       t_L1Decision[21]>0 ){ 
0475     if( t_L1CenJetPt->size()>0 && (*t_L1CenJetPt)[0]>maxL1Pt ) {
0476       maxL1Pt=(*t_L1CenJetPt)[0]; 
0477       maxL1Eta=(*t_L1CenJetEta)[0]; 
0478       maxL1Phi=(*t_L1CenJetPhi)[0];
0479     }
0480     if (t_L1FwdJetPt->size()>0 && (*t_L1FwdJetPt)[0]>maxL1Pt ) {
0481       maxL1Pt=(*t_L1FwdJetPt)[0]; 
0482       maxL1Eta=(*t_L1FwdJetEta)[0]; 
0483       maxL1Phi=(*t_L1FwdJetPhi)[0];
0484     }
0485     checkL1=true;
0486       } 
0487     }
0488 
0489     if( L1Seed=="L1Tau" || L1Seed=="L1JetL1Tau" || L1Seed=="L1JetL1TauL1EM") {
0490       l1Seed = true;
0491       if (t_L1Decision[30]>0 ||  t_L1Decision[31]>0  || t_L1Decision[32]>0 || 
0492       t_L1Decision[33]>0 ) {
0493     if (t_L1TauJetPt->size()>0 && (*t_L1TauJetPt)[0]>maxL1Pt ) {
0494       maxL1Pt=(*t_L1TauJetPt)[0]; 
0495       maxL1Eta=(*t_L1TauJetEta)[0]; 
0496       maxL1Phi=(*t_L1TauJetPhi)[0];
0497     }      
0498     checkL1=true;
0499       } 
0500     }
0501 
0502     if( L1Seed=="L1EM" || L1Seed=="L1JetL1TauL1EM") {      
0503       l1Seed = true;
0504       if (t_L1Decision[46]>0 || t_L1Decision[47]>0 || t_L1Decision[48]>0 ||  
0505       t_L1Decision[49]>0 || t_L1Decision[50]>0 || t_L1Decision[51]>0 ||  
0506       t_L1Decision[52]>0      ) {
0507     if(t_L1IsoEMPt->size()>0 && (*t_L1IsoEMPt)[0]>maxL1Pt ) {
0508       maxL1Pt=(*t_L1IsoEMPt)[0]; 
0509       maxL1Eta=(*t_L1IsoEMEta)[0]; 
0510       maxL1Phi=(*t_L1IsoEMPhi)[0];
0511     }      
0512     checkL1=true;
0513       } 
0514     }
0515     if (maxL1Pt<0 ) checkL1=false;
0516     if (!l1Seed)    checkL1=true;
0517 
0518     if (runNTrkList.find(t_RunNo) != runNTrkList.end() ) {
0519       runNTrkList[t_RunNo] += t_trackP->size(); 
0520     } else {
0521       runNTrkList.insert( std::pair<unsigned int, unsigned int>(t_RunNo, t_trackP->size()) );
0522     }
0523     
0524     unsigned int NIsoTrk=0;
0525     for(int itrk=0; itrk<t_trackP->size(); itrk++ ){
0526       
0527       // reject soft tracks
0528       if( (*t_trackPt)[itrk] < 1.0 ) continue;
0529 
0530       double p1              = (*t_trackP)[itrk];
0531       double pt1             = (*t_trackPt)[itrk];
0532       double eta1            = (*t_trackEta)[itrk];
0533       double phi1            = (*t_trackPhi)[itrk];
0534       double etaEcal1        = (*t_trackEcalEta)[itrk];
0535       double phiEcal1        = (*t_trackEcalPhi)[itrk];
0536       double etaHcal1        = (*t_trackHcalEta)[itrk];
0537       double phiHcal1        = (*t_trackHcalPhi)[itrk];
0538       int    trackNOuterHits = (*t_trackNOuterHits)[itrk];
0539       int    NLayersCrossed  = (*t_NLayersCrossed)[itrk];
0540       int    ecalSpike11x11  = (*t_ecalSpike11x11)[itrk];
0541 
0542       double maxNearP31x31   = (*t_maxNearP31x31)[itrk];
0543 
0544       double e7x7            = (*t_e7x7)[itrk]; 
0545       double e9x9            = (*t_e9x9)[itrk]; 
0546       double e11x11          = (*t_e11x11)[itrk]; 
0547       double e15x15          = (*t_e15x15)[itrk];
0548       
0549       double maxNearHcalP3x3 = (*t_maxNearHcalP3x3)[itrk];
0550       double maxNearHcalP5x5 = (*t_maxNearHcalP5x5)[itrk];
0551       double maxNearHcalP7x7 = (*t_maxNearHcalP7x7)[itrk];
0552       
0553       double h3x3            = (*t_h3x3)[itrk];
0554       double h5x5            = (*t_h5x5)[itrk];
0555       double h7x7            = (*t_h7x7)[itrk];
0556       
0557       if ( (*t_infoHcal)[itrk] < 1 ) {
0558     h3x3            = 0.0;
0559     h5x5            = 0.0;
0560     h7x7            = 0.0;
0561       }
0562 
0563       h_trackP_2       ->Fill((*t_trackP)[itrk]    );  
0564       h_trackPt_2      ->Fill((*t_trackPt)[itrk]   );  
0565       h_trackEta_2     ->Fill((*t_trackEta)[itrk]  );  
0566       h_trackPhi_2     ->Fill((*t_trackPhi)[itrk]  );  
0567       h_trackChisq_2   ->Fill((*t_trackChiSq)[itrk]);  
0568       h_trackDxyPV_2   ->Fill((*t_trackDxyPV)[itrk]);  
0569       h_trackDzPV_2    ->Fill((*t_trackDzPV )[itrk]);  
0570       
0571       int iTrkEtaBin=-1, iTrkMomBin=-1;
0572       for(int ieta=0; ieta<NEtaBins; ieta++)   {
0573     if(etaEcal1>genPartEtaBins[ieta] && etaEcal1<genPartEtaBins[ieta+1] ) iTrkEtaBin = ieta;
0574       }
0575       for(int ipt=0;  ipt<NPBins;   ipt++)  {
0576     if( p1>genPartPBins[ipt] &&  p1<genPartPBins[ipt+1] )  iTrkMomBin = ipt;
0577       }
0578       
0579       bool trackChargeIso = true, EcalChargeIso = true, HcalChargeIso = true;
0580       if(ecalCharIso=="maxNearP31X31"     && maxNearP31x31>0)     EcalChargeIso=false;
0581 
0582       if(hcalCharIso=="maxNearHcalP7X7"   && maxNearHcalP7x7>0)   HcalChargeIso=false;
0583       if(hcalCharIso=="maxNearHcalP5x5"   && maxNearHcalP5x5>0)   HcalChargeIso=false;
0584       if(hcalCharIso=="maxNearHcalP3x3"   && maxNearHcalP3x3>0)   HcalChargeIso=false;
0585       if(EcalChargeIso==false || HcalChargeIso==false)trackChargeIso=false;
0586       bool hcalNeutIso = true;
0587       
0588       if( h7x7-h5x5   > hhNeutIso )   hcalNeutIso=false;
0589       
0590       if( iTrkMomBin>=0 && iTrkEtaBin>=0) {  
0591     h_trackPhi_2_2[iTrkEtaBin] ->Fill((*t_trackPhi)[itrk]);
0592       }
0593 
0594       bool trackSel = true;
0595       if( ecalSpike11x11==0 ) trackSel=false ;
0596 
0597       if(dataType=="Data" && (std::abs((*t_trackDxyPV)[itrk])>0.2 || std::abs((*t_trackDzPV)[itrk])>0.2 || (*t_trackChiSq)[itrk]>5.0 ) ) 
0598     trackSel = false;
0599       if(dataType=="MC"   && (std::abs((*t_trackDxyPV)[itrk])>0.2 || std::abs((*t_trackDzPV)[itrk])>0.2 || (*t_trackChiSq)[itrk]>5.0 ) ) 
0600     trackSel = false;
0601       if(dataType=="DiPion" && (std::abs((*t_trackDxy)[itrk])>0.2 || std::abs((*t_trackDz)[itrk])>0.2 || (*t_trackChiSq)[itrk]>5.0 ) ) 
0602     trackSel = false;
0603 
0604       int trackPVid = (*t_trackPVIdx)[itrk];
0605       if(!VtxQlty[trackPVid])trackSel=false;
0606       if( trackSel)nTrk_trksel_1++;
0607       if( iTrkMomBin>=0 && iTrkEtaBin>=0 && trackSel)nTrk_trksel_2++;
0608 
0609       if( trackSel ) {
0610     h_trackP_3       ->Fill((*t_trackP)[itrk]    );  
0611     h_trackPt_3      ->Fill((*t_trackPt)[itrk]   );  
0612     h_trackEta_3     ->Fill((*t_trackEta)[itrk]  );  
0613     h_trackPhi_3     ->Fill((*t_trackPhi)[itrk]  );  
0614     h_trackChisq_3   ->Fill((*t_trackChiSq)[itrk]);  
0615     h_trackDxyPV_3   ->Fill((*t_trackDxyPV)[itrk]);  
0616     h_trackDzPV_3    ->Fill((*t_trackDzPV )[itrk]);  
0617     if( iTrkMomBin>=0 && iTrkEtaBin>=0) {    
0618       h_trackPhi_3_3[iTrkEtaBin] ->Fill((*t_trackPhi)[itrk]);
0619     }
0620       }
0621 
0622       if( iTrkMomBin>=0 && iTrkEtaBin>=0 && trackSel)nTrk_trksel_3++;
0623 
0624       if( NLayersCrossed<7 ) trackSel=false ;      
0625       if( iTrkMomBin>=0 && iTrkEtaBin>=0 && trackSel)nTrk_trksel_4++;
0626 
0627       // reject interactions in tracker
0628       if( (*t_trackHitInMissTOB )[itrk]>0 || (*t_trackHitInMissTEC )[itrk]>0 || 
0629       (*t_trackHitInMissTIB )[itrk]>0 || (*t_trackHitInMissTID )[itrk]>0 || 
0630       (*t_trackHitOutMissTOB)[itrk]>0 || (*t_trackHitOutMissTEC)[itrk]>0 ) {
0631     trackSel = false;
0632       }
0633       if( iTrkMomBin>=0 && iTrkEtaBin>=0 && trackSel)nTrk_trksel_5++;
0634       
0635       if( iTrkMomBin>=0 && iTrkEtaBin>=0) {  
0636     if((*t_trackHitInMissTIB )[itrk]<1 && (*t_trackHitInMissTID )[itrk]<1 )
0637       h_trackPhi_3_Inner[iTrkEtaBin] ->Fill((*t_trackPhi)[itrk]);
0638     if((*t_trackHitOutMissTOB)[itrk]<1 && (*t_trackHitOutMissTEC)[itrk]<1 )
0639       h_trackPhi_3_Outer[iTrkEtaBin] ->Fill((*t_trackPhi)[itrk]);
0640       }
0641       
0642       if( trackSel ) {
0643       h_trackP_4       ->Fill((*t_trackP)[itrk]    );  
0644       h_trackPt_4      ->Fill((*t_trackPt)[itrk]   );  
0645       h_trackEta_4     ->Fill((*t_trackEta)[itrk]  );  
0646       h_trackPhi_4     ->Fill((*t_trackPhi)[itrk]  );  
0647       h_trackChisq_4   ->Fill((*t_trackChiSq)[itrk]);  
0648       h_trackDxyPV_4   ->Fill((*t_trackDxyPV)[itrk]);  
0649       h_trackDzPV_4    ->Fill((*t_trackDzPV )[itrk]);  
0650       }
0651 
0652       double drTrackL1=0.0;
0653       if( checkL1 ){
0654     drTrackL1 = DeltaR(eta1, phi1, maxL1Eta, maxL1Phi);     
0655     if( drTrackL1<dRL1Jet) trackSel=false;
0656       }
0657       
0658       if( iTrkMomBin>=0 && iTrkEtaBin>=0 && trackSel) {  
0659      
0660     h_maxNearP31x31[iTrkMomBin][iTrkEtaBin]->Fill( maxNearP31x31 );
0661 
0662     if (EcalChargeIso) {
0663       nTrk_ecalcharIso++;
0664       if(HcalChargeIso)nTrk_hcalcharIso++;
0665     }
0666 
0667     if (trackChargeIso) {
0668       h_diff_e15x15e11x11      [iTrkMomBin][iTrkEtaBin]->Fill(e15x15-e11x11);
0669       h_diff_e15x15e11x11_20Sig[iTrkMomBin][iTrkEtaBin]->Fill((*t_e15x15_20Sig)[itrk]-(*t_e11x11_20Sig)[itrk]);
0670       h_diff_h7x7h5x5          [iTrkMomBin][iTrkEtaBin]->Fill( h7x7-h5x5 ) ;
0671     }
0672     
0673     if (trackChargeIso) { 
0674 
0675       bool ecalNeutIso = true;
0676       if( std::abs(eta1)<1.47 ) if( e15x15-e11x11 > ebNeutIso ) ecalNeutIso=false;
0677       if( std::abs(eta1)>1.47 ) if( e15x15-e11x11 > eeNeutIso ) ecalNeutIso=false;
0678       if(ecalNeutIso && hcalNeutIso) {
0679           h_trackP_5       ->Fill((*t_trackP)[itrk]    );  
0680           h_trackPt_5      ->Fill((*t_trackPt)[itrk]   );  
0681           h_trackEta_5     ->Fill((*t_trackEta)[itrk]  );  
0682           h_trackPhi_5     ->Fill((*t_trackPhi)[itrk]  );  
0683           h_trackChisq_5   ->Fill((*t_trackChiSq)[itrk]);  
0684           h_trackDxyPV_5   ->Fill((*t_trackDxyPV)[itrk]);  
0685           h_trackDzPV_5    ->Fill((*t_trackDzPV )[itrk]);  
0686           
0687           if( pt1>1.0 && std::abs(eta1)<2.3 ) {
0688         double etot11x11=0, etot9x9=0, etot7x7=0;
0689         if( std::abs(eta1)<1.7 ){
0690           etot11x11=e11x11+h3x3;
0691           etot9x9  =e9x9+h3x3;
0692           etot7x7  =e7x7+h3x3;
0693         } else {
0694           etot11x11=(*t_e11x11_20Sig)[itrk]+h3x3;
0695           etot9x9  =(*t_e9x9_20Sig)[itrk]+h3x3;
0696           etot7x7  =(*t_e7x7_20Sig)[itrk]+h3x3;
0697         }         
0698         h_trackPCaloE11x11H3x3_0->Fill( p1, etot11x11);
0699         h_trackPCaloE9x9H3x3_0  ->Fill( p1, etot9x9);
0700         h_trackPCaloE7x7H3x3_0  ->Fill( p1, etot7x7);
0701         
0702         if(std::abs(eta1)<1.1){
0703           h_trackPCaloE11x11H3x3_1->Fill( p1, etot11x11);
0704           h_trackPCaloE9x9H3x3_1  ->Fill( p1, etot9x9);
0705           h_trackPCaloE7x7H3x3_1  ->Fill( p1, etot7x7);        
0706         } else if(std::abs(eta1)<1.7) {
0707           h_trackPCaloE11x11H3x3_2->Fill( p1, etot11x11);
0708           h_trackPCaloE9x9H3x3_2  ->Fill( p1, etot9x9);
0709           h_trackPCaloE7x7H3x3_2  ->Fill( p1, etot7x7);        
0710         } else if( std::abs(eta1)<2.3 ){
0711           h_trackPCaloE11x11H3x3_3->Fill( p1, etot11x11);
0712           h_trackPCaloE9x9H3x3_3  ->Fill( p1, etot9x9);
0713           h_trackPCaloE7x7H3x3_3  ->Fill( p1, etot7x7);        
0714         }
0715           }
0716           
0717           if( p1>3.0 &&  p1<4.0)  h_eECAL11x11VsHCAL3x3[iTrkEtaBin]->Fill(e11x11, h3x3);
0718           // Ecal tranverse profile
0719           h_eECAL7x7_Frac  [iTrkMomBin][iTrkEtaBin]->Fill(e7x7/p1);
0720           h_eECAL9x9_Frac  [iTrkMomBin][iTrkEtaBin]->Fill(e9x9/p1);
0721           h_eECAL11x11_Frac[iTrkMomBin][iTrkEtaBin]->Fill(e11x11/p1);
0722           h_eECAL15x15_Frac[iTrkMomBin][iTrkEtaBin]->Fill(e15x15/p1);
0723           
0724           hh_eECAL7x7_Frac  [iTrkEtaBin]->Fill(e7x7/p1);
0725           hh_eECAL9x9_Frac  [iTrkEtaBin]->Fill(e9x9/p1);
0726           hh_eECAL11x11_Frac[iTrkEtaBin]->Fill(e11x11/p1);
0727           hh_eECAL15x15_Frac[iTrkEtaBin]->Fill(e15x15/p1);
0728           
0729           // Hcal transverse profile
0730           h_eHCAL3x3_Frac[iTrkMomBin][iTrkEtaBin]->Fill(h3x3/p1);
0731           h_eHCAL5x5_Frac[iTrkMomBin][iTrkEtaBin]->Fill(h5x5/p1);
0732           h_eHCAL7x7_Frac[iTrkMomBin][iTrkEtaBin]->Fill(h7x7/p1);
0733           if( e7x7<0.7) {
0734         h_eHCAL3x3MIP_Frac[iTrkMomBin][iTrkEtaBin]->Fill(h3x3/p1);
0735         h_eHCAL5x5MIP_Frac[iTrkMomBin][iTrkEtaBin]->Fill(h5x5/p1);
0736         h_eHCAL7x7MIP_Frac[iTrkMomBin][iTrkEtaBin]->Fill(h7x7/p1);
0737           }
0738           hh_eHCAL3x3_Frac[iTrkEtaBin]->Fill(h3x3/p1);
0739           hh_eHCAL5x5_Frac[iTrkEtaBin]->Fill(h5x5/p1);
0740           hh_eHCAL7x7_Frac[iTrkEtaBin]->Fill(h7x7/p1);
0741           
0742           // Response : Ecal+Hcal
0743           h_eHCAL3x3_eECAL11x11_response[iTrkMomBin][iTrkEtaBin]->Fill((h3x3+e11x11)/p1);
0744           h_eHCAL5x5_eECAL11x11_response[iTrkMomBin][iTrkEtaBin]->Fill((h5x5+e11x11)/p1);
0745           h_eHCAL7x7_eECAL11x11_response[iTrkMomBin][iTrkEtaBin]->Fill((h7x7+e11x11)/p1);
0746           if( e7x7<0.7) {
0747         h_eHCAL3x3_eECAL11x11_responseMIP[iTrkMomBin][iTrkEtaBin]->Fill((h3x3+e11x11)/p1);
0748         h_eHCAL5x5_eECAL11x11_responseMIP[iTrkMomBin][iTrkEtaBin]->Fill((h5x5+e11x11)/p1);
0749         h_eHCAL7x7_eECAL11x11_responseMIP[iTrkMomBin][iTrkEtaBin]->Fill((h7x7+e11x11)/p1);
0750           } else {
0751         h_eHCAL3x3_eECAL11x11_responseInteract[iTrkMomBin][iTrkEtaBin]->Fill((h3x3+e11x11)/p1);
0752         h_eHCAL5x5_eECAL11x11_responseInteract[iTrkMomBin][iTrkEtaBin]->Fill((h5x5+e11x11)/p1);
0753         h_eHCAL7x7_eECAL11x11_responseInteract[iTrkMomBin][iTrkEtaBin]->Fill((h7x7+e11x11)/p1);
0754           }
0755           
0756           h_eHCAL3x3_eECAL9x9_response[iTrkMomBin][iTrkEtaBin]->Fill((h3x3+e9x9)/p1);
0757           h_eHCAL5x5_eECAL9x9_response[iTrkMomBin][iTrkEtaBin]->Fill((h5x5+e9x9)/p1);
0758           h_eHCAL7x7_eECAL9x9_response[iTrkMomBin][iTrkEtaBin]->Fill((h7x7+e9x9)/p1);
0759           if( e7x7<0.700) {
0760         h_eHCAL3x3_eECAL9x9_responseMIP[iTrkMomBin][iTrkEtaBin]->Fill((h3x3+e9x9)/p1);
0761         h_eHCAL5x5_eECAL9x9_responseMIP[iTrkMomBin][iTrkEtaBin]->Fill((h5x5+e9x9)/p1);
0762         h_eHCAL7x7_eECAL9x9_responseMIP[iTrkMomBin][iTrkEtaBin]->Fill((h7x7+e9x9)/p1);
0763           } else {
0764         h_eHCAL3x3_eECAL9x9_responseInteract[iTrkMomBin][iTrkEtaBin]->Fill((h3x3+e9x9)/p1);
0765         h_eHCAL5x5_eECAL9x9_responseInteract[iTrkMomBin][iTrkEtaBin]->Fill((h5x5+e9x9)/p1);
0766         h_eHCAL7x7_eECAL9x9_responseInteract[iTrkMomBin][iTrkEtaBin]->Fill((h7x7+e9x9)/p1);
0767           }
0768           
0769           h_eHCAL3x3_eECAL7x7_response[iTrkMomBin][iTrkEtaBin]->Fill((h3x3+e7x7)/p1);
0770           h_eHCAL5x5_eECAL7x7_response[iTrkMomBin][iTrkEtaBin]->Fill((h5x5+e7x7)/p1);
0771           h_eHCAL7x7_eECAL7x7_response[iTrkMomBin][iTrkEtaBin]->Fill((h7x7+e7x7)/p1);
0772           if( e7x7<0.700) {
0773         h_eHCAL3x3_eECAL7x7_responseMIP[iTrkMomBin][iTrkEtaBin]->Fill((h3x3+e7x7)/p1);
0774         h_eHCAL5x5_eECAL7x7_responseMIP[iTrkMomBin][iTrkEtaBin]->Fill((h5x5+e7x7)/p1);
0775         h_eHCAL7x7_eECAL7x7_responseMIP[iTrkMomBin][iTrkEtaBin]->Fill((h7x7+e7x7)/p1);
0776           } else {
0777         h_eHCAL3x3_eECAL7x7_responseInteract[iTrkMomBin][iTrkEtaBin]->Fill((h3x3+e7x7)/p1);
0778         h_eHCAL5x5_eECAL7x7_responseInteract[iTrkMomBin][iTrkEtaBin]->Fill((h5x5+e7x7)/p1);
0779         h_eHCAL7x7_eECAL7x7_responseInteract[iTrkMomBin][iTrkEtaBin]->Fill((h7x7+e7x7)/p1);
0780           }     
0781       } // ecal neutral isolation
0782       
0783       bool ecalNeutIso20Sig = true;
0784       if( std::abs(eta1)<1.47 ) if( (*t_e15x15_20Sig)[itrk]-(*t_e11x11_20Sig)[itrk] > ebNeutIso ) ecalNeutIso20Sig=false;
0785       if( std::abs(eta1)>1.47 ) if( (*t_e15x15_20Sig)[itrk]-(*t_e11x11_20Sig)[itrk] > eeNeutIso ) ecalNeutIso20Sig=false;     
0786       if(ecalNeutIso20Sig) {
0787         nTrk_ecalNeutIso++;
0788         if(hcalNeutIso){
0789           nTrk_hcalNeutIso++;
0790           NIsoTrk++;  
0791           h_meanTrackP[iTrkMomBin][iTrkEtaBin]->Fill(p1);
0792           h_eECAL11x11_Frac_20Sig[iTrkMomBin][iTrkEtaBin]->Fill((*t_e11x11_20Sig)[itrk]/p1);
0793           h_eECAL9x9_Frac_20Sig  [iTrkMomBin][iTrkEtaBin]->Fill((*t_e9x9_20Sig)[itrk]/p1);
0794           h_eECAL7x7_Frac_20Sig  [iTrkMomBin][iTrkEtaBin]->Fill((*t_e7x7_20Sig)[itrk]/p1);
0795           
0796           hh_eECAL11x11_Frac_20Sig[iTrkEtaBin]->Fill((*t_e11x11_20Sig)[itrk]/p1);
0797           hh_eECAL9x9_Frac_20Sig  [iTrkEtaBin]->Fill((*t_e9x9_20Sig)[itrk]/p1);
0798           hh_eECAL7x7_Frac_20Sig  [iTrkEtaBin]->Fill((*t_e7x7_20Sig)[itrk]/p1);
0799         
0800           h_eHCAL3x3_Frac_20Sig[iTrkMomBin][iTrkEtaBin]->Fill(h3x3/p1);
0801           h_eHCAL5x5_Frac_20Sig[iTrkMomBin][iTrkEtaBin]->Fill(h5x5/p1);
0802           h_eHCAL7x7_Frac_20Sig[iTrkMomBin][iTrkEtaBin]->Fill(h7x7/p1);
0803           if((*t_e7x7_20Sig)[itrk]<0.7 ) {
0804         h_eHCAL3x3MIP_Frac_20Sig[iTrkMomBin][iTrkEtaBin]->Fill(h3x3/p1);
0805         h_eHCAL5x5MIP_Frac_20Sig[iTrkMomBin][iTrkEtaBin]->Fill(h5x5/p1);
0806         h_eHCAL7x7MIP_Frac_20Sig[iTrkMomBin][iTrkEtaBin]->Fill(h7x7/p1);
0807           }
0808           hh_eHCAL3x3_Frac_20Sig[iTrkEtaBin]->Fill(h3x3/p1);
0809           hh_eHCAL5x5_Frac_20Sig[iTrkEtaBin]->Fill(h5x5/p1);
0810           hh_eHCAL7x7_Frac_20Sig[iTrkEtaBin]->Fill(h7x7/p1);
0811           
0812           h_eHCAL3x3_eECAL11x11_response_20Sig[iTrkMomBin][iTrkEtaBin]->Fill((h3x3+(*t_e11x11_20Sig)[itrk])/p1);
0813           if((*t_e7x7_20Sig)[itrk]<0.7 ) {
0814         h_eHCAL3x3_eECAL11x11_responseMIP_20Sig[iTrkMomBin][iTrkEtaBin]->Fill((h3x3+(*t_e11x11_20Sig)[itrk])/p1);
0815           } else {
0816         h_eHCAL3x3_eECAL11x11_responseInteract_20Sig[iTrkMomBin][iTrkEtaBin]->Fill((h3x3+(*t_e11x11_20Sig)[itrk])/p1);
0817           }
0818           h_eHCAL3x3_eECAL9x9_response_20Sig  [iTrkMomBin][iTrkEtaBin]->Fill((h3x3+(*t_e9x9_20Sig)[itrk])/p1);
0819           if((*t_e7x7_20Sig)[itrk]<0.7 ) {
0820           h_eHCAL3x3_eECAL9x9_responseMIP_20Sig  [iTrkMomBin][iTrkEtaBin]->Fill((h3x3+(*t_e9x9_20Sig)[itrk])/p1);
0821           } else {
0822         h_eHCAL3x3_eECAL9x9_responseInteract_20Sig  [iTrkMomBin][iTrkEtaBin]->Fill((h3x3+(*t_e9x9_20Sig)[itrk])/p1);
0823           }
0824           h_eHCAL3x3_eECAL7x7_response_20Sig  [iTrkMomBin][iTrkEtaBin]->Fill((h3x3+(*t_e7x7_20Sig)[itrk])/p1);
0825           if((*t_e7x7_20Sig)[itrk]<0.7 ) {
0826         h_eHCAL3x3_eECAL7x7_responseMIP_20Sig  [iTrkMomBin][iTrkEtaBin]->Fill((h3x3+(*t_e7x7_20Sig)[itrk])/p1);
0827           } else {
0828         h_eHCAL3x3_eECAL7x7_responseInteract_20Sig  [iTrkMomBin][iTrkEtaBin]->Fill((h3x3+(*t_e7x7_20Sig)[itrk])/p1);
0829           }
0830           
0831           h_eHCAL5x5_eECAL11x11_response_20Sig[iTrkMomBin][iTrkEtaBin]->Fill((h5x5+(*t_e11x11_20Sig)[itrk])/p1);
0832           if((*t_e7x7_20Sig)[itrk]<0.7 ) {
0833         h_eHCAL5x5_eECAL11x11_responseMIP_20Sig[iTrkMomBin][iTrkEtaBin]->Fill((h5x5+(*t_e11x11_20Sig)[itrk])/p1);
0834           } else {
0835         h_eHCAL5x5_eECAL11x11_responseInteract_20Sig[iTrkMomBin][iTrkEtaBin]->Fill((h5x5+(*t_e11x11_20Sig)[itrk])/p1);
0836           }
0837           
0838           h_eHCAL5x5_eECAL7x7_response_20Sig  [iTrkMomBin][iTrkEtaBin]->Fill((h5x5+(*t_e7x7_20Sig)[itrk])/p1);
0839           if((*t_e7x7_20Sig)[itrk]<0.7 ) {
0840         h_eHCAL5x5_eECAL7x7_responseMIP_20Sig  [iTrkMomBin][iTrkEtaBin]->Fill((h5x5+(*t_e7x7_20Sig)[itrk])/p1);
0841           } else {
0842         h_eHCAL5x5_eECAL7x7_responseInteract_20Sig  [iTrkMomBin][iTrkEtaBin]->Fill((h5x5+(*t_e7x7_20Sig)[itrk])/p1);
0843           }
0844         }
0845       }
0846     } // if charged 
0847     
0848       } // momentum and eta bins 
0849       
0850     } // loop over tracks in the event 
0851 
0852     if( runNIsoTrkList.find(t_RunNo) != runNIsoTrkList.end() ) {
0853       runNIsoTrkList[t_RunNo] += NIsoTrk;
0854     } else {
0855       runNIsoTrkList.insert( std::pair<unsigned int, unsigned int>(t_RunNo,NIsoTrk) );
0856     }
0857   } //loop over tree entries
0858 
0859   std::cout << "Number of entries in tree "   << nentries
0860         << "\nnEventsGoodRuns           " << nEventsGoodRuns
0861         << "\nnEventsValidPV            " << nEventsValidPV  
0862         << " " << (double)(nEventsValidPV/nEventsGoodRuns)
0863         << "\nnEventsPVTracks           " << nEventsPVTracks 
0864         << " " << (double)(nEventsPVTracks/nEventsGoodRuns) << std::endl;
0865 
0866   std::cout << "saved runEvtList " << std::endl;
0867   std::map<unsigned int, unsigned int>::iterator runEvtListItr = runEvtList.begin();
0868   for(runEvtListItr=runEvtList.begin(); runEvtListItr != runEvtList.end(); runEvtListItr++) {
0869     std::cout<<runEvtListItr->first << " "<< runEvtListItr->second << std::endl;
0870   }
0871 
0872   std::cout << "Number of tracks in runs " << std::endl;
0873   std::map<unsigned int, unsigned int>::iterator runNTrkListItr = runNTrkList.begin();
0874   for(runNTrkListItr=runNTrkList.begin(); runNTrkListItr != runNTrkList.end(); runNTrkListItr++) {
0875     std::cout<<runNTrkListItr->first << " "<< runNTrkListItr->second << std::endl;
0876   }
0877 
0878   std::cout << "number of tracks_tracksel_1 " << nTrk_trksel_1 << std::endl
0879         << "number of tracks_tracksel_2 (etabin&mombin>0) " << nTrk_trksel_2 << std::endl
0880         << "number of tracks_tracksel_3 (TrackPVId=goodQlty) " << nTrk_trksel_3 << std::endl
0881         << "number of tracks_tracksel_4 (nlayers>6) " << nTrk_trksel_4 << std::endl
0882         << "number of tracks_tracksel_5 (noMissingHit) " << nTrk_trksel_5 << std::endl
0883         << "number of tracks_tracksel_ecalcharIso " << nTrk_ecalcharIso << std::endl
0884         << "number of tracks_tracksel_hcalcharIso " << nTrk_hcalcharIso << std::endl
0885         << "number of tracks_tracksel_ecalNeutIso " << nTrk_ecalNeutIso << std::endl
0886         << "number of tracks_tracksel_hcalNeutIso " << nTrk_hcalNeutIso << std::endl;
0887 
0888   std::cout << "Number of isolated tracks in runs " << std::endl;
0889   std::map<unsigned int, unsigned int>::iterator runNIsoTrkListItr = runNIsoTrkList.begin();
0890   for(runNIsoTrkListItr=runNIsoTrkList.begin(); runNIsoTrkListItr != runNIsoTrkList.end(); runNIsoTrkListItr++) {
0891     std::cout<<runNIsoTrkListItr->first << " "<< runNIsoTrkListItr->second << std::endl;
0892   }
0893   
0894 }
0895 
0896 Bool_t TreeAnalysisRecoXtalsTh::Notify() {
0897 
0898   // The Notify() function is called when a new file is opened. This
0899   // can be either for a new TTree in a TChain or when when a new TTree
0900   // is started when using PROOF. It is normally not necessary to make changes
0901   // to the generated code, but the routine can be extended by the
0902   // user if needed. The return value is currently not used.
0903 
0904   return kTRUE;
0905 }
0906 
0907 void TreeAnalysisRecoXtalsTh::Show(Long64_t entry) {
0908 
0909   // Print contents of entry.
0910   // If entry is not specified, print current entry
0911   if (!fChain) return;
0912   fChain->Show(entry);
0913 }
0914 
0915 double TreeAnalysisRecoXtalsTh::DeltaPhi(double v1, double v2) {
0916 
0917   // Computes the correctly normalized phi difference 
0918   // v1, v2 = phi of object 1 and 2
0919   double pi    = 3.141592654;
0920   double twopi = 6.283185307;
0921 
0922   double diff = std::abs(v2 - v1);
0923   double corr = twopi - diff;
0924   if (diff < pi) { return diff;} 
0925   else           { return corr;}
0926 }
0927 
0928 double TreeAnalysisRecoXtalsTh::DeltaR(double eta1, double phi1, double eta2, double phi2) {
0929   double deta = eta1 - eta2;
0930   double dphi = DeltaPhi(phi1, phi2);
0931   return std::sqrt(deta*deta + dphi*dphi);
0932 }
0933 
0934 void TreeAnalysisRecoXtalsTh::BookHistograms(const char *outFileName) {
0935 
0936   fout = new TFile(outFileName, "RECREATE");
0937 
0938   fout->cd();
0939 
0940   char hname[100], htit[100];
0941 
0942   TDirectory *d_caloTrkPProfile = fout->mkdir( "ProfileTrackPcaloEne" );
0943   d_caloTrkPProfile->cd();
0944   h_trackPCaloE11x11H3x3_0  = new TProfile("h_trackPCaloE11x11H3x3_0","CaloEne(E11x11H3x3) as trackP:|#eta|<2.4",   15, genPartPBins);
0945   h_trackPCaloE9x9H3x3_0    = new TProfile("h_trackPCaloE9x9H3x3_0",  "CaloEne(E9x9H3x3)   as trackP:|#eta|<2.4",   15, genPartPBins);
0946   h_trackPCaloE7x7H3x3_0    = new TProfile("h_trackPCaloE7x7H3x3_0",  "CaloEne(E7x7H3x3)   as trackP:|#eta|<2.4",   15, genPartPBins);
0947 
0948   h_trackPCaloE11x11H3x3_1  = new TProfile("h_trackPCaloE11x11H3x3_1","CaloEne(E11x11H3x3) as trackP:0.0<|#eta|<1.1",   15, genPartPBins);
0949   h_trackPCaloE9x9H3x3_1    = new TProfile("h_trackPCaloE9x9H3x3_1",  "CaloEne(E9x9H3x3)   as trackP:0.0<|#eta|<1.1",   15, genPartPBins);
0950   h_trackPCaloE7x7H3x3_1    = new TProfile("h_trackPCaloE7x7H3x3_1",  "CaloEne(E7x7H3x3)   as trackP:0.0<|#eta|<1.1",   15, genPartPBins);
0951 
0952   h_trackPCaloE11x11H3x3_2  = new TProfile("h_trackPCaloE11x11H3x3_2","CaloEne(E11x11H3x3) as trackP:1.1<|#eta|<1.7",   15, genPartPBins);
0953   h_trackPCaloE9x9H3x3_2    = new TProfile("h_trackPCaloE9x9H3x3_2",  "CaloEne(E9x9H3x3)   as trackP:1.1<|#eta|<1.7",   15, genPartPBins);
0954   h_trackPCaloE7x7H3x3_2    = new TProfile("h_trackPCaloE7x7H3x3_2",  "CaloEne(E7x7H3x3)   as trackP:1.1<|#eta|<1.7",   15, genPartPBins);
0955   
0956   h_trackPCaloE11x11H3x3_3  = new TProfile("h_trackPCaloE11x11H3x3_3","CaloEne(E11x11H3x3) as trackP:1.7<|#eta|<2.0",   15, genPartPBins);
0957   h_trackPCaloE9x9H3x3_3    = new TProfile("h_trackPCaloE9x9H3x3_3",  "CaloEne(E9x9H3x3)   as trackP:1.7<|#eta|<2.0",   15, genPartPBins);
0958   h_trackPCaloE7x7H3x3_3    = new TProfile("h_trackPCaloE7x7H3x3_3",  "CaloEne(E7x7H3x3)   as trackP:1.7<|#eta|<2.0",   15, genPartPBins);
0959 
0960   fout->cd();
0961   h_NPV_AnyGoodPV       = new TH1F("h_NPV_AnyGoodPV",       "h_NPV_AnyGoodPV",        10, -0.5,   9.5); 
0962   h_NPV_FirstGoodPV     = new TH1F("h_NPV_FirstGoodPV",     "h_NPV_FirstGoodPV",      10, -0.5,   9.5); 
0963   h_NPV_1               = new TH1F("h_NPV_1",               "h_NPV_1",                10, -0.5,   9.5); 
0964   h_nGoodPV             = new TH1F("h_nGoodPV",             "h_nGoodPV",              10, -0.5,   9.5); 
0965   h_nQltyVtx            = new TH1F("h_nQltyVtx",            "h_nQltyVtx",             10, -0.5,   9.5); 
0966   h_PVx_1               = new TH1F("h_PVx_1",               "h_PVx_1",                40, -2.0,   2.0); 
0967   h_PVy_1               = new TH1F("h_PVy_1",               "h_PVy_1",                40, -2.0,   2.0); 
0968   h_PVr_1               = new TH1F("h_PVr_1",               "h_PVr_1",               100,  0.0,   3.0); 
0969   h_PVz_1               = new TH1F("h_PVz_1",               "h_PVz_1",               100,-20.0,  20.0); 
0970   h_PVNDOF_1            = new TH1F("h_PVNDOF_1",            "h_PVNDOF_1",            300, -0.5, 299.5); 
0971   h_PVNTracks_1         = new TH1F("h_PVNTracks_1",         "h_PVNTracks_1",         300,  0.0, 300.0); 
0972   h_PVTracksSumPt_1     = new TH1F("h_PVTracksSumPt_1",     "h_PVTracksSumPt_1",     300,  0.0, 300.0);
0973   h_PVNTracksWt_1       = new TH1F("h_PVNTracksWt_1",       "h_PVNTracksWt_1",       300,  0.0, 300.0); 
0974   h_PVTracksSumPtWt_1   = new TH1F("h_PVTracksSumPtWt_1",   "h_PVTracksSumPtWt_1",   300,  0.0, 300.0);
0975   h_PVNTracksHP_1       = new TH1F("h_PVNTracksHP_1",       "h_PVNTracksHP_1",       300,  0.0, 300.0); 
0976   h_PVTracksSumPtHP_1   = new TH1F("h_PVTracksSumPtHP_1",   "h_PVTracksSumPtHP_1",   300,  0.0, 300.0);
0977   h_PVNTracksHPWt_1     = new TH1F("h_PVNTracksHPWt_1",     "h_PVNTracksHPWt_1",     300,  0.0, 300.0); 
0978   h_PVTracksSumPtHPWt_1 = new TH1F("h_PVTracksSumPtHPWt_1", "h_PVTracksSumPtHPWt_1", 300,  0.0, 300.0);
0979   h_NPV_1               ->Sumw2(); 
0980   h_PVx_1               ->Sumw2(); 
0981   h_PVy_1               ->Sumw2(); 
0982   h_PVr_1               ->Sumw2(); 
0983   h_PVz_1               ->Sumw2(); 
0984   h_PVNDOF_1            ->Sumw2();
0985   h_PVNTracks_1         ->Sumw2(); 
0986   h_PVTracksSumPt_1     ->Sumw2();
0987   h_PVNTracksWt_1       ->Sumw2(); 
0988   h_PVTracksSumPtWt_1   ->Sumw2();
0989   h_PVNTracksHP_1       ->Sumw2(); 
0990   h_PVTracksSumPtHP_1   ->Sumw2();
0991   h_PVNTracksHPWt_1     ->Sumw2(); 
0992   h_PVTracksSumPtHPWt_1 ->Sumw2();
0993 
0994   h_NPV_2             = new TH1F("h_NPV_2",             "h_NPV_2",              10, -0.5,   9.5); 
0995   h_PVx_2             = new TH1F("h_PVx_2",             "h_PVx_2",              80, -2.0,   2.0); 
0996   h_PVy_2             = new TH1F("h_PVy_2",             "h_PVy_2",              80, -2.0,   2.0); 
0997   h_PVr_2             = new TH1F("h_PVr_2",             "h_PVr_2",             100,  0.0,   3.0); 
0998   h_PVz_2             = new TH1F("h_PVz_2",             "h_PVz_2",             100,-20.0,  20.0); 
0999   h_PVNDOF_2          = new TH1F("h_PVNDOF_2",          "h_PVNDOF_2",          300, -0.5, 299.5); 
1000   h_PVNTracks_2       = new TH1F("h_PVNTracks_2",       "h_PVNTracks_2",       300,  0.0, 300.0); 
1001   h_PVTracksSumPt_2   = new TH1F("h_PVTracksSumPt_2",   "h_PVTracksSumPt_2",   300,  0.0, 300.0);
1002   h_PVNTracksWt_2     = new TH1F("h_PVNTracksWt_2",     "h_PVNTracksWt_2",     300,  0.0, 300.0); 
1003   h_PVTracksSumPtWt_2 = new TH1F("h_PVTracksSumPtWt_2", "h_PVTracksSumPtWt_2", 300,  0.0, 300.0);
1004   h_PVNTracksHP_2       = new TH1F("h_PVNTracksHP_2",       "h_PVNTracksHP_2",       300,  0.0, 300.0); 
1005   h_PVTracksSumPtHP_2   = new TH1F("h_PVTracksSumPtHP_2",   "h_PVTracksSumPtHP_2",   300,  0.0, 300.0);
1006   h_PVNTracksHPWt_2     = new TH1F("h_PVNTracksHPWt_2",     "h_PVNTracksHPWt_2",     300,  0.0, 300.0); 
1007   h_PVTracksSumPtHPWt_2 = new TH1F("h_PVTracksSumPtHPWt_2", "h_PVTracksSumPtHPWt_2", 300,  0.0, 300.0);
1008   h_NPV_2             ->Sumw2(); 
1009   h_PVx_2             ->Sumw2(); 
1010   h_PVy_2             ->Sumw2(); 
1011   h_PVr_2             ->Sumw2(); 
1012   h_PVz_2             ->Sumw2(); 
1013   h_PVNDOF_2          ->Sumw2();
1014   h_PVNTracks_2       ->Sumw2(); 
1015   h_PVTracksSumPt_2   ->Sumw2();
1016   h_PVNTracksWt_2     ->Sumw2(); 
1017   h_PVTracksSumPtWt_2 ->Sumw2();
1018   h_PVNTracksHP_2       ->Sumw2(); 
1019   h_PVTracksSumPtHP_2   ->Sumw2();
1020   h_PVNTracksHPWt_2     ->Sumw2(); 
1021   h_PVTracksSumPtHPWt_2 ->Sumw2();
1022 
1023   h_PVNTracksSumPt_1 = new TH2F("h_PVNTracksSumPt_1", "h_PVNTracksSumPt_1", 100,  0.0, 100.0, 100,  0.0, 100.0);
1024 
1025   h_trackPAll_1     = new TH1F("trackPAll_1",           "P:   All HighPirity Tracks",   15, genPartPBins);
1026   h_trackPtAll_1    = new TH1F("trackPtAll_1",          "Pt:  All HighPirity Tracks",   15, genPartPBins);
1027   h_trackEtaAll_1   = new TH1F("trackEtaAll_1",         "Eta: All HighPirity Tracks",   100, -3.0, 3.0);
1028   h_trackPhiAll_1   = new TH1F("trackPhiAll_1",         "Phi: All HighPirity Tracks",   100, -3.14159, 3.14159);
1029   h_trackChiSqAll_1 = new TH1F("trackChiSqAll_1",       "Chisq:All HighPirity Tracks",  200,  0.0,20.0);
1030   h_trackDxyAll_1   = new TH1F("h_trackDxyPVAll_1",     "DxyPV:All HighPirity Tracks",  200, -1.0, 1.0);
1031   h_trackDzAll_1    = new TH1F("h_trackDzPVAll_1",      "DzPV: All HighPirity Tracks",  200, -1.0, 1.0);
1032   h_trackPAll_1    ->Sumw2();
1033   h_trackPtAll_1    ->Sumw2();
1034   h_trackEtaAll_1  ->Sumw2();
1035   h_trackPhiAll_1  ->Sumw2();
1036   h_trackChiSqAll_1->Sumw2();  
1037   h_trackDxyAll_1  ->Sumw2();  
1038   h_trackDzAll_1   ->Sumw2();  
1039 
1040   h_trackP_1        = new TH1F("trackP_1",              "P:    good first PV & Iso in Ecal31x31",   15, genPartPBins);
1041   h_trackPt_1       = new TH1F("trackPt_1",             "Pt:   good first PV & Iso in Ecal31x31",   15, genPartPBins);
1042   h_trackEta_1      = new TH1F("trackEta_1",            "Eta:  good first PV & Iso in Ecal31x31",  100, -3.0, 3.0);
1043   h_trackPhi_1      = new TH1F("trackPhi_1",            "Phi:  good first PV & Iso in Ecal31x31",  100, -3.14159, 3.14159);
1044   h_trackChisq_1    = new TH1F("trackChisq_1",          "Chisq:good first PV & Iso in Ecal31x31",  200,  0.0,20.0);
1045   h_trackDxyPV_1    = new TH1F("h_trackDxyPV_1",        "DxyPV:good first PV & Iso in Ecal31x31",  200, -1.0, 1.0);
1046   h_trackDzPV_1     = new TH1F("h_trackDzPV_1",         "DzPV: good first PV & Iso in Ecal31x31",  200, -1.0, 1.0);
1047   h_trackP_1       ->Sumw2();  
1048   h_trackPt_1      ->Sumw2();  
1049   h_trackEta_1     ->Sumw2();  
1050   h_trackPhi_1     ->Sumw2();  
1051   h_trackChisq_1   ->Sumw2();  
1052   h_trackDxyPV_1   ->Sumw2();  
1053   h_trackDzPV_1    ->Sumw2();  
1054 
1055   h_trackP_2        = new TH1F("trackP_2",              "P:    good first PV & Iso in Ecal31x31 & Pt>1.0",   15, genPartPBins);
1056   h_trackPt_2       = new TH1F("trackPt_2",             "Pt:   good first PV & Iso in Ecal31x31 & Pt>1.0",   15, genPartPBins);
1057   h_trackEta_2      = new TH1F("trackEta_2",            "Eta:  good first PV & Iso in Ecal31x31 & Pt>1.0",  100, -3.0, 3.0);
1058   h_trackPhi_2      = new TH1F("trackPhi_2",            "Phi:  good first PV & Iso in Ecal31x31 & Pt>1.0",  100, -3.14159, 3.14159);
1059   h_trackChisq_2    = new TH1F("trackChisq_2",          "Chisq:good first PV & Iso in Ecal31x31 & Pt>1.0",  200,  0.0,20.0);
1060   h_trackDxyPV_2    = new TH1F("h_trackDxyPV_2",        "DxyPV:good first PV & Iso in Ecal31x31 & Pt>1.0",  200, -1.0, 1.0);
1061   h_trackDzPV_2     = new TH1F("h_trackDzPV_2",         "DzPV: good first PV & Iso in Ecal31x31 & Pt>1.0",  200, -1.0, 1.0);
1062   h_trackP_2       ->Sumw2();  
1063   h_trackPt_2      ->Sumw2();  
1064   h_trackEta_2     ->Sumw2();  
1065   h_trackPhi_2     ->Sumw2();  
1066   h_trackChisq_2   ->Sumw2();  
1067   h_trackDxyPV_2   ->Sumw2();  
1068   h_trackDzPV_2    ->Sumw2();  
1069 
1070   h_trackP_3        = new TH1F("trackP_3",              "P:    good first PV & Iso in Ecal31x31 & Pt>1.0 & NoIntInTracker",   15, genPartPBins);
1071   h_trackPt_3       = new TH1F("trackPt_3",             "Pt:   good first PV & Iso in Ecal31x31 & Pt>1.0 & NoIntInTracker",   15, genPartPBins);
1072   h_trackEta_3      = new TH1F("trackEta_3",            "Eta:  good first PV & Iso in Ecal31x31 & Pt>1.0 & NoIntInTracker",  100, -3.0, 3.0);
1073   h_trackPhi_3      = new TH1F("trackPhi_3",            "Phi:  good first PV & Iso in Ecal31x31 & Pt>1.0 & NoIntInTracker",  100, -3.14159, 3.14159);
1074   h_trackChisq_3    = new TH1F("trackChisq_3",          "Chisq:good first PV & Iso in Ecal31x31 & Pt>1.0 & NoIntInTracker",  200,  0.0,20.0);
1075   h_trackDxyPV_3    = new TH1F("h_trackDxyPV_3",        "DxyPV:good first PV & Iso in Ecal31x31 & Pt>1.0 & NoIntInTracker",  200, -1.0, 1.0);
1076   h_trackDzPV_3     = new TH1F("h_trackDzPV_3",         "DzPV: good first PV & Iso in Ecal31x31 & Pt>1.0 & NoIntInTracker",  200, -1.0, 1.0);
1077   h_trackP_3       ->Sumw2();  
1078   h_trackPt_3      ->Sumw2();  
1079   h_trackEta_3     ->Sumw2();  
1080   h_trackPhi_3     ->Sumw2();  
1081   h_trackChisq_3   ->Sumw2();  
1082   h_trackDxyPV_3   ->Sumw2();  
1083   h_trackDzPV_3    ->Sumw2();  
1084 
1085   h_trackP_4        = new TH1F("trackP_4",              "P:    good first PV & Iso in Ecal31x31 & Pt>1.0 & NoIntInTracker & |dxy|<0.2,|dz|<0.2,chisq<5.0",   15, genPartPBins);
1086   h_trackPt_4       = new TH1F("trackPt_4",             "Pt:   good first PV & Iso in Ecal31x31 & Pt>1.0 & NoIntInTracker & |dxy|<0.2,|dz|<0.2,chisq<5.0",   15, genPartPBins);
1087   h_trackEta_4      = new TH1F("trackEta_4",            "Eta:  good first PV & Iso in Ecal31x31 & Pt>1.0 & NoIntInTracker & |dxy|<0.2,|dz|<0.2,chisq<5.0",  100, -3.0, 3.0);
1088   h_trackPhi_4      = new TH1F("trackPhi_4",            "Phi:  good first PV & Iso in Ecal31x31 & Pt>1.0 & NoIntInTracker & |dxy|<0.2,|dz|<0.2,chisq<5.0",  100, -3.14159, 3.14159);
1089   h_trackChisq_4    = new TH1F("trackChisq_4",          "Chisq:good first PV & Iso in Ecal31x31 & Pt>1.0 & NoIntInTracker & |dxy|<0.2,|dz|<0.2,chisq<5.0",  200,  0.0,20.0);
1090   h_trackDxyPV_4    = new TH1F("h_trackDxyPV_4",        "DxyPV:good first PV & Iso in Ecal31x31 & Pt>1.0 & NoIntInTracker & |dxy|<0.2,|dz|<0.2,chisq<5.0",  200, -1.0, 1.0);
1091   h_trackDzPV_4     = new TH1F("h_trackDzPV_4",         "DzPV: good first PV & Iso in Ecal31x31 & Pt>1.0 & NoIntInTracker & |dxy|<0.2,|dz|<0.2,chisq<5.0",  200, -1.0, 1.0);
1092   h_trackP_4       ->Sumw2();  
1093   h_trackPt_4      ->Sumw2();  
1094   h_trackEta_4     ->Sumw2();  
1095   h_trackPhi_4     ->Sumw2();  
1096   h_trackChisq_4   ->Sumw2();  
1097   h_trackDxyPV_4   ->Sumw2();  
1098   h_trackDzPV_4    ->Sumw2();  
1099 
1100   h_trackP_5        = new TH1F("trackP_5",              "P:    Fully Isolated",   15, genPartPBins);
1101   h_trackPt_5       = new TH1F("trackPt_5",             "Pt:   Fully Isolated",   15, genPartPBins);
1102   h_trackEta_5      = new TH1F("trackEta_5",            "Eta:  Fully Isolated",  100, -3.0, 3.0);
1103   h_trackPhi_5      = new TH1F("trackPhi_5",            "Phi:  Fully Isolated",  100, -3.14159, 3.14159);
1104   h_trackChisq_5    = new TH1F("trackChisq_5",          "Chisq:Fully Isolated",  200,  0.0,20.0);
1105   h_trackDxyPV_5    = new TH1F("h_trackDxyPV_5",        "DxyPV:Fully Isolated",  200, -1.0, 1.0);
1106   h_trackDzPV_5     = new TH1F("h_trackDzPV_5",         "DzPV: Fully Isolated",  200, -1.0, 1.0);
1107   h_trackP_5       ->Sumw2();  
1108   h_trackPt_5      ->Sumw2();  
1109   h_trackEta_5     ->Sumw2();  
1110   h_trackPhi_5     ->Sumw2();  
1111   h_trackChisq_5   ->Sumw2();  
1112   h_trackDxyPV_5   ->Sumw2();  
1113   h_trackDzPV_5    ->Sumw2();  
1114 
1115 
1116   for(int ieta=0; ieta<NEtaBins; ieta++) {
1117     double lowEta=-5.0, highEta= 5.0;
1118     lowEta  = genPartEtaBins[ieta];
1119     highEta = genPartEtaBins[ieta+1];
1120     
1121      sprintf(hname, "h_trackPhi_2_2_etaBin%i",ieta);
1122      sprintf(htit,  "h_trackPhi(pt>1.0, iso in 31x31) : %3.2f<|#eta|<%3.2f)", lowEta, highEta);
1123     h_trackPhi_2_2[ieta]= new TH1F(hname, htit,  100, -3.14159, 3.14159);
1124      sprintf(hname, "h_trackPhi_3_3_etaBin%i",ieta);
1125      sprintf(htit,  "h_trackPhi(pt>1.0, iso in 31x31,no missing hit) : %3.2f<|#eta|<%3.2f)", lowEta, highEta);
1126     h_trackPhi_3_3[ieta]= new TH1F(hname, htit,  100, -3.14159, 3.14159);
1127 
1128      sprintf(hname, "h_trackPhi_3_Inner_etaBin%i",ieta);
1129      sprintf(htit,  "h_trackPhi(pt>1.0, iso in 31x31,no missing Inner hit) : %3.2f<|#eta|<%3.2f)", lowEta, highEta);
1130     h_trackPhi_3_Inner[ieta]= new TH1F(hname, htit,  100, -3.14159, 3.14159);
1131 
1132      sprintf(hname, "h_trackPhi_3_Outer_etaBin%i",ieta);
1133      sprintf(htit,  "h_trackPhi(pt>1.0, iso in 31x31,no missing Outer hit) : %3.2f<|#eta|<%3.2f)", lowEta, highEta);
1134     h_trackPhi_3_Outer[ieta]= new TH1F(hname, htit,  100, -3.14159, 3.14159);
1135 
1136   }
1137 
1138   TDirectory *d_meanTrackP   = fout->mkdir( "MeanTrackP" );
1139   TDirectory *d_maxNearP     = fout->mkdir( "MaxNearP" );
1140   TDirectory *d_neutralIso   = fout->mkdir( "NeutralIsolation" );
1141   TDirectory* d_trProf1      = fout->mkdir( "EcalTranverseProfile" );
1142   TDirectory* d_trProf2      = fout->mkdir( "HcalTranverseProfile" );
1143   TDirectory* d_response     = fout->mkdir( "Response" );
1144   TDirectory* d_response20Sig= fout->mkdir( "Response20Sig" );
1145 
1146   for(int ieta=0; ieta<NEtaBins; ieta++) {
1147     double lowEta=-5.0, highEta= 5.0;
1148     lowEta  = genPartEtaBins[ieta];
1149     highEta = genPartEtaBins[ieta+1];
1150 
1151     sprintf(hname, "h_eECAL11x11VsHCAL3x3_etaBin%i",ieta);
1152     sprintf(htit,  "eEvseH : %3.2f<|#eta|<%3.2f)", lowEta, highEta);
1153     h_eECAL11x11VsHCAL3x3[ieta] = new TH2F(hname, hname, 220, -1.0, 10.0, 220, -1.0, 10.0);
1154 
1155     d_trProf1->cd();
1156     sprintf(hname, "hh_eECAL3x3Frac_etaBin%i", ieta);
1157     sprintf(htit,  "eECAL3x3/trkP (%3.2f<|#eta|<%3.2f)", lowEta, highEta );
1158     hh_eECAL3x3_Frac[ieta] =  new TH1F(hname, htit, 500, -1.0, 4.0);
1159     hh_eECAL3x3_Frac[ieta] ->Sumw2();
1160     sprintf(hname, "hh_eECAL5x5Frac_etaBin%i", ieta);
1161     sprintf(htit,  "eECAL5x5/trkP (%3.2f<|#eta|<%3.2f)", lowEta, highEta );
1162     hh_eECAL5x5_Frac[ieta] =  new TH1F(hname, htit, 500, -1.0, 4.0);
1163     hh_eECAL5x5_Frac[ieta] ->Sumw2();
1164     sprintf(hname, "hh_eECAL7x7Frac_etaBin%i", ieta);
1165     sprintf(htit,  "eECAL7x7/trkP (%3.2f<|#eta|<%3.2f)", lowEta, highEta );
1166     hh_eECAL7x7_Frac[ieta] =  new TH1F(hname, htit, 500, -1.0, 4.0);
1167     hh_eECAL7x7_Frac[ieta] ->Sumw2();
1168     sprintf(hname, "hh_eECAL9x9Frac_etaBin%i", ieta);
1169     sprintf(htit,  "eECAL9x9/trkP (%3.2f<|#eta|<%3.2f)", lowEta, highEta );
1170     hh_eECAL9x9_Frac[ieta] =  new TH1F(hname, htit, 500, -1.0, 4.0);
1171     hh_eECAL9x9_Frac[ieta] ->Sumw2();
1172     sprintf(hname, "hh_eECAL11x11Frac_etaBin%i", ieta);
1173     sprintf(htit,  "eECAL11x11/trkP (%3.2f<|#eta|<%3.2f)", lowEta, highEta );
1174     hh_eECAL11x11_Frac[ieta] =  new TH1F(hname, htit, 500, -1.0, 4.0);
1175     hh_eECAL11x11_Frac[ieta] ->Sumw2();
1176     sprintf(hname, "hh_eECAL13x13Frac_etaBin%i", ieta);
1177     sprintf(htit,  "eECAL13x13/trkP (%3.2f<|#eta|<%3.2f)", lowEta, highEta );
1178     hh_eECAL13x13_Frac[ieta] =  new TH1F(hname, htit, 500, -1.0, 4.0);
1179     hh_eECAL13x13_Frac[ieta] ->Sumw2();
1180     sprintf(hname, "hh_eECAL15x15Frac_etaBin%i", ieta);
1181     sprintf(htit,  "eECAL15x15/trkP (%3.2f<|#eta|<%3.2f)", lowEta, highEta );
1182     hh_eECAL15x15_Frac[ieta] =  new TH1F(hname, htit, 500, -1.0, 4.0);
1183     hh_eECAL15x15_Frac[ieta] ->Sumw2();
1184     sprintf(hname, "hh_eECAL21x21Frac_etaBin%i", ieta);
1185     sprintf(htit,  "eECAL21x21/trkP (%3.2f<|#eta|<%3.2f)", lowEta, highEta );
1186     hh_eECAL21x21_Frac[ieta] =  new TH1F(hname, htit, 500, -1.0, 4.0);
1187     hh_eECAL21x21_Frac[ieta] ->Sumw2();
1188     sprintf(hname, "hh_eECAL25x25Frac_etaBin%i", ieta);
1189     sprintf(htit,  "eECAL25x25/trkP (%3.2f<|#eta|<%3.2f)", lowEta, highEta );
1190     hh_eECAL25x25_Frac[ieta] =  new TH1F(hname, htit, 500, -1.0, 4.0);
1191     hh_eECAL25x25_Frac[ieta] ->Sumw2();
1192     sprintf(hname, "hh_eECAL31x31Frac_etaBin%i", ieta);
1193     sprintf(htit,  "eECAL31x31/trkP (%3.2f<|#eta|<%3.2f)", lowEta, highEta );
1194     hh_eECAL31x31_Frac[ieta] =  new TH1F(hname, htit, 500, -1.0, 4.0);
1195     hh_eECAL31x31_Frac[ieta] ->Sumw2();
1196 
1197     sprintf(hname, "hh_eECAL7x7Frac20Sig_etaBin%i", ieta);
1198     sprintf(htit,  "eECAL7x7/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f)", lowEta, highEta );
1199     hh_eECAL7x7_Frac_20Sig[ieta] =  new TH1F(hname, htit, 500, -1.0, 4.0);
1200     hh_eECAL7x7_Frac_20Sig[ieta] ->Sumw2();
1201     sprintf(hname, "hh_eECAL9x9Frac20Sig_etaBin%i", ieta);
1202     sprintf(htit,  "eECAL9x9/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f)", lowEta, highEta );
1203     hh_eECAL9x9_Frac_20Sig[ieta] =  new TH1F(hname, htit, 500, -1.0, 4.0);
1204     hh_eECAL9x9_Frac_20Sig[ieta] ->Sumw2();
1205     sprintf(hname, "hh_eECAL11x11Frac20Sig_etaBin%i", ieta);
1206     sprintf(htit,  "eECAL11x11/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f)", lowEta, highEta );
1207     hh_eECAL11x11_Frac_20Sig[ieta] =  new TH1F(hname, htit, 500, -1.0, 4.0);
1208     hh_eECAL11x11_Frac_20Sig[ieta] ->Sumw2();
1209 
1210     d_trProf2->cd();
1211     sprintf(hname, "hh_eHCAL3x3Frac_etaBin%i", ieta);
1212     sprintf(htit,  "eHCAL3x3/trkP (%3.2f<|#eta|<%3.2f)", lowEta, highEta );
1213     hh_eHCAL3x3_Frac[ieta] = new TH1F(hname, htit, 500, -2.0, 4.0);
1214     hh_eHCAL3x3_Frac[ieta] ->Sumw2();
1215     sprintf(hname, "hh_eHCAL5x5Frac_etaBin%i", ieta);
1216     sprintf(htit,  "eHCAL5x5/trkP (%3.2f<|#eta|<%3.2f)", lowEta, highEta );
1217     hh_eHCAL5x5_Frac[ieta] = new TH1F(hname, htit, 500, -2.0, 4.0);
1218     hh_eHCAL5x5_Frac[ieta] ->Sumw2();
1219     sprintf(hname, "hh_eHCAL7x7Frac_etaBin%i", ieta);
1220     sprintf(htit,  "eHCAL7x7/trkP (%3.2f<|#eta|<%3.2f)", lowEta, highEta );
1221     hh_eHCAL7x7_Frac[ieta] = new TH1F(hname, htit, 500, -2.0, 4.0);
1222     hh_eHCAL7x7_Frac[ieta] ->Sumw2();
1223 
1224     sprintf(hname, "hh_eHCAL3x3Frac20Sig_etaBin%i", ieta);
1225     sprintf(htit,  "eHCAL3x3/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f)", lowEta, highEta );
1226     hh_eHCAL3x3_Frac_20Sig[ieta] = new TH1F(hname, htit, 500, -2.0, 4.0);
1227     hh_eHCAL3x3_Frac_20Sig[ieta] ->Sumw2();
1228     sprintf(hname, "hh_eHCAL5x5Frac20Sig_etaBin%i", ieta);
1229     sprintf(htit,  "eHCAL5x5/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f)", lowEta, highEta );
1230     hh_eHCAL5x5_Frac_20Sig[ieta] = new TH1F(hname, htit, 500, -2.0, 4.0);
1231     hh_eHCAL5x5_Frac_20Sig[ieta] ->Sumw2();
1232     sprintf(hname, "hh_eHCAL7x7Frac20Sig_etaBin%i", ieta);
1233     sprintf(htit,  "eHCAL7x7/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f)", lowEta, highEta );
1234     hh_eHCAL7x7_Frac_20Sig[ieta] = new TH1F(hname, htit, 500, -2.0, 4.0);
1235     hh_eHCAL7x7_Frac_20Sig[ieta] ->Sumw2();
1236 
1237 
1238     for(int ipt=0; ipt<NPBins; ipt++) {
1239       double lowP=0.0, highP=300.0;
1240       lowP    = genPartPBins[ipt];
1241       highP   = genPartPBins[ipt+1];
1242 
1243       d_neutralIso->cd();
1244       sprintf(hname, "h_diff_e15x15e11x11_ptBin%i_etaBin%i", ipt, ieta);
1245       sprintf(htit,  "h_diff_e15x15e11x11: (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1246       h_diff_e15x15e11x11      [ipt][ieta] = new TH1F(hname, htit, 600, -10.0, 50.0);
1247       sprintf(hname, "h_diff20Sig_e15x15e11x11_ptBin%i_etaBin%i", ipt, ieta);
1248       sprintf(htit,  "h_diff20Sig_e15x15e11x11: (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1249       h_diff_e15x15e11x11_20Sig[ipt][ieta] = new TH1F(hname, htit, 600, -10.0, 50.0);
1250 
1251       sprintf(hname, "h_diff_h7x7h5x5_ptBin%i_etaBin%i", ipt, ieta);
1252       sprintf(htit,  "h_diff_h7x7h5x5: (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1253       h_diff_h7x7h5x5[ipt][ieta] = new TH1F(hname, htit, 600, -10.0, 50.0);
1254 
1255       d_maxNearP->cd();
1256       sprintf(hname, "h_maxNearP31x31_ptBin%i_etaBin%i",ipt, ieta);
1257       sprintf(htit,  "maxNearP in 31x31 (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1258       h_maxNearP31x31[ipt][ieta] = new TH1F(hname, htit, 220, -2.0, 20.0);
1259       h_maxNearP31x31[ipt][ieta] ->Sumw2();
1260       sprintf(hname, "h_maxNearP25x25_ptBin%i_etaBin%i",ipt, ieta);
1261       sprintf(htit,  "maxNearP in 25x25 (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1262       h_maxNearP25x25[ipt][ieta] = new TH1F(hname, htit, 220, -2.0, 20.0);
1263       h_maxNearP25x25[ipt][ieta] ->Sumw2();
1264       sprintf(hname, "h_maxNearP21x21_ptBin%i_etaBin%i",ipt, ieta);
1265       sprintf(htit,  "maxNearP in 21x21 (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1266       h_maxNearP21x21[ipt][ieta] = new TH1F(hname, htit, 220, -2.0, 20.0);
1267       h_maxNearP21x21[ipt][ieta] ->Sumw2();
1268       sprintf(hname, "h_maxNearP15x15_ptBin%i_etaBin%i",ipt, ieta);
1269       sprintf(htit,  "maxNearP in 15x15 (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1270       h_maxNearP15x15[ipt][ieta] = new TH1F(hname, htit, 220, -2.0, 20.0);
1271       h_maxNearP15x15[ipt][ieta] ->Sumw2();
1272       sprintf(hname, "h_maxNearP11x11_ptBin%i_etaBin%i",ipt, ieta);
1273       sprintf(htit,  "maxNearP in 11x11 (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1274       h_maxNearP11x11[ipt][ieta] = new TH1F(hname, htit, 220, -2.0, 20.0);
1275       h_maxNearP11x11[ipt][ieta] ->Sumw2();
1276       sprintf(hname, "h_maxNearP9x9_ptBin%i_etaBin%i",ipt, ieta);
1277       sprintf(htit,  "maxNearP in 9x9 (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1278       h_maxNearP9x9[ipt][ieta] = new TH1F(hname, htit, 220, -2.0, 20.0);
1279       h_maxNearP9x9[ipt][ieta] ->Sumw2();
1280       sprintf(hname, "h_maxNearP7x7_ptBin%i_etaBin%i",ipt, ieta);
1281       sprintf(htit,  "maxNearP in 7x7 (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1282       h_maxNearP7x7[ipt][ieta] = new TH1F(hname, htit, 220, -2.0, 20.0);
1283       h_maxNearP7x7[ipt][ieta] ->Sumw2();
1284 
1285       //==============================
1286       // Ecal plots
1287       //==============================
1288       d_trProf1->cd();
1289       sprintf(hname, "h_eECAL3x3Frac_ptBin%i_etaBin%i", ipt, ieta);
1290       sprintf(htit,  "eECAL3x3/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1291       h_eECAL3x3_Frac[ipt][ieta] =  new TH1F(hname, htit, 1500, -1.0, 4.0);
1292       h_eECAL3x3_Frac[ipt][ieta] ->Sumw2();
1293       sprintf(hname, "h_eECAL5x5Frac_ptBin%i_etaBin%i", ipt, ieta);
1294       sprintf(htit,  "eECAL5x5/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1295       h_eECAL5x5_Frac[ipt][ieta] =  new TH1F(hname, htit, 1500, -1.0, 4.0);
1296       h_eECAL5x5_Frac[ipt][ieta] ->Sumw2();
1297       sprintf(hname, "h_eECAL7x7Frac_ptBin%i_etaBin%i", ipt, ieta);
1298       sprintf(htit,  "eECAL7x7/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1299       h_eECAL7x7_Frac[ipt][ieta] =  new TH1F(hname, htit, 1500, -1.0, 4.0);
1300       h_eECAL7x7_Frac[ipt][ieta] ->Sumw2();
1301       sprintf(hname, "h_eECAL9x9Frac_ptBin%i_etaBin%i", ipt, ieta);
1302       sprintf(htit,  "eECAL9x9/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1303       h_eECAL9x9_Frac[ipt][ieta] =  new TH1F(hname, htit, 1500, -1.0, 4.0);
1304       h_eECAL9x9_Frac[ipt][ieta] ->Sumw2();
1305       sprintf(hname, "h_eECAL11x11Frac_ptBin%i_etaBin%i", ipt, ieta);
1306       sprintf(htit,  "eECAL11x11/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1307       h_eECAL11x11_Frac[ipt][ieta] =  new TH1F(hname, htit, 1500, -1.0, 4.0);
1308       h_eECAL11x11_Frac[ipt][ieta] ->Sumw2();
1309       sprintf(hname, "h_eECAL13x13Frac_ptBin%i_etaBin%i", ipt, ieta);
1310       sprintf(htit,  "eECAL13x13/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1311       h_eECAL13x13_Frac[ipt][ieta] =  new TH1F(hname, htit, 1500, -1.0, 4.0);
1312       h_eECAL13x13_Frac[ipt][ieta] ->Sumw2();
1313       sprintf(hname, "h_eECAL15x15Frac_ptBin%i_etaBin%i", ipt, ieta);
1314       sprintf(htit,  "eECAL15x15/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1315       h_eECAL15x15_Frac[ipt][ieta] =  new TH1F(hname, htit, 1500, -1.0, 4.0);
1316       h_eECAL15x15_Frac[ipt][ieta] ->Sumw2();
1317       sprintf(hname, "h_eECAL21x21Frac_ptBin%i_etaBin%i", ipt, ieta);
1318       sprintf(htit,  "eECAL21x21/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1319       h_eECAL21x21_Frac[ipt][ieta] =  new TH1F(hname, htit, 1500, -1.0, 4.0);
1320       h_eECAL21x21_Frac[ipt][ieta] ->Sumw2();
1321       sprintf(hname, "h_eECAL25x25Frac_ptBin%i_etaBin%i", ipt, ieta);
1322       sprintf(htit,  "eECAL25x25/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1323       h_eECAL25x25_Frac[ipt][ieta] =  new TH1F(hname, htit, 1500, -1.0, 4.0);
1324       h_eECAL25x25_Frac[ipt][ieta] ->Sumw2();
1325       sprintf(hname, "h_eECAL31x31Frac_ptBin%i_etaBin%i", ipt, ieta);
1326       sprintf(htit,  "eECAL31x31/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1327       h_eECAL31x31_Frac[ipt][ieta] =  new TH1F(hname, htit, 1500, -1.0, 4.0);
1328       h_eECAL31x31_Frac[ipt][ieta] ->Sumw2();
1329 
1330       sprintf(hname, "h_eECAL7x7Frac20Sig_ptBin%i_etaBin%i", ipt, ieta);
1331       sprintf(htit,  "eECAL7x7/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1332       h_eECAL7x7_Frac_20Sig[ipt][ieta] =  new TH1F(hname, htit, 1500, -1.0, 4.0);
1333       h_eECAL7x7_Frac_20Sig[ipt][ieta] ->Sumw2();
1334       sprintf(hname, "h_eECAL9x9Frac20Sig_ptBin%i_etaBin%i", ipt, ieta);
1335       sprintf(htit,  "eECAL9x9/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1336       h_eECAL9x9_Frac_20Sig[ipt][ieta] =  new TH1F(hname, htit, 1500, -1.0, 4.0);
1337       h_eECAL9x9_Frac_20Sig[ipt][ieta] ->Sumw2();
1338       sprintf(hname, "h_eECAL11x11Frac20Sig_ptBin%i_etaBin%i", ipt, ieta);
1339       sprintf(htit,  "eECAL11x11/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1340       h_eECAL11x11_Frac_20Sig[ipt][ieta] =  new TH1F(hname, htit, 1500, -1.0, 4.0);
1341       h_eECAL11x11_Frac_20Sig[ipt][ieta] ->Sumw2();
1342 
1343       //==============================
1344       // Hcal plots
1345       //==============================
1346       d_trProf2->cd();
1347       sprintf(hname, "h_eHCAL3x3Frac_ptBin%i_etaBin%i", ipt, ieta);
1348       sprintf(htit,  "eHCAL3x3/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1349       h_eHCAL3x3_Frac[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1350       h_eHCAL3x3_Frac[ipt][ieta] ->Sumw2();
1351       sprintf(hname, "h_eHCAL5x5Frac_ptBin%i_etaBin%i", ipt, ieta);
1352       sprintf(htit,  "eHCAL5x5/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1353       h_eHCAL5x5_Frac[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1354       h_eHCAL5x5_Frac[ipt][ieta] ->Sumw2();
1355       sprintf(hname, "h_eHCAL7x7Frac_ptBin%i_etaBin%i", ipt, ieta);
1356       sprintf(htit,  "eHCAL7x7/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1357       h_eHCAL7x7_Frac[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1358       h_eHCAL7x7_Frac[ipt][ieta] ->Sumw2();
1359 
1360       sprintf(hname, "h_eHCAL3x3Frac20Sig_ptBin%i_etaBin%i", ipt, ieta);
1361       sprintf(htit,  "eHCAL3x3/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1362       h_eHCAL3x3_Frac_20Sig[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1363       h_eHCAL3x3_Frac_20Sig[ipt][ieta] ->Sumw2();
1364       sprintf(hname, "h_eHCAL5x5Frac20Sig_ptBin%i_etaBin%i", ipt, ieta);
1365       sprintf(htit,  "eHCAL5x5/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1366       h_eHCAL5x5_Frac_20Sig[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1367       h_eHCAL5x5_Frac_20Sig[ipt][ieta] ->Sumw2();
1368       sprintf(hname, "h_eHCAL7x7Frac20Sig_ptBin%i_etaBin%i", ipt, ieta);
1369       sprintf(htit,  "eHCAL7x7/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1370       h_eHCAL7x7_Frac_20Sig[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1371       h_eHCAL7x7_Frac_20Sig[ipt][ieta] ->Sumw2();
1372 
1373       sprintf(hname, "h_eHCAL3x3FracMIP_ptBin%i_etaBin%i", ipt, ieta);
1374       sprintf(htit,  "eHCAL3x3MIP/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1375       h_eHCAL3x3MIP_Frac[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1376       h_eHCAL3x3MIP_Frac[ipt][ieta] ->Sumw2();
1377       sprintf(hname, "h_eHCAL5x5FracMIP_ptBin%i_etaBin%i", ipt, ieta);
1378       sprintf(htit,  "eHCAL5x5MIP/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1379       h_eHCAL5x5MIP_Frac[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1380       h_eHCAL5x5MIP_Frac[ipt][ieta] ->Sumw2();
1381       sprintf(hname, "h_eHCAL7x7FracMIP_ptBin%i_etaBin%i", ipt, ieta);
1382       sprintf(htit,  "eHCAL7x7MIP/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1383       h_eHCAL7x7MIP_Frac[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1384       h_eHCAL7x7MIP_Frac[ipt][ieta] ->Sumw2();
1385 
1386       sprintf(hname, "h_eHCAL3x3FracMIP20Sig_ptBin%i_etaBin%i", ipt, ieta);
1387       sprintf(htit,  "eHCAL3x3MIP/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1388       h_eHCAL3x3MIP_Frac_20Sig[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1389       h_eHCAL3x3MIP_Frac_20Sig[ipt][ieta] ->Sumw2();
1390       sprintf(hname, "h_eHCAL5x5FracMIP20Sig_ptBin%i_etaBin%i", ipt, ieta);
1391       sprintf(htit,  "eHCAL5x5MIP/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1392       h_eHCAL5x5MIP_Frac_20Sig[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1393       h_eHCAL5x5MIP_Frac_20Sig[ipt][ieta] ->Sumw2();
1394       sprintf(hname, "h_eHCAL7x7FracMIP20Sig_ptBin%i_etaBin%i", ipt, ieta);
1395       sprintf(htit,  "eHCAL7x7MIP/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1396       h_eHCAL7x7MIP_Frac_20Sig[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1397       h_eHCAL7x7MIP_Frac_20Sig[ipt][ieta] ->Sumw2();
1398 
1399       //=================== Ecal+Hcal Response ====================
1400       d_response->cd();
1401       sprintf(hname, "h_Response_eHCAL3x3_eECAL11x11_ptBin%i_etaBin%i", ipt, ieta);
1402       sprintf(htit,  "(eHCAL3x3+eECAL11x11)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1403       h_eHCAL3x3_eECAL11x11_response[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1404       h_eHCAL3x3_eECAL11x11_response[ipt][ieta] ->Sumw2();
1405       sprintf(hname, "h_Response_eHCAL5x5_eECAL11x11_ptBin%i_etaBin%i", ipt, ieta);
1406       sprintf(htit,  "(eHCAL5x5+eECAL11x11)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1407       h_eHCAL5x5_eECAL11x11_response[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1408       h_eHCAL5x5_eECAL11x11_response[ipt][ieta] ->Sumw2();
1409       sprintf(hname, "h_Response_eHCAL7x7_eECAL11x11_ptBin%i_etaBin%i", ipt, ieta);
1410       sprintf(htit,  "(eHCAL7x7+eECAL11x11)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1411       h_eHCAL7x7_eECAL11x11_response[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1412       h_eHCAL7x7_eECAL11x11_response[ipt][ieta] ->Sumw2();
1413 
1414       sprintf(hname, "h_ResponseMIP_eHCAL3x3_eECAL11x11_ptBin%i_etaBin%i", ipt, ieta);
1415       sprintf(htit,  "(eHCAL3x3MIP+eECAL11x11)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1416       h_eHCAL3x3_eECAL11x11_responseMIP[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1417       h_eHCAL3x3_eECAL11x11_responseMIP[ipt][ieta] ->Sumw2();
1418       sprintf(hname, "h_ResponseMIP_eHCAL5x5_eECAL11x11_ptBin%i_etaBin%i", ipt, ieta);
1419       sprintf(htit,  "(eHCAL5x5MIP+eECAL11x11)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1420       h_eHCAL5x5_eECAL11x11_responseMIP[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1421       h_eHCAL5x5_eECAL11x11_responseMIP[ipt][ieta] ->Sumw2();
1422       sprintf(hname, "h_ResponseMIP_eHCAL7x7_eECAL11x11_ptBin%i_etaBin%i", ipt, ieta);
1423       sprintf(htit,  "(eHCAL7x7MIP+eECAL11x11)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1424       h_eHCAL7x7_eECAL11x11_responseMIP[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1425       h_eHCAL7x7_eECAL11x11_responseMIP[ipt][ieta] ->Sumw2();
1426       
1427       sprintf(hname, "h_ResponseInteract_eHCAL3x3_eECAL11x11_ptBin%i_etaBin%i", ipt, ieta);
1428       sprintf(htit,  "(eHCAL3x3Interact+eECAL11x11)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1429       h_eHCAL3x3_eECAL11x11_responseInteract[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1430       h_eHCAL3x3_eECAL11x11_responseInteract[ipt][ieta] ->Sumw2();
1431       sprintf(hname, "h_ResponseInteract_eHCAL5x5_eECAL11x11_ptBin%i_etaBin%i", ipt, ieta);
1432       sprintf(htit,  "(eHCAL5x5Interact+eECAL11x11)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1433       h_eHCAL5x5_eECAL11x11_responseInteract[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1434       h_eHCAL5x5_eECAL11x11_responseInteract[ipt][ieta] ->Sumw2();
1435       sprintf(hname, "h_ResponseInteract_eHCAL7x7_eECAL11x11_ptBin%i_etaBin%i", ipt, ieta);
1436       sprintf(htit,  "(eHCAL7x7Interact+eECAL11x11)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1437       h_eHCAL7x7_eECAL11x11_responseInteract[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1438       h_eHCAL7x7_eECAL11x11_responseInteract[ipt][ieta] ->Sumw2();
1439 
1440       sprintf(hname, "h_Response_eHCAL3x3_eECAL9x9_ptBin%i_etaBin%i", ipt, ieta);
1441       sprintf(htit,  "(eHCAL3x3+eECAL9x9)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1442       h_eHCAL3x3_eECAL9x9_response[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1443       h_eHCAL3x3_eECAL9x9_response[ipt][ieta] ->Sumw2();
1444       sprintf(hname, "h_Response_eHCAL5x5_eECAL9x9_ptBin%i_etaBin%i", ipt, ieta);
1445       sprintf(htit,  "(eHCAL5x5+eECAL9x9)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1446       h_eHCAL5x5_eECAL9x9_response[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1447       h_eHCAL5x5_eECAL9x9_response[ipt][ieta] ->Sumw2();
1448       sprintf(hname, "h_Response_eHCAL7x7_eECAL9x9_ptBin%i_etaBin%i", ipt, ieta);
1449       sprintf(htit,  "(eHCAL7x7+eECAL9x9)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1450       h_eHCAL7x7_eECAL9x9_response[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1451       h_eHCAL7x7_eECAL9x9_response[ipt][ieta] ->Sumw2();
1452 
1453       sprintf(hname, "h_ResponseMIP_eHCAL3x3_eECAL9x9_ptBin%i_etaBin%i", ipt, ieta);
1454       sprintf(htit,  "(eHCAL3x3MIP+eECAL9x9)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1455       h_eHCAL3x3_eECAL9x9_responseMIP[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1456       h_eHCAL3x3_eECAL9x9_responseMIP[ipt][ieta] ->Sumw2();
1457       sprintf(hname, "h_ResponseMIP_eHCAL5x5_eECAL9x9_ptBin%i_etaBin%i", ipt, ieta);
1458       sprintf(htit,  "(eHCAL5x5MIP+eECAL9x9)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1459       h_eHCAL5x5_eECAL9x9_responseMIP[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1460       h_eHCAL5x5_eECAL9x9_responseMIP[ipt][ieta] ->Sumw2();
1461       sprintf(hname, "h_ResponseMIP_eHCAL7x7_eECAL9x9_ptBin%i_etaBin%i", ipt, ieta);
1462       sprintf(htit,  "(eHCAL7x7MIP+eECAL9x9)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1463       h_eHCAL7x7_eECAL9x9_responseMIP[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1464       h_eHCAL7x7_eECAL9x9_responseMIP[ipt][ieta] ->Sumw2();
1465       
1466       sprintf(hname, "h_ResponseInteract_eHCAL3x3_eECAL9x9_ptBin%i_etaBin%i", ipt, ieta);
1467       sprintf(htit,  "(eHCAL3x3Interact+eECAL9x9)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1468       h_eHCAL3x3_eECAL9x9_responseInteract[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1469       h_eHCAL3x3_eECAL9x9_responseInteract[ipt][ieta] ->Sumw2();
1470       sprintf(hname, "h_ResponseInteract_eHCAL5x5_eECAL9x9_ptBin%i_etaBin%i", ipt, ieta);
1471       sprintf(htit,  "(eHCAL5x5Interact+eECAL9x9)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1472       h_eHCAL5x5_eECAL9x9_responseInteract[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1473       h_eHCAL5x5_eECAL9x9_responseInteract[ipt][ieta] ->Sumw2();
1474       sprintf(hname, "h_ResponseInteract_eHCAL7x7_eECAL9x9_ptBin%i_etaBin%i", ipt, ieta);
1475       sprintf(htit,  "(eHCAL7x7Interact+eECAL9x9)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1476       h_eHCAL7x7_eECAL9x9_responseInteract[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1477       h_eHCAL7x7_eECAL9x9_responseInteract[ipt][ieta] ->Sumw2();
1478 
1479       sprintf(hname, "h_Response_eHCAL3x3_eECAL7x7_ptBin%i_etaBin%i", ipt, ieta);
1480       sprintf(htit,  "(eHCAL3x3+eECAL7x7)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1481       h_eHCAL3x3_eECAL7x7_response[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1482       h_eHCAL3x3_eECAL7x7_response[ipt][ieta] ->Sumw2();
1483       sprintf(hname, "h_Response_eHCAL5x5_eECAL7x7_ptBin%i_etaBin%i", ipt, ieta);
1484       sprintf(htit,  "(eHCAL5x5+eECAL7x7)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1485       h_eHCAL5x5_eECAL7x7_response[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1486       h_eHCAL5x5_eECAL7x7_response[ipt][ieta] ->Sumw2();
1487       sprintf(hname, "h_Response_eHCAL7x7_eECAL7x7_ptBin%i_etaBin%i", ipt, ieta);
1488       sprintf(htit,  "(eHCAL7x7+eECAL7x7)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1489       h_eHCAL7x7_eECAL7x7_response[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1490       h_eHCAL7x7_eECAL7x7_response[ipt][ieta] ->Sumw2();
1491 
1492       sprintf(hname, "h_ResponseMIP_eHCAL3x3_eECAL7x7_ptBin%i_etaBin%i", ipt, ieta);
1493       sprintf(htit,  "(eHCAL3x3MIP+eECAL7x7)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1494       h_eHCAL3x3_eECAL7x7_responseMIP[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1495       h_eHCAL3x3_eECAL7x7_responseMIP[ipt][ieta] ->Sumw2();
1496       sprintf(hname, "h_ResponseMIP_eHCAL5x5_eECAL7x7_ptBin%i_etaBin%i", ipt, ieta);
1497       sprintf(htit,  "(eHCAL5x5MIP+eECAL7x7)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1498       h_eHCAL5x5_eECAL7x7_responseMIP[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1499       h_eHCAL5x5_eECAL7x7_responseMIP[ipt][ieta] ->Sumw2();
1500       sprintf(hname, "h_ResponseMIP_eHCAL7x7_eECAL7x7_ptBin%i_etaBin%i", ipt, ieta);
1501       sprintf(htit,  "(eHCAL7x7MIP+eECAL7x7)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1502       h_eHCAL7x7_eECAL7x7_responseMIP[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1503       h_eHCAL7x7_eECAL7x7_responseMIP[ipt][ieta] ->Sumw2();
1504       
1505       sprintf(hname, "h_ResponseInteract_eHCAL3x3_eECAL7x7_ptBin%i_etaBin%i", ipt, ieta);
1506       sprintf(htit,  "(eHCAL3x3Interact+eECAL7x7)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1507       h_eHCAL3x3_eECAL7x7_responseInteract[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1508       h_eHCAL3x3_eECAL7x7_responseInteract[ipt][ieta] ->Sumw2();
1509       sprintf(hname, "h_ResponseInteract_eHCAL5x5_eECAL7x7_ptBin%i_etaBin%i", ipt, ieta);
1510       sprintf(htit,  "(eHCAL5x5Interact+eECAL7x7)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1511       h_eHCAL5x5_eECAL7x7_responseInteract[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1512       h_eHCAL5x5_eECAL7x7_responseInteract[ipt][ieta] ->Sumw2();
1513       sprintf(hname, "h_ResponseInteract_eHCAL7x7_eECAL7x7_ptBin%i_etaBin%i", ipt, ieta);
1514       sprintf(htit,  "(eHCAL7x7Interact+eECAL7x7)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1515       h_eHCAL7x7_eECAL7x7_responseInteract[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1516       h_eHCAL7x7_eECAL7x7_responseInteract[ipt][ieta] ->Sumw2();
1517       
1518       d_response20Sig->cd();
1519       sprintf(hname, "h_Response20Sig_eHCAL3x3_eECAL11x11_ptBin%i_etaBin%i", ipt, ieta);
1520       sprintf(htit,  "(eHCAL3x3+eECAL11x11)(Xtal>2.0#sigma)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1521       h_eHCAL3x3_eECAL11x11_response_20Sig[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1522       h_eHCAL3x3_eECAL11x11_response_20Sig[ipt][ieta] ->Sumw2();
1523       sprintf(hname, "h_ResponseMIP20Sig_eHCAL3x3_eECAL11x11_ptBin%i_etaBin%i", ipt, ieta);
1524       sprintf(htit,  "(eHCAL3x3MIP+eECAL11x11)/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1525       h_eHCAL3x3_eECAL11x11_responseMIP_20Sig[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1526       h_eHCAL3x3_eECAL11x11_responseMIP_20Sig[ipt][ieta] ->Sumw2();
1527       sprintf(hname, "h_ResponseInteract20Sig_eHCAL3x3_eECAL11x11_ptBin%i_etaBin%i", ipt, ieta);
1528       sprintf(htit,  "(eHCAL3x3Interact+eECAL11x11)/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1529       h_eHCAL3x3_eECAL11x11_responseInteract_20Sig[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1530       h_eHCAL3x3_eECAL11x11_responseInteract_20Sig[ipt][ieta] ->Sumw2();
1531 
1532       sprintf(hname, "h_Response20Sig_eHCAL3x3_eECAL9x9_ptBin%i_etaBin%i", ipt, ieta);
1533       sprintf(htit,  "(eHCAL3x3+eECAL9x9)(Xtal>2.0#sigma)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1534       h_eHCAL3x3_eECAL9x9_response_20Sig[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1535       h_eHCAL3x3_eECAL9x9_response_20Sig[ipt][ieta] ->Sumw2();
1536       sprintf(hname, "h_ResponseMIP20Sig_eHCAL3x3_eECAL9x9_ptBin%i_etaBin%i", ipt, ieta);
1537       sprintf(htit,  "(eHCAL3x3MIP+eECAL9x9)/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1538       h_eHCAL3x3_eECAL9x9_responseMIP_20Sig[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1539       h_eHCAL3x3_eECAL9x9_responseMIP_20Sig[ipt][ieta] ->Sumw2();
1540       sprintf(hname, "h_ResponseInteract20Sig_eHCAL3x3_eECAL9x9_ptBin%i_etaBin%i", ipt, ieta);
1541       sprintf(htit,  "(eHCAL3x3Interact+eECAL9x9)/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1542       h_eHCAL3x3_eECAL9x9_responseInteract_20Sig[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1543       h_eHCAL3x3_eECAL9x9_responseInteract_20Sig[ipt][ieta] ->Sumw2();
1544 
1545       sprintf(hname, "h_Response20Sig_eHCAL3x3_eECAL7x7_ptBin%i_etaBin%i", ipt, ieta);
1546       sprintf(htit,  "(eHCAL3x3+eECAL7x7)(Xtal>2.0#sigma)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1547       h_eHCAL3x3_eECAL7x7_response_20Sig[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1548       h_eHCAL3x3_eECAL7x7_response_20Sig[ipt][ieta] ->Sumw2();
1549       sprintf(hname, "h_ResponseMIP20Sig_eHCAL3x3_eECAL7x7_ptBin%i_etaBin%i", ipt, ieta);
1550       sprintf(htit,  "(eHCAL3x3MIP+eECAL7x7)/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1551       h_eHCAL3x3_eECAL7x7_responseMIP_20Sig[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1552       h_eHCAL3x3_eECAL7x7_responseMIP_20Sig[ipt][ieta] ->Sumw2();
1553       sprintf(hname, "h_ResponseInteract20Sig_eHCAL3x3_eECAL7x7_ptBin%i_etaBin%i", ipt, ieta);
1554       sprintf(htit,  "(eHCAL3x3Interact+eECAL7x7)/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1555       h_eHCAL3x3_eECAL7x7_responseInteract_20Sig[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1556       h_eHCAL3x3_eECAL7x7_responseInteract_20Sig[ipt][ieta] ->Sumw2();
1557       
1558       sprintf(hname, "h_Response20Sig_eHCAL5x5_eECAL7x7_ptBin%i_etaBin%i", ipt, ieta);
1559       sprintf(htit,  "(eHCAL5x5+eECAL7x7)(Xtal>2.0#sigma)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1560       h_eHCAL5x5_eECAL7x7_response_20Sig[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1561       h_eHCAL5x5_eECAL7x7_response_20Sig[ipt][ieta] ->Sumw2();
1562       sprintf(hname, "h_ResponseMIP20Sig_eHCAL5x5_eECAL7x7_ptBin%i_etaBin%i", ipt, ieta);
1563       sprintf(htit,  "(eHCAL5x5MIP+eECAL7x7)/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1564       h_eHCAL5x5_eECAL7x7_responseMIP_20Sig[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1565       h_eHCAL5x5_eECAL7x7_responseMIP_20Sig[ipt][ieta] ->Sumw2();
1566       sprintf(hname, "h_ResponseInteract20Sig_eHCAL5x5_eECAL7x7_ptBin%i_etaBin%i", ipt, ieta);
1567       sprintf(htit,  "(eHCAL5x5Interact+eECAL7x7)/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1568       h_eHCAL5x5_eECAL7x7_responseInteract_20Sig[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1569       h_eHCAL5x5_eECAL7x7_responseInteract_20Sig[ipt][ieta] ->Sumw2();
1570  
1571       sprintf(hname, "h_Response20Sig_eHCAL5x5_eECAL11x11_ptBin%i_etaBin%i", ipt, ieta);
1572       sprintf(htit,  "(eHCAL5x5+eECAL11x11)(Xtal>2.0#sigma)/trkP (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1573       h_eHCAL5x5_eECAL11x11_response_20Sig[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1574       h_eHCAL5x5_eECAL11x11_response_20Sig[ipt][ieta] ->Sumw2();
1575       sprintf(hname, "h_ResponseMIP20Sig_eHCAL5x5_eECAL11x11_ptBin%i_etaBin%i", ipt, ieta);
1576       sprintf(htit,  "(eHCAL5x5MIP+eECAL11x11)/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1577       h_eHCAL5x5_eECAL11x11_responseMIP_20Sig[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1578       h_eHCAL5x5_eECAL11x11_responseMIP_20Sig[ipt][ieta] ->Sumw2();
1579       sprintf(hname, "h_ResponseInteract20Sig_eHCAL5x5_eECAL11x11_ptBin%i_etaBin%i", ipt, ieta);
1580       sprintf(htit,  "(eHCAL5x5Interact+eECAL11x11)/trkP(Xtal>2.0#sigma) (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );
1581       h_eHCAL5x5_eECAL11x11_responseInteract_20Sig[ipt][ieta] = new TH1F(hname, htit, 1500, -1.0, 4.0);
1582       h_eHCAL5x5_eECAL11x11_responseInteract_20Sig[ipt][ieta] ->Sumw2();
1583 
1584       d_maxNearP->cd();
1585       sprintf(hname, "h_meanTrackP_ptBin%i_etaBin%i", ipt, ieta);
1586       sprintf(htit,  "Track Momentum (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP );             
1587       h_meanTrackP[ipt][ieta] = new TH1F(hname, htit, 100, lowP*0.9, highP*1.1);
1588       h_meanTrackP[ipt][ieta] ->Sumw2();
1589     }
1590   }
1591  
1592   fout->cd();
1593   double hcalEtaBins[55] = {-2.650,-2.500,-2.322,-2.172,-2.043,-1.930,
1594                 -1.830,-1.740,-1.653,-1.566,-1.479,-1.392,-1.305,
1595                 -1.218,-1.131,-1.044,-0.957,-0.879,-0.783,-0.696,
1596                 -0.609,-0.522,-0.435,-0.348,-0.261,-0.174,-0.087,
1597                 0.000,
1598                 0.087, 0.174, 0.261, 0.348, 0.435, 0.522, 0.609,
1599                 0.696, 0.783, 0.879, 0.957, 1.044, 1.131, 1.218,
1600                 1.305, 1.392, 1.479, 1.566, 1.653, 1.740, 1.830,
1601                 1.930, 2.043, 2.172, 2.322, 2.500, 2.650};
1602   h_HcalMeanEneVsEta = new TProfile("h_HcalMeanEneVsEta", "h_HcalMeanEneVsEta", 54, hcalEtaBins);
1603   h_HcalMeanEneVsEta->Sumw2();
1604 
1605 }
1606