ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/Jeng/PVStudy/plugins/TrackSplitterProducer.cc
Revision: 1.3
Committed: Wed Nov 11 19:27:50 2009 UTC (15 years, 6 months ago) by yygao
Content type: text/plain
Branch: MAIN
CVS Tags: V00-01-01, V00-01-00, CMSSW_3_1_2_PileUp_v2
Changes since 1.2: +50 -3 lines
Log Message:
Add the option of cut the minzSep between pvtx[0] and the rest

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.2 #include "DataFormats/VertexReco/interface/Vertex.h"
45     #include <DataFormats/VertexReco/interface/VertexFwd.h>
46    
47 yygao 1.1 #include "UserCode/PVStudy/interface/TrackSorter.h"
48 yygao 1.2 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
49     #include "TrackingTools/Records/interface/TransientTrackRecord.h"
50     #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
51     #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
52 yygao 1.1
53 yygao 1.2 #include "PhysicsTools/UtilAlgos/interface/TFileService.h"
54 yygao 1.1
55     TrackSplitterProducer::TrackSplitterProducer(const edm::ParameterSet & iConfig)
56     : theTrackFilter(iConfig.getParameter<edm::ParameterSet>("TkFilterParameters"))
57     {
58     //now do what ever initialization is needed
59     debug_ = iConfig.getParameter<bool> ("Debug");
60     TrackCollectionTag_ = iConfig.getParameter<edm::InputTag>("TrackCollectionTag");
61     bsSrc = iConfig.getParameter< edm::InputTag >("beamSpot");
62 yygao 1.3 vertexCollectionTag_ = iConfig.getUntrackedParameter<edm::InputTag>("vertexCollection");
63     min_pvtx_zsep_ = iConfig.getParameter<double> ("minPvtxZsep");
64 yygao 1.2 pixelVertexCollectionTag_ = iConfig.getParameter<edm::InputTag>("pixelVertexCollectionTag");
65 yygao 1.3 trkFilterZcut_ = iConfig.getParameter<double> ("trkFilterZcut");
66 yygao 1.1
67     produces<TrackCollection>("SplittedTracks1");
68 yygao 1.2 produces<TrackCollection>("SplittedTracks2");
69 yygao 1.1 }
70    
71    
72     void TrackSplitterProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
73     {
74     using namespace edm;
75     using namespace reco;
76     using namespace std;
77    
78 yygao 1.2 if(debug_)
79     cout<<"Starting TrackSplitterProducer::produce()..."<<endl;
80 yygao 1.3
81    
82     //=======================================================
83     // PVTX accessors
84     //=======================================================
85     //vertexCollection
86     static const reco::VertexCollection s_empty_vertexColl;
87     const reco::VertexCollection *vertexColl = &s_empty_vertexColl;
88     edm::Handle<reco::VertexCollection> vertexCollectionHandle;
89     iEvent.getByLabel(vertexCollectionTag_, vertexCollectionHandle);
90     if( iEvent.getByLabel(vertexCollectionTag_, vertexCollectionHandle)) {
91     vertexColl = vertexCollectionHandle.product();
92     } else {
93     cout << "vertexCollection cannot be found -> using empty collection of same type." <<endl;
94     }
95    
96     if (vertexColl->size() == 0) return;
97     // The event must have a valid leading vertex with proper fit (!isFake())
98     if (vertexColl->begin()->isFake() || !(vertexColl->begin()->isValid())) return;
99    
100     // == Apply event pre-selection based on if the leading pvtx has nearby pvtx within max_zsep
101     // zsep is the mininum separation of the secondary pvtxs to the leading one
102     if(vertexColl->size()>1) {
103     double zsep = 9999.0;
104     for(reco::VertexCollection::const_iterator v=vertexColl->begin() + 1 ;
105     v!=vertexColl->end(); ++v) {
106     if(v->isValid() && ! v->isFake()) {
107     // zsep
108     if(fabs(v->z()- vertexColl->begin()->z())< zsep)
109     zsep = fabs(v->z()- vertexColl->begin()->z());
110     }
111     }
112     if( zsep < min_pvtx_zsep_ ) return;
113     }
114 yygao 1.2
115 yygao 1.3 //=======================================================
116     // Track accessors
117     //=======================================================
118    
119 yygao 1.2 edm::Handle<reco::TrackCollection> trackCollectionHandle;
120 yygao 1.1 iEvent.getByLabel(TrackCollectionTag_, trackCollectionHandle);
121    
122 yygao 1.3
123     //=======================================================
124     // BeamSpot accessors
125     //=======================================================
126    
127 yygao 1.1 edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
128     iEvent.getByLabel(bsSrc,recoBeamSpotHandle);
129     reco::BeamSpot bs = *recoBeamSpotHandle;
130    
131     if(trackCollectionHandle->size() == 0 ) return;
132    
133     // Sort the Tracks in the order of descending pT
134     TrackSorter trksorter;
135     TrackCollection sortedTrackColl = trksorter.sortedList(*(trackCollectionHandle.product()));
136    
137     // Filter the track for PV
138     std::vector<reco::Track> filtedTrackColl;
139     edm::ESHandle<TransientTrackBuilder> theB;
140     iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
141    
142 yygao 1.2 //Get pixelVertices
143     static const reco::VertexCollection s_empty_pixelVertexColl;
144     const reco::VertexCollection *pixelVertexColl = &s_empty_pixelVertexColl;
145     edm::Handle<reco::VertexCollection> pixelVertexCollectionHandle;
146     iEvent.getByLabel(pixelVertexCollectionTag_, pixelVertexCollectionHandle);
147     if( iEvent.getByLabel(pixelVertexCollectionTag_, pixelVertexCollectionHandle)) {
148     pixelVertexColl = pixelVertexCollectionHandle.product();
149     } else {
150     cout << "pixelVertexCollection cannot be found. -> using empty collection of same type." <<endl;
151     }
152    
153     if(debug_)
154     cout<<"Getting pixelVertices: size() = "<< pixelVertexColl->size() <<endl;
155    
156 yygao 1.3 if(pixelVertexColl->size() == 0 ) return;
157 yygao 1.1
158 yygao 1.2 //Get the z and sigma z of the first vertex of pixelVertexColl
159     double z_pixelpvtx = pixelVertexColl->begin()->z();
160     double zerr_pixelpvtx = pixelVertexColl->begin()->zError();
161     if(debug_) {
162     std::cout<<"TrackSplitterProducer pixelVertex z = "<<z_pixelpvtx<<" +/- "<<zerr_pixelpvtx<<std::endl;
163     // Fill the track->dz() in the PV difference from first pixelpvtx
164     std::cout<<"Looping through "<< pixelVertexColl->begin()->tracksSize() <<" pixelTracks"<<endl;
165     for(reco::Vertex::trackRef_iterator t = (pixelVertexColl->begin())->tracks_begin();
166     t!= (pixelVertexColl->begin())->tracks_end(); t++) {
167     // illegal charge
168     if ( (**t).charge() < -1 || (**t).charge() > 1 ) {
169     // h["trkPtPV"]->Fill(0.);
170     }
171     else
172     std::cout<<"pixelTrack dz() = "<<(**t).dz()<<" +/- "<<(**t).dzError()<<endl;
173     }
174     }
175    
176 yygao 1.3
177     // Loop over generalTracks and apply trkFilter
178 yygao 1.1 for(TrackCollection::const_iterator track = sortedTrackColl.begin(); track!= sortedTrackColl.end(); ++track)
179     {
180     TransientTrack t_track = (*theB).build(*track);
181     t_track.setBeamSpot(bs);
182 yygao 1.2 //Filter tracks based on default track filter by offlinePrimaryVertices
183     if (theTrackFilter(t_track)) {
184     //Filter tracks based on the z provided by pixelVertices
185     if (fabs((track->dz() - z_pixelpvtx)) < trkFilterZcut_ *track->dzError())
186     filtedTrackColl.push_back(Track (*track));
187     else {
188     if(debug_) {
189     cout<<"WARNING: track dz is not within "<< trkFilterZcut_ << " cm." <<endl;
190     cout<<" track->dz = "<<track->dz()<<" +/- "<<track->dzError()
191     <<" track->dz - z_pixelpvtx = " << track->dz()-z_pixelpvtx<<endl;
192     }
193     }
194     }
195 yygao 1.1 }
196    
197 yygao 1.3
198 yygao 1.1 // create the two splitted trackCollections
199     splittedTColl1_ = std::auto_ptr<TrackCollection> (new TrackCollection );
200     splittedTColl2_ = std::auto_ptr<TrackCollection> (new TrackCollection );
201    
202     unsigned int i = 0;
203     unsigned int irand = 0;
204    
205     if(debug_)
206     cout<<"TrackSplitterProducer looping over "<< filtedTrackColl.size()<< "tracks." << endl;
207    
208     for(TrackCollection::const_iterator track = filtedTrackColl.begin(); track!= filtedTrackColl.end(); ++track, ++i)
209     {
210 yygao 1.2 if(debug_) {
211     cout<< "Track "<<i<<" : pT = "<<track->pt()
212     <<"; dz = "<< track->dz() <<" +/- "<< track->dzError() <<endl;
213     }
214 yygao 1.1
215     //Iterate every other track
216     if (track != (filtedTrackColl.end() -1) ) {
217     // iterate every other track
218     irand = rand() % 2; // A random number from [0,1]
219     if( irand == 0 )
220     {
221     splittedTColl1_->push_back( Track( *track ) );
222     splittedTColl2_->push_back( Track( *(++track) ) );
223     ++i;
224     }
225     if( irand == 1 )
226     {
227     splittedTColl2_->push_back( Track( *track ) );
228     splittedTColl1_->push_back( Track( *(++track) ) );
229     ++i;
230     }
231     if(debug_)
232     cout<< "Track "<<i<<" : pT = "<<track->pt()<<endl;
233     }
234     }
235    
236     //Iterate over splittedTColl1
237     if(debug_)
238     {
239     unsigned int j = 0;
240     cout<<"TrackSplitterProducer Looping over splittedTColl1: " <<endl;
241     for(TrackCollection::const_iterator track = splittedTColl1_->begin(); track!= splittedTColl1_->end(); ++track, ++j){
242     cout<< "Track "<<j<<" : pT = "<<track->pt()<<endl;
243     }
244     j = 0;
245     cout<<"TrackSplitterProducer Looping over splittedTColl2: " <<endl;
246     for(TrackCollection::const_iterator track = splittedTColl2_->begin(); track!= splittedTColl2_->end(); ++track, ++j){
247     cout<< "Track "<<j<<" : pT = "<<track->pt()<<endl;
248     }
249     }// End if(debug_)
250    
251     // Put Splitted TrackCollection into Event
252     iEvent.put(splittedTColl1_, "SplittedTracks1");
253     iEvent.put(splittedTColl2_, "SplittedTracks2");
254     }
255    
256     //DEFINE_FWK_MODULE(TrackSplitterProducer);