ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/Jeng/PVStudy/plugins/TrackSplitterProducer.cc
Revision: 1.4
Committed: Tue Jan 26 12:07:20 2010 UTC (15 years, 3 months ago) by yygao
Content type: text/plain
Branch: MAIN
CVS Tags: CMSSW_3_6_1_patch3_May27thReReco, CMSSW_3_5_8_patch3_May6thReReco, PASPlot_27Feb2010, CMSSW_3_3_6_patch3_Dec19thReReco, HEAD
Changes since 1.3: +11 -113 lines
Log Message:
Remove the selection based on pixelVertices

File Contents

# User Rev Content
1 yygao 1.1
2     // -*- C++ -*-
3     //
4     // Package: TrackSplitterProducer
5     // Class: TrackSplitterProducer
6     //
7     /**\class TrackSplitterProducer TrackSplitterProducer.cc MySub/TrackSplitterProducer/src/TrackSplitterProducer.cc
8    
9     Description: <one line class summary>
10    
11     Implementation:
12     <Notes on implementation>
13     */
14     //
15    
16     // updated to 19/9/2009 8.30 pm
17    
18     //
19     //
20    
21     // system include files
22     #include <memory>
23    
24    
25     // user include files
26     #include "UserCode/PVStudy/interface/TrackSplitterProducer.h"
27     #include "FWCore/Framework/interface/Frameworkfwd.h"
28     #include "FWCore/Framework/interface/EDProducer.h"
29     #include "FWCore/Framework/interface/Event.h"
30     #include "FWCore/Framework/interface/ESHandle.h"
31     #include "FWCore/Framework/interface/EventSetup.h"
32     #include "FWCore/Framework/interface/MakerMacros.h"
33     #include "FWCore/ParameterSet/interface/ParameterSet.h"
34    
35     #include "TH1F.h"
36     #include "TH2F.h"
37     #include "TFile.h"
38     #include "TROOT.h"
39     #include "TChain.h"
40     #include "TNtuple.h"
41    
42     #include "DataFormats/TrackReco/interface/Track.h"
43     #include "DataFormats/TrackReco/interface/TrackFwd.h"
44 yygao 1.4 #include "UserCode/PVStudy/interface/TrackSorter.h"
45 yygao 1.2
46     #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
47 yygao 1.4 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
48     #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
49     #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
50 yygao 1.1
51     TrackSplitterProducer::TrackSplitterProducer(const edm::ParameterSet & iConfig)
52     : theTrackFilter(iConfig.getParameter<edm::ParameterSet>("TkFilterParameters"))
53     {
54     //now do what ever initialization is needed
55     debug_ = iConfig.getParameter<bool> ("Debug");
56     TrackCollectionTag_ = iConfig.getParameter<edm::InputTag>("TrackCollectionTag");
57     bsSrc = iConfig.getParameter< edm::InputTag >("beamSpot");
58    
59     produces<TrackCollection>("SplittedTracks1");
60 yygao 1.4 produces<TrackCollection>("SplittedTracks2");
61 yygao 1.1 }
62    
63    
64     void TrackSplitterProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
65     {
66     using namespace edm;
67     using namespace reco;
68     using namespace std;
69    
70 yygao 1.4 Handle< TrackCollection> trackCollectionHandle;
71 yygao 1.1 iEvent.getByLabel(TrackCollectionTag_, trackCollectionHandle);
72    
73     edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
74     iEvent.getByLabel(bsSrc,recoBeamSpotHandle);
75     reco::BeamSpot bs = *recoBeamSpotHandle;
76    
77     if(trackCollectionHandle->size() == 0 ) return;
78    
79     // Sort the Tracks in the order of descending pT
80     TrackSorter trksorter;
81     TrackCollection sortedTrackColl = trksorter.sortedList(*(trackCollectionHandle.product()));
82    
83 yygao 1.4
84 yygao 1.1 // Filter the track for PV
85     std::vector<reco::Track> filtedTrackColl;
86     edm::ESHandle<TransientTrackBuilder> theB;
87     iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
88    
89    
90     for(TrackCollection::const_iterator track = sortedTrackColl.begin(); track!= sortedTrackColl.end(); ++track)
91     {
92     TransientTrack t_track = (*theB).build(*track);
93     t_track.setBeamSpot(bs);
94 yygao 1.4 if (theTrackFilter(t_track)) filtedTrackColl.push_back(Track (*track));
95 yygao 1.1 }
96    
97     // create the two splitted trackCollections
98     splittedTColl1_ = std::auto_ptr<TrackCollection> (new TrackCollection );
99     splittedTColl2_ = std::auto_ptr<TrackCollection> (new TrackCollection );
100    
101     unsigned int i = 0;
102     unsigned int irand = 0;
103    
104     if(debug_)
105     cout<<"TrackSplitterProducer looping over "<< filtedTrackColl.size()<< "tracks." << endl;
106    
107 yygao 1.4 //for(TrackCollection::const_iterator track = sortedTrackColl.begin(); track!= sortedTrackColl.end(); ++track, ++i)
108 yygao 1.1 for(TrackCollection::const_iterator track = filtedTrackColl.begin(); track!= filtedTrackColl.end(); ++track, ++i)
109     {
110 yygao 1.4 if(debug_)
111     cout<< "Track "<<i<<" : pT = "<<track->pt()<<endl;
112 yygao 1.1
113     //Iterate every other track
114     if (track != (filtedTrackColl.end() -1) ) {
115     // iterate every other track
116     irand = rand() % 2; // A random number from [0,1]
117     if( irand == 0 )
118     {
119     splittedTColl1_->push_back( Track( *track ) );
120     splittedTColl2_->push_back( Track( *(++track) ) );
121     ++i;
122     }
123     if( irand == 1 )
124     {
125     splittedTColl2_->push_back( Track( *track ) );
126     splittedTColl1_->push_back( Track( *(++track) ) );
127     ++i;
128     }
129     if(debug_)
130     cout<< "Track "<<i<<" : pT = "<<track->pt()<<endl;
131     }
132     }
133    
134     //Iterate over splittedTColl1
135     if(debug_)
136     {
137     unsigned int j = 0;
138     cout<<"TrackSplitterProducer Looping over splittedTColl1: " <<endl;
139     for(TrackCollection::const_iterator track = splittedTColl1_->begin(); track!= splittedTColl1_->end(); ++track, ++j){
140     cout<< "Track "<<j<<" : pT = "<<track->pt()<<endl;
141     }
142     j = 0;
143     cout<<"TrackSplitterProducer Looping over splittedTColl2: " <<endl;
144     for(TrackCollection::const_iterator track = splittedTColl2_->begin(); track!= splittedTColl2_->end(); ++track, ++j){
145     cout<< "Track "<<j<<" : pT = "<<track->pt()<<endl;
146     }
147     }// End if(debug_)
148    
149     // Put Splitted TrackCollection into Event
150     iEvent.put(splittedTColl1_, "SplittedTracks1");
151     iEvent.put(splittedTColl2_, "SplittedTracks2");
152     }
153    
154     //DEFINE_FWK_MODULE(TrackSplitterProducer);