CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/DataFormats/PatCandidates/src/Jet.cc

Go to the documentation of this file.
00001 //
00002 // $Id: Jet.cc,v 1.44 2011/06/08 20:40:19 rwolf Exp $
00003 //
00004 
00005 #include "DataFormats/PatCandidates/interface/Jet.h"
00006 #include "DataFormats/RecoCandidate/interface/RecoCaloTowerCandidate.h"
00007 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00009 
00010 using namespace pat;
00011 
00013 Jet::Jet() :
00014   PATObject<reco::Jet>(reco::Jet()),
00015   embeddedCaloTowers_(false),
00016   embeddedPFCandidates_(false),
00017   partonFlavour_(0),
00018   jetCharge_(0.),
00019   isCaloTowerCached_(false),
00020   isPFCandidateCached_(false)
00021 {
00022 }
00023 
00025 Jet::Jet(const reco::Jet & aJet) :
00026   PATObject<reco::Jet>(aJet),
00027   embeddedCaloTowers_(false),
00028   embeddedPFCandidates_(false),
00029   partonFlavour_(0),
00030   jetCharge_(0.0),
00031   isCaloTowerCached_(false),
00032   isPFCandidateCached_(false)
00033 {
00034   tryImportSpecific(aJet);
00035 }
00036 
00038 Jet::Jet(const edm::Ptr<reco::Jet> & aJetRef) :
00039   PATObject<reco::Jet>(aJetRef),
00040   embeddedCaloTowers_(false),
00041   embeddedPFCandidates_(false),
00042   partonFlavour_(0),
00043   jetCharge_(0.0),
00044   isCaloTowerCached_(false),
00045   isPFCandidateCached_(false)
00046 {
00047   tryImportSpecific(*aJetRef);
00048 }
00049 
00051 Jet::Jet(const edm::RefToBase<reco::Jet> & aJetRef) :
00052   PATObject<reco::Jet>(aJetRef),
00053   embeddedCaloTowers_(false),
00054   embeddedPFCandidates_(false),
00055   partonFlavour_(0),
00056   jetCharge_(0.0),
00057   isCaloTowerCached_(false),
00058   isPFCandidateCached_(false)
00059 {
00060   tryImportSpecific(*aJetRef);
00061 }
00062 
00063 std::ostream& 
00064 reco::operator<<(std::ostream& out, const pat::Jet& obj) 
00065 {
00066   if(!out) return out;
00067   
00068   out << "\tpat::Jet: ";
00069   out << std::setiosflags(std::ios::right);
00070   out << std::setiosflags(std::ios::fixed);
00071   out << std::setprecision(3);
00072   out << " E/pT/eta/phi " 
00073       << obj.energy()<<"/"
00074       << obj.pt()<<"/"
00075       << obj.eta()<<"/"
00076       << obj.phi();
00077   return out; 
00078 }
00079 
00081 void Jet::tryImportSpecific(const reco::Jet& source)
00082 {
00083   const std::type_info & type = typeid(source);
00084   if( type == typeid(reco::CaloJet) ){
00085     specificCalo_.push_back( (static_cast<const reco::CaloJet&>(source)).getSpecific() );
00086   } else if( type == typeid(reco::JPTJet) ){
00087     reco::JPTJet const & jptJet = static_cast<reco::JPTJet const &>(source);
00088     specificJPT_.push_back( jptJet.getSpecific() );
00089     reco::CaloJet const * caloJet = 0;
00090     if ( jptJet.getCaloJetRef().isNonnull() && jptJet.getCaloJetRef().isAvailable() ) {
00091       caloJet = dynamic_cast<reco::CaloJet const *>( jptJet.getCaloJetRef().get() );
00092     }
00093     if ( caloJet != 0 ) {
00094       specificCalo_.push_back( caloJet->getSpecific() );
00095     }
00096     else {
00097       edm::LogWarning("OptionalProductNotFound") << " in pat::Jet, Attempted to add Calo Specifics to JPT Jets, but failed."
00098                                                  << " Jet ID for JPT Jets will not be available for you." << std::endl;
00099     }
00100   } else if( type == typeid(reco::PFJet) ){
00101     specificPF_.push_back( (static_cast<const reco::PFJet&>(source)).getSpecific() );
00102   }
00103 }
00104 
00106 Jet::~Jet() {
00107 }
00108 
00110 
00111 CaloTowerPtr Jet::getCaloConstituent (unsigned fIndex) const {
00112     if (embeddedCaloTowers_) {
00113       // Refactorized PAT access
00114       if ( caloTowersFwdPtr_.size() > 0 ) {
00115         return (fIndex < caloTowersFwdPtr_.size() ?
00116                 caloTowersFwdPtr_[fIndex].ptr() : CaloTowerPtr());
00117       }
00118       // Compatibility PAT access
00119       else {
00120         if ( caloTowers_.size() > 0 ) {
00121           return (fIndex < caloTowers_.size() ?
00122                   CaloTowerPtr(&caloTowers_, fIndex) : CaloTowerPtr());
00123 
00124         }
00125       }
00126     }
00127     // Non-embedded access
00128     else {
00129       Constituent dau = daughterPtr (fIndex);
00130       const CaloTower* caloTower = dynamic_cast <const CaloTower*> (dau.get());
00131       if (caloTower != 0) {
00132         return CaloTowerPtr(dau.id(), caloTower, dau.key() );
00133       }
00134       else {
00135         throw cms::Exception("Invalid Constituent") << "CaloJet constituent is not of CaloTower type";
00136       }
00137 
00138     }
00139 
00140     return CaloTowerPtr ();
00141 }
00142 
00143 
00144 
00145 std::vector<CaloTowerPtr> const & Jet::getCaloConstituents () const {
00146   if ( !isCaloTowerCached_ || caloTowers_.size() > 0 ) cacheCaloTowers();
00147   return caloTowersTemp_;
00148 }
00149 
00150 
00152 
00153 reco::PFCandidatePtr Jet::getPFConstituent (unsigned fIndex) const {
00154     if (embeddedPFCandidates_) {
00155       // Refactorized PAT access
00156       if ( pfCandidatesFwdPtr_.size() > 0 ) {
00157         return (fIndex < pfCandidatesFwdPtr_.size() ?
00158                 pfCandidatesFwdPtr_[fIndex].ptr() : reco::PFCandidatePtr());
00159       }
00160       // Compatibility PAT access
00161       else {
00162         if ( pfCandidates_.size() > 0 ) {
00163           return (fIndex < pfCandidates_.size() ?
00164                   reco::PFCandidatePtr(&pfCandidates_, fIndex) : reco::PFCandidatePtr());
00165 
00166         }
00167       }
00168     }
00169     // Non-embedded access
00170     else {
00171       Constituent dau = daughterPtr (fIndex);
00172       const reco::PFCandidate* pfCandidate = dynamic_cast <const reco::PFCandidate*> (dau.get());
00173       if (pfCandidate) {
00174         return reco::PFCandidatePtr(dau.id(), pfCandidate, dau.key() );
00175       }
00176       else {
00177         throw cms::Exception("Invalid Constituent") << "PFJet constituent is not of PFCandidate type";
00178       }
00179 
00180     }
00181 
00182     return reco::PFCandidatePtr ();
00183 }
00184 
00185 std::vector<reco::PFCandidatePtr> const & Jet::getPFConstituents () const {
00186   if ( !isPFCandidateCached_ || pfCandidates_.size() > 0 ) cachePFCandidates();
00187   return pfCandidatesTemp_;
00188 }
00189 
00190 
00191 
00193 const reco::GenJet * Jet::genJet() const {
00194   if (genJet_.size()) return  &(genJet_.front());
00195   else if ( genJetRef_.size() ) return genJetRef_[0].get();
00196   else return genJetFwdRef_.get();
00197 }
00198 
00200 int Jet::partonFlavour() const {
00201   return partonFlavour_;
00202 }
00203 
00205 
00206 // initialize the jet to a given JEC level during creation starting from Uncorrected
00207 void Jet::initializeJEC(unsigned int level, const JetCorrFactors::Flavor& flavor, unsigned int set)
00208 {
00209   currentJECSet(set);
00210   currentJECLevel(level);
00211   currentJECFlavor(flavor);
00212   setP4(jec_[set].correction(level, flavor)*p4());
00213 }
00214 
00216 int Jet::jecSet(const std::string& set) const
00217 {
00218   for(std::vector<pat::JetCorrFactors>::const_iterator corrFactor=jec_.begin(); corrFactor!=jec_.end(); ++corrFactor)
00219     if( corrFactor->jecSet()==set ){ return (corrFactor-jec_.begin()); }
00220   return -1;
00221 }
00222 
00224 const std::vector<std::string> Jet::availableJECSets() const
00225 {
00226   std::vector<std::string> sets;
00227   for(std::vector<pat::JetCorrFactors>::const_iterator corrFactor=jec_.begin(); corrFactor!=jec_.end(); ++corrFactor)
00228     sets.push_back(corrFactor->jecSet());
00229   return sets;
00230 }
00231 
00232 const std::vector<std::string> Jet::availableJECLevels(const int& set) const
00233 {
00234   return set>=0 ? jec_.at(set).correctionLabels() : std::vector<std::string>();
00235 }
00236 
00239 float Jet::jecFactor(const std::string& level, const std::string& flavor, const std::string& set) const
00240 {
00241   for(unsigned int idx=0; idx<jec_.size(); ++idx){
00242     if(set.empty() || jec_.at(idx).jecSet()==set){
00243       if(jec_[idx].jecLevel(level)>=0){
00244         return jecFactor(jec_[idx].jecLevel(level), jec_[idx].jecFlavor(flavor), idx);
00245       }
00246       else{
00247         throw cms::Exception("InvalidRequest") << "This JEC level " << level << " does not exist. \n";
00248       }
00249     }
00250   }
00251   throw cms::Exception("InvalidRequest") << "This jet does not carry any jet energy correction factor information \n"
00252                                          << "for a jet energy correction set with label " << set << "\n";
00253 }
00254 
00257 float Jet::jecFactor(const unsigned int& level, const JetCorrFactors::Flavor& flavor, const unsigned int& set) const
00258 {
00259   if(!jecSetsAvailable()){
00260     throw cms::Exception("InvalidRequest") << "This jet does not carry any jet energy correction factor information \n";
00261   }
00262   if(!jecSetAvailable(set)){
00263     throw cms::Exception("InvalidRequest") << "This jet does not carry any jet energy correction factor information \n"
00264                                            << "for a jet energy correction set with index " << set << "\n";
00265   }
00266   return jec_.at(set).correction(level, flavor)/jec_.at(currentJECSet_).correction(currentJECLevel_, currentJECFlavor_);
00267 }
00268 
00271 Jet Jet::correctedJet(const std::string& level, const std::string& flavor, const std::string& set) const
00272 {
00273   // rescale p4 of the jet; the update of current values is
00274   // done within the called jecFactor function
00275   for(unsigned int idx=0; idx<jec_.size(); ++idx){
00276     if(set.empty() || jec_.at(idx).jecSet()==set){
00277       if(jec_[idx].jecLevel(level)>=0){
00278         return correctedJet(jec_[idx].jecLevel(level), jec_[idx].jecFlavor(flavor), idx);
00279       }
00280       else{
00281         throw cms::Exception("InvalidRequest") << "This JEC level " << level << " does not exist. \n";
00282       }
00283     }
00284   }
00285   throw cms::Exception("InvalidRequest") << "This JEC set " << set << " does not exist. \n";
00286 }
00287 
00290 Jet Jet::correctedJet(const unsigned int& level, const JetCorrFactors::Flavor& flavor, const unsigned int& set) const
00291 {
00292   Jet correctedJet(*this);
00293   //rescale p4 of the jet
00294   correctedJet.setP4(jecFactor(level, flavor, set)*p4());
00295   // update current level, flavor and set
00296   correctedJet.currentJECSet(set); correctedJet.currentJECLevel(level); correctedJet.currentJECFlavor(flavor);
00297   return correctedJet;
00298 }
00299 
00300 
00302 
00303 const std::vector<std::pair<std::string, float> > & Jet::getPairDiscri() const {
00304    return pairDiscriVector_;
00305 }
00306 
00308 float Jet::bDiscriminator(const std::string & aLabel) const {
00309   float discriminator = -1000.;
00310   const std::string & theLabel = ((aLabel == "" || aLabel == "default")) ? "trackCountingHighEffBJetTags" : aLabel;
00311   for(unsigned int i=0; i!=pairDiscriVector_.size(); i++){
00312     if(pairDiscriVector_[i].first == theLabel){
00313       discriminator = pairDiscriVector_[i].second;
00314     }
00315   }
00316   return discriminator;
00317 }
00318 
00319 const reco::BaseTagInfo * Jet::tagInfo(const std::string &label) const {
00320     std::vector<std::string>::const_iterator it = std::find(tagInfoLabels_.begin(), tagInfoLabels_.end(), label);
00321     if (it != tagInfoLabels_.end()) {
00322       if ( tagInfosFwdPtr_.size() > 0 ) return tagInfosFwdPtr_[it - tagInfoLabels_.begin()].get();
00323       else if ( tagInfos_.size() > 0 )  return & tagInfos_[it - tagInfoLabels_.begin()];
00324       return 0;
00325     }
00326     return 0;
00327 }
00328 
00329 
00330 template<typename T>
00331 const T *  Jet::tagInfoByType() const {
00332   // First check the factorized PAT version
00333     for (size_t i = 0, n = tagInfosFwdPtr_.size(); i < n; ++i) {
00334       TagInfoFwdPtrCollection::value_type const & val = tagInfosFwdPtr_[i];
00335       reco::BaseTagInfo const * baseTagInfo = val.get();
00336       if ( typeid(*baseTagInfo) == typeid(T) ) {
00337         return static_cast<const T *>( baseTagInfo );
00338       }
00339     }
00340     // Then check compatibility version
00341     for (size_t i = 0, n = tagInfos_.size(); i < n; ++i) {
00342       edm::OwnVector<reco::BaseTagInfo>::value_type const & val = tagInfos_[i];
00343       reco::BaseTagInfo const * baseTagInfo = &val;
00344       if ( typeid(*baseTagInfo) == typeid(T) ) {
00345         return static_cast<const T *>( baseTagInfo );
00346       }
00347     }
00348     return 0;
00349 }
00350 
00351 
00352 
00353 const reco::TrackIPTagInfo *
00354 Jet::tagInfoTrackIP(const std::string &label) const {
00355     return (label.empty() ? tagInfoByType<reco::TrackIPTagInfo>()
00356                           : dynamic_cast<const reco::TrackIPTagInfo *>(tagInfo(label)) );
00357 }
00358 
00359 const reco::SoftLeptonTagInfo *
00360 Jet::tagInfoSoftLepton(const std::string &label) const {
00361     return (label.empty() ? tagInfoByType<reco::SoftLeptonTagInfo>()
00362                           : dynamic_cast<const reco::SoftLeptonTagInfo *>(tagInfo(label)) );
00363 }
00364 
00365 const reco::SecondaryVertexTagInfo *
00366 Jet::tagInfoSecondaryVertex(const std::string &label) const {
00367     return (label.empty() ? tagInfoByType<reco::SecondaryVertexTagInfo>()
00368                           : dynamic_cast<const reco::SecondaryVertexTagInfo *>(tagInfo(label)) );
00369 }
00370 
00371 void
00372 Jet::addTagInfo(const std::string &label,
00373                 const TagInfoFwdPtrCollection::value_type &info) {
00374     std::string::size_type idx = label.find("TagInfos");
00375     if (idx == std::string::npos) {
00376       tagInfoLabels_.push_back(label);
00377     } else {
00378         tagInfoLabels_.push_back(label.substr(0,idx));
00379     }
00380     tagInfosFwdPtr_.push_back(info);
00381 }
00382 
00383 
00384 
00386 float Jet::jetCharge() const {
00387   return jetCharge_;
00388 }
00389 
00391 const reco::TrackRefVector & Jet::associatedTracks() const {
00392   return associatedTracks_;
00393 }
00394 
00396 void Jet::setAssociatedTracks(const reco::TrackRefVector &tracks) {
00397     associatedTracks_ = tracks;
00398 }
00399 
00401 void Jet::setCaloTowers(const CaloTowerFwdPtrCollection & caloTowers) {
00402   for(unsigned int i = 0; i < caloTowers.size(); ++i) {
00403     caloTowersFwdPtr_.push_back( caloTowers.at(i) );
00404   }
00405   embeddedCaloTowers_ = true;
00406   isCaloTowerCached_ = false;
00407 }
00408 
00409 
00411 void Jet::setPFCandidates(const PFCandidateFwdPtrCollection & pfCandidates) {
00412   for(unsigned int i = 0; i < pfCandidates.size(); ++i) {
00413     pfCandidatesFwdPtr_.push_back(pfCandidates.at(i));
00414   }
00415   embeddedPFCandidates_ = true;
00416   isPFCandidateCached_ = false;
00417 }
00418 
00419 
00421 void Jet::setGenJetRef(const edm::FwdRef<reco::GenJetCollection> & gj)
00422 {
00423   genJetFwdRef_ = gj;
00424 }
00425 
00426 
00427 
00429 void Jet::setPartonFlavour(int partonFl) {
00430   partonFlavour_ = partonFl;
00431 }
00432 
00434 void Jet::addBDiscriminatorPair(const std::pair<std::string, float> & thePair) {
00435   pairDiscriVector_.push_back(thePair);
00436 }
00437 
00439 void Jet::setJetCharge(float jetCharge) {
00440   jetCharge_ = jetCharge;
00441 }
00442 
00443 
00444 
00446 void Jet::cacheCaloTowers() const {
00447   // Clear the cache
00448   caloTowersTemp_.clear();
00449   // Here is where we've embedded constituents
00450   if ( embeddedCaloTowers_ ) {
00451     // Refactorized PAT access
00452     if ( caloTowersFwdPtr_.size() > 0 ) {
00453       for ( CaloTowerFwdPtrVector::const_iterator ibegin=caloTowersFwdPtr_.begin(),
00454               iend = caloTowersFwdPtr_.end(),
00455               icalo = ibegin;
00456             icalo != iend; ++icalo ) {
00457         caloTowersTemp_.push_back( CaloTowerPtr(icalo->ptr() ) );
00458       }
00459     }
00460     // Compatibility access
00461     else if ( caloTowers_.size() > 0 ) {
00462       for ( CaloTowerCollection::const_iterator ibegin=caloTowers_.begin(),
00463               iend = caloTowers_.end(),
00464               icalo = ibegin;
00465             icalo != iend; ++icalo ) {
00466         caloTowersTemp_.push_back( CaloTowerPtr(&caloTowers_, icalo - ibegin ) );
00467       }
00468     }
00469   }
00470   // Non-embedded access
00471   else {
00472     for ( unsigned fIndex = 0; fIndex < numberOfDaughters(); ++fIndex ) {
00473       Constituent const & dau = daughterPtr (fIndex);
00474       const CaloTower* caloTower = dynamic_cast <const CaloTower*> (dau.get());
00475       if (caloTower) {
00476         caloTowersTemp_.push_back( CaloTowerPtr(dau.id(), caloTower,dau.key() ) );
00477       }
00478       else {
00479         throw cms::Exception("Invalid Constituent") << "CaloJet constituent is not of CaloTower type";
00480       }
00481     }
00482   }
00483   // Set the cache flag
00484   isCaloTowerCached_=true;
00485 }
00486 
00488 void Jet::cachePFCandidates() const {
00489   // Clear the cache
00490   pfCandidatesTemp_.clear();
00491   // Here is where we've embedded constituents
00492   if ( embeddedPFCandidates_ ) {
00493     // Refactorized PAT access
00494     if ( pfCandidatesFwdPtr_.size() > 0 ) {
00495       for ( PFCandidateFwdPtrCollection::const_iterator ibegin=pfCandidatesFwdPtr_.begin(),
00496               iend = pfCandidatesFwdPtr_.end(),
00497               ipf = ibegin;
00498             ipf != iend; ++ipf ) {
00499         pfCandidatesTemp_.push_back( reco::PFCandidatePtr(ipf->ptr() ) );
00500       }
00501     }
00502     // Compatibility access
00503     else if ( pfCandidates_.size() > 0 ) {
00504       for ( reco::PFCandidateCollection::const_iterator ibegin=pfCandidates_.begin(),
00505               iend = pfCandidates_.end(),
00506               ipf = ibegin;
00507             ipf != iend; ++ipf ) {
00508         pfCandidatesTemp_.push_back( reco::PFCandidatePtr(&pfCandidates_, ipf - ibegin ) );
00509       }
00510     }
00511   }
00512   // Non-embedded access
00513   else {
00514     for ( unsigned fIndex = 0; fIndex < numberOfDaughters(); ++fIndex ) {
00515       Constituent const & dau = daughterPtr (fIndex);
00516       const reco::PFCandidate* pfCandidate = dynamic_cast <const reco::PFCandidate*> (dau.get());
00517       if (pfCandidate) {
00518         pfCandidatesTemp_.push_back( reco::PFCandidatePtr(dau.id(), pfCandidate,dau.key() ) );
00519       }
00520       else {
00521         throw cms::Exception("Invalid Constituent") << "PFJet constituent is not of PFCandidate type";
00522       }
00523     }
00524   }
00525   // Set the cache flag
00526   isPFCandidateCached_=true;
00527 }