00001 /* 00002 * "Unembed" the pizeros in a reco::PFTau. 00003 * 00004 * This converts a collection of PFTaus which have their PiZeros stored 00005 * as std::vector<RecoTauPiZeros>s to an output collection which has 00006 * the PiZeros stored in a separate product, with the PiZeros stored as Refs 00007 * within the tau. This will improve the de-serialization speed of the taus. 00008 * 00009 * Author: Evan K. Friis, UW Madison 00010 * 00011 */ 00012 00013 #include "FWCore/Framework/interface/Frameworkfwd.h" 00014 #include "FWCore/Framework/interface/Event.h" 00015 #include "FWCore/Framework/interface/EventSetup.h" 00016 #include "FWCore/ParameterSet/interface/ParameterSet.h" 00017 #include "FWCore/Framework/interface/EDProducer.h" 00018 00019 #include "RecoTauTag/RecoTau/interface/RecoTauCommonUtilities.h" 00020 00021 #include "DataFormats/TauReco/interface/PFTau.h" 00022 #include "DataFormats/TauReco/interface/PFTauFwd.h" 00023 #include "DataFormats/TauReco/interface/RecoTauPiZero.h" 00024 #include "DataFormats/TauReco/interface/RecoTauPiZeroFwd.h" 00025 00026 class RecoTauPiZeroUnembedder : public edm::EDProducer { 00027 public: 00028 RecoTauPiZeroUnembedder(const edm::ParameterSet& pset); 00029 virtual ~RecoTauPiZeroUnembedder(){} 00030 void produce(edm::Event& evt, const edm::EventSetup& es); 00031 private: 00032 edm::InputTag src_; 00033 }; 00034 00035 RecoTauPiZeroUnembedder::RecoTauPiZeroUnembedder(const edm::ParameterSet& pset) { 00036 src_ = pset.getParameter<edm::InputTag>("src"); 00037 produces<reco::RecoTauPiZeroCollection>("pizeros"); 00038 produces<reco::PFTauCollection>(); 00039 } 00040 void RecoTauPiZeroUnembedder::produce(edm::Event& evt, const edm::EventSetup& es) { 00041 std::auto_ptr<reco::RecoTauPiZeroCollection> piZerosOut( 00042 new reco::RecoTauPiZeroCollection); 00043 std::auto_ptr<reco::PFTauCollection> tausOut(new reco::PFTauCollection); 00044 00045 edm::Handle<reco::CandidateView> tauView; 00046 evt.getByLabel(src_, tauView); 00047 00048 reco::PFTauRefVector taus = 00049 reco::tau::castView<reco::PFTauRefVector>(tauView); 00050 00051 // Get the reference to the product of where the final pizeros will end up 00052 reco::RecoTauPiZeroRefProd piZeroProd = 00053 evt.getRefBeforePut<reco::RecoTauPiZeroCollection>("pizeros"); 00054 00055 for (size_t iTau = 0; iTau < taus.size(); ++iTau) { 00056 // Make a copy 00057 reco::PFTau myTau = *taus[iTau]; 00058 // The ref vectors that will be filled 00059 reco::RecoTauPiZeroRefVector signalPiZeroRefs; 00060 reco::RecoTauPiZeroRefVector isolationPiZeroRefs; 00061 00062 // Copy the PiZeros into the new vector, while updating what refs they will 00063 // have 00064 const reco::RecoTauPiZeroCollection& signalPiZeros = 00065 myTau.signalPiZeroCandidates(); 00066 00067 for (size_t iPiZero = 0; iPiZero < signalPiZeros.size(); ++iPiZero) { 00068 piZerosOut->push_back(signalPiZeros[iPiZero]); 00069 // Figure out what the ref for this pizero will be in the new coll. 00070 signalPiZeroRefs.push_back( 00071 reco::RecoTauPiZeroRef(piZeroProd, piZerosOut->size()-1)); 00072 } 00073 00074 const reco::RecoTauPiZeroCollection& isolationPiZeroCandidates = 00075 myTau.isolationPiZeroCandidates(); 00076 for (size_t iPiZero = 0; iPiZero < isolationPiZeroCandidates.size(); ++iPiZero) { 00077 piZerosOut->push_back(isolationPiZeroCandidates[iPiZero]); 00078 // Figure out what the ref for this pizero will be in the new coll. 00079 isolationPiZeroRefs.push_back( 00080 reco::RecoTauPiZeroRef(piZeroProd, piZerosOut->size()-1)); 00081 } 00082 00083 myTau.setSignalPiZeroCandidatesRefs(signalPiZeroRefs); 00084 myTau.setIsolationPiZeroCandidatesRefs(isolationPiZeroRefs); 00085 00086 tausOut->push_back(myTau); 00087 } 00088 00089 evt.put(piZerosOut, "pizeros"); 00090 evt.put(tausOut); 00091 } 00092 00093 #include "FWCore/Framework/interface/MakerMacros.h" 00094 DEFINE_FWK_MODULE(RecoTauPiZeroUnembedder);