CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2_patch1/src/PhysicsTools/PatAlgos/plugins/PATConversionProducer.cc

Go to the documentation of this file.
00001 //
00002 // $Id: PATConversionProducer.cc,v 1.1 2012/04/14 02:12:39 tjkim Exp $
00003 //
00004 #include "PhysicsTools/PatAlgos/plugins/PATConversionProducer.h"
00005 
00006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00007 #include "FWCore/ParameterSet/interface/FileInPath.h"
00008 
00009 #include "DataFormats/Common/interface/Association.h"
00010 #include "DataFormats/Common/interface/ValueMap.h"
00011 #include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h"
00012 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00013 
00014 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
00015 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00016 
00017 #include "PhysicsTools/PatUtils/interface/TrackerIsolationPt.h"
00018 #include "PhysicsTools/PatUtils/interface/CaloIsolationEnergy.h"
00019 
00020 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00021 #include "DataFormats/VertexReco/interface/Vertex.h"
00022 
00023 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
00024 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
00025 
00026 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
00027 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00028 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00029 #include "TrackingTools/IPTools/interface/IPTools.h"
00030 
00031 #include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h"
00032 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
00033 #include "DataFormats/PatCandidates/interface/Conversion.h"
00034 
00035 #include "RecoEcal/EgammaCoreTools/interface/EcalClusterLazyTools.h"
00036 #include "RecoEgamma/EgammaTools/interface/ConversionTools.h"
00037 #include "DataFormats/Math/interface/LorentzVector.h"
00038 
00039 #include <vector>
00040 #include <memory>
00041 #include "TMath.h"
00042 
00043 using namespace pat;
00044 using namespace std;
00045 
00046 
00047 PATConversionProducer::PATConversionProducer(const edm::ParameterSet & iConfig){
00048 
00049   // general configurables
00050   electronSrc_      = iConfig.getParameter<edm::InputTag>( "electronSource" );
00051 
00052   // produces vector of muons
00053   produces<std::vector<Conversion> >();
00054 
00055 }
00056 
00057 
00058 PATConversionProducer::~PATConversionProducer() {
00059 }
00060 
00061 
00062 void PATConversionProducer::produce(edm::Event & iEvent, const edm::EventSetup & iSetup) {
00063 
00064   // Get the collection of electrons from the event
00065   edm::Handle<edm::View<reco::GsfElectron> > electrons;
00066   iEvent.getByLabel(electronSrc_, electrons);
00067 
00068   edm::Handle<reco::BeamSpot> bsHandle;
00069   iEvent.getByLabel("offlineBeamSpot", bsHandle);
00070   const reco::BeamSpot &beamspot = *bsHandle.product();
00071 
00072   // for conversion veto selection  
00073   edm::Handle<reco::ConversionCollection> hConversions;
00074   iEvent.getByLabel("allConversions", hConversions);
00075 
00076   std::vector<Conversion> * patConversions = new std::vector<Conversion>();
00077 
00078   for (reco::ConversionCollection::const_iterator conv = hConversions->begin(); conv!= hConversions->end(); ++conv) {
00079 
00080     reco::Vertex vtx = conv->conversionVertex();
00081 
00082     int index = 0; 
00083     for (edm::View<reco::GsfElectron>::const_iterator itElectron = electrons->begin(); itElectron != electrons->end(); ++itElectron) {
00084 
00085       //find matched conversions with electron and save those conversions with matched electron index
00086       if (ConversionTools::matchesConversion(*itElectron, *conv)) {
00087 
00088         double vtxProb = TMath::Prob( vtx.chi2(), vtx.ndof());
00089         math::XYZVector mom(conv->refittedPairMomentum());
00090         double dbsx = vtx.x() - beamspot.position().x();   
00091         double dbsy = vtx.y() - beamspot.position().y();
00092         double lxy = (mom.x()*dbsx + mom.y()*dbsy)/mom.rho();
00093         int nHitsMax = 0;
00094 
00095         for (std::vector<uint8_t>::const_iterator it = conv->nHitsBeforeVtx().begin(); it!=conv->nHitsBeforeVtx().end(); ++it) {
00096           if ((*it) > nHitsMax) nHitsMax = (*it);
00097         }
00098 
00099         pat::Conversion anConversion( index ); 
00100         anConversion.setVtxProb( vtxProb );
00101         anConversion.setLxy( lxy );
00102         anConversion.setNHitsMax(  nHitsMax );
00103 
00104         patConversions->push_back(anConversion);
00105         break;
00106       }
00107       index++;
00108     }
00109     
00110   }
00111 
00112   // add the electrons to the event output
00113   std::auto_ptr<std::vector<Conversion> > ptr(patConversions);
00114   iEvent.put(ptr);
00115 
00116 }
00117 
00118 #include "FWCore/Framework/interface/MakerMacros.h"
00119 
00120 DEFINE_FWK_MODULE(PATConversionProducer);