CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FastjetJetProducer.cc
Go to the documentation of this file.
1 //
3 // FastjetJetProducer
4 // ------------------
5 //
6 // 04/21/2009 Philipp Schieferdecker <philipp.schieferdecker@cern.ch>
8 
10 
12 
18 
29 
30 
33 
34 #include "fastjet/SISConePlugin.hh"
35 #include "fastjet/CMSIterativeConePlugin.hh"
36 #include "fastjet/ATLASConePlugin.hh"
37 #include "fastjet/CDFMidPointPlugin.hh"
38 
39 
40 #include <iostream>
41 #include <memory>
42 #include <algorithm>
43 #include <limits>
44 #include <cmath>
45 
46 using namespace std;
47 
48 
49 
51 // construction / destruction
53 
54 //______________________________________________________________________________
56  : VirtualJetProducer( iConfig )
57 {
58 
59  if ( iConfig.exists("UseOnlyVertexTracks") )
60  useOnlyVertexTracks_ = iConfig.getParameter<bool>("UseOnlyVertexTracks");
61  else
62  useOnlyVertexTracks_ = false;
63 
64  if ( iConfig.exists("UseOnlyOnePV") )
65  useOnlyOnePV_ = iConfig.getParameter<bool>("UseOnlyOnePV");
66  else
67  useOnlyOnePV_ = false;
68 
69  if ( iConfig.exists("DzTrVtxMax") )
70  dzTrVtxMax_ = iConfig.getParameter<double>("DzTrVtxMax");
71  else
72  dzTrVtxMax_ = 999999.;
73  if ( iConfig.exists("DxyTrVtxMax") )
74  dxyTrVtxMax_ = iConfig.getParameter<double>("DxyTrVtxMax");
75  else
76  dxyTrVtxMax_ = 999999.;
77  if ( iConfig.exists("MinVtxNdof") )
78  minVtxNdof_ = iConfig.getParameter<int>("MinVtxNdof");
79  else
80  minVtxNdof_ = 5;
81  if ( iConfig.exists("MaxVtxZ") )
82  maxVtxZ_ = iConfig.getParameter<double>("MaxVtxZ");
83  else
84  maxVtxZ_ = 15;
85 
86 }
87 
88 
89 //______________________________________________________________________________
91 {
92 }
93 
94 
96 // implementation of member functions
98 
100 {
101 
102  // for everything but track jets
103  if (!makeTrackJet(jetTypeE)) {
104 
105  // use the default production from one collection
106  VirtualJetProducer::produce( iEvent, iSetup );
107 
108  } else { // produce trackjets from tracks grouped per primary vertex
109 
110  produceTrackJets(iEvent, iSetup);
111 
112  }
113 
114 }
115 
116 
118 {
119 
120  // read in the track candidates
122  iEvent.getByLabel(src_, inputsHandle);
123  // make collection with pointers so we can play around with it
124  std::vector<edm::Ptr<reco::RecoChargedRefCandidate> > allInputs;
125  std::vector<edm::Ptr<reco::Candidate> > origInputs;
126  for (size_t i = 0; i < inputsHandle->size(); ++i) {
127  allInputs.push_back(inputsHandle->ptrAt(i));
128  origInputs.push_back(inputsHandle->ptrAt(i));
129  }
130 
131  // read in the PV collection
133  iEvent.getByLabel(srcPVs_, pvCollection);
134  // define the overall output jet container
135  std::auto_ptr<std::vector<reco::TrackJet> > jets(new std::vector<reco::TrackJet>() );
136 
137  // loop over the good vertices, clustering for each vertex separately
138  for (reco::VertexCollection::const_iterator itVtx = pvCollection->begin(); itVtx != pvCollection->end(); ++itVtx) {
139  if (itVtx->isFake() || itVtx->ndof() < minVtxNdof_ || fabs(itVtx->z()) > maxVtxZ_) continue;
140 
141  // clear the intermediate containers
142  inputs_.clear();
143  fjInputs_.clear();
144  fjJets_.clear();
145 
146  // if only vertex-associated tracks should be used
147  if (useOnlyVertexTracks_) {
148  // loop over the tracks associated to the vertex
149  for (reco::Vertex::trackRef_iterator itTr = itVtx->tracks_begin(); itTr != itVtx->tracks_end(); ++itTr) {
150  // whether a match was found in the track candidate input
151  bool found = false;
152  // loop over input track candidates
153  for (std::vector<edm::Ptr<reco::RecoChargedRefCandidate> >::iterator itIn = allInputs.begin(); itIn != allInputs.end(); ++itIn) {
154  // match the input track candidate to the track from the vertex
155  reco::TrackRef trref(itTr->castTo<reco::TrackRef>());
156  // check if the tracks match
157  if ((*itIn)->track() == trref) {
158  found = true;
159  // add this track candidate to the input for clustering
160  inputs_.push_back(*itIn);
161  // erase the track candidate from the total list of input, so we don't reuse it later
162  allInputs.erase(itIn);
163  // found the candidate track corresponding to the vertex track, so stop the loop
164  break;
165  } // end if match found
166  } // end loop over input tracks
167  // give an info message in case no match is found (can happen if candidates are subset of tracks used for clustering)
168  if (!found) edm::LogInfo("FastjetTrackJetProducer") << "Ignoring a track at vertex which is not in input track collection!";
169  } // end loop over tracks associated to vertex
170  // if all inpt track candidates should be used
171  } else {
172  // loop over input track candidates
173  for (std::vector<edm::Ptr<reco::RecoChargedRefCandidate> >::iterator itIn = allInputs.begin(); itIn != allInputs.end(); ++itIn) {
174  // check if the track is close enough to the vertex
175  float dz = (*itIn)->track()->dz(itVtx->position());
176  float dxy = (*itIn)->track()->dxy(itVtx->position());
177  if (fabs(dz) > dzTrVtxMax_) continue;
178  if (fabs(dxy) > dxyTrVtxMax_) continue;
179  bool closervtx = false;
180  // now loop over the good vertices a second time
181  for (reco::VertexCollection::const_iterator itVtx2 = pvCollection->begin(); itVtx2 != pvCollection->end(); ++itVtx2) {
182  if (itVtx->isFake() || itVtx->ndof() < minVtxNdof_ || fabs(itVtx->z()) > maxVtxZ_) continue;
183  // and check this track is closer to any other vertex (if more than 1 vertex considered)
184  if (!useOnlyOnePV_ &&
185  itVtx != itVtx2 &&
186  fabs((*itIn)->track()->dz(itVtx2->position())) < fabs(dz)) {
187  closervtx = true;
188  break; // 1 closer vertex makes the track already not matched, so break
189  }
190  }
191  // don't add this track if another vertex is found closer
192  if (closervtx) continue;
193  // add this track candidate to the input for clustering
194  inputs_.push_back(*itIn);
195  // erase the track candidate from the total list of input, so we don't reuse it later
196  allInputs.erase(itIn);
197  // take a step back in the loop since we just erased
198  --itIn;
199  }
200  }
201 
202  // convert candidates in inputs_ to fastjet::PseudoJets in fjInputs_
203  fjInputs_.reserve(inputs_.size());
204  inputTowers();
205  LogDebug("FastjetTrackJetProducer") << "Inputted towers\n";
206 
207  // run algorithm, using fjInputs_, modifying fjJets_ and allocating fjClusterSeq_
208  runAlgorithm(iEvent, iSetup);
209  LogDebug("FastjetTrackJetProducer") << "Ran algorithm\n";
210 
211  // convert our jets and add to the overall jet vector
212  for (unsigned int ijet=0;ijet<fjJets_.size();++ijet) {
213  // get the constituents from fastjet
214  std::vector<fastjet::PseudoJet> fjConstituents = sorted_by_pt(fjClusterSeq_->constituents(fjJets_[ijet]));
215  // convert them to CandidatePtr vector
216  std::vector<reco::CandidatePtr> constituents = getConstituents(fjConstituents);
217  // fill the trackjet
219  // write the specifics to the jet (simultaneously sets 4-vector, vertex).
220  writeSpecific( jet,
221  reco::Particle::LorentzVector(fjJets_[ijet].px(), fjJets_[ijet].py(), fjJets_[ijet].pz(), fjJets_[ijet].E()),
222  vertex_, constituents, iSetup);
223  jet.setJetArea(0);
224  jet.setPileup(0);
225  jet.setPrimaryVertex(edm::Ref<reco::VertexCollection>(pvCollection, (int) (itVtx-pvCollection->begin())));
226  jet.setVertex(itVtx->position());
227  jets->push_back(jet);
228  }
229 
230  if (useOnlyOnePV_) break; // stop vertex loop if only one vertex asked for
231  } // end loop over vertices
232 
233  // put the jets in the collection
234  LogDebug("FastjetTrackJetProducer") << "Put " << jets->size() << " jets in the event.\n";
235  iEvent.put(jets);
236 
237 }
238 
239 
240 //______________________________________________________________________________
242 {
243  // run algorithm
244  if ( !doAreaFastjet_ && !doRhoFastjet_) {
245  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequence( fjInputs_, *fjJetDefinition_ ) );
246  } else if (voronoiRfact_ <= 0) {
247  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceArea( fjInputs_, *fjJetDefinition_ , *fjActiveArea_ ) );
248  } else {
249  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceVoronoiArea( fjInputs_, *fjJetDefinition_ , fastjet::VoronoiAreaSpec(voronoiRfact_) ) );
250  }
251 
252  fjJets_ = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
253 
254 }
255 
256 
257 
259 // define as cmssw plugin
261 
263 
#define LogDebug(id)
JetType::Type jetTypeE
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
reco::Particle::Point vertex_
virtual std::vector< reco::CandidatePtr > getConstituents(const std::vector< fastjet::PseudoJet > &fjConstituents)
std::vector< fastjet::PseudoJet > fjJets_
virtual void setPileup(float fEnergy)
Set pileup energy contribution as calculated by algorithm.
Definition: Jet.h:96
virtual void inputTowers()
virtual void runAlgorithm(edm::Event &iEvent, const edm::EventSetup &iSetup)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool exists(std::string const &parameterName) const
checks if a parameter exists
FastjetJetProducer(const edm::ParameterSet &iConfig)
bool makeTrackJet(const JetType::Type &fTag)
virtual void setJetArea(float fArea)
set jet area
Definition: Jet.h:91
std::vector< fastjet::PseudoJet > fjInputs_
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)
int iEvent
Definition: GenABIO.cc:243
std::vector< edm::Ptr< reco::Candidate > > inputs_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:84
ClusterSequencePtr fjClusterSeq_
virtual void produceTrackJets(edm::Event &iEvent, const edm::EventSetup &iSetup)
void setPrimaryVertex(const reco::VertexRef &vtx)
set associated primary vertex
Definition: TrackJet.cc:81
ActiveAreaSpecPtr fjActiveArea_
virtual void setVertex(const Point &vertex)
set vertex
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:355
Jets made out of tracks.
Definition: TrackJet.h:28
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:38
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:26
boost::shared_ptr< fastjet::ClusterSequence > ClusterSequencePtr
void writeSpecific(reco::CaloJet &jet, reco::Particle::LorentzVector const &p4, reco::Particle::Point const &point, std::vector< reco::CandidatePtr > const &constituents, edm::EventSetup const &c)
Definition: JetSpecific.cc:41