CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2_patch1/src/RecoVZero/VZeroFinding/plugins/VZeroProducer.cc

Go to the documentation of this file.
00001 #include "VZeroProducer.h"
00002 
00003 #include "RecoVZero/VZeroFinding/interface/VZeroFinder.h"
00004 
00005 #include "FWCore/Framework/interface/MakerMacros.h"
00006 #include "FWCore/Framework/interface/ESHandle.h"
00007 #include "FWCore/Framework/interface/Frameworkfwd.h"
00008 #include "FWCore/ServiceRegistry/interface/Service.h"
00009 
00010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00012 #include "FWCore/Utilities/interface/InputTag.h"
00013 
00014 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00015 
00016 using namespace std;
00017 using namespace edm;
00018 
00019 /*****************************************************************************/
00020 VZeroProducer::VZeroProducer(const ParameterSet& pset)
00021   : pset_(pset)
00022 {
00023   LogInfo("VZeroProducer") << " constructor";
00024   produces<reco::VZeroCollection>();
00025 
00026   // Get track level cuts
00027   minImpactPositiveDaughter =
00028     pset.getParameter<double>("minImpactPositiveDaughter");
00029   minImpactNegativeDaughter =
00030     pset.getParameter<double>("minImpactNegativeDaughter");
00031 }
00032 
00033 /*****************************************************************************/
00034 VZeroProducer::~VZeroProducer()
00035 {
00036   LogInfo("VZeroProducer") << " destructor";
00037 }
00038 
00039 /*****************************************************************************/
00040 void VZeroProducer::produce(Event& ev, const EventSetup& es)
00041 {
00042   // Get tracks
00043   Handle<reco::TrackCollection> trackCollection;
00044   ev.getByLabel(pset_.getParameter<InputTag>("trackCollection"),
00045                                               trackCollection);
00046   const reco::TrackCollection tracks = *(trackCollection.product());
00047 
00048   // Get primary vertices
00049   Handle<reco::VertexCollection> vertexCollection;
00050   ev.getByLabel(pset_.getParameter<InputTag>("vertexCollection"),
00051                                               vertexCollection);
00052   const reco::VertexCollection* vertices = vertexCollection.product();
00053 
00054   // Find vzeros
00055   VZeroFinder theFinder(es,pset_);
00056 
00057   // Selection based on track impact parameter
00058   reco::TrackRefVector positives;
00059   reco::TrackRefVector negatives;
00060 
00061   for(unsigned int i=0; i<tracks.size(); i++)
00062   {
00063     if(tracks[i].charge() > 0 &&
00064        fabs(tracks[i].d0()) > minImpactPositiveDaughter)
00065       positives.push_back(reco::TrackRef(trackCollection, i));
00066 
00067     if(tracks[i].charge() < 0 &&
00068        fabs(tracks[i].d0()) > minImpactNegativeDaughter)
00069       negatives.push_back(reco::TrackRef(trackCollection, i));
00070   }
00071 
00072   LogTrace("MinBiasTracking") << "[VZeroProducer] using tracks :"
00073        << " +" << positives.size()
00074        << " -" << negatives.size();
00075 
00076   auto_ptr<reco::VZeroCollection> result(new reco::VZeroCollection);
00077 
00078   // Check all combination of positives and negatives
00079   if(positives.size() > 0 && negatives.size() > 0)
00080     for(reco::track_iterator ipos = positives.begin();
00081                              ipos!= positives.end(); ipos++)
00082     for(reco::track_iterator ineg = negatives.begin();
00083                              ineg!= negatives.end(); ineg++)
00084     {
00085       reco::VZeroData data;
00086 
00087       if(theFinder.checkTrackPair(**ipos,**ineg, vertices, data) == true)
00088       {
00089         // Create vertex (creation point)
00090         reco::Vertex vertex(reco::Vertex::Point(data.crossingPoint.x(),
00091                                                 data.crossingPoint.y(),
00092                                                 data.crossingPoint.z()),
00093                             reco::Vertex::Error(), 0.,0.,0);
00094 
00095         // Add references to daughters
00096         vertex.add(reco::TrackBaseRef(*ipos));
00097         vertex.add(reco::TrackBaseRef(*ineg));
00098 
00099         // Store vzero
00100         result->push_back(reco::VZero(vertex,data));
00101       }
00102     }
00103 
00104   LogTrace("MinBiasTracking")
00105     << "[VZeroProducer] found candidates : " << result->size();
00106 
00107   // Put result back to the event
00108   ev.put(result);
00109 }
00110