Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-10-25 09:34:25

0001 ////////////////////////////////////////////////////////////////
0002 // This class skeleton has been automatically generated on
0003 // from TTree hcalCalibTree/Tree for IsoTrack Calibration
0004 // with ROOT version 5.17/02
0005 //
0006 //  TSelector-based code for getting the HCAL resp. correction
0007 //  from physics events. Works for DiJet and IsoTrack calibration.
0008 //
0009 //  Anton Anastassov (Northwestern)
0010 //  Email: aa@fnal.gov
0011 //
0012 //
0013 ///////////////////////////////////////////////////////////////
0014 
0015 #ifndef hcalCalib_h
0016 #define hcalCalib_h
0017 
0018 #include <vector>
0019 #include <map>
0020 
0021 #include <TROOT.h>
0022 #include <TChain.h>
0023 #include <TFile.h>
0024 #include <TSelector.h>
0025 #include <TROOT.h>
0026 #include <TChain.h>
0027 #include <TFile.h>
0028 #include <TSelector.h>
0029 
0030 #include <TLorentzVector.h>
0031 #include <TClonesArray.h>
0032 #include <TRefArray.h>
0033 
0034 #include <TH2.h>
0035 #include <TStyle.h>
0036 #include <TFile.h>
0037 #include <TMatrixF.h>
0038 #include <TMatrixD.h>
0039 #include <TDecompSVD.h>
0040 #include <TDecompQRH.h>
0041 
0042 // needed to get cell coordinates
0043 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
0044 #include "Geometry/CaloTopology/interface/HcalTopology.h"
0045 
0046 class hcalCalib : public TSelector {
0047 public:
0048   TTree *fChain;  //!pointer to the analyzed TTree or TChain
0049 
0050   UInt_t eventNumber;
0051   UInt_t runNumber;
0052   Int_t iEtaHit;
0053   UInt_t iPhiHit;
0054   TClonesArray *cells;
0055   Float_t emEnergy;
0056   Float_t targetE;
0057   Float_t etVetoJet;
0058 
0059   Float_t xTrkHcal;
0060   Float_t yTrkHcal;
0061   Float_t zTrkHcal;
0062   Float_t xTrkEcal;
0063   Float_t yTrkEcal;
0064   Float_t zTrkEcal;
0065 
0066   TLorentzVector *tagJetP4;
0067   TLorentzVector *probeJetP4;
0068 
0069   Float_t tagJetEmFrac;
0070   Float_t probeJetEmFrac;
0071 
0072   // List of branches
0073   TBranch *b_eventNumber;  //!
0074   TBranch *b_runNumber;    //!
0075   TBranch *b_iEtaHit;      //!
0076   TBranch *b_iPhiHit;      //!
0077   TBranch *b_cells;        //!
0078   TBranch *b_emEnergy;     //!
0079   TBranch *b_targetE;      //!
0080   TBranch *b_etVetoJet;    //!
0081 
0082   TBranch *b_xTrkHcal;
0083   TBranch *b_yTrkHcal;
0084   TBranch *b_zTrkHcal;
0085   TBranch *b_xTrkEcal;
0086   TBranch *b_yTrkEcal;
0087   TBranch *b_zTrkEcal;
0088 
0089   TBranch *b_tagJetEmFrac;    //!
0090   TBranch *b_probeJetEmFrac;  //!
0091 
0092   TBranch *b_tagJetP4;    //!
0093   TBranch *b_probeJetP4;  //!
0094 
0095   UInt_t nEvents;
0096 
0097   TFile *histoFile;
0098 
0099   // sanity check histograms
0100   TH1F *h1_trkP;
0101   TH1F *h1_allTrkP;
0102 
0103   TH1F *h1_selTrkP_iEta10;
0104 
0105   TH1F *h1_rawSumE;
0106   TH1F *h1_rawResp;
0107   TH1F *h1_corResp;
0108   TH1F *h1_rawRespBarrel;
0109   TH1F *h1_corRespBarrel;
0110   TH1F *h1_rawRespEndcap;
0111   TH1F *h1_corRespEndcap;
0112   TH1F *h1_numEventsTwrIEta;
0113 
0114   TH2F *h2_dHitRefBarrel;
0115   TH2F *h2_dHitRefEndcap;
0116 
0117   // histograms based on iEta, iPhi of refPosition forthe cluster (at the moment: hottest tower)
0118   // expect range |iEta|<=24 (to do: add flexibility for arbitrary range)
0119   TH1F *h1_corRespIEta[48];
0120 
0121   hcalCalib(TTree * /*tree*/ = nullptr) {}
0122   ~hcalCalib() override {}
0123   Int_t Version() const override { return 2; }
0124   void Begin(TTree *tree) override;
0125   //   virtual void    SlaveBegin(TTree *tree);
0126   void Init(TTree *tree) override;
0127   Bool_t Notify() override;
0128   Bool_t Process(Long64_t entry) override;
0129   Int_t GetEntry(Long64_t entry, Int_t getall = 0) override {
0130     return fChain ? fChain->GetTree()->GetEntry(entry, getall) : 0;
0131   }
0132   void SetOption(const char *option) override { fOption = option; }
0133   void SetObject(TObject *obj) override { fObject = obj; }
0134   void SetInputList(TList *input) override { fInput = input; }
0135   TList *GetOutputList() const override { return fOutput; }
0136   //   virtual void    SlaveTerminate();
0137   void Terminate() override;
0138 
0139   //------------ CUTS ---------------
0140   Float_t MIN_TARGET_E;
0141   Float_t MAX_TARGET_E;
0142 
0143   Float_t MIN_CELL_E;
0144   Float_t MIN_EOVERP;
0145   Float_t MAX_EOVERP;
0146   Float_t MAX_TRK_EME;
0147 
0148   Float_t MAX_ET_THIRD_JET;
0149   Float_t MIN_DPHI_DIJETS;
0150 
0151   Bool_t SUM_DEPTHS;
0152   Bool_t SUM_SMALL_DEPTHS;
0153   Bool_t COMBINE_PHI;
0154 
0155   Int_t HB_CLUSTER_SIZE;
0156   Int_t HE_CLUSTER_SIZE;
0157 
0158   Bool_t USE_CONE_CLUSTERING;
0159   Float_t MAX_CONE_DIST;
0160 
0161   Int_t CALIB_ABS_IETA_MAX;
0162   Int_t CALIB_ABS_IETA_MIN;
0163 
0164   Float_t MAX_PROBEJET_EMFRAC;
0165   Float_t MAX_TAGJET_EMFRAC;
0166   Float_t MAX_TAGJET_ABSETA;
0167   Float_t MIN_TAGJET_ET;
0168 
0169   Float_t MIN_PROBEJET_ABSETA;
0170 
0171   TString CALIB_TYPE;    // "ISO_TRACK" or "DI_JET"
0172   TString CALIB_METHOD;  // L3, matrix inversion, everage then matrix inversion,...
0173 
0174   TString PHI_SYM_COR_FILENAME;
0175   Bool_t APPLY_PHI_SYM_COR_FLAG;
0176 
0177   TString OUTPUT_COR_COEF_FILENAME;
0178   TString HISTO_FILENAME;
0179 
0180   const CaloGeometry *theCaloGeometry;
0181   const HcalTopology *topo_;
0182 
0183   void SetMinTargetE(Float_t e) { MIN_TARGET_E = e; }
0184   void SetMaxTargetE(Float_t e) { MAX_TARGET_E = e; }
0185   void SetSumDepthsFlag(Bool_t b) { SUM_DEPTHS = b; }
0186   void SetSumSmallDepthsFlag(Bool_t b) { SUM_SMALL_DEPTHS = b; }
0187   void SetCombinePhiFlag(Bool_t b) { COMBINE_PHI = b; }
0188   void SetMinCellE(Float_t e) { MIN_CELL_E = e; }
0189   void SetMinEOverP(Float_t e) { MIN_EOVERP = e; }
0190   void SetMaxEOverP(Float_t e) { MAX_EOVERP = e; }
0191   void SetMaxTrkEmE(Float_t e) { MAX_TRK_EME = e; }
0192   void SetCalibType(const TString &s) { CALIB_TYPE = s; }
0193   void SetCalibMethod(const TString &s) { CALIB_METHOD = s; }
0194   void SetHbClusterSize(Int_t i) { HB_CLUSTER_SIZE = i; }
0195   void SetHeClusterSize(Int_t i) { HE_CLUSTER_SIZE = i; }
0196 
0197   void SetUseConeClustering(Bool_t b) { USE_CONE_CLUSTERING = b; }
0198   void SetConeMaxDist(Float_t d) { MAX_CONE_DIST = d; }
0199 
0200   void SetCalibAbsIEtaMax(Int_t i) { CALIB_ABS_IETA_MAX = i; }
0201   void SetCalibAbsIEtaMin(Int_t i) { CALIB_ABS_IETA_MIN = i; }
0202   void SetMaxEtThirdJet(Float_t et) { MAX_ET_THIRD_JET = et; }
0203   void SetMinDPhiDiJets(Float_t dphi) { MIN_DPHI_DIJETS = dphi; }
0204   void SetApplyPhiSymCorFlag(Bool_t b) { APPLY_PHI_SYM_COR_FLAG = b; }
0205   void SetPhiSymCorFileName(const TString &filename) { PHI_SYM_COR_FILENAME = filename; }
0206   void SetMaxProbeJetEmFrac(Float_t f) { MAX_PROBEJET_EMFRAC = f; }
0207   void SetMaxTagJetEmFrac(Float_t f) { MAX_TAGJET_EMFRAC = f; }
0208   void SetMaxTagJetAbsEta(Float_t e) { MAX_TAGJET_ABSETA = e; }
0209   void SetMinTagJetEt(Float_t e) { MIN_TAGJET_ET = e; }
0210   void SetMinProbeJetAbsEta(Float_t e) { MIN_PROBEJET_ABSETA = e; }
0211   void SetOutputCorCoefFileName(const TString &filename) { OUTPUT_COR_COEF_FILENAME = filename; }
0212   void SetHistoFileName(const TString &filename) { HISTO_FILENAME = filename; }
0213 
0214   void SetCaloGeometry(const CaloGeometry *g, const HcalTopology *topo) {
0215     theCaloGeometry = g;
0216     topo_ = topo;
0217   }
0218 
0219   void GetCoefFromMtrxInvOfAve();
0220 
0221   Bool_t ReadPhiSymCor();
0222 
0223   void makeTextFile();
0224 
0225   // --------- containers passed to minimizers ----------------
0226   std::vector<std::vector<Float_t> > cellEnergies;
0227   std::vector<std::vector<UInt_t> > cellIds;
0228   std::vector<std::pair<Int_t, UInt_t> > refIEtaIPhi;  // centroid of jet or hottest tower iEta, iPhi
0229   std::vector<Float_t> targetEnergies;
0230 
0231   std::map<UInt_t, Float_t> phiSymCor;  // holds the phi symmetry corrections read from the file
0232 
0233   std::map<UInt_t, Float_t> solution;  // correction coef: solution from L3, holds final coef for hybrid methods as well
0234   std::map<Int_t, Float_t>
0235       iEtaCoefMap;  // correction coef: from matrix inversion AFTER averaging, also intermediate results for hybrid methods
0236 
0237   //   ClassDef(hcalCalib,0);
0238 };
0239 
0240 #endif