CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/RecoJets/FFTJetProducers/plugins/FFTJetPFPileupCleaner.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    FFTJetProducers
00004 // Class:      FFTJetPFPileupCleaner
00005 // 
00013 //
00014 // Original Author:  Igor Volobouev
00015 //         Created:  Thu Jul 14 17:50:33 CDT 2011
00016 // $Id: FFTJetPFPileupCleaner.cc,v 1.1 2011/07/15 04:26:34 igv Exp $
00017 //
00018 //
00019 
00020 // framework include files
00021 #include "FWCore/Framework/interface/Frameworkfwd.h"
00022 #include "FWCore/Framework/interface/EDProducer.h"
00023 #include "FWCore/Framework/interface/Event.h"
00024 #include "FWCore/Framework/interface/MakerMacros.h"
00025 
00026 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00027 
00028 #include "DataFormats/Common/interface/View.h"
00029 #include "DataFormats/Common/interface/Handle.h"
00030 
00031 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00032 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
00033 
00034 #include "DataFormats/VertexReco/interface/Vertex.h"
00035 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00036 
00037 #define init_param(type, varname) varname (ps.getParameter< type >( #varname ))
00038 
00039 //
00040 // class declaration
00041 //
00042 class FFTJetPFPileupCleaner : public edm::EDProducer
00043 {
00044 public:
00045     explicit FFTJetPFPileupCleaner(const edm::ParameterSet&);
00046     ~FFTJetPFPileupCleaner();
00047 
00048 protected:
00049     // methods
00050     void beginJob();
00051     void produce(edm::Event&, const edm::EventSetup&);
00052     void endJob();
00053 
00054 private:
00055     FFTJetPFPileupCleaner();
00056     FFTJetPFPileupCleaner(const FFTJetPFPileupCleaner&);
00057     FFTJetPFPileupCleaner& operator=(const FFTJetPFPileupCleaner&);
00058 
00059     bool isRemovable(reco::PFCandidate::ParticleType ptype) const;
00060     void setRemovalBit(reco::PFCandidate::ParticleType ptype, bool onOff);
00061     void buildRemovalMask();
00062 
00063     reco::VertexRef findSomeVertex(
00064         const edm::Handle<reco::VertexCollection>& vertices,
00065         const reco::PFCandidate& pfcand) const;
00066 
00067     edm::InputTag PFCandidates;
00068     edm::InputTag Vertices;
00069 
00070     // The following, if true, will cause association of a candidate
00071     // with some vertex no matter what
00072     bool checkClosestZVertex;
00073 
00074     // The following, if true, will cause removal of candidates
00075     // associated with the main vertex
00076     bool removeMainVertex;
00077 
00078     // The following will tell us to remove candidates not associated
00079     // with any vertex
00080     bool removeUnassociated;
00081 
00082     // Do we want to reverse the decision?
00083     bool reverseRemovalDecision;
00084 
00085     // Flags for removing things. See PFCandidate header
00086     // for particle types.
00087     bool remove_X;         // undefined
00088     bool remove_h;         // charged hadron
00089     bool remove_e;         // electron
00090     bool remove_mu;        // muon
00091     bool remove_gamma;     // photon
00092     bool remove_h0;        // neutral hadron
00093     bool remove_h_HF;      // HF tower identified as a hadron
00094     bool remove_egamma_HF; // HF tower identified as an EM particle
00095 
00096     // Mask for removing things
00097     unsigned removalMask;
00098 
00099     // Min and max eta for keeping things
00100     double etaMin;
00101     double etaMax;
00102 
00103     // Cut for the vertex Ndof
00104     double vertexNdofCut;
00105 };
00106 
00107 //
00108 // constructors and destructor
00109 //
00110 FFTJetPFPileupCleaner::FFTJetPFPileupCleaner(const edm::ParameterSet& ps)
00111     : init_param(edm::InputTag, PFCandidates),
00112       init_param(edm::InputTag, Vertices),
00113       init_param(bool, checkClosestZVertex),
00114       init_param(bool, removeMainVertex),
00115       init_param(bool, removeUnassociated),
00116       init_param(bool, reverseRemovalDecision),
00117       init_param(bool, remove_X        ),
00118       init_param(bool, remove_h        ),
00119       init_param(bool, remove_e        ),
00120       init_param(bool, remove_mu       ),
00121       init_param(bool, remove_gamma    ),
00122       init_param(bool, remove_h0       ),
00123       init_param(bool, remove_h_HF     ),
00124       init_param(bool, remove_egamma_HF),
00125       removalMask(0),
00126       init_param(double, etaMin),
00127       init_param(double, etaMax),
00128       init_param(double, vertexNdofCut)
00129 {
00130     buildRemovalMask();
00131     produces<reco::PFCandidateCollection>();
00132 }
00133 
00134 
00135 FFTJetPFPileupCleaner::~FFTJetPFPileupCleaner()
00136 {
00137 }
00138 
00139 
00140 // ------------ method called to produce the data  ------------
00141 void FFTJetPFPileupCleaner::produce(
00142     edm::Event& iEvent, const edm::EventSetup& iSetup)
00143 {
00144     // get PFCandidates
00145     std::auto_ptr<reco::PFCandidateCollection> 
00146         pOutput(new reco::PFCandidateCollection);
00147 
00148     edm::Handle<reco::PFCandidateCollection> pfCandidates;
00149     iEvent.getByLabel(PFCandidates, pfCandidates);
00150 
00151     // get vertices
00152     edm::Handle<reco::VertexCollection> vertices;
00153     iEvent.getByLabel(Vertices, vertices);
00154 
00155     const unsigned ncand = pfCandidates->size();
00156     for (unsigned i=0; i<ncand; ++i)
00157     {
00158         reco::PFCandidatePtr candptr(pfCandidates, i);
00159         bool remove = false;
00160 
00161         if (isRemovable(candptr->particleId()))
00162         {
00163             reco::VertexRef vertexref(findSomeVertex(vertices, *candptr));
00164             if (vertexref.isNull())
00165                 remove = removeUnassociated;
00166             else if (vertexref.key() == 0)
00167                 remove = removeMainVertex;
00168             else
00169                 remove = true;
00170         }
00171 
00172         // Check the eta range
00173         if (!remove)
00174         {
00175             const double eta = candptr->p4().Eta();
00176             remove = eta < etaMin || eta > etaMax;
00177         }
00178 
00179         // Should we remember this candidate?
00180         if (reverseRemovalDecision)
00181             remove = !remove;
00182         if (!remove)
00183         {
00184             const reco::PFCandidate& cand = (*pfCandidates)[i];
00185             pOutput->push_back(cand);
00186             pOutput->back().setSourceCandidatePtr(candptr);
00187         }
00188     }
00189 
00190     iEvent.put(pOutput);
00191 }
00192 
00193 
00194 bool FFTJetPFPileupCleaner::isRemovable(
00195     const reco::PFCandidate::ParticleType ptype) const
00196 {
00197     const unsigned shift = static_cast<unsigned>(ptype);
00198     assert(shift < 32U);
00199     return removalMask & (1U << shift);
00200 }
00201 
00202 
00203 void FFTJetPFPileupCleaner::setRemovalBit(
00204     const reco::PFCandidate::ParticleType ptype, const bool value)
00205 {
00206     const unsigned shift = static_cast<unsigned>(ptype);
00207     assert(shift < 32U);
00208     const unsigned mask = (1U << shift);
00209     if (value)
00210         removalMask |= mask;
00211     else
00212         removalMask &= ~mask;
00213 }
00214 
00215 
00216 // The following essentially duplicates the code in PFPileUp.cc,
00217 // with added cut on ndof.
00218 reco::VertexRef FFTJetPFPileupCleaner::findSomeVertex(
00219     const edm::Handle<reco::VertexCollection>& vertices,
00220     const reco::PFCandidate& pfcand) const
00221 {  
00222     typedef reco::VertexCollection::const_iterator IV;
00223     typedef reco::Vertex::trackRef_iterator IT;
00224 
00225     reco::TrackBaseRef trackBaseRef(pfcand.trackRef());
00226 
00227     size_t iVertex = 0;
00228     unsigned index = 0;
00229     unsigned nFoundVertex = 0;
00230     double bestweight = 0;
00231 
00232     const IV vertend(vertices->end());
00233     for (IV iv=vertices->begin(); iv!=vertend; ++iv, ++index)
00234     {
00235         const double ndof = iv->ndof();
00236         if (!iv->isFake() && ndof > vertexNdofCut)
00237         {
00238             const reco::Vertex& vtx = *iv;
00239 
00240             // loop on tracks in vertices
00241             IT trackend(vtx.tracks_end());
00242             for (IT iTrack=vtx.tracks_begin(); iTrack!=trackend; ++iTrack)
00243             {
00244                 const reco::TrackBaseRef& baseRef = *iTrack;
00245 
00246                 // one of the tracks in the vertex is the same as 
00247                 // the track considered in the function
00248                 if (baseRef == trackBaseRef)
00249                 {
00250                     // select the vertex for which the track has the highest weight
00251                     const double w = vtx.trackWeight(baseRef);
00252                     if (w > bestweight)
00253                     {
00254                         bestweight=w;
00255                         iVertex=index;
00256                         nFoundVertex++;
00257                     }   
00258                 }
00259             }
00260         }
00261     }
00262 
00263     if (nFoundVertex > 0)
00264     {
00265         if (nFoundVertex != 1)
00266             edm::LogWarning("TrackOnTwoVertex")
00267                 << "a track is shared by at least two vertices. "
00268                 << "Used to be an assert";
00269         return reco::VertexRef(vertices, iVertex);
00270     }
00271 
00272     // optional: as a secondary solution, associate the closest vertex in z
00273     if (checkClosestZVertex) 
00274     {
00275         double dzmin = 10000;
00276         double ztrack = pfcand.vertex().z();
00277         bool foundVertex = false;
00278         index = 0;
00279         for (IV iv=vertices->begin(); iv!=vertend; ++iv, ++index)
00280         {
00281             const double ndof = iv->ndof();
00282             if (!iv->isFake() && ndof > vertexNdofCut)
00283             {
00284                 const double dz = fabs(ztrack - iv->z());
00285                 if (dz < dzmin)
00286                 {
00287                     dzmin = dz; 
00288                     iVertex = index;
00289                     foundVertex = true;
00290                 }
00291             }
00292         }
00293 
00294         if (foundVertex) 
00295             return reco::VertexRef(vertices, iVertex);
00296     }
00297 
00298     return reco::VertexRef();
00299 }
00300 
00301 
00302 void FFTJetPFPileupCleaner::buildRemovalMask()
00303 {
00304     setRemovalBit(reco::PFCandidate::X,         remove_X        );
00305     setRemovalBit(reco::PFCandidate::h,         remove_h        );
00306     setRemovalBit(reco::PFCandidate::e,         remove_e        );
00307     setRemovalBit(reco::PFCandidate::mu,        remove_mu       );
00308     setRemovalBit(reco::PFCandidate::gamma,     remove_gamma    );
00309     setRemovalBit(reco::PFCandidate::h0,        remove_h0       );
00310     setRemovalBit(reco::PFCandidate::h_HF,      remove_h_HF     );
00311     setRemovalBit(reco::PFCandidate::egamma_HF, remove_egamma_HF);
00312 }
00313 
00314 
00315 // ------------ method called once each job just before starting event loop
00316 void FFTJetPFPileupCleaner::beginJob()
00317 {
00318 }
00319 
00320 
00321 // ------------ method called once each job just after ending the event loop
00322 void FFTJetPFPileupCleaner::endJob()
00323 {
00324 }
00325 
00326 
00327 //define this as a plug-in
00328 DEFINE_FWK_MODULE(FFTJetPFPileupCleaner);