CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SimTrackIdProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: FastSimulation/SimTrackIdProducer
4 // Class: SimTrackIdProducer
5 //
11 //
12 // Original Author: Vilius Kripas
13 // Created: Fri, 24 Oct 2014 09:47:25 GMT
14 //
15 //
16 #include <memory>
27 #include <vector>
28 #include <stdio.h>
29 
31 {
32  //Main products
33  produces<std::vector<unsigned int> >();
34 
35  // Input Tag
36  edm::InputTag trackCollectionTag = conf.getParameter<edm::InputTag>("trackCollection");
37 
38  maxChi2_ = conf.getParameter<double>("maxChi2");
39 
40 
41  if (conf.exists("overrideTrkQuals")) {
42  edm::InputTag overrideTrkQuals = conf.getParameter<edm::InputTag>("overrideTrkQuals");
43  if ( !(overrideTrkQuals==edm::InputTag("")) )
44  overrideTrkQuals_.push_back( consumes<edm::ValueMap<int> >(overrideTrkQuals) );
45  }
47  filterTracks_=false;
48  if (conf.exists("TrackQuality")){
49  filterTracks_=true;
50  std::string trackQuality = conf.getParameter<std::string>("TrackQuality");
51  if ( !trackQuality.empty() ) {
53  }
54  }
55 
56  // consumes
57  trackToken = consumes<reco::TrackCollection>(trackCollectionTag);
58 }
59 
60 void
62 {
63  // The produced object
64  std::auto_ptr<std::vector<unsigned int> > SimTrackIds(new std::vector<unsigned int>());
65 
66  // The input track collection handle
68  e.getByToken(trackToken,trackCollection);
69 
70  std::vector<edm::Handle<edm::ValueMap<int> > > quals;
71  if ( overrideTrkQuals_.size() > 0) {
72  quals.resize(1);
73  e.getByToken(overrideTrkQuals_[0],quals[0]);
74  }
75 
76  for (size_t i = 0 ; i!=trackCollection->size();++i)
77  {
78  const reco::Track & track = trackCollection->at(i);
79  reco::TrackRef trackRef(trackCollection,i);
80  if (filterTracks_) {
81  bool goodTk = true;
82 
83  if ( quals.size()!=0) {
84  int qual=(*(quals[0]))[trackRef];
85  //std::cout << qual << std::endl;
86  if ( qual < 0 ) {goodTk=false;}
87  //note that this does not work for some trackquals (goodIterative or undefQuality)
88  else
89  goodTk = ( qual & (1<<trackQuality_))>>trackQuality_;
90  }
91  else
92  goodTk=(track.quality(trackQuality_));
93  if ( !goodTk) continue;
94  }
95  if(track.chi2()>maxChi2_) continue ;
96 
97  const TrackingRecHit * hit = *(track.recHitsBegin());
98  if (hit)
99  {
100  const SiTrackerGSMatchedRecHit2D* fsimhit = dynamic_cast<const SiTrackerGSMatchedRecHit2D*>(hit);
101  if (fsimhit)
102  {
103  SimTrackIds->push_back(fsimhit->simtrackId());
104  }
105  }
106 
107  }
108  e.put(SimTrackIds);
109 }
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::TrackCollection > trackToken
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
std::vector< edm::EDGetTokenT< edm::ValueMap< int > > > overrideTrkQuals_
reco::TrackBase::TrackQuality trackQuality_
bool exists(std::string const &parameterName) const
checks if a parameter exists
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
virtual void produce(edm::Event &e, const edm::EventSetup &es) override
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:512
SimTrackIdProducer(const edm::ParameterSet &conf)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:104
tuple conf
Definition: dbtoconf.py:185
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:123
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:473