Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:31

0001 //
0002 // Jet Analyzer class for heavy ion jets. for DQM jet analysis monitoring
0003 // For CMSSW_7_5_X, especially reading background subtracted jets - jet matching 1 - 1
0004 // author: Raghav Kunnawalkam Elayavalli,
0005 //         April 6th 2015
0006 //         Rutgers University, email: raghav.k.e at CERN dot CH
0007 //
0008 
0009 #include "DQMOffline/JetMET/interface/JetAnalyzer_HeavyIons_matching.h"
0010 
0011 using namespace edm;
0012 using namespace reco;
0013 using namespace std;
0014 
0015 // declare the constructors:
0016 
0017 JetAnalyzer_HeavyIons_matching::JetAnalyzer_HeavyIons_matching(const edm::ParameterSet& iConfig)
0018     : mInputJet1Collection(iConfig.getParameter<edm::InputTag>("src_Jet1")),
0019       mInputJet2Collection(iConfig.getParameter<edm::InputTag>("src_Jet2")),
0020       JetType1(iConfig.getUntrackedParameter<std::string>("Jet1")),
0021       JetType2(iConfig.getUntrackedParameter<std::string>("Jet2")),
0022       mRecoJetPtThreshold(iConfig.getParameter<double>("recoJetPtThreshold")),
0023       mRecoDelRMatch(iConfig.getParameter<double>("recoDelRMatch")),
0024       mRecoJetEtaCut(iConfig.getParameter<double>("recoJetEtaCut")) {
0025   std::string inputCollectionLabelJet1(mInputJet1Collection.label());
0026   std::string inputCollectionLabelJet2(mInputJet2Collection.label());
0027 
0028   //consumes
0029 
0030   if (std::string("VsCalo") == JetType1)
0031     caloJet1Token_ = consumes<reco::CaloJetCollection>(mInputJet1Collection);
0032   if (std::string("VsPF") == JetType1)
0033     pfJetsToken_ = consumes<reco::PFJetCollection>(mInputJet1Collection);
0034   if (std::string("PuCalo") == JetType1)
0035     caloJet2Token_ = consumes<reco::CaloJetCollection>(mInputJet1Collection);
0036   if (std::string("PuPF") == JetType1)
0037     basicJetsToken_ = consumes<reco::BasicJetCollection>(mInputJet1Collection);
0038 
0039   if (std::string("VsCalo") == JetType2)
0040     caloJet1Token_ = consumes<reco::CaloJetCollection>(mInputJet2Collection);
0041   if (std::string("VsPF") == JetType2)
0042     pfJetsToken_ = consumes<reco::PFJetCollection>(mInputJet2Collection);
0043   if (std::string("PuCalo") == JetType2)
0044     caloJet2Token_ = consumes<reco::CaloJetCollection>(mInputJet2Collection);
0045   if (std::string("PuPF") == JetType2)
0046     basicJetsToken_ = consumes<reco::BasicJetCollection>(mInputJet2Collection);
0047 
0048   // initialize the Jet matching histograms
0049 
0050   mpT_ratio_Jet1Jet2 = nullptr;
0051   mpT_Jet1_matched = nullptr;
0052   mpT_Jet2_matched = nullptr;
0053   mpT_Jet1_unmatched = nullptr;
0054   mpT_Jet2_unmatched = nullptr;
0055 
0056   // we need to add histograms which will hold the hadronic and electromagnetic energy content for the unmatched Jets.
0057   if (std::string("VsCalo") == JetType1 || std::string("PuCalo") == JetType1) {
0058     mHadEnergy_Jet1_unmatched = nullptr;
0059     mEmEnergy_Jet1_unmatched = nullptr;
0060   }
0061   if (std::string("VsCalo") == JetType2 || std::string("PuCalo") == JetType2) {
0062     mHadEnergy_Jet2_unmatched = nullptr;
0063     mEmEnergy_Jet2_unmatched = nullptr;
0064   }
0065 
0066   if (std::string("VsPF") == JetType1) {
0067     mChargedHadronEnergy_Jet1_unmatched = nullptr;
0068     mNeutralHadronEnergy_Jet1_unmatched = nullptr;
0069     mChargedEmEnergy_Jet1_unmatched = nullptr;
0070     mNeutralEmEnergy_Jet1_unmatched = nullptr;
0071     mChargedMuEnergy_Jet1_unmatched = nullptr;
0072     mChargedHadEnergyFraction_Jet1_unmatched = nullptr;
0073     mNeutralHadEnergyFraction_Jet1_unmatched = nullptr;
0074     mPhotonEnergyFraction_Jet1_unmatched = nullptr;
0075     mElectronEnergyFraction_Jet1_unmatched = nullptr;
0076     mMuonEnergyFraction_Jet1_unmatched = nullptr;
0077   }
0078 
0079   if (std::string("VsPF") == JetType2) {
0080     mChargedHadronEnergy_Jet2_unmatched = nullptr;
0081     mNeutralHadronEnergy_Jet2_unmatched = nullptr;
0082     mChargedEmEnergy_Jet2_unmatched = nullptr;
0083     mNeutralEmEnergy_Jet2_unmatched = nullptr;
0084     mChargedMuEnergy_Jet2_unmatched = nullptr;
0085     mChargedHadEnergyFraction_Jet2_unmatched = nullptr;
0086     mNeutralHadEnergyFraction_Jet2_unmatched = nullptr;
0087     mPhotonEnergyFraction_Jet2_unmatched = nullptr;
0088     mElectronEnergyFraction_Jet2_unmatched = nullptr;
0089     mMuonEnergyFraction_Jet2_unmatched = nullptr;
0090   }
0091 }
0092 
0093 void JetAnalyzer_HeavyIons_matching::bookHistograms(DQMStore::IBooker& ibooker,
0094                                                     edm::Run const& iRun,
0095                                                     edm::EventSetup const&) {
0096   ibooker.setCurrentFolder("JetMET/HIJetValidation/" + mInputJet1Collection.label() + "_DeltaRMatched_" +
0097                            mInputJet2Collection.label());
0098 
0099   mpT_ratio_Jet1Jet2 = ibooker.book1D(
0100       "Ratio_Jet1pT_vs_Jet2pT", Form(";Matched %s Jet pT / %s Jet pT;Counts", "JetType1", "JetType2"), 100, 0, 10);
0101   mpT_Jet1_matched =
0102       ibooker.book1D("Jet1_matched_jet_Spectra", Form(";Matched %s Spectra; counts", "JetType1"), 100, 0, 1000);
0103   mpT_Jet2_matched =
0104       ibooker.book1D("Jet2_matched_jet_Spectra", Form(";Matched %s Spectra; counts", "JetType2"), 100, 0, 1000);
0105   mpT_Jet1_unmatched =
0106       ibooker.book1D("Jet1_unmatched_jet_Spectra", Form(";Unmatched %s spectra;counts", "JetType1"), 100, 0, 1000);
0107   mpT_Jet2_unmatched =
0108       ibooker.book1D("Jet2_unmatched_jet_Spectra", Form(";Unmatched %s spectra;counts", "JetType2"), 100, 0, 1000);
0109 
0110   if (std::string("VsCalo") == JetType1 || std::string("PuCalo") == JetType1) {
0111     mHadEnergy_Jet1_unmatched =
0112         ibooker.book1D("HadEnergy_Jet1_unmatched",
0113                        Form("HadEnergy_Jet1_unmatched;HadEnergy unmatched %s;counts", "JetType1"),
0114                        50,
0115                        0,
0116                        200);
0117     mEmEnergy_Jet1_unmatched = ibooker.book1D(
0118         "EmEnergy_Jet1_unmatched", Form("EmEnergy_Jet1_unmatched;EMEnergy unmatched %s;counts", "JetType1"), 50, 0, 200);
0119   }
0120 
0121   if (std::string("VsCalo") == JetType2 || std::string("PuCalo") == JetType2) {
0122     mHadEnergy_Jet2_unmatched =
0123         ibooker.book1D("HadEnergy_Jet2_unmatched",
0124                        Form("HadEnergy_Jet2_unmatched;HadEnergy unmatched %s;counts", "JetType2"),
0125                        50,
0126                        0,
0127                        200);
0128     mEmEnergy_Jet2_unmatched = ibooker.book1D(
0129         "EmEnergy_Jet2_unmatched", Form("EmEnergy_Jet2_unmatched;EMEnergy unmatched %s;counts", "JetType2"), 50, 0, 200);
0130   }
0131 
0132   if (std::string("VsPF") == JetType1) {
0133     mChargedHadronEnergy_Jet1_unmatched = ibooker.book1D(
0134         "ChargedHadronEnergy_Jet1_unmatched", Form(";charged HAD energy unmatched %s;counts", "JetType1"), 100, 0, 300);
0135     mNeutralHadronEnergy_Jet1_unmatched = ibooker.book1D(
0136         "neutralHadronEnergy_Jet1_unmatched", Form(";neutral HAD energy unmatched %s;counts", "JetType1"), 100, 0, 300);
0137     mChargedEmEnergy_Jet1_unmatched = ibooker.book1D(
0138         "ChargedEmEnergy_Jet1_unmatched", Form(";charged EM energy unmatched %s;counts", "JetType1"), 100, 0, 300);
0139     mNeutralEmEnergy_Jet1_unmatched = ibooker.book1D(
0140         "neutralEmEnergy_Jet1_unmatched", Form(";neutral EM energy unmatched %s;counts", "JetType1"), 100, 0, 300);
0141     mChargedMuEnergy_Jet1_unmatched = ibooker.book1D(
0142         "ChargedMuEnergy_Jet1_unmatched", Form(";charged Mu energy unmatched %s;counts", "JetType1"), 100, 0, 300);
0143 
0144     mChargedHadEnergyFraction_Jet1_unmatched = ibooker.book1D(
0145         "ChargedHadEnergyFraction_Jet1_unmatched", Form(";h^{+/-} Energy Fraction %s;counts", "JetType1"), 50, 0, 1);
0146     mNeutralHadEnergyFraction_Jet1_unmatched = ibooker.book1D(
0147         "NeutralHadEnergyFraction_Jet1_unmatched", Form(";h^{0} Energy Fraction %s;counts", "JetType1"), 50, 0, 1);
0148     mPhotonEnergyFraction_Jet1_unmatched = ibooker.book1D(
0149         "PhotonEnergyFraction_Jet1_unmatched", Form(";#gamma Energy Fraction %s;counts", "JetType1"), 50, 0, 1);
0150     mElectronEnergyFraction_Jet1_unmatched = ibooker.book1D(
0151         "ElectronEnergyFraction_Jet1_unmatched", Form(";e Energy Fraction %s;counts", "JetType1"), 50, 0, 1);
0152     mMuonEnergyFraction_Jet1_unmatched = ibooker.book1D(
0153         "MuonoEnergyFraction_Jet1_unmatched", Form(";#mu Energy Fraction %s;counts", "JetType1"), 50, 0, 1);
0154   }
0155 
0156   if (std::string("VsPF") == JetType2) {
0157     mChargedHadronEnergy_Jet2_unmatched = ibooker.book1D(
0158         "ChargedHadronEnergy_Jet2_unmatched", Form(";charged HAD energy unmatched %s;counts", "JetType2"), 100, 0, 300);
0159     mNeutralHadronEnergy_Jet2_unmatched = ibooker.book1D(
0160         "neutralHadronEnergy_Jet2_unmatched", Form(";neutral HAD energy unmatched %s;counts", "JetType2"), 100, 0, 300);
0161     mChargedEmEnergy_Jet2_unmatched = ibooker.book1D(
0162         "ChargedEmEnergy_Jet2_unmatched", Form(";charged EM energy unmatched %s;counts", "JetType2"), 100, 0, 300);
0163     mNeutralEmEnergy_Jet2_unmatched = ibooker.book1D(
0164         "neutralEmEnergy_Jet2_unmatched", Form(";neutral EM energy unmatched %s;counts", "JetType2"), 100, 0, 300);
0165     mChargedMuEnergy_Jet2_unmatched = ibooker.book1D(
0166         "ChargedMuEnergy_Jet2_unmatched", Form(";charged Mu energy unmatched %s;counts", "JetType2"), 100, 0, 300);
0167 
0168     mChargedHadEnergyFraction_Jet2_unmatched = ibooker.book1D(
0169         "ChargedHadEnergyFraction_Jet2_unmatched", Form(";h^{+/-} Energy Fraction %s;counts", "JetType2"), 50, 0, 1);
0170     mNeutralHadEnergyFraction_Jet2_unmatched = ibooker.book1D(
0171         "NeutralHadEnergyFraction_Jet2_unmatched", Form(";h^{0} Energy Fraction %s;counts", "JetType2"), 50, 0, 1);
0172     mPhotonEnergyFraction_Jet2_unmatched = ibooker.book1D(
0173         "PhotonEnergyFraction_Jet2_unmatched", Form(";#gamma Energy Fraction %s;counts", "JetType2"), 50, 0, 1);
0174     mElectronEnergyFraction_Jet2_unmatched = ibooker.book1D(
0175         "ElectronEnergyFraction_Jet2_unmatched", Form(";e Energy Fraction %s;counts", "JetType2"), 50, 0, 1);
0176     mMuonEnergyFraction_Jet2_unmatched = ibooker.book1D(
0177         "MuonoEnergyFraction_Jet2_unmatched", Form(";#mu Energy Fraction %s;counts", "JetType2"), 50, 0, 1);
0178   }
0179 
0180   if (mOutputFile.empty())
0181     LogInfo("OutputInfo") << " Histograms will NOT be saved";
0182   else
0183     LogInfo("OutputInfo") << " Histograms will be saved to file:" << mOutputFile;
0184 }
0185 
0186 //------------------------------------------------------------------------------
0187 // ~JetAnalyzer_HeavyIons
0188 //------------------------------------------------------------------------------
0189 JetAnalyzer_HeavyIons_matching::~JetAnalyzer_HeavyIons_matching() {}
0190 
0191 //------------------------------------------------------------------------------
0192 // analyze
0193 //------------------------------------------------------------------------------
0194 void JetAnalyzer_HeavyIons_matching::analyze(const edm::Event& mEvent, const edm::EventSetup& mSetup) {
0195   // Get the Jet collection
0196   //----------------------------------------------------------------------------
0197 
0198   std::vector<const Jet*> recoJet1;
0199   recoJet1.clear();
0200   std::vector<const Jet*> recoJet2;
0201   recoJet2.clear();
0202 
0203   edm::Handle<CaloJetCollection> caloJet1;
0204   edm::Handle<CaloJetCollection> caloJet2;
0205   edm::Handle<JPTJetCollection> jptJets;
0206   edm::Handle<PFJetCollection> pfJets;
0207   edm::Handle<BasicJetCollection> basicJets;
0208 
0209   if (std::string("VsCalo") == JetType1) {
0210     mEvent.getByToken(caloJet1Token_, caloJet1);
0211     for (unsigned ijet = 0; ijet < caloJet1->size(); ++ijet)
0212       recoJet1.push_back(&(*caloJet1)[ijet]);
0213   }
0214   if (std::string("PuCalo") == JetType1) {
0215     mEvent.getByToken(caloJet2Token_, caloJet1);
0216     for (unsigned ijet = 0; ijet < caloJet1->size(); ++ijet)
0217       recoJet1.push_back(&(*caloJet1)[ijet]);
0218   }
0219   if (std::string("VsPF") == JetType1) {
0220     mEvent.getByToken(pfJetsToken_, pfJets);
0221     for (unsigned ijet = 0; ijet < pfJets->size(); ++ijet)
0222       recoJet1.push_back(&(*pfJets)[ijet]);
0223   }
0224   if (std::string("PuPF") == JetType1) {
0225     mEvent.getByToken(basicJetsToken_, basicJets);
0226     for (unsigned ijet = 0; ijet < basicJets->size(); ++ijet)
0227       recoJet1.push_back(&(*basicJets)[ijet]);
0228   }
0229 
0230   if (std::string("VsCalo") == JetType2) {
0231     mEvent.getByToken(caloJet1Token_, caloJet2);
0232     for (unsigned ijet = 0; ijet < caloJet2->size(); ++ijet)
0233       recoJet2.push_back(&(*caloJet2)[ijet]);
0234   }
0235   if (std::string("PuCalo") == JetType2) {
0236     mEvent.getByToken(caloJet2Token_, caloJet2);
0237     for (unsigned ijet = 0; ijet < caloJet2->size(); ++ijet)
0238       recoJet2.push_back(&(*caloJet2)[ijet]);
0239   }
0240   if (std::string("VsPF") == JetType2) {
0241     mEvent.getByToken(pfJetsToken_, pfJets);
0242     for (unsigned ijet = 0; ijet < pfJets->size(); ++ijet)
0243       recoJet2.push_back(&(*pfJets)[ijet]);
0244   }
0245   if (std::string("PuPF") == JetType2) {
0246     mEvent.getByToken(basicJetsToken_, basicJets);
0247     for (unsigned ijet = 0; ijet < basicJets->size(); ++ijet)
0248       recoJet2.push_back(&(*basicJets)[ijet]);
0249   }
0250 
0251   // start to perform the matching - between recoJet1 and recoJet2.
0252 
0253   Int_t Jet1_nref = recoJet1.size();
0254   Int_t Jet2_nref = recoJet2.size();
0255 
0256   int jet1 = 0;
0257   int jet2 = 0;
0258 
0259   std::vector<MyJet> vJet1, vJet2;
0260   std::vector<int> Jet1_ID(Jet1_nref), Jet2_ID(Jet2_nref);
0261 
0262   if (Jet1_nref == 0 || Jet2_nref == 0)
0263     return;
0264 
0265   for (unsigned ijet1 = 0; ijet1 < recoJet1.size(); ++ijet1) {
0266     if (recoJet1[ijet1]->pt() < mRecoJetPtThreshold)
0267       continue;
0268     if (fabs(recoJet1[ijet1]->eta()) < mRecoJetEtaCut)
0269       continue;
0270 
0271     MyJet JET1;
0272     JET1.eta = recoJet1[ijet1]->eta();
0273     JET1.phi = recoJet1[ijet1]->phi();
0274     JET1.pt = recoJet1[ijet1]->pt();
0275     JET1.id = ijet1;
0276 
0277     vJet1.push_back(JET1);
0278     jet1++;
0279 
0280   }  // first jet loop
0281 
0282   for (unsigned ijet2 = 0; ijet2 < recoJet2.size(); ++ijet2) {
0283     if (recoJet2[ijet2]->pt() < mRecoJetPtThreshold)
0284       continue;
0285     if (fabs(recoJet2[ijet2]->eta()) < mRecoJetEtaCut)
0286       continue;
0287 
0288     MyJet JET2;
0289     JET2.eta = recoJet2[ijet2]->eta();
0290     JET2.phi = recoJet2[ijet2]->phi();
0291     JET2.pt = recoJet2[ijet2]->pt();
0292     JET2.id = ijet2;
0293 
0294     vJet2.push_back(JET2);
0295     jet2++;
0296 
0297   }  // second jet loop
0298 
0299   bool onlyJet1 = (jet1 > 0 && jet2 == 0) ? true : false;
0300   bool onlyJet2 = (jet1 == 0 && jet2 > 0) ? true : false;
0301   bool bothJet1Jet2 = (jet1 > 0 && jet2 > 0) ? true : false;
0302 
0303   std::vector<MyJet>::const_iterator iJet, jJet;
0304 
0305   if (onlyJet1) {
0306     for (iJet = vJet1.begin(); iJet != vJet1.end(); ++iJet) {
0307       int pj = (*iJet).id;
0308 
0309       mpT_Jet1_unmatched->Fill(recoJet1[pj]->pt());
0310     }
0311 
0312   } else if (onlyJet2) {
0313     for (iJet = vJet2.begin(); iJet != vJet2.end(); ++iJet) {
0314       int cj = (*iJet).id;
0315 
0316       mpT_Jet2_unmatched->Fill(recoJet2[cj]->pt());
0317     }
0318 
0319   } else if (bothJet1Jet2) {
0320     ABMatchedJets mABMatchedJets;
0321 
0322     for (iJet = vJet1.begin(); iJet != vJet1.end(); ++iJet) {
0323       for (jJet = vJet2.begin(); jJet != vJet2.end(); ++jJet) {
0324         mABMatchedJets.insert(std::make_pair(*iJet, *jJet));
0325       }
0326     }
0327 
0328     ABItr itr;
0329     // matched Jets matching Jet 1 to Jet 2
0330     for (itr = mABMatchedJets.begin(); itr != mABMatchedJets.end(); ++itr) {
0331       ABJetPair jetpair = (*itr);
0332       MyJet Aj = jetpair.first;
0333       MyJet Bj = jetpair.second;
0334 
0335       float delr = JetAnalyzer_HeavyIons_matching::deltaRR(Bj.eta, Bj.phi, Aj.eta, Aj.phi);
0336 
0337       if (delr < mRecoDelRMatch && Jet1_ID[Aj.id] == 0) {
0338         mpT_ratio_Jet1Jet2->Fill((Float_t)recoJet2[Bj.id]->pt() / recoJet1[Aj.id]->pt());
0339 
0340         mpT_Jet1_matched->Fill(recoJet1[Aj.id]->pt());
0341         mpT_Jet2_matched->Fill(recoJet2[Bj.id]->pt());
0342 
0343         Jet1_ID[Aj.id] = 1;
0344         Jet2_ID[Bj.id] = 1;
0345       }
0346     }
0347 
0348     // for unmatched Jets
0349     for (itr = mABMatchedJets.begin(); itr != mABMatchedJets.end(); ++itr) {
0350       ABJetPair jetpair = (*itr);
0351 
0352       MyJet Aj = jetpair.first;
0353       MyJet Bj = jetpair.second;
0354 
0355       if (Jet1_ID[Aj.id] == 0) {
0356         mpT_Jet1_unmatched->Fill(recoJet1[Aj.id]->pt());
0357         Jet1_ID[Aj.id] = 1;
0358 
0359         if (std::string("VsCalo") == JetType1 || std::string("PuCalo") == JetType1) {
0360           mHadEnergy_Jet1_unmatched->Fill((*caloJet1)[Aj.id].hadEnergyInHO() + (*caloJet1)[Aj.id].hadEnergyInHB() +
0361                                           (*caloJet1)[Aj.id].hadEnergyInHF() + (*caloJet1)[Aj.id].hadEnergyInHE());
0362           mEmEnergy_Jet1_unmatched->Fill((*caloJet1)[Aj.id].emEnergyInEB() + (*caloJet1)[Aj.id].emEnergyInEE() +
0363                                          (*caloJet1)[Aj.id].emEnergyInHF());
0364         }
0365 
0366         if (std::string("VsPF") == JetType1) {
0367           mChargedHadronEnergy_Jet1_unmatched->Fill((*pfJets)[Aj.id].chargedHadronEnergy());
0368           mNeutralHadronEnergy_Jet1_unmatched->Fill((*pfJets)[Aj.id].neutralHadronEnergy());
0369           mChargedEmEnergy_Jet1_unmatched->Fill((*pfJets)[Aj.id].chargedEmEnergy());
0370           mNeutralEmEnergy_Jet1_unmatched->Fill((*pfJets)[Aj.id].neutralEmEnergy());
0371           mChargedMuEnergy_Jet1_unmatched->Fill((*pfJets)[Aj.id].chargedMuEnergy());
0372 
0373           mChargedHadEnergyFraction_Jet1_unmatched->Fill((*pfJets)[Aj.id].chargedHadronEnergyFraction());
0374           mNeutralHadEnergyFraction_Jet1_unmatched->Fill((*pfJets)[Aj.id].neutralHadronEnergyFraction());
0375           mPhotonEnergyFraction_Jet1_unmatched->Fill((*pfJets)[Aj.id].photonEnergyFraction());
0376           mElectronEnergyFraction_Jet1_unmatched->Fill((*pfJets)[Aj.id].electronEnergyFraction());
0377           mMuonEnergyFraction_Jet1_unmatched->Fill((*pfJets)[Aj.id].muonEnergyFraction());
0378         }
0379       }
0380 
0381       if (Jet2_ID[Bj.id] == 0) {
0382         mpT_Jet2_unmatched->Fill(recoJet2[Bj.id]->pt());
0383         Jet2_ID[Bj.id] = 2;
0384         if (std::string("VsCalo") == JetType2 || std::string("PuCalo") == JetType2) {
0385           mHadEnergy_Jet2_unmatched->Fill((*caloJet2)[Bj.id].hadEnergyInHO() + (*caloJet2)[Bj.id].hadEnergyInHB() +
0386                                           (*caloJet2)[Bj.id].hadEnergyInHF() + (*caloJet2)[Bj.id].hadEnergyInHE());
0387           mEmEnergy_Jet2_unmatched->Fill((*caloJet2)[Bj.id].emEnergyInEB() + (*caloJet2)[Bj.id].emEnergyInEE() +
0388                                          (*caloJet2)[Bj.id].emEnergyInHF());
0389         }
0390 
0391         if (std::string("VsPF") == JetType2) {
0392           mChargedHadronEnergy_Jet2_unmatched->Fill((*pfJets)[Bj.id].chargedHadronEnergy());
0393           mNeutralHadronEnergy_Jet2_unmatched->Fill((*pfJets)[Bj.id].neutralHadronEnergy());
0394           mChargedEmEnergy_Jet2_unmatched->Fill((*pfJets)[Bj.id].chargedEmEnergy());
0395           mNeutralEmEnergy_Jet2_unmatched->Fill((*pfJets)[Bj.id].neutralEmEnergy());
0396           mChargedMuEnergy_Jet2_unmatched->Fill((*pfJets)[Bj.id].chargedMuEnergy());
0397 
0398           mChargedHadEnergyFraction_Jet2_unmatched->Fill((*pfJets)[Bj.id].chargedHadronEnergyFraction());
0399           mNeutralHadEnergyFraction_Jet2_unmatched->Fill((*pfJets)[Bj.id].neutralHadronEnergyFraction());
0400           mPhotonEnergyFraction_Jet2_unmatched->Fill((*pfJets)[Bj.id].photonEnergyFraction());
0401           mElectronEnergyFraction_Jet2_unmatched->Fill((*pfJets)[Bj.id].electronEnergyFraction());
0402           mMuonEnergyFraction_Jet2_unmatched->Fill((*pfJets)[Bj.id].muonEnergyFraction());
0403         }
0404       }
0405     }
0406 
0407   }  // both Jet1 and Jet2 in the event.
0408 }