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 #include "fastjet/tools/Filter.hh"
39 #include "fastjet/tools/Pruner.hh"
40 #include "fastjet/tools/MassDropTagger.hh"
41 
42 #include <iostream>
43 #include <memory>
44 #include <algorithm>
45 #include <limits>
46 #include <cmath>
47 //#include <fstream>
48 
49 using namespace std;
50 
51 
52 
54 // construction / destruction
56 
57 //______________________________________________________________________________
59  : VirtualJetProducer( iConfig ),
60  useMassDropTagger_(false),
61  useFiltering_(false),
62  useTrimming_(false),
63  usePruning_(false),
64  muCut_(-1.0),
65  yCut_(-1.0),
66  rFilt_(-1.0),
67  nFilt_(-1),
68  trimPtFracMin_(-1.0),
69  zCut_(-1.0),
70  RcutFactor_(-1.0)
71 {
72 
73  if ( iConfig.exists("UseOnlyVertexTracks") )
74  useOnlyVertexTracks_ = iConfig.getParameter<bool>("UseOnlyVertexTracks");
75  else
76  useOnlyVertexTracks_ = false;
77 
78  if ( iConfig.exists("UseOnlyOnePV") )
79  useOnlyOnePV_ = iConfig.getParameter<bool>("UseOnlyOnePV");
80  else
81  useOnlyOnePV_ = false;
82 
83  if ( iConfig.exists("DzTrVtxMax") )
84  dzTrVtxMax_ = iConfig.getParameter<double>("DzTrVtxMax");
85  else
86  dzTrVtxMax_ = 999999.;
87  if ( iConfig.exists("DxyTrVtxMax") )
88  dxyTrVtxMax_ = iConfig.getParameter<double>("DxyTrVtxMax");
89  else
90  dxyTrVtxMax_ = 999999.;
91  if ( iConfig.exists("MinVtxNdof") )
92  minVtxNdof_ = iConfig.getParameter<int>("MinVtxNdof");
93  else
94  minVtxNdof_ = 5;
95  if ( iConfig.exists("MaxVtxZ") )
96  maxVtxZ_ = iConfig.getParameter<double>("MaxVtxZ");
97  else
98  maxVtxZ_ = 15;
99 
100 
101  if ( iConfig.exists("useFiltering") ||
102  iConfig.exists("useTrimming") ||
103  iConfig.exists("usePruning") ||
104  iConfig.exists("useMassDropTagger") ) {
105  useMassDropTagger_=false;
106  useFiltering_=false;
107  useTrimming_=false;
108  usePruning_=false;
109  rFilt_=-1.0;
110  nFilt_=-1;
111  trimPtFracMin_=-1.0;
112  zCut_=-1.0;
113  RcutFactor_=-1.0;
114  muCut_=-1.0;
115  yCut_=-1.0;
116  useExplicitGhosts_ = true;
117 
118 
119  if ( iConfig.exists("useMassDropTagger") ) {
120  useMassDropTagger_ = true;
121  muCut_ = iConfig.getParameter<double>("muCut");
122  yCut_ = iConfig.getParameter<double>("yCut");
123  }
124 
125  if ( iConfig.exists("useFiltering") ) {
126  useFiltering_ = true;
127  rFilt_ = iConfig.getParameter<double>("rFilt");
128  nFilt_ = iConfig.getParameter<int>("nFilt");
129  }
130 
131  if ( iConfig.exists("useTrimming") ) {
132  useTrimming_ = true;
133  rFilt_ = iConfig.getParameter<double>("rFilt");
134  trimPtFracMin_ = iConfig.getParameter<double>("trimPtFracMin");
135  }
136 
137  if ( iConfig.exists("usePruning") ) {
138  usePruning_ = true;
139  zCut_ = iConfig.getParameter<double>("zcut");
140  RcutFactor_ = iConfig.getParameter<double>("rcut_factor");
141  nFilt_ = iConfig.getParameter<int>("nFilt");
142  }
143 
144  }
145 
146 }
147 
148 
149 //______________________________________________________________________________
151 {
152 }
153 
154 
156 // implementation of member functions
158 
160 {
161 
162  // for everything but track jets
163  if (!makeTrackJet(jetTypeE)) {
164 
165  // use the default production from one collection
166  VirtualJetProducer::produce( iEvent, iSetup );
167 
168  } else { // produce trackjets from tracks grouped per primary vertex
169 
170  produceTrackJets(iEvent, iSetup);
171 
172  }
173 
174 }
175 
176 
178 {
179 
180  // read in the track candidates
182  iEvent.getByLabel(src_, inputsHandle);
183  // make collection with pointers so we can play around with it
184  std::vector<edm::Ptr<reco::RecoChargedRefCandidate> > allInputs;
185  std::vector<edm::Ptr<reco::Candidate> > origInputs;
186  for (size_t i = 0; i < inputsHandle->size(); ++i) {
187  allInputs.push_back(inputsHandle->ptrAt(i));
188  origInputs.push_back(inputsHandle->ptrAt(i));
189  }
190 
191  // read in the PV collection
193  iEvent.getByLabel(srcPVs_, pvCollection);
194  // define the overall output jet container
195  std::auto_ptr<std::vector<reco::TrackJet> > jets(new std::vector<reco::TrackJet>() );
196 
197  // loop over the good vertices, clustering for each vertex separately
198  for (reco::VertexCollection::const_iterator itVtx = pvCollection->begin(); itVtx != pvCollection->end(); ++itVtx) {
199  if (itVtx->isFake() || itVtx->ndof() < minVtxNdof_ || fabs(itVtx->z()) > maxVtxZ_) continue;
200 
201  // clear the intermediate containers
202  inputs_.clear();
203  fjInputs_.clear();
204  fjJets_.clear();
205 
206  // if only vertex-associated tracks should be used
207  if (useOnlyVertexTracks_) {
208  // loop over the tracks associated to the vertex
209  for (reco::Vertex::trackRef_iterator itTr = itVtx->tracks_begin(); itTr != itVtx->tracks_end(); ++itTr) {
210  // whether a match was found in the track candidate input
211  bool found = false;
212  // loop over input track candidates
213  for (std::vector<edm::Ptr<reco::RecoChargedRefCandidate> >::iterator itIn = allInputs.begin(); itIn != allInputs.end(); ++itIn) {
214  // match the input track candidate to the track from the vertex
215  reco::TrackRef trref(itTr->castTo<reco::TrackRef>());
216  // check if the tracks match
217  if ((*itIn)->track() == trref) {
218  found = true;
219  // add this track candidate to the input for clustering
220  inputs_.push_back(*itIn);
221  // erase the track candidate from the total list of input, so we don't reuse it later
222  allInputs.erase(itIn);
223  // found the candidate track corresponding to the vertex track, so stop the loop
224  break;
225  } // end if match found
226  } // end loop over input tracks
227  // give an info message in case no match is found (can happen if candidates are subset of tracks used for clustering)
228  if (!found) edm::LogInfo("FastjetTrackJetProducer") << "Ignoring a track at vertex which is not in input track collection!";
229  } // end loop over tracks associated to vertex
230  // if all inpt track candidates should be used
231  } else {
232  // loop over input track candidates
233  for (std::vector<edm::Ptr<reco::RecoChargedRefCandidate> >::iterator itIn = allInputs.begin(); itIn != allInputs.end(); ++itIn) {
234  // check if the track is close enough to the vertex
235  float dz = (*itIn)->track()->dz(itVtx->position());
236  float dxy = (*itIn)->track()->dxy(itVtx->position());
237  if (fabs(dz) > dzTrVtxMax_) continue;
238  if (fabs(dxy) > dxyTrVtxMax_) continue;
239  bool closervtx = false;
240  // now loop over the good vertices a second time
241  for (reco::VertexCollection::const_iterator itVtx2 = pvCollection->begin(); itVtx2 != pvCollection->end(); ++itVtx2) {
242  if (itVtx->isFake() || itVtx->ndof() < minVtxNdof_ || fabs(itVtx->z()) > maxVtxZ_) continue;
243  // and check this track is closer to any other vertex (if more than 1 vertex considered)
244  if (!useOnlyOnePV_ &&
245  itVtx != itVtx2 &&
246  fabs((*itIn)->track()->dz(itVtx2->position())) < fabs(dz)) {
247  closervtx = true;
248  break; // 1 closer vertex makes the track already not matched, so break
249  }
250  }
251  // don't add this track if another vertex is found closer
252  if (closervtx) continue;
253  // add this track candidate to the input for clustering
254  inputs_.push_back(*itIn);
255  // erase the track candidate from the total list of input, so we don't reuse it later
256  allInputs.erase(itIn);
257  // take a step back in the loop since we just erased
258  --itIn;
259  }
260  }
261 
262  // convert candidates in inputs_ to fastjet::PseudoJets in fjInputs_
263  fjInputs_.reserve(inputs_.size());
264  inputTowers();
265  LogDebug("FastjetTrackJetProducer") << "Inputted towers\n";
266 
267  // run algorithm, using fjInputs_, modifying fjJets_ and allocating fjClusterSeq_
268  runAlgorithm(iEvent, iSetup);
269  LogDebug("FastjetTrackJetProducer") << "Ran algorithm\n";
270 
271  // convert our jets and add to the overall jet vector
272  for (unsigned int ijet=0;ijet<fjJets_.size();++ijet) {
273  // get the constituents from fastjet
274  std::vector<fastjet::PseudoJet> fjConstituents = sorted_by_pt(fjClusterSeq_->constituents(fjJets_[ijet]));
275  // convert them to CandidatePtr vector
276  std::vector<reco::CandidatePtr> constituents = getConstituents(fjConstituents);
277  // fill the trackjet
279  // write the specifics to the jet (simultaneously sets 4-vector, vertex).
280  writeSpecific( jet,
281  reco::Particle::LorentzVector(fjJets_[ijet].px(), fjJets_[ijet].py(), fjJets_[ijet].pz(), fjJets_[ijet].E()),
282  vertex_, constituents, iSetup);
283  jet.setJetArea(0);
284  jet.setPileup(0);
285  jet.setPrimaryVertex(edm::Ref<reco::VertexCollection>(pvCollection, (int) (itVtx-pvCollection->begin())));
286  jet.setVertex(itVtx->position());
287  jets->push_back(jet);
288  }
289 
290  if (useOnlyOnePV_) break; // stop vertex loop if only one vertex asked for
291  } // end loop over vertices
292 
293  // put the jets in the collection
294  LogDebug("FastjetTrackJetProducer") << "Put " << jets->size() << " jets in the event.\n";
295  iEvent.put(jets);
296 
297 }
298 
299 
300 //______________________________________________________________________________
302 {
303  // run algorithm
304  /*
305  fjInputs_.clear();
306  double px, py , pz, E;
307  string line;
308  std::ifstream fin("dump3.txt");
309  while (getline(fin, line)){
310  if (line == "#END") break;
311  if (line.substr(0,1) == "#") {continue;}
312  istringstream istr(line);
313  istr >> px >> py >> pz >> E;
314  // create a fastjet::PseudoJet with these components and put it onto
315  // back of the input_particles vector
316  fastjet::PseudoJet j(px,py,pz,E);
317  //if ( fabs(j.rap()) < inputEtaMax )
318  fjInputs_.push_back(fastjet::PseudoJet(px,py,pz,E));
319  }
320  fin.close();
321  */
322 
323  if ( !doAreaFastjet_ && !doRhoFastjet_) {
324  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequence( fjInputs_, *fjJetDefinition_ ) );
325  } else if (voronoiRfact_ <= 0) {
326  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceArea( fjInputs_, *fjJetDefinition_ , *fjAreaDefinition_ ) );
327  } else {
328  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceVoronoiArea( fjInputs_, *fjJetDefinition_ , fastjet::VoronoiAreaSpec(voronoiRfact_) ) );
329  }
330 
331  if ( !useTrimming_ && !useFiltering_ && !usePruning_ ) {
332  fjJets_ = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
333  }
334  else {
335  fjJets_.clear();
336  std::vector<fastjet::PseudoJet> tempJets = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
337 
338  fastjet::MassDropTagger md_tagger( muCut_, yCut_ );
339  fastjet::Filter trimmer( fastjet::Filter(fastjet::JetDefinition(fastjet::kt_algorithm, rFilt_), fastjet::SelectorPtFractionMin(trimPtFracMin_)));
340  fastjet::Filter filter( fastjet::Filter(fastjet::JetDefinition(fastjet::cambridge_algorithm, rFilt_), fastjet::SelectorNHardest(nFilt_)));
341  fastjet::Pruner pruner(fastjet::cambridge_algorithm, zCut_, RcutFactor_);
342 
343  std::vector<fastjet::Transformer const *> transformers;
344 
345  if ( useMassDropTagger_ ) {
346  transformers.push_back(&md_tagger);
347  }
348  if ( useTrimming_ ) {
349  transformers.push_back(&trimmer);
350  }
351  if ( useFiltering_ ) {
352  transformers.push_back(&filter);
353  }
354  if ( usePruning_ ) {
355  transformers.push_back(&pruner);
356  }
357 
358 
359 
360  for ( std::vector<fastjet::PseudoJet>::const_iterator ijet = tempJets.begin(),
361  ijetEnd = tempJets.end(); ijet != ijetEnd; ++ijet ) {
362 
363  fastjet::PseudoJet transformedJet = *ijet;
364  bool passed = true;
365  for ( std::vector<fastjet::Transformer const *>::const_iterator itransf = transformers.begin(),
366  itransfEnd = transformers.end(); itransf != itransfEnd; ++itransf ) {
367  if ( transformedJet != 0 ) {
368  transformedJet = (**itransf)(transformedJet);
369  }
370  else {
371  passed=false;
372  }
373  }
374  if ( passed )
375  fjJets_.push_back( transformedJet );
376  }
377  }
378 
379 }
380 
381 
382 
384 // define as cmssw plugin
386 
388 
#define LogDebug(id)
JetType::Type jetTypeE
T getParameter(std::string const &) const
bool useFiltering_
Mass-drop tagging for boosted Higgs.
int i
Definition: DBlmapReader.cc:9
double rFilt_
for mass-drop tagging, symmetry cut: min(pt1^2,pt2^2) * dR(1,2) / mjet &gt; ycut
tuple Filter
Definition: Filter_cff.py:5
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:109
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
double RcutFactor_
for pruning: constituent minimum pt fraction of parent cluster
FastjetJetProducer(const edm::ParameterSet &iConfig)
bool useTrimming_
Jet filtering technique.
bool makeTrackJet(const JetType::Type &fTag)
virtual void setJetArea(float fArea)
set jet area
Definition: Jet.h:104
std::vector< fastjet::PseudoJet > fjInputs_
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)
double trimPtFracMin_
for filtering, pruning: number of subjets expected
int iEvent
Definition: GenABIO.cc:243
bool usePruning_
Jet trimming technique.
std::vector< edm::Ptr< reco::Candidate > > inputs_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
ClusterSequencePtr fjClusterSeq_
virtual void produceTrackJets(edm::Event &iEvent, const edm::EventSetup &iSetup)
vector< PseudoJet > jets
double muCut_
Jet pruning technique.
void setPrimaryVertex(const reco::VertexRef &vtx)
set associated primary vertex
Definition: TrackJet.cc:81
double zCut_
for trimming: constituent minimum pt fraction of full jet
virtual void setVertex(const Point &vertex)
set vertex
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
Jets made out of tracks.
Definition: TrackJet.h:28
int nFilt_
for filtering, trimming: dR scale of sub-clustering
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)
double yCut_
for mass-drop tagging, m0/mjet (m0 = mass of highest mass subjet)
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:38
AreaDefinitionPtr fjAreaDefinition_
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:25
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