00001
00002
00003
00004
00005
00007
00008
00009 #include "RecoJets/JetProducers/interface/JetSpecific.h"
00010
00011
00012 #include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
00013 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
00014 #include "DataFormats/CaloTowers/interface/CaloTowerDetId.h"
00015 #include "DataFormats/RecoCandidate/interface/RecoCaloTowerCandidate.h"
00016 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00017 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00018 #include "Geometry/CaloTopology/interface/HcalTopology.h"
00019 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00020 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
00021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00022 #include "FWCore/Framework/interface/ESHandle.h"
00023 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00024 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00025
00026 #include <cmath>
00027
00028
00029 using namespace std;
00030
00031
00033
00035
00036
00037
00038
00039
00040
00041 void reco::writeSpecific(reco::CaloJet & jet,
00042 reco::Particle::LorentzVector const & p4,
00043 reco::Particle::Point const & point,
00044 std::vector<reco::CandidatePtr> const & constituents,
00045 edm::EventSetup const & c )
00046 {
00047
00048 edm::ESHandle<CaloGeometry> geometry;
00049 c.get<CaloGeometryRecord>().get(geometry);
00050 const CaloSubdetectorGeometry* towerGeometry =
00051 geometry->getSubdetectorGeometry(DetId::Calo, CaloTowerDetId::SubdetId);
00052
00053
00054 reco::CaloJet::Specific specific;
00055 makeSpecific (constituents, *towerGeometry, &specific);
00056
00057 jet = reco::CaloJet( p4, point, specific, constituents);
00058 }
00059
00060
00061
00062
00063 void reco::writeSpecific(reco::BasicJet & jet,
00064 reco::Particle::LorentzVector const & p4,
00065 reco::Particle::Point const & point,
00066 std::vector<reco::CandidatePtr> const & constituents,
00067 edm::EventSetup const & c )
00068 {
00069 jet = reco::BasicJet( p4, point, constituents);
00070 }
00071
00072
00073 void reco::writeSpecific(reco::GenJet & jet,
00074 reco::Particle::LorentzVector const & p4,
00075 reco::Particle::Point const & point,
00076 std::vector<reco::CandidatePtr> const & constituents,
00077 edm::EventSetup const & c )
00078 {
00079
00080
00081 reco::GenJet::Specific specific;
00082 makeSpecific (constituents, &specific);
00083
00084 jet = reco::GenJet( p4, point, specific, constituents);
00085 }
00086
00087
00088 void reco::writeSpecific(reco::PFJet & jet,
00089 reco::Particle::LorentzVector const & p4,
00090 reco::Particle::Point const & point,
00091 std::vector<reco::CandidatePtr> const & constituents,
00092 edm::EventSetup const & c )
00093 {
00094
00095 reco::PFJet::Specific specific;
00096 makeSpecific (constituents, &specific);
00097
00098 int charge = 0.;
00099 for ( std::vector<reco::CandidatePtr>::const_iterator ic = constituents.begin(),
00100 icbegin = constituents.begin(), icend = constituents.end();
00101 ic != icend; ++ic ) {
00102 charge += (*ic)->charge();
00103 }
00104 jet = reco::PFJet( p4, point, specific, constituents);
00105 jet.setCharge( charge );
00106 }
00107
00108
00109
00110 void reco::writeSpecific(reco::TrackJet & jet,
00111 reco::Particle::LorentzVector const & p4,
00112 reco::Particle::Point const & point,
00113 std::vector<reco::CandidatePtr> const & constituents,
00114 edm::EventSetup const & c )
00115 {
00116 jet = reco::TrackJet(p4, point, constituents);
00117 }
00118
00119
00120
00121
00122
00123 bool reco::makeSpecific(vector<reco::CandidatePtr> const & towers,
00124 const CaloSubdetectorGeometry& towerGeometry,
00125 CaloJet::Specific* caloJetSpecific)
00126 {
00127 if (0==caloJetSpecific) return false;
00128
00129
00130
00131
00132 vector<double> eECal_i;
00133 vector<double> eHCal_i;
00134 double eInHad = 0.;
00135 double eInEm = 0.;
00136 double eInHO = 0.;
00137 double eInHB = 0.;
00138 double eInHE = 0.;
00139 double eHadInHF = 0.;
00140 double eEmInHF = 0.;
00141 double eInEB = 0.;
00142 double eInEE = 0.;
00143 double jetArea = 0.;
00144
00145 vector<reco::CandidatePtr>::const_iterator itTower;
00146 for (itTower=towers.begin();itTower!=towers.end();++itTower) {
00147 if ( itTower->isNull() || !itTower->isAvailable() ) {
00148 edm::LogWarning("DataNotFound") << " JetSpecific: Tower is invalid\n";
00149 continue;
00150 }
00151 const CaloTower* tower = dynamic_cast<const CaloTower*>(itTower->get());
00152 if (tower) {
00153
00154 eECal_i.push_back(tower->emEnergy());
00155 eInEm += tower->emEnergy();
00156
00157 eHCal_i.push_back(tower->hadEnergy());
00158 eInHad += tower->hadEnergy();
00159
00160
00161 switch (reco::hcalSubdetector(tower->id().ieta())) {
00162 case HcalBarrel:
00163 eInHB += tower->hadEnergy();
00164 eInHO += tower->outerEnergy();
00165 eInEB += tower->emEnergy();
00166 break;
00167 case HcalEndcap:
00168 eInHE += tower->hadEnergy();
00169 eInEE += tower->emEnergy();
00170 break;
00171 case HcalForward:
00172 eHadInHF += tower->hadEnergy();
00173 eEmInHF += tower->emEnergy();
00174 break;
00175 default:
00176 break;
00177 }
00178
00179 const CaloCellGeometry* geometry = towerGeometry.getGeometry(tower->id());
00180 if (geometry) {
00181 float dEta = fabs(geometry->getCorners()[0].eta()-
00182 geometry->getCorners()[2].eta());
00183 float dPhi = fabs(geometry->getCorners()[0].phi() -
00184 geometry->getCorners()[2].phi());
00185 jetArea += dEta * dPhi;
00186 }
00187 else {
00188 edm::LogWarning("DataNotFound") <<"reco::makeCaloJetSpecific: Geometry for cell "
00189 <<tower->id()<<" can not be found. Ignoring cell\n";
00190 }
00191 }
00192 else {
00193 edm::LogWarning("DataNotFound")<<"reco::makeCaloJetSpecific: Constituent is not of "
00194 <<"CaloTower type\n";
00195 }
00196 }
00197
00198 double towerEnergy = eInHad + eInEm;
00199 caloJetSpecific->mHadEnergyInHO = eInHO;
00200 caloJetSpecific->mHadEnergyInHB = eInHB;
00201 caloJetSpecific->mHadEnergyInHE = eInHE;
00202 caloJetSpecific->mHadEnergyInHF = eHadInHF;
00203 caloJetSpecific->mEmEnergyInHF = eEmInHF;
00204 caloJetSpecific->mEmEnergyInEB = eInEB;
00205 caloJetSpecific->mEmEnergyInEE = eInEE;
00206 if (towerEnergy > 0) {
00207 caloJetSpecific->mEnergyFractionHadronic = eInHad / towerEnergy;
00208 caloJetSpecific->mEnergyFractionEm = eInEm / towerEnergy;
00209 }
00210 else {
00211 caloJetSpecific->mEnergyFractionHadronic = 1.;
00212 caloJetSpecific->mEnergyFractionEm = 0.;
00213 }
00214 caloJetSpecific->mTowersArea = jetArea;
00215 caloJetSpecific->mMaxEInEmTowers = 0;
00216 caloJetSpecific->mMaxEInHadTowers = 0;
00217
00218
00219 sort(eECal_i.begin(), eECal_i.end(), greater<double>());
00220 sort(eHCal_i.begin(), eHCal_i.end(), greater<double>());
00221
00222 if (!towers.empty()) {
00223
00224 caloJetSpecific->mMaxEInEmTowers = eECal_i.front();
00225 caloJetSpecific->mMaxEInHadTowers = eHCal_i.front();
00226 }
00227
00228 return true;
00229 }
00230
00231
00232
00233 bool reco::makeSpecific(vector<reco::CandidatePtr> const & particles,
00234 PFJet::Specific* pfJetSpecific)
00235 {
00236 if (0==pfJetSpecific) return false;
00237
00238
00239
00240
00241
00242 float chargedHadronEnergy=0.;
00243 float neutralHadronEnergy=0.;
00244 float photonEnergy=0.;
00245 float electronEnergy=0.;
00246 float muonEnergy=0.;
00247 float HFHadronEnergy=0.;
00248 float HFEMEnergy=0.;
00249 int chargedHadronMultiplicity=0;
00250 int neutralHadronMultiplicity=0;
00251 int photonMultiplicity=0;
00252 int electronMultiplicity=0;
00253 int muonMultiplicity=0;
00254 int HFHadronMultiplicity=0;
00255 int HFEMMultiplicity=0;
00256
00257 float chargedEmEnergy=0.;
00258 float neutralEmEnergy=0.;
00259 float chargedMuEnergy=0.;
00260 int chargedMultiplicity=0;
00261 int neutralMultiplicity=0;
00262
00263 vector<reco::CandidatePtr>::const_iterator itParticle;
00264 for (itParticle=particles.begin();itParticle!=particles.end();++itParticle){
00265 if ( itParticle->isNull() || !itParticle->isAvailable() ) {
00266 edm::LogWarning("DataNotFound") << " JetSpecific: PF Particle is invalid\n";
00267 continue;
00268 }
00269 const PFCandidate* pfCand = dynamic_cast<const PFCandidate*> (itParticle->get());
00270 if (pfCand) {
00271 switch (PFCandidate::ParticleType(pfCand->particleId())) {
00272 case PFCandidate::h:
00273 chargedHadronEnergy += pfCand->energy();
00274 chargedHadronMultiplicity++;
00275 chargedMultiplicity++;
00276 break;
00277
00278 case PFCandidate::h0 :
00279 neutralHadronEnergy += pfCand->energy();
00280 neutralHadronMultiplicity++;
00281 neutralMultiplicity++;
00282 break;
00283
00284 case PFCandidate::gamma:
00285 photonEnergy += pfCand->energy();
00286 photonMultiplicity++;
00287 neutralEmEnergy += pfCand->energy();
00288 neutralMultiplicity++;
00289 break;
00290
00291 case PFCandidate::e:
00292 electronEnergy += pfCand->energy();
00293 electronMultiplicity++;
00294 chargedEmEnergy += pfCand->energy();
00295 chargedMultiplicity++;
00296 break;
00297
00298 case PFCandidate::mu:
00299 muonEnergy += pfCand->energy();
00300 muonMultiplicity++;
00301 chargedMuEnergy += pfCand->energy();
00302 chargedMultiplicity++;
00303 break;
00304
00305 case PFCandidate::h_HF :
00306 HFHadronEnergy += pfCand->energy();
00307 HFHadronMultiplicity++;
00308 neutralMultiplicity++;
00309 break;
00310
00311 case PFCandidate::egamma_HF :
00312 HFEMEnergy += pfCand->energy();
00313 HFEMMultiplicity++;
00314 neutralEmEnergy += pfCand->energy();
00315 neutralMultiplicity++;
00316 break;
00317
00318
00319 default:
00320 edm::LogWarning("DataNotFound") <<"reco::makePFJetSpecific: Unknown PFCandidate::ParticleType: "
00321 <<pfCand->particleId()<<" is ignored\n";
00322 break;
00323 }
00324 }
00325 else {
00326 edm::LogWarning("DataNotFound") <<"reco::makePFJetSpecific: Referred constituent is not "
00327 <<"a PFCandidate\n";
00328 }
00329 }
00330
00331 pfJetSpecific->mChargedHadronEnergy=chargedHadronEnergy;
00332 pfJetSpecific->mNeutralHadronEnergy= neutralHadronEnergy;
00333 pfJetSpecific->mPhotonEnergy= photonEnergy;
00334 pfJetSpecific->mElectronEnergy= electronEnergy;
00335 pfJetSpecific->mMuonEnergy= muonEnergy;
00336 pfJetSpecific->mHFHadronEnergy= HFHadronEnergy;
00337 pfJetSpecific->mHFEMEnergy= HFEMEnergy;
00338
00339 pfJetSpecific->mChargedHadronMultiplicity=chargedHadronMultiplicity;
00340 pfJetSpecific->mNeutralHadronMultiplicity= neutralHadronMultiplicity;
00341 pfJetSpecific->mPhotonMultiplicity= photonMultiplicity;
00342 pfJetSpecific->mElectronMultiplicity= electronMultiplicity;
00343 pfJetSpecific->mMuonMultiplicity= muonMultiplicity;
00344 pfJetSpecific->mHFHadronMultiplicity= HFHadronMultiplicity;
00345 pfJetSpecific->mHFEMMultiplicity= HFEMMultiplicity;
00346
00347 pfJetSpecific->mChargedEmEnergy=chargedEmEnergy;
00348 pfJetSpecific->mChargedMuEnergy=chargedMuEnergy;
00349 pfJetSpecific->mNeutralEmEnergy=neutralEmEnergy;
00350 pfJetSpecific->mChargedMultiplicity=chargedMultiplicity;
00351 pfJetSpecific->mNeutralMultiplicity=neutralMultiplicity;
00352
00353 return true;
00354 }
00355
00356
00357
00358 bool reco::makeSpecific(vector<reco::CandidatePtr> const & mcparticles,
00359 GenJet::Specific* genJetSpecific)
00360 {
00361 if (0==genJetSpecific) return false;
00362
00363 vector<reco::CandidatePtr>::const_iterator itMcParticle=mcparticles.begin();
00364 for (;itMcParticle!=mcparticles.end();++itMcParticle) {
00365 if ( itMcParticle->isNull() || !itMcParticle->isAvailable() ) {
00366 edm::LogWarning("DataNotFound") << " JetSpecific: MC Particle is invalid\n";
00367 continue;
00368 }
00369 const Candidate* candidate = itMcParticle->get();
00370 if (candidate->hasMasterClone()) candidate = candidate->masterClone().get();
00371 const GenParticle* genParticle = GenJet::genParticle(candidate);
00372 if (genParticle) {
00373 double e = genParticle->energy();
00374 switch (abs (genParticle->pdgId ())) {
00375 case 22:
00376 case 11:
00377 genJetSpecific->m_EmEnergy += e;
00378 break;
00379 case 211:
00380 case 321:
00381 case 130:
00382 case 2212:
00383 case 2112:
00384 genJetSpecific->m_HadEnergy += e;
00385 break;
00386 case 13:
00387 case 12:
00388 case 14:
00389 case 16:
00390
00391 genJetSpecific->m_InvisibleEnergy += e;
00392 break;
00393 default:
00394 genJetSpecific->m_AuxiliaryEnergy += e;
00395 }
00396 }
00397 else {
00398 edm::LogWarning("DataNotFound") <<"reco::makeGenJetSpecific: Referred GenParticleCandidate "
00399 <<"is not available in the event\n";
00400 }
00401 }
00402
00403 return true;
00404 }
00405
00406
00407
00408 HcalSubdetector reco::hcalSubdetector(int iEta)
00409 {
00410 static const HcalTopology topology;
00411 int eta = std::abs(iEta);
00412 if (eta <= topology.lastHBRing()) return HcalBarrel;
00413 else if (eta <= topology.lastHERing()) return HcalEndcap;
00414 else if (eta <= topology.lastHFRing()) return HcalForward;
00415 return HcalEmpty;
00416 }
00417
00418
00419