Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-07 04:35:00

0001 //
0002 // DiJetAnalyzer.cc
0003 //
0004 // description: Create an ntuple necessary for dijet balance calibration for the HCAL
0005 //
0006 
0007 // system include files
0008 #include <memory>
0009 #include <string>
0010 #include <vector>
0011 #include <set>
0012 #include <map>
0013 
0014 #include "TTree.h"
0015 #include "TFile.h"
0016 #include "TH1D.h"
0017 #include "TH2D.h"
0018 #include "TClonesArray.h"
0019 
0020 // user include files
0021 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0022 #include "CondFormats/JetMETObjects/interface/JetCorrectorParameters.h"
0023 
0024 #include "FWCore/Framework/interface/Frameworkfwd.h"
0025 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0026 #include "FWCore/Framework/interface/Event.h"
0027 #include "FWCore/Framework/interface/EventSetup.h"
0028 #include "FWCore/Framework/interface/MakerMacros.h"
0029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0030 #include "FWCore/ServiceRegistry/interface/Service.h"
0031 #include "FWCore/Utilities/interface/EDGetToken.h"
0032 #include "FWCore/Utilities/interface/EDMException.h"
0033 
0034 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
0035 #include "DataFormats/JetReco/interface/PFJetCollection.h"
0036 #include "DataFormats/JetReco/interface/GenJetCollection.h"
0037 #include "DataFormats/HcalRecHit/interface/HBHERecHit.h"
0038 #include "DataFormats/HcalRecHit/interface/HFRecHit.h"
0039 #include "DataFormats/HcalRecHit/interface/HORecHit.h"
0040 #include "DataFormats/ParticleFlowReco/interface/PFBlockFwd.h"
0041 #include "DataFormats/ParticleFlowReco/interface/PFBlock.h"
0042 #include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
0043 #include "DataFormats/ParticleFlowReco/interface/PFRecHit.h"
0044 #include "DataFormats/ParticleFlowReco/interface/PFRecHitFwd.h"
0045 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0046 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0047 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0048 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0049 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0050 
0051 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
0052 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
0053 #include "Geometry/HcalTowerAlgo/interface/HcalGeometry.h"
0054 #include "Geometry/Records/interface/CaloGeometryRecord.h"
0055 
0056 #include "JetMETCorrections/JetCorrector/interface/JetCorrector.h"
0057 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
0058 
0059 //
0060 // class declarations
0061 //
0062 
0063 class JetCorretPair : protected std::pair<const reco::PFJet*, double> {
0064 public:
0065   JetCorretPair() {
0066     first = 0;
0067     second = 1.0;
0068   }
0069   JetCorretPair(const reco::PFJet* j, double s) {
0070     first = j;
0071     second = s;
0072   }
0073   ~JetCorretPair() {}
0074 
0075   inline const reco::PFJet* jet(void) const { return first; }
0076   inline void jet(const reco::PFJet* j) {
0077     first = j;
0078     return;
0079   }
0080   inline double scale(void) const { return second; }
0081   inline void scale(double d) {
0082     second = d;
0083     return;
0084   }
0085 
0086 private:
0087 };
0088 
0089 class DiJetAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0090 public:
0091   explicit DiJetAnalyzer(const edm::ParameterSet&);
0092   ~DiJetAnalyzer() override = default;
0093 
0094 private:
0095   void beginJob() override;  //(const edm::EventSetup&);
0096   void analyze(const edm::Event&, const edm::EventSetup&) override;
0097   void endJob() override;
0098 
0099   // parameters
0100   const std::string pfJetCollName_;     // label for the PF jet collection
0101   const std::string hbheRecHitName_;    // label for HBHERecHits collection
0102   const std::string hfRecHitName_;      // label for HFRecHit collection
0103   const std::string hoRecHitName_;      // label for HORecHit collection
0104   const std::string pvCollName_;        // label for primary vertex collection
0105   const std::string rootHistFilename_;  // name of the histogram file
0106   const double maxDeltaEta_;            // maximum delta-|Eta| between Jets
0107   const double minTagJetEta_;           // minimum |eta| of the tag jet
0108   const double maxTagJetEta_;           // maximum |eta| of the tag jet
0109   const double minSumJetEt_;            // minimum Sum of the tag and probe jet Et
0110   const double minJetEt_;               // minimum Jet Et
0111   const double maxThirdJetEt_;          // maximum 3rd jet Et
0112   const bool debug_;                    // print debug statements
0113 
0114   const edm::EDGetTokenT<reco::PFJetCollection> tok_PFJet_;
0115   const edm::EDGetTokenT<edm::SortedCollection<HBHERecHit, edm::StrictWeakOrdering<HBHERecHit>>> tok_HBHE_;
0116   const edm::EDGetTokenT<edm::SortedCollection<HFRecHit, edm::StrictWeakOrdering<HFRecHit>>> tok_HF_;
0117   const edm::EDGetTokenT<edm::SortedCollection<HORecHit, edm::StrictWeakOrdering<HORecHit>>> tok_HO_;
0118   const edm::EDGetTokenT<reco::VertexCollection> tok_Vertex_;
0119   const edm::EDGetTokenT<reco::JetCorrector> jetCorrectorToken_;
0120 
0121   const edm::ESGetToken<CaloGeometry, CaloGeometryRecord> tok_geom_;
0122 
0123   // root file/histograms
0124   TH1D* h_PassSelPF_;
0125   TTree* tree_;
0126 
0127   float tpfjet_pt_, tpfjet_p_, tpfjet_E_, tpfjet_eta_, tpfjet_phi_, tpfjet_EMfrac_, tpfjet_hadEcalEfrac_, tpfjet_scale_,
0128       tpfjet_area_;
0129   int tpfjet_jetID_;
0130   float tpfjet_EBE_, tpfjet_EEE_, tpfjet_HBE_, tpfjet_HEE_, tpfjet_HFE_;
0131   float tpfjet_unkown_E_, tpfjet_unkown_px_, tpfjet_unkown_py_, tpfjet_unkown_pz_, tpfjet_unkown_EcalE_;
0132   float tpfjet_electron_E_, tpfjet_electron_px_, tpfjet_electron_py_, tpfjet_electron_pz_, tpfjet_electron_EcalE_;
0133   float tpfjet_muon_E_, tpfjet_muon_px_, tpfjet_muon_py_, tpfjet_muon_pz_, tpfjet_muon_EcalE_;
0134   float tpfjet_photon_E_, tpfjet_photon_px_, tpfjet_photon_py_, tpfjet_photon_pz_, tpfjet_photon_EcalE_;
0135   int tpfjet_unkown_n_, tpfjet_electron_n_, tpfjet_muon_n_, tpfjet_photon_n_;
0136   int tpfjet_had_n_;
0137   std::vector<float> tpfjet_had_E_, tpfjet_had_px_, tpfjet_had_py_, tpfjet_had_pz_, tpfjet_had_EcalE_,
0138       tpfjet_had_rawHcalE_, tpfjet_had_emf_;
0139   std::vector<int> tpfjet_had_id_, tpfjet_had_candtrackind_, tpfjet_had_ntwrs_;
0140   int tpfjet_ntwrs_;
0141   std::vector<int> tpfjet_twr_ieta_, tpfjet_twr_iphi_, tpfjet_twr_depth_, tpfjet_twr_subdet_, tpfjet_twr_candtrackind_,
0142       tpfjet_twr_hadind_, tpfjet_twr_elmttype_, tpfjet_twr_clusterind_;
0143   std::vector<float> tpfjet_twr_hade_, tpfjet_twr_frac_, tpfjet_twr_dR_;
0144   int tpfjet_cluster_n_;
0145   std::vector<float> tpfjet_cluster_eta_, tpfjet_cluster_phi_, tpfjet_cluster_dR_;
0146   int tpfjet_ncandtracks_;
0147   std::vector<float> tpfjet_candtrack_px_, tpfjet_candtrack_py_, tpfjet_candtrack_pz_, tpfjet_candtrack_EcalE_;
0148   float ppfjet_pt_, ppfjet_p_, ppfjet_E_, ppfjet_eta_, ppfjet_phi_, ppfjet_EMfrac_, ppfjet_hadEcalEfrac_, ppfjet_scale_,
0149       ppfjet_area_;
0150   int ppfjet_jetID_;
0151   float ppfjet_EBE_, ppfjet_EEE_, ppfjet_HBE_, ppfjet_HEE_, ppfjet_HFE_;
0152   float ppfjet_unkown_E_, ppfjet_unkown_px_, ppfjet_unkown_py_, ppfjet_unkown_pz_, ppfjet_unkown_EcalE_;
0153   float ppfjet_electron_E_, ppfjet_electron_px_, ppfjet_electron_py_, ppfjet_electron_pz_, ppfjet_electron_EcalE_;
0154   float ppfjet_muon_E_, ppfjet_muon_px_, ppfjet_muon_py_, ppfjet_muon_pz_, ppfjet_muon_EcalE_;
0155   float ppfjet_photon_E_, ppfjet_photon_px_, ppfjet_photon_py_, ppfjet_photon_pz_, ppfjet_photon_EcalE_;
0156   int ppfjet_unkown_n_, ppfjet_electron_n_, ppfjet_muon_n_, ppfjet_photon_n_;
0157   int ppfjet_had_n_;
0158   std::vector<float> ppfjet_had_E_, ppfjet_had_px_, ppfjet_had_py_, ppfjet_had_pz_, ppfjet_had_EcalE_,
0159       ppfjet_had_rawHcalE_, ppfjet_had_emf_;
0160   std::vector<int> ppfjet_had_id_, ppfjet_had_candtrackind_, ppfjet_had_ntwrs_;
0161   int ppfjet_ntwrs_;
0162   std::vector<int> ppfjet_twr_ieta_, ppfjet_twr_iphi_, ppfjet_twr_depth_, ppfjet_twr_subdet_, ppfjet_twr_candtrackind_,
0163       ppfjet_twr_hadind_, ppfjet_twr_elmttype_, ppfjet_twr_clusterind_;
0164   std::vector<float> ppfjet_twr_hade_, ppfjet_twr_frac_, ppfjet_twr_dR_;
0165   int ppfjet_cluster_n_;
0166   std::vector<float> ppfjet_cluster_eta_, ppfjet_cluster_phi_, ppfjet_cluster_dR_;
0167   int ppfjet_ncandtracks_;
0168   std::vector<float> ppfjet_candtrack_px_, ppfjet_candtrack_py_, ppfjet_candtrack_pz_, ppfjet_candtrack_EcalE_;
0169   float pf_dijet_deta_, pf_dijet_dphi_, pf_dijet_balance_;
0170   float pf_thirdjet_px_, pf_thirdjet_py_, pf_realthirdjet_px_, pf_realthirdjet_py_, pf_realthirdjet_scale_;
0171   int pf_Run_, pf_Lumi_, pf_Event_;
0172   int pf_NPV_;
0173 
0174   // helper functions
0175   double deltaR(const reco::Jet* j1, const reco::Jet* j2);
0176   double deltaR(const double eta1, const double phi1, const double eta2, const double phi2);
0177   int getEtaPhi(const DetId id);
0178   int getEtaPhi(const HcalDetId id);
0179 
0180   struct JetCorretPairComp {
0181     inline bool operator()(const JetCorretPair& a, const JetCorretPair& b) const {
0182       return (a.jet()->pt() * a.scale()) > (b.jet()->pt() * b.scale());
0183     }
0184   };
0185 };
0186 
0187 DiJetAnalyzer::DiJetAnalyzer(const edm::ParameterSet& iConfig)
0188     : pfJetCollName_(iConfig.getParameter<std::string>("pfJetCollName")),
0189       hbheRecHitName_(iConfig.getParameter<std::string>("hbheRecHitName")),
0190       hfRecHitName_(iConfig.getParameter<std::string>("hfRecHitName")),
0191       hoRecHitName_(iConfig.getParameter<std::string>("hoRecHitName")),
0192       pvCollName_(iConfig.getParameter<std::string>("pvCollName")),
0193       rootHistFilename_(iConfig.getParameter<std::string>("rootHistFilename")),
0194       maxDeltaEta_(iConfig.getParameter<double>("maxDeltaEta")),
0195       minTagJetEta_(iConfig.getParameter<double>("minTagJetEta")),
0196       maxTagJetEta_(iConfig.getParameter<double>("maxTagJetEta")),
0197       minSumJetEt_(iConfig.getParameter<double>("minSumJetEt")),
0198       minJetEt_(iConfig.getParameter<double>("minJetEt")),
0199       maxThirdJetEt_(iConfig.getParameter<double>("maxThirdJetEt")),
0200       debug_(iConfig.getUntrackedParameter<bool>("debug", false)),
0201       tok_PFJet_(consumes<reco::PFJetCollection>(pfJetCollName_)),
0202       tok_HBHE_(consumes<edm::SortedCollection<HBHERecHit, edm::StrictWeakOrdering<HBHERecHit>>>(hbheRecHitName_)),
0203       tok_HF_(consumes<edm::SortedCollection<HFRecHit, edm::StrictWeakOrdering<HFRecHit>>>(hfRecHitName_)),
0204       tok_HO_(consumes<edm::SortedCollection<HORecHit, edm::StrictWeakOrdering<HORecHit>>>(hoRecHitName_)),
0205       tok_Vertex_(consumes<reco::VertexCollection>(pvCollName_)),
0206       jetCorrectorToken_(consumes<reco::JetCorrector>(iConfig.getParameter<edm::InputTag>("JetCorrections"))),
0207       tok_geom_(esConsumes<CaloGeometry, CaloGeometryRecord>()) {
0208   usesResource(TFileService::kSharedResource);
0209 }
0210 
0211 //
0212 // member functions
0213 //
0214 
0215 // ------------ method called to for each event  ------------
0216 void DiJetAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& evSetup) {
0217   pf_Run_ = iEvent.id().run();
0218   pf_Lumi_ = iEvent.id().luminosityBlock();
0219   pf_Event_ = iEvent.id().event();
0220 
0221   // Get PFJets
0222   const edm::Handle<reco::PFJetCollection> pfjets = iEvent.getHandle(tok_PFJet_);
0223   if (!pfjets.isValid()) {
0224     throw edm::Exception(edm::errors::ProductNotFound)
0225         << " could not find PFJetCollection named " << pfJetCollName_ << ".\n";
0226     return;
0227   }
0228 
0229   // Get RecHits in HB and HE
0230   const edm::Handle<edm::SortedCollection<HBHERecHit, edm::StrictWeakOrdering<HBHERecHit>>> hbhereco =
0231       iEvent.getHandle(tok_HBHE_);
0232   if (!hbhereco.isValid()) {
0233     throw edm::Exception(edm::errors::ProductNotFound)
0234         << " could not find HBHERecHit named " << hbheRecHitName_ << ".\n";
0235     return;
0236   }
0237 
0238   // Get RecHits in HF
0239   const edm::Handle<edm::SortedCollection<HFRecHit, edm::StrictWeakOrdering<HFRecHit>>> hfreco =
0240       iEvent.getHandle(tok_HF_);
0241   if (!hfreco.isValid()) {
0242     throw edm::Exception(edm::errors::ProductNotFound) << " could not find HFRecHit named " << hfRecHitName_ << ".\n";
0243     return;
0244   }
0245 
0246   // Get RecHits in HO
0247   const edm::Handle<edm::SortedCollection<HORecHit, edm::StrictWeakOrdering<HORecHit>>> horeco =
0248       iEvent.getHandle(tok_HO_);
0249   if (!horeco.isValid()) {
0250     throw edm::Exception(edm::errors::ProductNotFound) << " could not find HORecHit named " << hoRecHitName_ << ".\n";
0251     return;
0252   }
0253 
0254   // Get geometry
0255   const CaloGeometry* geo = &evSetup.getData(tok_geom_);
0256   const HcalGeometry* HBGeom = static_cast<const HcalGeometry*>(geo->getSubdetectorGeometry(DetId::Hcal, 1));
0257   const HcalGeometry* HEGeom = static_cast<const HcalGeometry*>(geo->getSubdetectorGeometry(DetId::Hcal, 2));
0258   const CaloSubdetectorGeometry* HOGeom = geo->getSubdetectorGeometry(DetId::Hcal, 3);
0259   const CaloSubdetectorGeometry* HFGeom = geo->getSubdetectorGeometry(DetId::Hcal, 4);
0260 
0261   // Get primary vertices
0262   const edm::Handle<std::vector<reco::Vertex>> pv = iEvent.getHandle(tok_Vertex_);
0263   if (!pv.isValid()) {
0264     throw edm::Exception(edm::errors::ProductNotFound) << " could not find Vertex named " << pvCollName_ << ".\n";
0265     return;
0266   }
0267   pf_NPV_ = 0;
0268   for (std::vector<reco::Vertex>::const_iterator it = pv->begin(); it != pv->end(); ++it) {
0269     if (!it->isFake() && it->ndof() > 4)
0270       ++pf_NPV_;
0271   }
0272 
0273   reco::JetCorrector const& correctorPF = iEvent.get(jetCorrectorToken_);
0274 
0275   //////////////////////////////
0276   // Event Selection
0277   //////////////////////////////
0278 
0279   // determine which cut results in failure
0280   int passSelPF = 0;
0281 
0282   // sort jets by corrected et
0283   std::set<JetCorretPair, JetCorretPairComp> pfjetcorretpairset;
0284   for (reco::PFJetCollection::const_iterator it = pfjets->begin(); it != pfjets->end(); ++it) {
0285     const reco::PFJet* jet = &(*it);
0286     double jec = correctorPF.correction(*it);
0287     pfjetcorretpairset.insert(JetCorretPair(jet, jec));
0288   }
0289 
0290   JetCorretPair pf_tag, pf_probe;
0291   pf_thirdjet_px_ = pf_thirdjet_py_ = 0.0;
0292   pf_realthirdjet_px_ = pf_realthirdjet_py_ = 0.0;
0293   pf_realthirdjet_px_ = 1;
0294   int cntr = 0;
0295   for (std::set<JetCorretPair, JetCorretPairComp>::const_iterator it = pfjetcorretpairset.begin();
0296        it != pfjetcorretpairset.end();
0297        ++it) {
0298     JetCorretPair jet = (*it);
0299     ++cntr;
0300     if (cntr == 1)
0301       pf_tag = jet;
0302     else if (cntr == 2)
0303       pf_probe = jet;
0304     else {
0305       pf_thirdjet_px_ += jet.scale() * jet.jet()->px();
0306       pf_thirdjet_py_ += jet.scale() * jet.jet()->py();
0307       if (cntr == 3) {
0308         pf_realthirdjet_px_ = jet.jet()->px();
0309         pf_realthirdjet_py_ = jet.jet()->py();
0310         pf_realthirdjet_scale_ = jet.scale();
0311       }
0312     }
0313   }
0314 
0315   if (pf_tag.jet() && pf_probe.jet()) {
0316     // require that the first two jets are above some minimum,
0317     // and the rest are below some maximum
0318     if ((pf_tag.jet()->et() + pf_probe.jet()->et()) < minSumJetEt_)
0319       passSelPF |= 0x1;
0320     if (pf_tag.jet()->et() < minJetEt_ || pf_probe.jet()->et() < minJetEt_)
0321       passSelPF |= 0x2;
0322     if (sqrt(pf_thirdjet_px_ * pf_thirdjet_px_ + pf_thirdjet_py_ * pf_thirdjet_py_) > maxThirdJetEt_)
0323       passSelPF |= 0x4;
0324 
0325     // force the tag jet to have the smaller |eta|
0326     if (std::fabs(pf_tag.jet()->eta()) > std::fabs(pf_probe.jet()->eta())) {
0327       JetCorretPair temp = pf_tag;
0328       pf_tag = pf_probe;
0329       pf_probe = temp;
0330     }
0331 
0332     // eta cuts
0333     double dAbsEta = std::fabs(std::fabs(pf_tag.jet()->eta()) - std::fabs(pf_probe.jet()->eta()));
0334     if (dAbsEta > maxDeltaEta_)
0335       passSelPF |= 0x8;
0336     if (fabs(pf_tag.jet()->eta()) < minTagJetEta_)
0337       passSelPF |= 0x10;
0338     if (fabs(pf_tag.jet()->eta()) > maxTagJetEta_)
0339       passSelPF |= 0x10;
0340   } else {
0341     passSelPF = 0x40;
0342   }
0343 
0344   h_PassSelPF_->Fill(passSelPF);
0345   if (passSelPF)
0346     return;
0347   // dump
0348   if (debug_) {
0349     edm::LogVerbatim("HcalCalib") << "Run: " << iEvent.id().run() << "; Event: " << iEvent.id().event();
0350     for (reco::PFJetCollection::const_iterator it = pfjets->begin(); it != pfjets->end(); ++it) {
0351       const reco::PFJet* jet = &(*it);
0352       edm::LogVerbatim("HcalCalib") << "istag=" << (jet == pf_tag.jet()) << "; isprobe=" << (jet == pf_probe.jet())
0353                                     << "; et=" << jet->et() << "; eta=" << jet->eta();
0354     }
0355   }
0356 
0357   // Reset particle variables
0358   tpfjet_unkown_E_ = tpfjet_unkown_px_ = tpfjet_unkown_py_ = tpfjet_unkown_pz_ = tpfjet_unkown_EcalE_ = 0.0;
0359   tpfjet_electron_E_ = tpfjet_electron_px_ = tpfjet_electron_py_ = tpfjet_electron_pz_ = tpfjet_electron_EcalE_ = 0.0;
0360   tpfjet_muon_E_ = tpfjet_muon_px_ = tpfjet_muon_py_ = tpfjet_muon_pz_ = tpfjet_muon_EcalE_ = 0.0;
0361   tpfjet_photon_E_ = tpfjet_photon_px_ = tpfjet_photon_py_ = tpfjet_photon_pz_ = tpfjet_photon_EcalE_ = 0.0;
0362   tpfjet_unkown_n_ = tpfjet_electron_n_ = tpfjet_muon_n_ = tpfjet_photon_n_ = 0;
0363   tpfjet_had_n_ = 0;
0364   tpfjet_cluster_n_ = 0;
0365   ppfjet_unkown_E_ = ppfjet_unkown_px_ = ppfjet_unkown_py_ = ppfjet_unkown_pz_ = ppfjet_unkown_EcalE_ = 0.0;
0366   ppfjet_electron_E_ = ppfjet_electron_px_ = ppfjet_electron_py_ = ppfjet_electron_pz_ = ppfjet_electron_EcalE_ = 0.0;
0367   ppfjet_muon_E_ = ppfjet_muon_px_ = ppfjet_muon_py_ = ppfjet_muon_pz_ = ppfjet_muon_EcalE_ = 0.0;
0368   ppfjet_photon_E_ = ppfjet_photon_px_ = ppfjet_photon_py_ = ppfjet_photon_pz_ = ppfjet_photon_EcalE_ = 0.0;
0369   ppfjet_unkown_n_ = ppfjet_electron_n_ = ppfjet_muon_n_ = ppfjet_photon_n_ = 0;
0370   ppfjet_had_n_ = 0;
0371   ppfjet_cluster_n_ = 0;
0372 
0373   tpfjet_had_E_.clear();
0374   tpfjet_had_px_.clear();
0375   tpfjet_had_py_.clear();
0376   tpfjet_had_pz_.clear();
0377   tpfjet_had_EcalE_.clear();
0378   tpfjet_had_rawHcalE_.clear();
0379   tpfjet_had_emf_.clear();
0380   tpfjet_had_id_.clear();
0381   tpfjet_had_candtrackind_.clear();
0382   tpfjet_had_ntwrs_.clear();
0383   tpfjet_twr_ieta_.clear();
0384   tpfjet_twr_iphi_.clear();
0385   tpfjet_twr_depth_.clear();
0386   tpfjet_twr_subdet_.clear();
0387   tpfjet_twr_candtrackind_.clear();
0388   tpfjet_twr_hadind_.clear();
0389   tpfjet_twr_elmttype_.clear();
0390   tpfjet_twr_hade_.clear();
0391   tpfjet_twr_frac_.clear();
0392   tpfjet_twr_dR_.clear();
0393   tpfjet_twr_clusterind_.clear();
0394   tpfjet_cluster_eta_.clear();
0395   tpfjet_cluster_phi_.clear();
0396   tpfjet_cluster_dR_.clear();
0397   tpfjet_candtrack_px_.clear();
0398   tpfjet_candtrack_py_.clear();
0399   tpfjet_candtrack_pz_.clear();
0400   tpfjet_candtrack_EcalE_.clear();
0401   ppfjet_had_E_.clear();
0402   ppfjet_had_px_.clear();
0403   ppfjet_had_py_.clear();
0404   ppfjet_had_pz_.clear();
0405   ppfjet_had_EcalE_.clear();
0406   ppfjet_had_rawHcalE_.clear();
0407   ppfjet_had_emf_.clear();
0408   ppfjet_had_id_.clear();
0409   ppfjet_had_candtrackind_.clear();
0410   ppfjet_had_ntwrs_.clear();
0411   ppfjet_twr_ieta_.clear();
0412   ppfjet_twr_iphi_.clear();
0413   ppfjet_twr_depth_.clear();
0414   ppfjet_twr_subdet_.clear();
0415   ppfjet_twr_candtrackind_.clear();
0416   ppfjet_twr_hadind_.clear();
0417   ppfjet_twr_elmttype_.clear();
0418   ppfjet_twr_hade_.clear();
0419   ppfjet_twr_frac_.clear();
0420   ppfjet_twr_dR_.clear();
0421   ppfjet_twr_clusterind_.clear();
0422   ppfjet_cluster_eta_.clear();
0423   ppfjet_cluster_phi_.clear();
0424   ppfjet_cluster_dR_.clear();
0425   ppfjet_candtrack_px_.clear();
0426   ppfjet_candtrack_py_.clear();
0427   ppfjet_candtrack_pz_.clear();
0428   ppfjet_candtrack_EcalE_.clear();
0429 
0430   std::map<int, std::pair<int, std::set<float>>> tpfjet_rechits;
0431   std::map<int, std::pair<int, std::set<float>>> ppfjet_rechits;
0432   std::map<float, int> tpfjet_clusters;
0433   std::map<float, int> ppfjet_clusters;
0434 
0435   // fill tag jet variables
0436   tpfjet_pt_ = pf_tag.jet()->pt();
0437   tpfjet_p_ = pf_tag.jet()->p();
0438   tpfjet_E_ = pf_tag.jet()->energy();
0439   tpfjet_eta_ = pf_tag.jet()->eta();
0440   tpfjet_phi_ = pf_tag.jet()->phi();
0441   tpfjet_scale_ = pf_tag.scale();
0442   tpfjet_area_ = pf_tag.jet()->jetArea();
0443   tpfjet_ntwrs_ = 0;
0444   tpfjet_ncandtracks_ = 0;
0445 
0446   tpfjet_jetID_ = 0;  // Not a loose, medium, or tight jet
0447   if (fabs(pf_tag.jet()->eta()) < 2.4) {
0448     if (pf_tag.jet()->chargedHadronEnergyFraction() > 0 && pf_tag.jet()->chargedMultiplicity() > 0 &&
0449         pf_tag.jet()->chargedEmEnergyFraction() < 0.99 &&
0450         (pf_tag.jet()->chargedMultiplicity() + pf_tag.jet()->neutralMultiplicity()) > 1) {
0451       if (pf_tag.jet()->neutralHadronEnergyFraction() < 0.9 && pf_tag.jet()->neutralEmEnergyFraction() < 0.9) {
0452         tpfjet_jetID_ = 3;  // Tight jet
0453       } else if (pf_tag.jet()->neutralHadronEnergyFraction() < 0.95 && pf_tag.jet()->neutralEmEnergyFraction() < 0.95) {
0454         tpfjet_jetID_ = 2;  // Medium jet
0455       } else if (pf_tag.jet()->neutralHadronEnergyFraction() < 0.99 && pf_tag.jet()->neutralEmEnergyFraction() < 0.99) {
0456         tpfjet_jetID_ = 1;  // Loose jet
0457       }
0458     }
0459   } else if ((pf_tag.jet()->chargedMultiplicity() + pf_tag.jet()->neutralMultiplicity()) > 1) {
0460     if (pf_tag.jet()->neutralHadronEnergyFraction() < 0.9 && pf_tag.jet()->neutralEmEnergyFraction() < 0.9) {
0461       tpfjet_jetID_ = 3;  // Tight jet
0462     } else if (pf_tag.jet()->neutralHadronEnergyFraction() < 0.95 && pf_tag.jet()->neutralEmEnergyFraction() < 0.95) {
0463       tpfjet_jetID_ = 2;  // Medium jet
0464     } else if (pf_tag.jet()->neutralHadronEnergyFraction() < 0.99 && pf_tag.jet()->neutralEmEnergyFraction() < 0.99) {
0465       tpfjet_jetID_ = 1;  // Loose jet
0466     }
0467   }
0468 
0469   /////////////////////////////////////////////
0470   // Get PF constituents and fill HCAL towers
0471   /////////////////////////////////////////////
0472 
0473   // Get tag PFCandidates
0474   std::vector<reco::PFCandidatePtr> tagconst = pf_tag.jet()->getPFConstituents();
0475   for (std::vector<reco::PFCandidatePtr>::const_iterator it = tagconst.begin(); it != tagconst.end(); ++it) {
0476     bool hasTrack = false;
0477     // Test PFCandidate type
0478     reco::PFCandidate::ParticleType candidateType = (*it)->particleId();
0479     switch (candidateType) {
0480       case reco::PFCandidate::X:
0481         tpfjet_unkown_E_ += (*it)->energy();
0482         tpfjet_unkown_px_ += (*it)->px();
0483         tpfjet_unkown_py_ += (*it)->py();
0484         tpfjet_unkown_pz_ += (*it)->pz();
0485         tpfjet_unkown_EcalE_ += (*it)->ecalEnergy();
0486         tpfjet_unkown_n_++;
0487         continue;
0488       case reco::PFCandidate::h: {
0489         tpfjet_had_E_.push_back((*it)->energy());
0490         tpfjet_had_px_.push_back((*it)->px());
0491         tpfjet_had_py_.push_back((*it)->py());
0492         tpfjet_had_pz_.push_back((*it)->pz());
0493         tpfjet_had_EcalE_.push_back((*it)->ecalEnergy());
0494         tpfjet_had_rawHcalE_.push_back((*it)->rawHcalEnergy());
0495         tpfjet_had_id_.push_back(0);
0496         tpfjet_had_ntwrs_.push_back(0);
0497         tpfjet_had_n_++;
0498 
0499         reco::TrackRef trackRef = (*it)->trackRef();
0500         if (trackRef.isNonnull()) {
0501           reco::Track track = *trackRef;
0502           tpfjet_candtrack_px_.push_back(track.px());
0503           tpfjet_candtrack_py_.push_back(track.py());
0504           tpfjet_candtrack_pz_.push_back(track.pz());
0505           tpfjet_candtrack_EcalE_.push_back((*it)->ecalEnergy());
0506           tpfjet_had_candtrackind_.push_back(tpfjet_ncandtracks_);
0507           hasTrack = true;
0508           tpfjet_ncandtracks_++;
0509         } else {
0510           tpfjet_had_candtrackind_.push_back(-2);
0511         }
0512       } break;
0513       case reco::PFCandidate::e:
0514         tpfjet_electron_E_ += (*it)->energy();
0515         tpfjet_electron_px_ += (*it)->px();
0516         tpfjet_electron_py_ += (*it)->py();
0517         tpfjet_electron_pz_ += (*it)->pz();
0518         tpfjet_electron_EcalE_ += (*it)->ecalEnergy();
0519         tpfjet_electron_n_++;
0520         continue;
0521       case reco::PFCandidate::mu:
0522         tpfjet_muon_E_ += (*it)->energy();
0523         tpfjet_muon_px_ += (*it)->px();
0524         tpfjet_muon_py_ += (*it)->py();
0525         tpfjet_muon_pz_ += (*it)->pz();
0526         tpfjet_muon_EcalE_ += (*it)->ecalEnergy();
0527         tpfjet_muon_n_++;
0528         continue;
0529       case reco::PFCandidate::gamma:
0530         tpfjet_photon_E_ += (*it)->energy();
0531         tpfjet_photon_px_ += (*it)->px();
0532         tpfjet_photon_py_ += (*it)->py();
0533         tpfjet_photon_pz_ += (*it)->pz();
0534         tpfjet_photon_EcalE_ += (*it)->ecalEnergy();
0535         tpfjet_photon_n_++;
0536         continue;
0537       case reco::PFCandidate::h0: {
0538         tpfjet_had_E_.push_back((*it)->energy());
0539         tpfjet_had_px_.push_back((*it)->px());
0540         tpfjet_had_py_.push_back((*it)->py());
0541         tpfjet_had_pz_.push_back((*it)->pz());
0542         tpfjet_had_EcalE_.push_back((*it)->ecalEnergy());
0543         tpfjet_had_rawHcalE_.push_back((*it)->rawHcalEnergy());
0544         tpfjet_had_id_.push_back(1);
0545         tpfjet_had_candtrackind_.push_back(-1);
0546         tpfjet_had_ntwrs_.push_back(0);
0547         tpfjet_had_n_++;
0548         break;
0549       }
0550       case reco::PFCandidate::h_HF: {
0551         tpfjet_had_E_.push_back((*it)->energy());
0552         tpfjet_had_px_.push_back((*it)->px());
0553         tpfjet_had_py_.push_back((*it)->py());
0554         tpfjet_had_pz_.push_back((*it)->pz());
0555         tpfjet_had_EcalE_.push_back((*it)->ecalEnergy());
0556         tpfjet_had_rawHcalE_.push_back((*it)->rawHcalEnergy());
0557         tpfjet_had_id_.push_back(2);
0558         tpfjet_had_candtrackind_.push_back(-1);
0559         tpfjet_had_ntwrs_.push_back(0);
0560         tpfjet_had_n_++;
0561         break;
0562       }
0563       case reco::PFCandidate::egamma_HF: {
0564         tpfjet_had_E_.push_back((*it)->energy());
0565         tpfjet_had_px_.push_back((*it)->px());
0566         tpfjet_had_py_.push_back((*it)->py());
0567         tpfjet_had_pz_.push_back((*it)->pz());
0568         tpfjet_had_EcalE_.push_back((*it)->ecalEnergy());
0569         tpfjet_had_rawHcalE_.push_back((*it)->rawHcalEnergy());
0570         tpfjet_had_id_.push_back(3);
0571         tpfjet_had_candtrackind_.push_back(-1);
0572         tpfjet_had_ntwrs_.push_back(0);
0573         tpfjet_had_n_++;
0574         break;
0575       }
0576     }
0577 
0578     std::map<int, int> twrietas;
0579     float HFHAD_E = 0;
0580     float HFEM_E = 0;
0581     int maxElement = (*it)->elementsInBlocks().size();
0582     for (int e = 0; e < maxElement; ++e) {
0583       // Get elements from block
0584       reco::PFBlockRef blockRef = (*it)->elementsInBlocks()[e].first;
0585       const edm::OwnVector<reco::PFBlockElement>& elements = blockRef->elements();
0586       for (unsigned iEle = 0; iEle < elements.size(); iEle++) {
0587         if (elements[iEle].index() == (*it)->elementsInBlocks()[e].second) {
0588           if (elements[iEle].type() == reco::PFBlockElement::HCAL) {  // Element is HB or HE
0589             // Get cluster and hits
0590             reco::PFClusterRef clusterref = elements[iEle].clusterRef();
0591             reco::PFCluster cluster = *clusterref;
0592             double cluster_dR = deltaR(tpfjet_eta_, tpfjet_phi_, cluster.eta(), cluster.phi());
0593             if (tpfjet_clusters.count(cluster_dR) == 0) {
0594               tpfjet_clusters[cluster_dR] = tpfjet_cluster_n_;
0595               tpfjet_cluster_eta_.push_back(cluster.eta());
0596               tpfjet_cluster_phi_.push_back(cluster.phi());
0597               tpfjet_cluster_dR_.push_back(cluster_dR);
0598               tpfjet_cluster_n_++;
0599             }
0600             int cluster_ind = tpfjet_clusters[cluster_dR];
0601 
0602             std::vector<std::pair<DetId, float>> hitsAndFracs = cluster.hitsAndFractions();
0603 
0604             // Run over hits and match
0605             int nHits = hitsAndFracs.size();
0606             for (int iHit = 0; iHit < nHits; iHit++) {
0607               int etaPhiPF = getEtaPhi(hitsAndFracs[iHit].first);
0608 
0609               for (edm::SortedCollection<HBHERecHit, edm::StrictWeakOrdering<HBHERecHit>>::const_iterator ith =
0610                        hbhereco->begin();
0611                    ith != hbhereco->end();
0612                    ++ith) {
0613                 int etaPhiRecHit = getEtaPhi((*ith).id());
0614                 if (etaPhiPF == etaPhiRecHit) {
0615                   tpfjet_had_ntwrs_.at(tpfjet_had_n_ - 1)++;
0616                   if (tpfjet_rechits.count((*ith).id()) == 0) {
0617                     tpfjet_twr_ieta_.push_back((*ith).id().ieta());
0618                     tpfjet_twr_iphi_.push_back((*ith).id().iphi());
0619                     tpfjet_twr_depth_.push_back((*ith).id().depth());
0620                     tpfjet_twr_subdet_.push_back((*ith).id().subdet());
0621                     if (hitsAndFracs[iHit].second > 0.05 && (*ith).energy() > 0.0)
0622                       twrietas[(*ith).id().ieta()]++;
0623                     tpfjet_twr_hade_.push_back((*ith).energy());
0624                     tpfjet_twr_frac_.push_back(hitsAndFracs[iHit].second);
0625                     tpfjet_rechits[(*ith).id()].second.insert(hitsAndFracs[iHit].second);
0626                     tpfjet_twr_hadind_.push_back(tpfjet_had_n_ - 1);
0627                     tpfjet_twr_elmttype_.push_back(0);
0628                     tpfjet_twr_clusterind_.push_back(cluster_ind);
0629                     if (hasTrack) {
0630                       tpfjet_twr_candtrackind_.push_back(tpfjet_ncandtracks_ - 1);
0631                     } else {
0632                       tpfjet_twr_candtrackind_.push_back(-1);
0633                     }
0634                     switch ((*ith).id().subdet()) {
0635                       case HcalSubdetector::HcalBarrel: {
0636                         CaloCellGeometry::CornersVec cv = HBGeom->getCorners((*ith).id());
0637                         float avgeta = (cv[0].eta() + cv[2].eta()) / 2.0;
0638                         float avgphi = (static_cast<double>(cv[0].phi()) + static_cast<double>(cv[2].phi())) / 2.0;
0639                         if (cv[0].phi() < cv[2].phi())
0640                           avgphi = (2.0 * 3.141592653 + static_cast<double>(cv[0].phi()) +
0641                                     static_cast<double>(cv[2].phi())) /
0642                                    2.0;
0643                         tpfjet_twr_dR_.push_back(deltaR(tpfjet_eta_, tpfjet_phi_, avgeta, avgphi));
0644                         break;
0645                       }
0646                       case HcalSubdetector::HcalEndcap: {
0647                         CaloCellGeometry::CornersVec cv = HEGeom->getCorners((*ith).id());
0648                         float avgeta = (cv[0].eta() + cv[2].eta()) / 2.0;
0649                         float avgphi = (static_cast<double>(cv[0].phi()) + static_cast<double>(cv[2].phi())) / 2.0;
0650                         if (cv[0].phi() < cv[2].phi())
0651                           avgphi = (2.0 * 3.141592653 + static_cast<double>(cv[0].phi()) +
0652                                     static_cast<double>(cv[2].phi())) /
0653                                    2.0;
0654                         tpfjet_twr_dR_.push_back(deltaR(tpfjet_eta_, tpfjet_phi_, avgeta, avgphi));
0655                         break;
0656                       }
0657                       default:
0658                         tpfjet_twr_dR_.push_back(-1);
0659                         break;
0660                     }
0661                     tpfjet_rechits[(*ith).id()].first = tpfjet_ntwrs_;
0662                     ++tpfjet_ntwrs_;
0663                   } else if (tpfjet_rechits[(*ith).id()].second.count(hitsAndFracs[iHit].second) == 0) {
0664                     tpfjet_twr_frac_.at(tpfjet_rechits[(*ith).id()].first) += hitsAndFracs[iHit].second;
0665                     if (cluster_dR <
0666                         tpfjet_cluster_dR_.at(tpfjet_twr_clusterind_.at(tpfjet_rechits[(*ith).id()].first))) {
0667                       tpfjet_twr_clusterind_.at(tpfjet_rechits[(*ith).id()].first) = cluster_ind;
0668                     }
0669                     tpfjet_rechits[(*ith).id()].second.insert(hitsAndFracs[iHit].second);
0670                   }
0671                 }  // Test if ieta,iphi matches
0672               }  // Loop over rechits
0673             }  // Loop over hits
0674           }  // Test if element is from HCAL
0675           else if (elements[iEle].type() == reco::PFBlockElement::HFHAD) {  // Element is HF
0676             for (edm::SortedCollection<HFRecHit, edm::StrictWeakOrdering<HFRecHit>>::const_iterator ith =
0677                      hfreco->begin();
0678                  ith != hfreco->end();
0679                  ++ith) {
0680               if ((*ith).id().depth() == 1)
0681                 continue;  // Remove long fibers
0682               auto thisCell = HFGeom->getGeometry((*ith).id().rawId());
0683               const CaloCellGeometry::CornersVec& cv = thisCell->getCorners();
0684 
0685               bool passMatch = false;
0686               if ((*it)->eta() < cv[0].eta() && (*it)->eta() > cv[2].eta()) {
0687                 if ((*it)->phi() < cv[0].phi() && (*it)->phi() > cv[2].phi())
0688                   passMatch = true;
0689                 else if (cv[0].phi() < cv[2].phi()) {
0690                   if ((*it)->phi() < cv[0].phi())
0691                     passMatch = true;
0692                   else if ((*it)->phi() > cv[2].phi())
0693                     passMatch = true;
0694                 }
0695               }
0696 
0697               if (passMatch) {
0698                 tpfjet_had_ntwrs_.at(tpfjet_had_n_ - 1)++;
0699                 tpfjet_twr_ieta_.push_back((*ith).id().ieta());
0700                 tpfjet_twr_iphi_.push_back((*ith).id().iphi());
0701                 tpfjet_twr_depth_.push_back((*ith).id().depth());
0702                 tpfjet_twr_subdet_.push_back((*ith).id().subdet());
0703                 tpfjet_twr_hade_.push_back((*ith).energy());
0704                 tpfjet_twr_frac_.push_back(1.0);
0705                 tpfjet_twr_hadind_.push_back(tpfjet_had_n_ - 1);
0706                 tpfjet_twr_elmttype_.push_back(1);
0707                 tpfjet_twr_clusterind_.push_back(-1);
0708                 tpfjet_twr_candtrackind_.push_back(-1);
0709                 float avgeta = (cv[0].eta() + cv[2].eta()) / 2.0;
0710                 float avgphi = (static_cast<double>(cv[0].phi()) + static_cast<double>(cv[2].phi())) / 2.0;
0711                 if (cv[0].phi() < cv[2].phi())
0712                   avgphi =
0713                       (2.0 * 3.141592653 + static_cast<double>(cv[0].phi()) + static_cast<double>(cv[2].phi())) / 2.0;
0714                 tpfjet_twr_dR_.push_back(deltaR(tpfjet_eta_, tpfjet_phi_, avgeta, avgphi));
0715                 ++tpfjet_ntwrs_;
0716                 HFHAD_E += (*ith).energy();
0717               }
0718             }
0719           } else if (elements[iEle].type() == reco::PFBlockElement::HFEM) {  // Element is HF
0720             for (edm::SortedCollection<HFRecHit, edm::StrictWeakOrdering<HFRecHit>>::const_iterator ith =
0721                      hfreco->begin();
0722                  ith != hfreco->end();
0723                  ++ith) {
0724               if ((*ith).id().depth() == 2)
0725                 continue;  // Remove short fibers
0726               auto thisCell = HFGeom->getGeometry((*ith).id().rawId());
0727               const CaloCellGeometry::CornersVec& cv = thisCell->getCorners();
0728 
0729               bool passMatch = false;
0730               if ((*it)->eta() < cv[0].eta() && (*it)->eta() > cv[2].eta()) {
0731                 if ((*it)->phi() < cv[0].phi() && (*it)->phi() > cv[2].phi())
0732                   passMatch = true;
0733                 else if (cv[0].phi() < cv[2].phi()) {
0734                   if ((*it)->phi() < cv[0].phi())
0735                     passMatch = true;
0736                   else if ((*it)->phi() > cv[2].phi())
0737                     passMatch = true;
0738                 }
0739               }
0740 
0741               if (passMatch) {
0742                 tpfjet_had_ntwrs_.at(tpfjet_had_n_ - 1)++;
0743                 tpfjet_twr_ieta_.push_back((*ith).id().ieta());
0744                 tpfjet_twr_iphi_.push_back((*ith).id().iphi());
0745                 tpfjet_twr_depth_.push_back((*ith).id().depth());
0746                 tpfjet_twr_subdet_.push_back((*ith).id().subdet());
0747                 tpfjet_twr_hade_.push_back((*ith).energy());
0748                 tpfjet_twr_frac_.push_back(1.0);
0749                 tpfjet_twr_hadind_.push_back(tpfjet_had_n_ - 1);
0750                 tpfjet_twr_elmttype_.push_back(2);
0751                 tpfjet_twr_clusterind_.push_back(-1);
0752                 tpfjet_twr_candtrackind_.push_back(-1);
0753                 float avgeta = (cv[0].eta() + cv[2].eta()) / 2.0;
0754                 float avgphi = (static_cast<double>(cv[0].phi()) + static_cast<double>(cv[2].phi())) / 2.0;
0755                 if (cv[0].phi() < cv[2].phi())
0756                   avgphi =
0757                       (2.0 * 3.141592653 + static_cast<double>(cv[0].phi()) + static_cast<double>(cv[2].phi())) / 2.0;
0758                 tpfjet_twr_dR_.push_back(deltaR(tpfjet_eta_, tpfjet_phi_, avgeta, avgphi));
0759                 ++tpfjet_ntwrs_;
0760                 HFEM_E += (*ith).energy();
0761               }
0762             }
0763           } else if (elements[iEle].type() == reco::PFBlockElement::HO) {  // Element is HO
0764             reco::PFClusterRef clusterref = elements[iEle].clusterRef();
0765             reco::PFCluster cluster = *clusterref;
0766             double cluster_dR = deltaR(tpfjet_eta_, tpfjet_phi_, cluster.eta(), cluster.phi());
0767             if (tpfjet_clusters.count(cluster_dR) == 0) {
0768               tpfjet_clusters[cluster_dR] = tpfjet_cluster_n_;
0769               tpfjet_cluster_eta_.push_back(cluster.eta());
0770               tpfjet_cluster_phi_.push_back(cluster.phi());
0771               tpfjet_cluster_dR_.push_back(cluster_dR);
0772               tpfjet_cluster_n_++;
0773             }
0774             int cluster_ind = tpfjet_clusters[cluster_dR];
0775 
0776             std::vector<std::pair<DetId, float>> hitsAndFracs = cluster.hitsAndFractions();
0777             int nHits = hitsAndFracs.size();
0778             for (int iHit = 0; iHit < nHits; iHit++) {
0779               int etaPhiPF = getEtaPhi(hitsAndFracs[iHit].first);
0780 
0781               for (edm::SortedCollection<HORecHit, edm::StrictWeakOrdering<HORecHit>>::const_iterator ith =
0782                        horeco->begin();
0783                    ith != horeco->end();
0784                    ++ith) {
0785                 int etaPhiRecHit = getEtaPhi((*ith).id());
0786                 if (etaPhiPF == etaPhiRecHit) {
0787                   tpfjet_had_ntwrs_.at(tpfjet_had_n_ - 1)++;
0788                   if (tpfjet_rechits.count((*ith).id()) == 0) {
0789                     tpfjet_twr_ieta_.push_back((*ith).id().ieta());
0790                     tpfjet_twr_iphi_.push_back((*ith).id().iphi());
0791                     tpfjet_twr_depth_.push_back((*ith).id().depth());
0792                     tpfjet_twr_subdet_.push_back((*ith).id().subdet());
0793                     if (hitsAndFracs[iHit].second > 0.05 && (*ith).energy() > 0.0)
0794                       twrietas[(*ith).id().ieta()]++;
0795                     tpfjet_twr_hade_.push_back((*ith).energy());
0796                     tpfjet_twr_frac_.push_back(hitsAndFracs[iHit].second);
0797                     tpfjet_rechits[(*ith).id()].second.insert(hitsAndFracs[iHit].second);
0798                     tpfjet_twr_hadind_.push_back(tpfjet_had_n_ - 1);
0799                     tpfjet_twr_elmttype_.push_back(3);
0800                     tpfjet_twr_clusterind_.push_back(cluster_ind);
0801                     if (hasTrack) {
0802                       tpfjet_twr_candtrackind_.push_back(tpfjet_ncandtracks_ - 1);
0803                     } else {
0804                       tpfjet_twr_candtrackind_.push_back(-1);
0805                     }
0806                     auto thisCell = HOGeom->getGeometry((*ith).id().rawId());
0807                     const CaloCellGeometry::CornersVec& cv = thisCell->getCorners();
0808                     float avgeta = (cv[0].eta() + cv[2].eta()) / 2.0;
0809                     float avgphi = (static_cast<double>(cv[0].phi()) + static_cast<double>(cv[2].phi())) / 2.0;
0810                     if (cv[0].phi() < cv[2].phi())
0811                       avgphi =
0812                           (2.0 * 3.141592653 + static_cast<double>(cv[0].phi()) + static_cast<double>(cv[2].phi())) /
0813                           2.0;
0814                     tpfjet_twr_dR_.push_back(deltaR(tpfjet_eta_, tpfjet_phi_, avgeta, avgphi));
0815                     tpfjet_rechits[(*ith).id()].first = tpfjet_ntwrs_;
0816                     ++tpfjet_ntwrs_;
0817                   } else if (tpfjet_rechits[(*ith).id()].second.count(hitsAndFracs[iHit].second) == 0) {
0818                     tpfjet_twr_frac_.at(tpfjet_rechits[(*ith).id()].first) += hitsAndFracs[iHit].second;
0819                     if (cluster_dR <
0820                         tpfjet_cluster_dR_.at(tpfjet_twr_clusterind_.at(tpfjet_rechits[(*ith).id()].first))) {
0821                       tpfjet_twr_clusterind_.at(tpfjet_rechits[(*ith).id()].first) = cluster_ind;
0822                     }
0823                     tpfjet_rechits[(*ith).id()].second.insert(hitsAndFracs[iHit].second);
0824                   }
0825                 }  // Test if ieta,iphi match
0826               }  // Loop over rechits
0827             }  // Loop over hits
0828           }  // Test if element is HO
0829         }  // Test for right element index
0830       }  // Loop over elements
0831     }  // Loop over elements in blocks
0832 
0833     switch (candidateType) {
0834       case reco::PFCandidate::h_HF:
0835         tpfjet_had_emf_.push_back(HFEM_E / (HFEM_E + HFHAD_E));
0836         break;
0837       case reco::PFCandidate::egamma_HF:
0838         tpfjet_had_emf_.push_back(-1);
0839         break;
0840       default:
0841         tpfjet_had_emf_.push_back(-1);
0842         break;
0843     }
0844   }  // Loop over PF constitutents
0845 
0846   int tag_had_EcalE = 0;
0847   int tag_had_rawHcalE = 0;
0848   for (int i = 0; i < tpfjet_had_n_; i++) {
0849     tag_had_EcalE += tpfjet_had_EcalE_[i];
0850     tag_had_rawHcalE += tpfjet_had_rawHcalE_[i];
0851   }
0852   tpfjet_EMfrac_ = 1.0 - tag_had_rawHcalE / (tag_had_rawHcalE + tag_had_EcalE + tpfjet_unkown_E_ + tpfjet_electron_E_ +
0853                                              tpfjet_muon_E_ + tpfjet_photon_E_);
0854   tpfjet_hadEcalEfrac_ = tag_had_EcalE / (tag_had_rawHcalE + tag_had_EcalE + tpfjet_unkown_E_ + tpfjet_electron_E_ +
0855                                           tpfjet_muon_E_ + tpfjet_photon_E_);
0856 
0857   // fill probe jet variables
0858   ppfjet_pt_ = pf_probe.jet()->pt();
0859   ppfjet_p_ = pf_probe.jet()->p();
0860   ppfjet_E_ = pf_probe.jet()->energy();
0861   ppfjet_eta_ = pf_probe.jet()->eta();
0862   ppfjet_phi_ = pf_probe.jet()->phi();
0863   ppfjet_scale_ = pf_probe.scale();
0864   ppfjet_area_ = pf_probe.jet()->jetArea();
0865   ppfjet_ntwrs_ = 0;
0866   ppfjet_ncandtracks_ = 0;
0867 
0868   ppfjet_jetID_ = 0;  // Not a loose, medium, or tight jet
0869   if (fabs(pf_probe.jet()->eta()) < 2.4) {
0870     if (pf_probe.jet()->chargedHadronEnergyFraction() > 0 && pf_probe.jet()->chargedMultiplicity() > 0 &&
0871         pf_probe.jet()->chargedEmEnergyFraction() < 0.99 &&
0872         (pf_probe.jet()->chargedMultiplicity() + pf_probe.jet()->neutralMultiplicity()) > 1) {
0873       if (pf_probe.jet()->neutralHadronEnergyFraction() < 0.9 && pf_probe.jet()->neutralEmEnergyFraction() < 0.9) {
0874         ppfjet_jetID_ = 3;  // Tight jet
0875       } else if (pf_probe.jet()->neutralHadronEnergyFraction() < 0.95 &&
0876                  pf_probe.jet()->neutralEmEnergyFraction() < 0.95) {
0877         ppfjet_jetID_ = 2;  // Medium jet
0878       } else if (pf_probe.jet()->neutralHadronEnergyFraction() < 0.99 &&
0879                  pf_probe.jet()->neutralEmEnergyFraction() < 0.99) {
0880         ppfjet_jetID_ = 1;  // Loose jet
0881       }
0882     }
0883   } else if ((pf_probe.jet()->chargedMultiplicity() + pf_probe.jet()->neutralMultiplicity()) > 1) {
0884     if (pf_probe.jet()->neutralHadronEnergyFraction() < 0.9 && pf_probe.jet()->neutralEmEnergyFraction() < 0.9) {
0885       ppfjet_jetID_ = 3;  // Tight jet
0886     } else if (pf_probe.jet()->neutralHadronEnergyFraction() < 0.95 &&
0887                pf_probe.jet()->neutralEmEnergyFraction() < 0.95) {
0888       ppfjet_jetID_ = 2;  // Medium jet
0889     } else if (pf_probe.jet()->neutralHadronEnergyFraction() < 0.99 &&
0890                pf_probe.jet()->neutralEmEnergyFraction() < 0.99) {
0891       ppfjet_jetID_ = 1;  // Loose jet
0892     }
0893   }
0894 
0895   // Get PF constituents and fill HCAL towers
0896   std::vector<reco::PFCandidatePtr> probeconst = pf_probe.jet()->getPFConstituents();
0897   for (std::vector<reco::PFCandidatePtr>::const_iterator it = probeconst.begin(); it != probeconst.end(); ++it) {
0898     bool hasTrack = false;
0899     reco::PFCandidate::ParticleType candidateType = (*it)->particleId();
0900     switch (candidateType) {
0901       case reco::PFCandidate::X:
0902         ppfjet_unkown_E_ += (*it)->energy();
0903         ppfjet_unkown_px_ += (*it)->px();
0904         ppfjet_unkown_py_ += (*it)->py();
0905         ppfjet_unkown_pz_ += (*it)->pz();
0906         ppfjet_unkown_EcalE_ += (*it)->ecalEnergy();
0907         ppfjet_unkown_n_++;
0908         continue;
0909       case reco::PFCandidate::h: {
0910         ppfjet_had_E_.push_back((*it)->energy());
0911         ppfjet_had_px_.push_back((*it)->px());
0912         ppfjet_had_py_.push_back((*it)->py());
0913         ppfjet_had_pz_.push_back((*it)->pz());
0914         ppfjet_had_EcalE_.push_back((*it)->ecalEnergy());
0915         ppfjet_had_rawHcalE_.push_back((*it)->rawHcalEnergy());
0916         ppfjet_had_id_.push_back(0);
0917         ppfjet_had_ntwrs_.push_back(0);
0918         ppfjet_had_n_++;
0919 
0920         reco::TrackRef trackRef = (*it)->trackRef();
0921         if (trackRef.isNonnull()) {
0922           reco::Track track = *trackRef;
0923           ppfjet_candtrack_px_.push_back(track.px());
0924           ppfjet_candtrack_py_.push_back(track.py());
0925           ppfjet_candtrack_pz_.push_back(track.pz());
0926           ppfjet_candtrack_EcalE_.push_back((*it)->ecalEnergy());
0927           ppfjet_had_candtrackind_.push_back(ppfjet_ncandtracks_);
0928           hasTrack = true;
0929           ppfjet_ncandtracks_++;
0930         } else {
0931           ppfjet_had_candtrackind_.push_back(-2);
0932         }
0933       } break;
0934       case reco::PFCandidate::e:
0935         ppfjet_electron_E_ += (*it)->energy();
0936         ppfjet_electron_px_ += (*it)->px();
0937         ppfjet_electron_py_ += (*it)->py();
0938         ppfjet_electron_pz_ += (*it)->pz();
0939         ppfjet_electron_EcalE_ += (*it)->ecalEnergy();
0940         ppfjet_electron_n_++;
0941         continue;
0942       case reco::PFCandidate::mu:
0943         ppfjet_muon_E_ += (*it)->energy();
0944         ppfjet_muon_px_ += (*it)->px();
0945         ppfjet_muon_py_ += (*it)->py();
0946         ppfjet_muon_pz_ += (*it)->pz();
0947         ppfjet_muon_EcalE_ += (*it)->ecalEnergy();
0948         ppfjet_muon_n_++;
0949         continue;
0950       case reco::PFCandidate::gamma:
0951         ppfjet_photon_E_ += (*it)->energy();
0952         ppfjet_photon_px_ += (*it)->px();
0953         ppfjet_photon_py_ += (*it)->py();
0954         ppfjet_photon_pz_ += (*it)->pz();
0955         ppfjet_photon_EcalE_ += (*it)->ecalEnergy();
0956         ppfjet_photon_n_++;
0957         continue;
0958       case reco::PFCandidate::h0: {
0959         ppfjet_had_E_.push_back((*it)->energy());
0960         ppfjet_had_px_.push_back((*it)->px());
0961         ppfjet_had_py_.push_back((*it)->py());
0962         ppfjet_had_pz_.push_back((*it)->pz());
0963         ppfjet_had_EcalE_.push_back((*it)->ecalEnergy());
0964         ppfjet_had_rawHcalE_.push_back((*it)->rawHcalEnergy());
0965         ppfjet_had_id_.push_back(1);
0966         ppfjet_had_candtrackind_.push_back(-1);
0967         ppfjet_had_ntwrs_.push_back(0);
0968         ppfjet_had_n_++;
0969         break;
0970       }
0971       case reco::PFCandidate::h_HF: {
0972         ppfjet_had_E_.push_back((*it)->energy());
0973         ppfjet_had_px_.push_back((*it)->px());
0974         ppfjet_had_py_.push_back((*it)->py());
0975         ppfjet_had_pz_.push_back((*it)->pz());
0976         ppfjet_had_EcalE_.push_back((*it)->ecalEnergy());
0977         ppfjet_had_rawHcalE_.push_back((*it)->rawHcalEnergy());
0978         ppfjet_had_id_.push_back(2);
0979         ppfjet_had_candtrackind_.push_back(-1);
0980         ppfjet_had_ntwrs_.push_back(0);
0981         ppfjet_had_n_++;
0982         break;
0983       }
0984       case reco::PFCandidate::egamma_HF: {
0985         ppfjet_had_E_.push_back((*it)->energy());
0986         ppfjet_had_px_.push_back((*it)->px());
0987         ppfjet_had_py_.push_back((*it)->py());
0988         ppfjet_had_pz_.push_back((*it)->pz());
0989         ppfjet_had_EcalE_.push_back((*it)->ecalEnergy());
0990         ppfjet_had_rawHcalE_.push_back((*it)->rawHcalEnergy());
0991         ppfjet_had_id_.push_back(3);
0992         ppfjet_had_candtrackind_.push_back(-1);
0993         ppfjet_had_ntwrs_.push_back(0);
0994         ppfjet_had_n_++;
0995         break;
0996       }
0997     }
0998 
0999     float HFHAD_E = 0;
1000     float HFEM_E = 0;
1001     int maxElement = (*it)->elementsInBlocks().size();
1002     for (int e = 0; e < maxElement; ++e) {
1003       // Get elements from block
1004       reco::PFBlockRef blockRef = (*it)->elementsInBlocks()[e].first;
1005       const edm::OwnVector<reco::PFBlockElement>& elements = blockRef->elements();
1006       for (unsigned iEle = 0; iEle < elements.size(); iEle++) {
1007         if (elements[iEle].index() == (*it)->elementsInBlocks()[e].second) {
1008           if (elements[iEle].type() == reco::PFBlockElement::HCAL) {  // Element is HB or HE
1009             // Get cluster and hits
1010             reco::PFClusterRef clusterref = elements[iEle].clusterRef();
1011             reco::PFCluster cluster = *clusterref;
1012             double cluster_dR = deltaR(ppfjet_eta_, ppfjet_phi_, cluster.eta(), cluster.phi());
1013             if (ppfjet_clusters.count(cluster_dR) == 0) {
1014               ppfjet_clusters[cluster_dR] = ppfjet_cluster_n_;
1015               ppfjet_cluster_eta_.push_back(cluster.eta());
1016               ppfjet_cluster_phi_.push_back(cluster.phi());
1017               ppfjet_cluster_dR_.push_back(cluster_dR);
1018               ppfjet_cluster_n_++;
1019             }
1020             int cluster_ind = ppfjet_clusters[cluster_dR];
1021             std::vector<std::pair<DetId, float>> hitsAndFracs = cluster.hitsAndFractions();
1022 
1023             // Run over hits and match
1024             int nHits = hitsAndFracs.size();
1025             for (int iHit = 0; iHit < nHits; iHit++) {
1026               int etaPhiPF = getEtaPhi(hitsAndFracs[iHit].first);
1027 
1028               for (edm::SortedCollection<HBHERecHit, edm::StrictWeakOrdering<HBHERecHit>>::const_iterator ith =
1029                        hbhereco->begin();
1030                    ith != hbhereco->end();
1031                    ++ith) {
1032                 int etaPhiRecHit = getEtaPhi((*ith).id());
1033                 if (etaPhiPF == etaPhiRecHit) {
1034                   ppfjet_had_ntwrs_.at(ppfjet_had_n_ - 1)++;
1035                   if (ppfjet_rechits.count((*ith).id()) == 0) {
1036                     ppfjet_twr_ieta_.push_back((*ith).id().ieta());
1037                     ppfjet_twr_iphi_.push_back((*ith).id().iphi());
1038                     ppfjet_twr_depth_.push_back((*ith).id().depth());
1039                     ppfjet_twr_subdet_.push_back((*ith).id().subdet());
1040                     ppfjet_twr_hade_.push_back((*ith).energy());
1041                     ppfjet_twr_frac_.push_back(hitsAndFracs[iHit].second);
1042                     ppfjet_rechits[(*ith).id()].second.insert(hitsAndFracs[iHit].second);
1043                     ppfjet_twr_hadind_.push_back(ppfjet_had_n_ - 1);
1044                     ppfjet_twr_elmttype_.push_back(0);
1045                     ppfjet_twr_clusterind_.push_back(cluster_ind);
1046                     if (hasTrack) {
1047                       ppfjet_twr_candtrackind_.push_back(ppfjet_ncandtracks_ - 1);
1048                     } else {
1049                       ppfjet_twr_candtrackind_.push_back(-1);
1050                     }
1051                     switch ((*ith).id().subdet()) {
1052                       case HcalSubdetector::HcalBarrel: {
1053                         CaloCellGeometry::CornersVec cv = HBGeom->getCorners((*ith).id());
1054                         float avgeta = (cv[0].eta() + cv[2].eta()) / 2.0;
1055                         float avgphi = (static_cast<double>(cv[0].phi()) + static_cast<double>(cv[2].phi())) / 2.0;
1056                         if (cv[0].phi() < cv[2].phi())
1057                           avgphi = (2.0 * 3.141592653 + static_cast<double>(cv[0].phi()) +
1058                                     static_cast<double>(cv[2].phi())) /
1059                                    2.0;
1060                         ppfjet_twr_dR_.push_back(deltaR(ppfjet_eta_, ppfjet_phi_, avgeta, avgphi));
1061                         break;
1062                       }
1063                       case HcalSubdetector::HcalEndcap: {
1064                         CaloCellGeometry::CornersVec cv = HEGeom->getCorners((*ith).id());
1065                         float avgeta = (cv[0].eta() + cv[2].eta()) / 2.0;
1066                         float avgphi = (static_cast<double>(cv[0].phi()) + static_cast<double>(cv[2].phi())) / 2.0;
1067                         if (cv[0].phi() < cv[2].phi())
1068                           avgphi = (2.0 * 3.141592653 + static_cast<double>(cv[0].phi()) +
1069                                     static_cast<double>(cv[2].phi())) /
1070                                    2.0;
1071                         ppfjet_twr_dR_.push_back(deltaR(ppfjet_eta_, ppfjet_phi_, avgeta, avgphi));
1072                         break;
1073                       }
1074                       default:
1075                         ppfjet_twr_dR_.push_back(-1);
1076                         break;
1077                     }
1078                     ppfjet_rechits[(*ith).id()].first = ppfjet_ntwrs_;
1079                     ++ppfjet_ntwrs_;
1080                   } else if (ppfjet_rechits[(*ith).id()].second.count(hitsAndFracs[iHit].second) == 0) {
1081                     ppfjet_twr_frac_.at(ppfjet_rechits[(*ith).id()].first) += hitsAndFracs[iHit].second;
1082                     if (cluster_dR <
1083                         ppfjet_cluster_dR_.at(ppfjet_twr_clusterind_.at(ppfjet_rechits[(*ith).id()].first))) {
1084                       ppfjet_twr_clusterind_.at(ppfjet_rechits[(*ith).id()].first) = cluster_ind;
1085                     }
1086                     ppfjet_rechits[(*ith).id()].second.insert(hitsAndFracs[iHit].second);
1087                   }
1088                 }  // Test if ieta,iphi matches
1089               }  // Loop over rechits
1090             }  // Loop over hits
1091           }  // Test if element is from HCAL
1092           else if (elements[iEle].type() == reco::PFBlockElement::HFHAD) {  // Element is HF
1093             for (edm::SortedCollection<HFRecHit, edm::StrictWeakOrdering<HFRecHit>>::const_iterator ith =
1094                      hfreco->begin();
1095                  ith != hfreco->end();
1096                  ++ith) {
1097               if ((*ith).id().depth() == 1)
1098                 continue;  // Remove long fibers
1099               auto thisCell = HFGeom->getGeometry((*ith).id().rawId());
1100               const CaloCellGeometry::CornersVec& cv = thisCell->getCorners();
1101 
1102               bool passMatch = false;
1103               if ((*it)->eta() < cv[0].eta() && (*it)->eta() > cv[2].eta()) {
1104                 if ((*it)->phi() < cv[0].phi() && (*it)->phi() > cv[2].phi())
1105                   passMatch = true;
1106                 else if (cv[0].phi() < cv[2].phi()) {
1107                   if ((*it)->phi() < cv[0].phi())
1108                     passMatch = true;
1109                   else if ((*it)->phi() > cv[2].phi())
1110                     passMatch = true;
1111                 }
1112               }
1113 
1114               if (passMatch) {
1115                 ppfjet_had_ntwrs_.at(ppfjet_had_n_ - 1)++;
1116                 ppfjet_twr_ieta_.push_back((*ith).id().ieta());
1117                 ppfjet_twr_iphi_.push_back((*ith).id().iphi());
1118                 ppfjet_twr_depth_.push_back((*ith).id().depth());
1119                 ppfjet_twr_subdet_.push_back((*ith).id().subdet());
1120                 ppfjet_twr_hade_.push_back((*ith).energy());
1121                 ppfjet_twr_frac_.push_back(1.0);
1122                 ppfjet_twr_hadind_.push_back(ppfjet_had_n_ - 1);
1123                 ppfjet_twr_elmttype_.push_back(1);
1124                 ppfjet_twr_clusterind_.push_back(-1);
1125                 ppfjet_twr_candtrackind_.push_back(-1);
1126                 float avgeta = (cv[0].eta() + cv[2].eta()) / 2.0;
1127                 float avgphi = (static_cast<double>(cv[0].phi()) + static_cast<double>(cv[2].phi())) / 2.0;
1128                 if (cv[0].phi() < cv[2].phi())
1129                   avgphi =
1130                       (2.0 * 3.141592653 + static_cast<double>(cv[0].phi()) + static_cast<double>(cv[2].phi())) / 2.0;
1131                 ppfjet_twr_dR_.push_back(deltaR(ppfjet_eta_, ppfjet_phi_, avgeta, avgphi));
1132                 ++ppfjet_ntwrs_;
1133                 HFHAD_E += (*ith).energy();
1134               }
1135             }
1136           } else if (elements[iEle].type() == reco::PFBlockElement::HFEM) {  // Element is HF
1137             for (edm::SortedCollection<HFRecHit, edm::StrictWeakOrdering<HFRecHit>>::const_iterator ith =
1138                      hfreco->begin();
1139                  ith != hfreco->end();
1140                  ++ith) {
1141               if ((*ith).id().depth() == 2)
1142                 continue;  // Remove short fibers
1143               auto thisCell = HFGeom->getGeometry((*ith).id().rawId());
1144               const CaloCellGeometry::CornersVec& cv = thisCell->getCorners();
1145 
1146               bool passMatch = false;
1147               if ((*it)->eta() < cv[0].eta() && (*it)->eta() > cv[2].eta()) {
1148                 if ((*it)->phi() < cv[0].phi() && (*it)->phi() > cv[2].phi())
1149                   passMatch = true;
1150                 else if (cv[0].phi() < cv[2].phi()) {
1151                   if ((*it)->phi() < cv[0].phi())
1152                     passMatch = true;
1153                   else if ((*it)->phi() > cv[2].phi())
1154                     passMatch = true;
1155                 }
1156               }
1157 
1158               if (passMatch) {
1159                 ppfjet_had_ntwrs_.at(ppfjet_had_n_ - 1)++;
1160                 ppfjet_twr_ieta_.push_back((*ith).id().ieta());
1161                 ppfjet_twr_iphi_.push_back((*ith).id().iphi());
1162                 ppfjet_twr_depth_.push_back((*ith).id().depth());
1163                 ppfjet_twr_subdet_.push_back((*ith).id().subdet());
1164                 ppfjet_twr_hade_.push_back((*ith).energy());
1165                 ppfjet_twr_frac_.push_back(1.0);
1166                 ppfjet_twr_hadind_.push_back(ppfjet_had_n_ - 1);
1167                 ppfjet_twr_elmttype_.push_back(2);
1168                 ppfjet_twr_clusterind_.push_back(-1);
1169                 ppfjet_twr_candtrackind_.push_back(-1);
1170                 float avgeta = (cv[0].eta() + cv[2].eta()) / 2.0;
1171                 float avgphi = (static_cast<double>(cv[0].phi()) + static_cast<double>(cv[2].phi())) / 2.0;
1172                 if (cv[0].phi() < cv[2].phi())
1173                   avgphi =
1174                       (2.0 * 3.141592653 + static_cast<double>(cv[0].phi()) + static_cast<double>(cv[2].phi())) / 2.0;
1175                 ppfjet_twr_dR_.push_back(deltaR(ppfjet_eta_, ppfjet_phi_, avgeta, avgphi));
1176                 ++ppfjet_ntwrs_;
1177                 HFEM_E += (*ith).energy();
1178               }
1179             }
1180           } else if (elements[iEle].type() == reco::PFBlockElement::HO) {  // Element is HO
1181             reco::PFClusterRef clusterref = elements[iEle].clusterRef();
1182             reco::PFCluster cluster = *clusterref;
1183             double cluster_dR = deltaR(ppfjet_eta_, ppfjet_phi_, cluster.eta(), cluster.phi());
1184             if (ppfjet_clusters.count(cluster_dR) == 0) {
1185               ppfjet_clusters[cluster_dR] = ppfjet_cluster_n_;
1186               ppfjet_cluster_eta_.push_back(cluster.eta());
1187               ppfjet_cluster_phi_.push_back(cluster.phi());
1188               ppfjet_cluster_dR_.push_back(cluster_dR);
1189               ppfjet_cluster_n_++;
1190             }
1191             int cluster_ind = ppfjet_clusters[cluster_dR];
1192 
1193             std::vector<std::pair<DetId, float>> hitsAndFracs = cluster.hitsAndFractions();
1194             int nHits = hitsAndFracs.size();
1195             for (int iHit = 0; iHit < nHits; iHit++) {
1196               int etaPhiPF = getEtaPhi(hitsAndFracs[iHit].first);
1197 
1198               for (edm::SortedCollection<HORecHit, edm::StrictWeakOrdering<HORecHit>>::const_iterator ith =
1199                        horeco->begin();
1200                    ith != horeco->end();
1201                    ++ith) {
1202                 int etaPhiRecHit = getEtaPhi((*ith).id());
1203                 if (etaPhiPF == etaPhiRecHit) {
1204                   ppfjet_had_ntwrs_.at(ppfjet_had_n_ - 1)++;
1205                   if (ppfjet_rechits.count((*ith).id()) == 0) {
1206                     ppfjet_twr_ieta_.push_back((*ith).id().ieta());
1207                     ppfjet_twr_iphi_.push_back((*ith).id().iphi());
1208                     ppfjet_twr_depth_.push_back((*ith).id().depth());
1209                     ppfjet_twr_subdet_.push_back((*ith).id().subdet());
1210                     ppfjet_twr_hade_.push_back((*ith).energy());
1211                     ppfjet_twr_frac_.push_back(hitsAndFracs[iHit].second);
1212                     ppfjet_rechits[(*ith).id()].second.insert(hitsAndFracs[iHit].second);
1213                     ppfjet_twr_hadind_.push_back(ppfjet_had_n_ - 1);
1214                     ppfjet_twr_elmttype_.push_back(3);
1215                     ppfjet_twr_clusterind_.push_back(cluster_ind);
1216                     if (hasTrack) {
1217                       ppfjet_twr_candtrackind_.push_back(ppfjet_ncandtracks_ - 1);
1218                     } else {
1219                       ppfjet_twr_candtrackind_.push_back(-1);
1220                     }
1221                     auto thisCell = HOGeom->getGeometry((*ith).id().rawId());
1222                     const CaloCellGeometry::CornersVec& cv = thisCell->getCorners();
1223                     float avgeta = (cv[0].eta() + cv[2].eta()) / 2.0;
1224                     float avgphi = (static_cast<double>(cv[0].phi()) + static_cast<double>(cv[2].phi())) / 2.0;
1225                     if (cv[0].phi() < cv[2].phi())
1226                       avgphi =
1227                           (2.0 * 3.141592653 + static_cast<double>(cv[0].phi()) + static_cast<double>(cv[2].phi())) /
1228                           2.0;
1229                     ppfjet_twr_dR_.push_back(deltaR(ppfjet_eta_, ppfjet_phi_, avgeta, avgphi));
1230                     ppfjet_rechits[(*ith).id()].first = ppfjet_ntwrs_;
1231                     ++ppfjet_ntwrs_;
1232                   } else if (ppfjet_rechits[(*ith).id()].second.count(hitsAndFracs[iHit].second) == 0) {
1233                     ppfjet_twr_frac_.at(ppfjet_rechits[(*ith).id()].first) += hitsAndFracs[iHit].second;
1234                     if (cluster_dR <
1235                         ppfjet_cluster_dR_.at(ppfjet_twr_clusterind_.at(ppfjet_rechits[(*ith).id()].first))) {
1236                       ppfjet_twr_clusterind_.at(ppfjet_rechits[(*ith).id()].first) = cluster_ind;
1237                     }
1238                     ppfjet_rechits[(*ith).id()].second.insert(hitsAndFracs[iHit].second);
1239                   }
1240                 }  // Test if ieta,iphi match
1241               }  // Loop over rechits
1242             }  // Loop over hits
1243           }  // Test if element is from HO
1244         }  // Test for right element index
1245       }  // Loop over elements
1246     }  // Loop over elements in blocks
1247     switch (candidateType) {
1248       case reco::PFCandidate::h_HF:
1249         ppfjet_had_emf_.push_back(HFEM_E / (HFEM_E + HFHAD_E));
1250         break;
1251       case reco::PFCandidate::egamma_HF:
1252         ppfjet_had_emf_.push_back(-1);
1253         break;
1254       default:
1255         ppfjet_had_emf_.push_back(-1);
1256         break;
1257     }
1258   }  // Loop over PF constitutents
1259 
1260   int probe_had_EcalE = 0;
1261   int probe_had_rawHcalE = 0;
1262   for (int i = 0; i < ppfjet_had_n_; i++) {
1263     probe_had_EcalE += ppfjet_had_EcalE_[i];
1264     probe_had_rawHcalE += ppfjet_had_rawHcalE_[i];
1265   }
1266   ppfjet_EMfrac_ = 1.0 - probe_had_rawHcalE / (probe_had_rawHcalE + probe_had_EcalE + ppfjet_unkown_E_ +
1267                                                ppfjet_electron_E_ + ppfjet_muon_E_ + ppfjet_photon_E_);
1268   ppfjet_hadEcalEfrac_ = probe_had_EcalE / (probe_had_rawHcalE + probe_had_EcalE + ppfjet_unkown_E_ +
1269                                             ppfjet_electron_E_ + ppfjet_muon_E_ + ppfjet_photon_E_);
1270 
1271   // fill dijet variables
1272   pf_dijet_deta_ = std::fabs(std::fabs(pf_tag.jet()->eta()) - std::fabs(pf_probe.jet()->eta()));
1273   pf_dijet_dphi_ = pf_tag.jet()->phi() - pf_probe.jet()->phi();
1274   if (pf_dijet_dphi_ > 3.1415)
1275     pf_dijet_dphi_ = 6.2832 - pf_dijet_dphi_;
1276   pf_dijet_balance_ = (tpfjet_pt_ - ppfjet_pt_) / (tpfjet_pt_ + ppfjet_pt_);
1277 
1278   tree_->Fill();
1279 
1280   return;
1281 }
1282 
1283 // ------------ method called once each job just before starting event loop  ------------
1284 void DiJetAnalyzer::beginJob() {
1285   // book histograms
1286 
1287   edm::Service<TFileService> fs;
1288   h_PassSelPF_ = fs->make<TH1D>("h_PassSelectionPF", "Selection Pass Failures PFJets", 200, -0.5, 199.5);
1289 
1290   tree_ = fs->make<TTree>("dijettree", "tree for dijet balancing");
1291 
1292   tree_->Branch("tpfjet_pt", &tpfjet_pt_, "tpfjet_pt/F");
1293   tree_->Branch("tpfjet_p", &tpfjet_p_, "tpfjet_p/F");
1294   tree_->Branch("tpfjet_E", &tpfjet_E_, "tpfjet_E/F");
1295   tree_->Branch("tpfjet_eta", &tpfjet_eta_, "tpfjet_eta/F");
1296   tree_->Branch("tpfjet_phi", &tpfjet_phi_, "tpfjet_phi/F");
1297   tree_->Branch("tpfjet_EMfrac", &tpfjet_EMfrac_, "tpfjet_EMfrac/F");
1298   tree_->Branch("tpfjet_hadEcalEfrac", &tpfjet_hadEcalEfrac_, "tpfjet_hadEcalEfrac/F");
1299   tree_->Branch("tpfjet_scale", &tpfjet_scale_, "tpfjet_scale/F");
1300   tree_->Branch("tpfjet_area", &tpfjet_area_, "tpfjet_area/F");
1301   tree_->Branch("tpfjet_jetID", &tpfjet_jetID_, "tpfjet_jetID/I");
1302   tree_->Branch("tpfjet_unkown_E", &tpfjet_unkown_E_, "tpfjet_unkown_E/F");
1303   tree_->Branch("tpfjet_electron_E", &tpfjet_electron_E_, "tpfjet_electron_E/F");
1304   tree_->Branch("tpfjet_muon_E", &tpfjet_muon_E_, "tpfjet_muon_E/F");
1305   tree_->Branch("tpfjet_photon_E", &tpfjet_photon_E_, "tpfjet_photon_E/F");
1306   tree_->Branch("tpfjet_unkown_px", &tpfjet_unkown_px_, "tpfjet_unkown_px/F");
1307   tree_->Branch("tpfjet_electron_px", &tpfjet_electron_px_, "tpfjet_electron_px/F");
1308   tree_->Branch("tpfjet_muon_px", &tpfjet_muon_px_, "tpfjet_muon_px/F");
1309   tree_->Branch("tpfjet_photon_px", &tpfjet_photon_px_, "tpfjet_photon_px/F");
1310   tree_->Branch("tpfjet_unkown_py", &tpfjet_unkown_py_, "tpfjet_unkown_py/F");
1311   tree_->Branch("tpfjet_electron_py", &tpfjet_electron_py_, "tpfjet_electron_py/F");
1312   tree_->Branch("tpfjet_muon_py", &tpfjet_muon_py_, "tpfjet_muon_py/F");
1313   tree_->Branch("tpfjet_photon_py", &tpfjet_photon_py_, "tpfjet_photon_py/F");
1314   tree_->Branch("tpfjet_unkown_pz", &tpfjet_unkown_pz_, "tpfjet_unkown_pz/F");
1315   tree_->Branch("tpfjet_electron_pz", &tpfjet_electron_pz_, "tpfjet_electron_pz/F");
1316   tree_->Branch("tpfjet_muon_pz", &tpfjet_muon_pz_, "tpfjet_muon_pz/F");
1317   tree_->Branch("tpfjet_photon_pz", &tpfjet_photon_pz_, "tpfjet_photon_pz/F");
1318   tree_->Branch("tpfjet_unkown_EcalE", &tpfjet_unkown_EcalE_, "tpfjet_unkown_EcalE/F");
1319   tree_->Branch("tpfjet_electron_EcalE", &tpfjet_electron_EcalE_, "tpfjet_electron_EcalE/F");
1320   tree_->Branch("tpfjet_muon_EcalE", &tpfjet_muon_EcalE_, "tpfjet_muon_EcalE/F");
1321   tree_->Branch("tpfjet_photon_EcalE", &tpfjet_photon_EcalE_, "tpfjet_photon_EcalE/F");
1322   tree_->Branch("tpfjet_unkown_n", &tpfjet_unkown_n_, "tpfjet_unkown_n/I");
1323   tree_->Branch("tpfjet_electron_n", &tpfjet_electron_n_, "tpfjet_electron_n/I");
1324   tree_->Branch("tpfjet_muon_n", &tpfjet_muon_n_, "tpfjet_muon_n/I");
1325   tree_->Branch("tpfjet_photon_n", &tpfjet_photon_n_, "tpfjet_photon_n/I");
1326   tree_->Branch("tpfjet_had_n", &tpfjet_had_n_, "tpfjet_had_n/I");
1327   tree_->Branch("tpfjet_had_E", &tpfjet_had_E_);
1328   tree_->Branch("tpfjet_had_px", &tpfjet_had_px_);
1329   tree_->Branch("tpfjet_had_py", &tpfjet_had_py_);
1330   tree_->Branch("tpfjet_had_pz", &tpfjet_had_pz_);
1331   tree_->Branch("tpfjet_had_EcalE", &tpfjet_had_EcalE_);
1332   tree_->Branch("tpfjet_had_rawHcalE", &tpfjet_had_rawHcalE_);
1333   tree_->Branch("tpfjet_had_emf", &tpfjet_had_emf_);
1334   tree_->Branch("tpfjet_had_id", &tpfjet_had_id_);
1335   tree_->Branch("tpfjet_had_candtrackind", &tpfjet_had_candtrackind_);
1336   tree_->Branch("tpfjet_had_ntwrs", &tpfjet_had_ntwrs_);
1337   tree_->Branch("tpfjet_ntwrs", &tpfjet_ntwrs_, "tpfjet_ntwrs/I");
1338   tree_->Branch("tpfjet_twr_ieta", &tpfjet_twr_ieta_);
1339   tree_->Branch("tpfjet_twr_iphi", &tpfjet_twr_iphi_);
1340   tree_->Branch("tpfjet_twr_depth", &tpfjet_twr_depth_);
1341   tree_->Branch("tpfjet_twr_subdet", &tpfjet_twr_subdet_);
1342   tree_->Branch("tpfjet_twr_hade", &tpfjet_twr_hade_);
1343   tree_->Branch("tpfjet_twr_frac", &tpfjet_twr_frac_);
1344   tree_->Branch("tpfjet_twr_candtrackind", &tpfjet_twr_candtrackind_);
1345   tree_->Branch("tpfjet_twr_hadind", &tpfjet_twr_hadind_);
1346   tree_->Branch("tpfjet_twr_elmttype", &tpfjet_twr_elmttype_);
1347   tree_->Branch("tpfjet_twr_dR", &tpfjet_twr_dR_);
1348   tree_->Branch("tpfjet_twr_clusterind", &tpfjet_twr_clusterind_);
1349   tree_->Branch("tpfjet_cluster_n", &tpfjet_cluster_n_, "tpfjet_cluster_n/I");
1350   tree_->Branch("tpfjet_cluster_eta", &tpfjet_cluster_eta_);
1351   tree_->Branch("tpfjet_cluster_phi", &tpfjet_cluster_phi_);
1352   tree_->Branch("tpfjet_cluster_dR", &tpfjet_cluster_dR_);
1353   tree_->Branch("tpfjet_ncandtracks", &tpfjet_ncandtracks_, "tpfjet_ncandtracks/I");
1354   tree_->Branch("tpfjet_candtrack_px", &tpfjet_candtrack_px_);
1355   tree_->Branch("tpfjet_candtrack_py", &tpfjet_candtrack_py_);
1356   tree_->Branch("tpfjet_candtrack_pz", &tpfjet_candtrack_pz_);
1357   tree_->Branch("tpfjet_candtrack_EcalE", &tpfjet_candtrack_EcalE_);
1358   tree_->Branch("ppfjet_pt", &ppfjet_pt_, "ppfjet_pt/F");
1359   tree_->Branch("ppfjet_p", &ppfjet_p_, "ppfjet_p/F");
1360   tree_->Branch("ppfjet_E", &ppfjet_E_, "ppfjet_E/F");
1361   tree_->Branch("ppfjet_eta", &ppfjet_eta_, "ppfjet_eta/F");
1362   tree_->Branch("ppfjet_phi", &ppfjet_phi_, "ppfjet_phi/F");
1363   tree_->Branch("ppfjet_EMfrac", &ppfjet_EMfrac_, "ppfjet_EMfrac/F");
1364   tree_->Branch("ppfjet_hadEcalEfrac", &ppfjet_hadEcalEfrac_, "ppfjet_hadEcalEfrac/F");
1365   tree_->Branch("ppfjet_scale", &ppfjet_scale_, "ppfjet_scale/F");
1366   tree_->Branch("ppfjet_area", &ppfjet_area_, "ppfjet_area/F");
1367   tree_->Branch("ppfjet_jetID", &ppfjet_jetID_, "ppfjet_jetID/I");
1368   tree_->Branch("ppfjet_unkown_E", &ppfjet_unkown_E_, "ppfjet_unkown_E/F");
1369   tree_->Branch("ppfjet_electron_E", &ppfjet_electron_E_, "ppfjet_electron_E/F");
1370   tree_->Branch("ppfjet_muon_E", &ppfjet_muon_E_, "ppfjet_muon_E/F");
1371   tree_->Branch("ppfjet_photon_E", &ppfjet_photon_E_, "ppfjet_photon_E/F");
1372   tree_->Branch("ppfjet_unkown_px", &ppfjet_unkown_px_, "ppfjet_unkown_px/F");
1373   tree_->Branch("ppfjet_electron_px", &ppfjet_electron_px_, "ppfjet_electron_px/F");
1374   tree_->Branch("ppfjet_muon_px", &ppfjet_muon_px_, "ppfjet_muon_px/F");
1375   tree_->Branch("ppfjet_photon_px", &ppfjet_photon_px_, "ppfjet_photon_px/F");
1376   tree_->Branch("ppfjet_unkown_py", &ppfjet_unkown_py_, "ppfjet_unkown_py/F");
1377   tree_->Branch("ppfjet_electron_py", &ppfjet_electron_py_, "ppfjet_electron_py/F");
1378   tree_->Branch("ppfjet_muon_py", &ppfjet_muon_py_, "ppfjet_muon_py/F");
1379   tree_->Branch("ppfjet_photon_py", &ppfjet_photon_py_, "ppfjet_photon_py/F");
1380   tree_->Branch("ppfjet_unkown_pz", &ppfjet_unkown_pz_, "ppfjet_unkown_pz/F");
1381   tree_->Branch("ppfjet_electron_pz", &ppfjet_electron_pz_, "ppfjet_electron_pz/F");
1382   tree_->Branch("ppfjet_muon_pz", &ppfjet_muon_pz_, "ppfjet_muon_pz/F");
1383   tree_->Branch("ppfjet_photon_pz", &ppfjet_photon_pz_, "ppfjet_photon_pz/F");
1384   tree_->Branch("ppfjet_unkown_EcalE", &ppfjet_unkown_EcalE_, "ppfjet_unkown_EcalE/F");
1385   tree_->Branch("ppfjet_electron_EcalE", &ppfjet_electron_EcalE_, "ppfjet_electron_EcalE/F");
1386   tree_->Branch("ppfjet_muon_EcalE", &ppfjet_muon_EcalE_, "ppfjet_muon_EcalE/F");
1387   tree_->Branch("ppfjet_photon_EcalE", &ppfjet_photon_EcalE_, "ppfjet_photon_EcalE/F");
1388   tree_->Branch("ppfjet_unkown_n", &ppfjet_unkown_n_, "ppfjet_unkown_n/I");
1389   tree_->Branch("ppfjet_electron_n", &ppfjet_electron_n_, "ppfjet_electron_n/I");
1390   tree_->Branch("ppfjet_muon_n", &ppfjet_muon_n_, "ppfjet_muon_n/I");
1391   tree_->Branch("ppfjet_photon_n", &ppfjet_photon_n_, "ppfjet_photon_n/I");
1392   tree_->Branch("ppfjet_had_n", &ppfjet_had_n_, "ppfjet_had_n/I");
1393   tree_->Branch("ppfjet_had_E", &ppfjet_had_E_);
1394   tree_->Branch("ppfjet_had_px", &ppfjet_had_px_);
1395   tree_->Branch("ppfjet_had_py", &ppfjet_had_py_);
1396   tree_->Branch("ppfjet_had_pz", &ppfjet_had_pz_);
1397   tree_->Branch("ppfjet_had_EcalE", &ppfjet_had_EcalE_);
1398   tree_->Branch("ppfjet_had_rawHcalE", &ppfjet_had_rawHcalE_);
1399   tree_->Branch("ppfjet_had_emf", &ppfjet_had_emf_);
1400   tree_->Branch("ppfjet_had_id", &ppfjet_had_id_);
1401   tree_->Branch("ppfjet_had_candtrackind", &ppfjet_had_candtrackind_);
1402   tree_->Branch("ppfjet_had_ntwrs", &ppfjet_had_ntwrs_);
1403   tree_->Branch("ppfjet_ntwrs", &ppfjet_ntwrs_, "ppfjet_ntwrs/I");
1404   tree_->Branch("ppfjet_twr_ieta", &ppfjet_twr_ieta_);
1405   tree_->Branch("ppfjet_twr_iphi", &ppfjet_twr_iphi_);
1406   tree_->Branch("ppfjet_twr_depth", &ppfjet_twr_depth_);
1407   tree_->Branch("ppfjet_twr_subdet", &ppfjet_twr_subdet_);
1408   tree_->Branch("ppfjet_twr_hade", &ppfjet_twr_hade_);
1409   tree_->Branch("ppfjet_twr_frac", &ppfjet_twr_frac_);
1410   tree_->Branch("ppfjet_twr_candtrackind", &ppfjet_twr_candtrackind_);
1411   tree_->Branch("ppfjet_twr_hadind", &ppfjet_twr_hadind_);
1412   tree_->Branch("ppfjet_twr_elmttype", &ppfjet_twr_elmttype_);
1413   tree_->Branch("ppfjet_twr_dR", &ppfjet_twr_dR_);
1414   tree_->Branch("ppfjet_twr_clusterind", &ppfjet_twr_clusterind_);
1415   tree_->Branch("ppfjet_cluster_n", &ppfjet_cluster_n_, "ppfjet_cluster_n/I");
1416   tree_->Branch("ppfjet_cluster_eta", &ppfjet_cluster_eta_);
1417   tree_->Branch("ppfjet_cluster_phi", &ppfjet_cluster_phi_);
1418   tree_->Branch("ppfjet_cluster_dR", &ppfjet_cluster_dR_);
1419   tree_->Branch("ppfjet_ncandtracks", &ppfjet_ncandtracks_, "ppfjet_ncandtracks/I");
1420   tree_->Branch("ppfjet_candtrack_px", &ppfjet_candtrack_px_);
1421   tree_->Branch("ppfjet_candtrack_py", &ppfjet_candtrack_py_);
1422   tree_->Branch("ppfjet_candtrack_pz", &ppfjet_candtrack_pz_);
1423   tree_->Branch("ppfjet_candtrack_EcalE", &ppfjet_candtrack_EcalE_);
1424   tree_->Branch("pf_dijet_deta", &pf_dijet_deta_, "pf_dijet_deta/F");
1425   tree_->Branch("pf_dijet_dphi", &pf_dijet_dphi_, "pf_dijet_dphi/F");
1426   tree_->Branch("pf_dijet_balance", &pf_dijet_balance_, "pf_dijet_balance/F");
1427   tree_->Branch("pf_thirdjet_px", &pf_thirdjet_px_, "pf_thirdjet_px/F");
1428   tree_->Branch("pf_thirdjet_py", &pf_thirdjet_py_, "pf_thirdjet_py/F");
1429   tree_->Branch("pf_realthirdjet_px", &pf_realthirdjet_px_, "pf_realthirdjet_px/F");
1430   tree_->Branch("pf_realthirdjet_py", &pf_realthirdjet_py_, "pf_realthirdjet_py/F");
1431   tree_->Branch("pf_realthirdjet_scale", &pf_realthirdjet_scale_, "pf_realthirdjet_scale/F");
1432   tree_->Branch("pf_Run", &pf_Run_, "pf_Run/I");
1433   tree_->Branch("pf_Lumi", &pf_Lumi_, "pf_Lumi/I");
1434   tree_->Branch("pf_Event", &pf_Event_, "pf_Event/I");
1435   tree_->Branch("pf_NPV", &pf_NPV_, "pf_NPV/I");
1436 
1437   return;
1438 }
1439 
1440 // ------------ method called once each job just after ending the event loop  ------------
1441 void DiJetAnalyzer::endJob() {}
1442 
1443 // helper function
1444 
1445 double DiJetAnalyzer::deltaR(const reco::Jet* j1, const reco::Jet* j2) {
1446   double deta = j1->eta() - j2->eta();
1447   double dphi = std::fabs(j1->phi() - j2->phi());
1448   if (dphi > 3.1415927)
1449     dphi = 2 * 3.1415927 - dphi;
1450   return std::sqrt(deta * deta + dphi * dphi);
1451 }
1452 
1453 double DiJetAnalyzer::deltaR(const double eta1, const double phi1, const double eta2, const double phi2) {
1454   double deta = eta1 - eta2;
1455   double dphi = std::fabs(phi1 - phi2);
1456   if (dphi > 3.1415927)
1457     dphi = 2 * 3.1415927 - dphi;
1458   return std::sqrt(deta * deta + dphi * dphi);
1459 }
1460 
1461 int DiJetAnalyzer::getEtaPhi(const DetId id) {
1462   return id.rawId() & 0x3FFF;  // Get 14 least-significant digits
1463 }
1464 
1465 int DiJetAnalyzer::getEtaPhi(const HcalDetId id) {
1466   return id.rawId() & 0x3FFF;  // Get 14 least-significant digits
1467 }
1468 
1469 //define this as a plug-in
1470 DEFINE_FWK_MODULE(DiJetAnalyzer);