CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFTauPrimaryVertexProducer.cc
Go to the documentation of this file.
1 /* class PFTauPrimaryVertexProducer
2  * EDProducer of the
3  * authors: Ian M. Nugent
4  * This work is based on the impact parameter work by Rosamaria Venditti and reconstructing the 3 prong taus.
5  * The idea of the fully reconstructing the tau using a kinematic fit comes from
6  * Lars Perchalla and Philip Sauerland Theses under Achim Stahl supervision. This
7  * work was continued by Ian M. Nugent and Vladimir Cherepanov.
8  * Thanks goes to Christian Veelken and Evan Klose Friis for their help and suggestions.
9  */
10 
11 
21 
27 
42 
46 
49 #include <memory>
50 #include <boost/foreach.hpp>
51 #include <TFormula.h>
52 
53 #include <memory>
54 
55 using namespace reco;
56 using namespace edm;
57 using namespace std;
58 
60  public:
61  enum Alg{useInputPV=0, useFontPV};
62 
63  struct DiscCutPair{
64  DiscCutPair():discr_(0),cutFormula_(0){}
65  ~DiscCutPair(){delete cutFormula_;}
68  double cut_;
69  TFormula* cutFormula_;
70  };
71  typedef std::vector<DiscCutPair*> DiscCutPairVec;
72 
73  explicit PFTauPrimaryVertexProducer(const edm::ParameterSet& iConfig);
75  virtual void produce(edm::Event&,const edm::EventSetup&);
76  private:
97  std::auto_ptr<StringCutObjectSelector<reco::PFTau> > cut_;
98  std::auto_ptr<tau::RecoTauVertexAssociator> vertexAssociator_;
99 };
100 
102  PFTauTag_(iConfig.getParameter<edm::InputTag>("PFTauTag")),
103  PFTauToken_(consumes<std::vector<reco::PFTau> >(PFTauTag_)),
104  ElectronTag_(iConfig.getParameter<edm::InputTag>("ElectronTag")),
105  ElectronToken_(consumes<std::vector<reco::Electron> >(ElectronTag_)),
106  MuonTag_(iConfig.getParameter<edm::InputTag>("MuonTag")),
107  MuonToken_(consumes<std::vector<reco::Muon> >(MuonTag_)),
108  PVTag_(iConfig.getParameter<edm::InputTag>("PVTag")),
109  PVToken_(consumes<reco::VertexCollection>(PVTag_)),
110  beamSpotTag_(iConfig.getParameter<edm::InputTag>("beamSpot")),
111  beamSpotToken_(consumes<reco::BeamSpot>(beamSpotTag_)),
112  TrackCollectionTag_(iConfig.getParameter<edm::InputTag>("TrackCollectionTag")),
113  TrackCollectionToken_(consumes<reco::TrackCollection>(TrackCollectionTag_)),
114  Algorithm_(iConfig.getParameter<int>("Algorithm")),
115  qualityCutsPSet_(iConfig.getParameter<edm::ParameterSet>("qualityCuts")),
116  useBeamSpot_(iConfig.getParameter<bool>("useBeamSpot")),
117  useSelectedTaus_(iConfig.getParameter<bool>("useSelectedTaus")),
118  RemoveMuonTracks_(iConfig.getParameter<bool>("RemoveMuonTracks")),
119  RemoveElectronTracks_(iConfig.getParameter<bool>("RemoveElectronTracks"))
120 {
122  std::vector<edm::ParameterSet> discriminators =iConfig.getParameter<std::vector<edm::ParameterSet> >("discriminators");
123  // Build each of our cuts
124  BOOST_FOREACH(const edm::ParameterSet &pset, discriminators) {
125  DiscCutPair* newCut = new DiscCutPair();
126  newCut->inputToken_ =consumes<reco::PFTauDiscriminator>(pset.getParameter<edm::InputTag>("discriminator"));
127 
128  if ( pset.existsAs<std::string>("selectionCut") ) newCut->cutFormula_ = new TFormula("selectionCut", pset.getParameter<std::string>("selectionCut").data());
129  else newCut->cut_ = pset.getParameter<double>("selectionCut");
130  discriminators_.push_back(newCut);
131  }
132  // Build a string cut if desired
133  if (iConfig.exists("cut")) cut_.reset(new StringCutObjectSelector<reco::PFTau>(iConfig.getParameter<std::string>( "cut" )));
135  produces<edm::AssociationVector<PFTauRefProd, std::vector<reco::VertexRef> > >();
136  produces<VertexCollection>("PFTauPrimaryVertices");
137 
139 }
140 
142 
143 }
144 
146  // Obtain Collections
147  edm::ESHandle<TransientTrackBuilder> transTrackBuilder;
148  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",transTrackBuilder);
149 
151  iEvent.getByToken(PFTauToken_,Tau);
152 
154  iEvent.getByToken(ElectronToken_,Electron);
155 
157  iEvent.getByToken(MuonToken_,Mu);
158 
160  iEvent.getByToken(PVToken_,PV);
161 
163  iEvent.getByToken(beamSpotToken_,beamSpot);
164 
165  edm::Handle<reco::TrackCollection> trackCollection;
166  iEvent.getByToken(TrackCollectionToken_,trackCollection);
167 
168  // Set Association Map
169  auto_ptr<edm::AssociationVector<PFTauRefProd, std::vector<reco::VertexRef> > > AVPFTauPV(new edm::AssociationVector<PFTauRefProd, std::vector<reco::VertexRef> >(PFTauRefProd(Tau)));
170  std::auto_ptr<VertexCollection> VertexCollection_out= std::auto_ptr<VertexCollection>(new VertexCollection);
171  reco::VertexRefProd VertexRefProd_out = iEvent.getRefBeforePut<reco::VertexCollection>("PFTauPrimaryVertices");
172 
173  // Load each discriminator
174  BOOST_FOREACH(DiscCutPair *disc, discriminators_) {
176  iEvent.getByToken(disc->inputToken_, discr);
177  disc->discr_ = &(*discr);
178  }
179 
180  // For each Tau Run Algorithim
181  if(Tau.isValid()){
182  for(reco::PFTauCollection::size_type iPFTau = 0; iPFTau < Tau->size(); iPFTau++) {
183  reco::PFTauRef tau(Tau, iPFTau);
184  reco::Vertex thePV;
185  if(useInputPV==Algorithm_){
186  vertexAssociator_->setEvent(iEvent);
187  thePV =(*( vertexAssociator_->associatedVertex(*tau)));
188  }
189  else if(useFontPV==Algorithm_){
190  thePV=PV->front();
191  }
193  // Check if it passed all the discrimiantors
194  bool passed(true);
195  BOOST_FOREACH(const DiscCutPair* disc, discriminators_) {
196  // Check this discriminator passes
197  bool passedDisc = true;
198  if ( disc->cutFormula_ )passedDisc = (disc->cutFormula_->Eval((*disc->discr_)[tau]) > 0.5);
199  else passedDisc = ((*disc->discr_)[tau] > disc->cut_);
200  if ( !passedDisc ){passed = false; break;}
201  }
202  if (passed && cut_.get()){passed = (*cut_)(*tau);}
203  if (passed){
204  std::vector<reco::TrackBaseRef> SignalTracks;
205  for(reco::PFTauCollection::size_type jPFTau = 0; jPFTau < Tau->size(); jPFTau++) {
206  if(useSelectedTaus_ || iPFTau==jPFTau){
207  reco::PFTauRef RefPFTau(Tau, jPFTau);
209  // Get tracks from PFTau daugthers
210  const std::vector<edm::Ptr<reco::PFCandidate> > cands = RefPFTau->signalPFChargedHadrCands();
211  for (std::vector<edm::Ptr<reco::PFCandidate> >::const_iterator iter = cands.begin(); iter!=cands.end(); iter++){
212  if(iter->get()->trackRef().isNonnull()) SignalTracks.push_back(reco::TrackBaseRef(iter->get()->trackRef()));
213  else if(iter->get()->gsfTrackRef().isNonnull()){SignalTracks.push_back(reco::TrackBaseRef(((iter)->get()->gsfTrackRef())));}
214  }
215  }
216  }
217  // Get Muon tracks
218  if(RemoveMuonTracks_){
219 
220  if(Mu.isValid()) {
221  for(reco::MuonCollection::size_type iMuon = 0; iMuon< Mu->size(); iMuon++){
222  reco::MuonRef RefMuon(Mu, iMuon);
223  if(RefMuon->track().isNonnull()) SignalTracks.push_back(reco::TrackBaseRef(RefMuon->track()));
224  }
225  }
226  }
227  // Get Electron Tracks
229  if(Electron.isValid()) {
230  for(reco::ElectronCollection::size_type iElectron = 0; iElectron<Electron->size(); iElectron++){
231  reco::ElectronRef RefElectron(Electron, iElectron);
232  if(RefElectron->track().isNonnull()) SignalTracks.push_back(reco::TrackBaseRef(RefElectron->track()));
233  }
234  }
235  }
237  // Get Non-Tau tracks
238  reco::TrackCollection nonTauTracks;
239 // if (trackCollection.isValid()) {
240 // // remove tau tracks and only tracks associated with the vertex
241 // unsigned int idx = 0;
242 // for (reco::TrackCollection::const_iterator iTrk = trackCollection->begin(); iTrk != trackCollection->end(); ++iTrk, idx++) {
243 // reco::TrackRef tmpRef(trackCollection, idx);
244 // reco::TrackRef tmpRefForBase=tmpRef;
245 // bool isSigTrk = false;
246 // bool fromVertex=false;
247 // for (unsigned int sigTrk = 0; sigTrk < SignalTracks.size(); sigTrk++) {
248 // if (reco::TrackBaseRef(tmpRefForBase)==SignalTracks.at(sigTrk)){isSigTrk = true; break;}
249 // }
250 // if (!isSigTrk) nonSigTracks.push_back(*iTrk);
251 // }
252 // }
253 
254  for(std::vector<reco::TrackBaseRef>::const_iterator vtxTrkRef=thePV.tracks_begin();vtxTrkRef<thePV.tracks_end();vtxTrkRef++){
255  for (unsigned int sigTrk = 0; sigTrk < SignalTracks.size(); sigTrk++) {
256  if((*vtxTrkRef)!=SignalTracks[sigTrk] ){
257  nonTauTracks.push_back(**vtxTrkRef);
258  }
259  }
260  }
262  // Refit the vertex
263  TransientVertex transVtx;
264  std::vector<reco::TransientTrack> transTracks;
265  for (reco::TrackCollection::iterator iter=nonTauTracks.begin(); iter!=nonTauTracks.end(); ++iter){
266  transTracks.push_back(transTrackBuilder->build(*iter));
267  }
268  bool FitOk(true);
269  if ( transTracks.size() >= 3 ) {
271  avf.setWeightThreshold(0.1); //weight per track. allow almost every fit, else --> exception
272  try {
273  if ( !useBeamSpot_ ){
274  transVtx = avf.vertex(transTracks);
275  } else {
276  transVtx = avf.vertex(transTracks, *beamSpot);
277  }
278  } catch (...) {
279  FitOk = false;
280  }
281  } else FitOk = false;
282  if ( FitOk ) thePV = transVtx;
283  }
284  VertexRef VRef = reco::VertexRef(VertexRefProd_out, VertexCollection_out->size());
285  VertexCollection_out->push_back(thePV);
286  AVPFTauPV->setValue(iPFTau, VRef);
287  }
288  }
289  iEvent.put(VertexCollection_out,"PFTauPrimaryVertices");
290  iEvent.put(AVPFTauPV);
291 }
292 
T getParameter(std::string const &) const
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:184
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:44
edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< reco::PFTauDiscriminator > inputToken_
edm::EDGetTokenT< std::vector< reco::PFTau > > PFTauToken_
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:10
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::RefProd< PFTauCollection > PFTauRefProd
references to PFTau collection
Definition: PFTauFwd.h:15
uint16_t size_type
edm::EDGetTokenT< reco::VertexCollection > PVToken_
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
edm::EDGetTokenT< std::vector< reco::Muon > > MuonToken_
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
std::auto_ptr< tau::RecoTauVertexAssociator > vertexAssociator_
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const
edm::EDGetTokenT< reco::TrackCollection > TrackCollectionToken_
bool isValid() const
Definition: HandleBase.h:76
edm::Ref< VertexCollection > VertexRef
persistent reference to a Vertex
Definition: VertexFwd.h:13
RefProd< PROD > getRefBeforePut()
Definition: Event.h:128
std::auto_ptr< StringCutObjectSelector< reco::PFTau > > cut_
Definition: L1GtObject.h:30
const T & get() const
Definition: EventSetup.h:55
std::vector< DiscCutPair * > DiscCutPairVec
virtual void produce(edm::Event &, const edm::EventSetup &)
PFTauPrimaryVertexProducer(const edm::ParameterSet &iConfig)
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:39
tuple discr
Definition: listHistos.py:51
edm::EDGetTokenT< std::vector< reco::Electron > > ElectronToken_