CMS 3D CMS Logo

CATopJetAlgorithm.cc
Go to the documentation of this file.
1 // Original author: Brock Tweedie (JHU)
2 // Ported to CMSSW by: Sal Rappoccio (JHU)
3 
11 
12 using namespace std;
13 using namespace reco;
14 using namespace edm;
15 
16 
17 // Run the algorithm
18 // ------------------
19 void CATopJetAlgorithm::run( const vector<fastjet::PseudoJet> & cell_particles,
20  vector<fastjet::PseudoJet> & hardjetsOutput,
21  boost::shared_ptr<fastjet::ClusterSequence> & fjClusterSeq
22  )
23 {
24  if ( verbose_ ) cout << "Welcome to CATopSubJetAlgorithm::run" << endl;
25 
26  // Sum Et of the event
27  double sumEt = 0.;
28 
29  //make a list of input objects ordered by ET and calculate sum et
30  // list of fastjet pseudojet constituents
31  for (unsigned i = 0; i < cell_particles.size(); ++i) {
32  sumEt += cell_particles[i].perp();
33  }
34 
35  // Determine which bin we are in for et clustering
36  int sumEtBinId = -1;
37  for ( unsigned int i = 0; i < sumEtBins_.size(); ++i ) {
38  if ( sumEt > sumEtBins_[i] ) sumEtBinId = i;
39  }
40  if ( verbose_ ) cout << "Using sumEt = " << sumEt << ", bin = " << sumEtBinId << endl;
41 
42  // If the sum et is too low, exit
43  if ( sumEtBinId < 0 ) {
44  return;
45  }
46 
47  // empty 4-vector
48  fastjet::PseudoJet blankJetA(0,0,0,0);
49  blankJetA.set_user_index(-1);
50  const fastjet::PseudoJet blankJet = blankJetA;
51 
52  // Define adjacency variables which depend on which sumEtBin we are in
53  double deltarcut = deltarBins_[sumEtBinId];
54  double nCellMin = nCellBins_[sumEtBinId];
55 
56  if ( verbose_ )cout<<"useAdjacency_ = "<<useAdjacency_<<endl;
57  if ( verbose_ && useAdjacency_==0)cout<<"No Adjacency"<<endl;
58  if ( verbose_ && useAdjacency_==1)cout<<"using deltar adjacency"<<endl;
59  if ( verbose_ && useAdjacency_==2)cout<<"using modified adjacency"<<endl;
60  if ( verbose_ && useAdjacency_==3)cout<<"using calorimeter nearest neighbor based adjacency"<<endl;
61  if ( verbose_ && useAdjacency_==1)cout << "Using deltarcut = " << deltarcut << endl;
62  if ( verbose_ && useAdjacency_==3)cout << "Using nCellMin = " << nCellMin << endl;
63 
64  if ( verbose_ ) cout << "About to do jet clustering in CA" << endl;
65  // run the jet clustering
66 
67  //cluster the jets with the jet definition jetDef:
68  // run algorithm
69  // boost::shared_ptr<fastjet::ClusterSequence> fjClusterSeq;
70  // if ( !doAreaFastjet_ ) {
71  // fjClusterSeq = boost::shared_ptr<fastjet::ClusterSequence>( new fastjet::ClusterSequence( cell_particles, jetDef ) );
72  // } else if (voronoiRfact_ <= 0) {
73  // fjClusterSeq = boost::shared_ptr<fastjet::ClusterSequence>( new fastjet::ClusterSequenceArea( cell_particles, jetDef , *fjActiveArea_ ) );
74  // } else {
75  // fjClusterSeq = boost::shared_ptr<fastjet::ClusterSequence>( new fastjet::ClusterSequenceVoronoiArea( cell_particles, jetDef , fastjet::VoronoiAreaSpec(voronoiRfact_) ) );
76  // }
77 
78  if ( verbose_ ) cout << "Getting inclusive jets" << endl;
79  // Get the transient inclusive jets
80  vector<fastjet::PseudoJet> inclusiveJets = fjClusterSeq->inclusive_jets(ptMin_);
81 
82  if ( verbose_ ) cout << "Getting central jets" << endl;
83  // Find the transient central jets
84  vector<fastjet::PseudoJet> centralJets;
85  for (unsigned int i = 0; i < inclusiveJets.size(); i++) {
86 
87  if (inclusiveJets[i].perp() > ptMin_ && fabs(inclusiveJets[i].rapidity()) < centralEtaCut_) {
88  centralJets.push_back(inclusiveJets[i]);
89  }
90  }
91  // Sort the transient central jets in Et
92  sort( centralJets.begin(), centralJets.end(), greaterByEtPseudoJet );
93 
94  // These will store the 4-vectors of each hard jet
95  vector<math::XYZTLorentzVector> p4_hardJets;
96 
97  // These will store the indices of each subjet that
98  // are present in each jet
99  vector<vector<int> > indices( centralJets.size() );
100 
101  // Loop over central jets, attempt to find substructure
102  vector<fastjet::PseudoJet>::iterator jetIt = centralJets.begin(),
103  centralJetsEnd = centralJets.end();
104  if ( verbose_ )cout<<"Loop over jets"<<endl;
105  int i=0;
106  for ( ; jetIt != centralJetsEnd; ++jetIt ) {
107  if ( verbose_ )cout<<"\nJet "<<i<<endl;
108  i++;
109  fastjet::PseudoJet localJet = *jetIt;
110 
111  // Get the 4-vector for this jet
112  p4_hardJets.push_back( math::XYZTLorentzVector(localJet.px(), localJet.py(), localJet.pz(), localJet.e() ));
113 
114  // jet decomposition. try to find 3 or 4 hard, well-localized subjets, characteristic of a boosted top.
115  if ( verbose_ )cout<<"local jet pt = "<<localJet.perp()<<endl;
116  if ( verbose_ )cout<<"deltap = "<<ptFracBins_[sumEtBinId]<<endl;
117 
118  double ptHard = ptFracBins_[sumEtBinId]*localJet.perp();
119  vector<fastjet::PseudoJet> leftoversAll;
120 
121  // stage 1: primary decomposition. look for when the jet declusters into two hard subjets
122  if ( verbose_ ) cout << "Doing decomposition 1" << endl;
123  fastjet::PseudoJet ja, jb;
124  vector<fastjet::PseudoJet> leftovers1;
125  bool hardBreak1 = decomposeJet(localJet,*fjClusterSeq,cell_particles,ptHard,nCellMin,deltarcut,ja,jb,leftovers1);
126  leftoversAll.insert(leftoversAll.end(),leftovers1.begin(),leftovers1.end());
127 
128  // stage 2: secondary decomposition. look for when the hard subjets found above further decluster into two hard sub-subjets
129  //
130  // ja -> jaa+jab ?
131  if ( verbose_ ) cout << "Doing decomposition 2. ja->jaa+jab?" << endl;
132  fastjet::PseudoJet jaa, jab;
133  vector<fastjet::PseudoJet> leftovers2a;
134  bool hardBreak2a = false;
135  if (hardBreak1) hardBreak2a = decomposeJet(ja,*fjClusterSeq,cell_particles,ptHard,nCellMin,deltarcut,jaa,jab,leftovers2a);
136  leftoversAll.insert(leftoversAll.end(),leftovers2a.begin(),leftovers2a.end());
137  // jb -> jba+jbb ?
138  if ( verbose_ ) cout << "Doing decomposition 2. ja->jba+jbb?" << endl;
139  fastjet::PseudoJet jba, jbb;
140  vector<fastjet::PseudoJet> leftovers2b;
141  bool hardBreak2b = false;
142  if (hardBreak1) hardBreak2b = decomposeJet(jb,*fjClusterSeq,cell_particles,ptHard,nCellMin,deltarcut,jba,jbb,leftovers2b);
143  leftoversAll.insert(leftoversAll.end(),leftovers2b.begin(),leftovers2b.end());
144 
145  // NOTE: it might be good to consider some checks for whether these subjets can be further decomposed. e.g., the above procedure leaves
146  // open the possibility of "subjets" that actually consist of two or more distinct hard clusters. however, this kind of thing
147  // is a rarity for the simulations so far considered.
148 
149  // proceed if one or both of the above hard subjets successfully decomposed
150  if ( verbose_ ) cout << "Done with decomposition" << endl;
151 
152  if ( verbose_ ) cout<<"hardBreak1 = "<<hardBreak1<<endl;
153  if ( verbose_ ) cout<<"hardBreak2a = "<<hardBreak2a<<endl;
154  if ( verbose_ ) cout<<"hardBreak2b = "<<hardBreak2b<<endl;
155 
156  fastjet::PseudoJet hardA = blankJet, hardB = blankJet, hardC = blankJet, hardD = blankJet;
157  if (!hardBreak1) {
158  hardA = localJet;
159  hardB = blankJet;
160  hardC = blankJet;
161  hardD = blankJet;
162  if(verbose_)cout<<"Hardbreak failed. Save subjet1=localJet"<<endl;
163  }
164  if (hardBreak1 && !hardBreak2a && !hardBreak2b) {
165  hardA = ja;
166  hardB = jb;
167  hardC = blankJet;
168  hardD = blankJet;
169  if(verbose_)cout<<"First decomposition succeeded, both second decompositions failed. Save subjet1=ja subjet2=jb"<<endl;
170  }
171  if (hardBreak1 && hardBreak2a && !hardBreak2b) {
172  hardA = jaa;
173  hardB = jab;
174  hardC = jb;
175  hardD = blankJet;
176  if(verbose_)cout<<"First decomposition succeeded, ja split succesfully, jb did not split. Save subjet1=jaa subjet2=jab subjet3=jb"<<endl;
177  }
178  if (hardBreak1 && !hardBreak2a && hardBreak2b) {
179  hardA = jba;
180  hardB = jbb;
181  hardC = ja;
182  hardD = blankJet;
183  if(verbose_)cout<<"First decomposition succeeded, jb split succesfully, ja did not split. Save subjet1=jba subjet2=jbb subjet3=ja"<<endl;
184  }
185  if (hardBreak1 && hardBreak2a && hardBreak2b) {
186  hardA = jaa;
187  hardB = jab;
188  hardC = jba;
189  hardD = jbb;
190  if(verbose_)cout<<"First decomposition and both secondary decompositions succeeded. Save subjet1=jaa subjet2=jab subjet3=jba subjet4=jbb"<<endl;
191  }
192 
193  // check if we are left with >= 3 hard subjets
194  fastjet::PseudoJet subjet1 = blankJet;
195  fastjet::PseudoJet subjet2 = blankJet;
196  fastjet::PseudoJet subjet3 = blankJet;
197  fastjet::PseudoJet subjet4 = blankJet;
198  subjet1 = hardA; subjet2 = hardB; subjet3 = hardC; subjet4 = hardD;
199 
200  // record the hard subjets
201  vector<fastjet::PseudoJet> hardSubjets;
202 
203  if ( verbose_ ) {
204  std::cout << "HardA : user_index = " << hardA.user_index() << ", (Pt,Y,Phi,M) = ("
205  << hardA.pt() << ", " << hardA.rapidity() << ", "
206  << hardA.phi() << ", " << hardA.m() << ")" << std::endl;
207 
208  std::cout << "HardB : user_index = " << hardB.user_index() << ", (Pt,Y,Phi,M) = ("
209  << hardB.pt() << ", " << hardB.rapidity() << ", "
210  << hardB.phi() << ", " << hardB.m() << ")" << std::endl;
211 
212  std::cout << "HardC : user_index = " << hardC.user_index() << ", (Pt,Y,Phi,M) = ("
213  << hardC.pt() << ", " << hardC.rapidity() << ", "
214  << hardC.phi() << ", " << hardC.m() << ")" << std::endl;
215 
216  std::cout << "HardD : user_index = " << hardD.user_index() << ", (Pt,Y,Phi,M) = ("
217  << hardD.pt() << ", " << hardD.rapidity() << ", "
218  << hardD.phi() << ", " << hardD.m() << ")" << std::endl;
219  }
220 
221  // Check to see if any subjects are counted amongst the "hard" subjets from previous
222  // lines. NOTE: In Fastjet 3.0, the default "user_index" changed from 0 to -1, so
223  // this can no longer be used as a designator for the veto of "blankJet" subjets,
224  // and now switch to pt > some small value.
225  if ( subjet1.pt() > 0.0001 )
226  hardSubjets.push_back(subjet1);
227  if ( subjet2.pt() > 0.0001 )
228  hardSubjets.push_back(subjet2);
229  if ( subjet3.pt() > 0.0001 )
230  hardSubjets.push_back(subjet3);
231  if ( subjet4.pt() > 0.0001 )
232  hardSubjets.push_back(subjet4);
233  sort(hardSubjets.begin(), hardSubjets.end(), greaterByEtPseudoJet );
234 
235  // Use new fastjet functionality to create a Pseudojet from constituents
236  fastjet::PseudoJet candidate = join(hardSubjets);
237  // Reset the jet's 4-vector to the "ungroomed" value
238  candidate.reset_momentum( jetIt->px(), jetIt->py(), jetIt->pz(), jetIt->e() );
239 
240  if ( verbose_ ) {
241  std::cout << "Final top-jet candidate: (Pt,Y,Phi,M) = ("
242  << candidate.pt() << ", " << candidate.rapidity() << ", "
243  << candidate.phi() << ", " << candidate.m() << ")" << std::endl;
244  std::vector<fastjet::PseudoJet> pieces = candidate.pieces();
245  std::cout << "Number of pieces = " << pieces.size() << std::endl;
246  for ( std::vector<fastjet::PseudoJet>::const_iterator ibegin = pieces.begin(), iend = pieces.end(), i = ibegin;
247  i != iend; ++i ) {
248  std::cout << " Piece : " << i - ibegin << ", (Pt,Y,Phi,M) = ("
249  << i->pt() << ", " << i->rapidity() << ", "
250  << i->phi() << ", " << i->m() << ")" << std::endl;
251  }
252  }
253  // Add to the list
254  hardjetsOutput.push_back( candidate );
255  }
256 }
257 
258 
259 
260 
261 //-----------------------------------------------------------------------
262 // determine whether two clusters (made of calorimeter towers) are living on "adjacent" cells. if they are, then
263 // we probably shouldn't consider them to be independent objects!
264 //
265 // From Sal: Ignoring genjet case
266 //
267 bool CATopJetAlgorithm::adjacentCells(const fastjet::PseudoJet & jet1, const fastjet::PseudoJet & jet2,
268  const vector<fastjet::PseudoJet> & cell_particles,
269  const fastjet::ClusterSequence & theClusterSequence,
270  double nCellMin ) const {
271 
272 
273  double eta1 = jet1.rapidity();
274  double phi1 = jet1.phi();
275  double eta2 = jet2.rapidity();
276  double phi2 = jet2.phi();
277 
278  double deta = abs(eta2 - eta1) / 0.087;
279  double dphi = fabs( reco::deltaPhi(phi2,phi1) ) / 0.087;
280 
281  return ( ( deta + dphi ) <= nCellMin );
282 }
283 
284 
285 
286 //-------------------------------------------------------------------------
287 // attempt to decompose a jet into "hard" subjets, where hardness is set by ptHard
288 //
289 bool CATopJetAlgorithm::decomposeJet(const fastjet::PseudoJet & theJet,
290  const fastjet::ClusterSequence & theClusterSequence,
291  const vector<fastjet::PseudoJet> & cell_particles,
292  double ptHard, double nCellMin, double deltarcut,
293  fastjet::PseudoJet & ja, fastjet::PseudoJet & jb,
294  vector<fastjet::PseudoJet> & leftovers) const {
295 
296  bool goodBreak;
297  fastjet::PseudoJet j = theJet;
298  double InputObjectPt = j.perp();
299  if ( verbose_ )cout<<"Input Object Pt = "<<InputObjectPt<<endl;
300  if ( verbose_ )cout<<"ptHard = "<<ptHard<<endl;
301  leftovers.clear();
302  if ( verbose_ )cout<<"start while loop"<<endl;
303 
304  while (true) { // watch out for infinite loop!
305  goodBreak = theClusterSequence.has_parents(j,ja,jb);
306  if (!goodBreak){
307  if ( verbose_ )cout<<"bad break. this is one cell. can't decluster anymore."<<endl;
308  break; // this is one cell, can't decluster anymore
309  }
310 
311  if ( verbose_ )cout<<"good break. ja Pt = "<<ja.perp()<<" jb Pt = "<<jb.perp()<<endl;
312 
314 
315  // check if clusters are adjacent using a constant deltar adjacency.
316  double clusters_deltar=fabs(ja.eta()-jb.eta())+fabs(deltaPhi(ja.phi(),jb.phi()));
317 
318  if ( verbose_ && useAdjacency_ ==1)cout<<"clusters_deltar = "<<clusters_deltar<<endl;
319  if ( verbose_ && useAdjacency_ ==1)cout<<"deltar cut = "<<deltarcut<<endl;
320 
321  if ( useAdjacency_==1 && clusters_deltar < deltarcut){
322  if ( verbose_ )cout<<"clusters too close. consant adj. break."<<endl;
323  break;
324  }
325 
326  // Check if clusters are adjacent using a DeltaR adjacency which is a function of pT.
327  double clusters_deltaR=deltaR( ja.rapidity(), ja.phi(), jb.rapidity(), jb.phi() );
328 
329  if ( verbose_ && useAdjacency_ ==2)cout<<"clusters_deltaR = "<<clusters_deltaR<<endl;
330  if ( verbose_ && useAdjacency_ ==2)cout<<"0.4-0.0004*InputObjectPt = "<<0.4-0.0004*InputObjectPt<<endl;
331 
332  if ( useAdjacency_==2 && clusters_deltaR < 0.4-0.0004*InputObjectPt)
333  {
334  if ( verbose_ )cout<<"clusters too close. modified adj. break."<<endl;
335  break;
336  }
337 
338  // Check if clusters are adjacent in the calorimeter.
339  if ( useAdjacency_==3 && adjacentCells(ja,jb,cell_particles,theClusterSequence,nCellMin) ){
340  if ( verbose_ )cout<<"clusters too close in the calorimeter. calorimeter adj. break."<<endl;
341  break; // the clusters are "adjacent" in the calorimeter => shouldn't have decomposed
342  }
343 
344  if ( verbose_ )cout<<"clusters pass distance cut"<<endl;
345 
347 
348  if ( verbose_ )cout<<"ptHard = "<<ptHard<<endl;
349 
350  if (ja.perp() < ptHard && jb.perp() < ptHard){
351  if ( verbose_ )cout<<"two soft clusters. dead end"<<endl;
352  break; // broke into two soft clusters, dead end
353  }
354 
355  if (ja.perp() > ptHard && jb.perp() > ptHard){
356  if ( verbose_ )cout<<"two hard clusters. done"<<endl;
357  return true; // broke into two hard clusters, we're done!
358  }
359 
360  else if (ja.perp() > jb.perp()) { // broke into one hard and one soft, ditch the soft one and try again
361  if ( verbose_ )cout<<"ja hard jb soft. try to split hard. j = ja"<<endl;
362  j = ja;
363  vector<fastjet::PseudoJet> particles = theClusterSequence.constituents(jb);
364  leftovers.insert(leftovers.end(),particles.begin(),particles.end());
365  }
366  else {
367  if ( verbose_ )cout<<"ja hard jb soft. try to split hard. j = jb"<<endl;
368  j = jb;
369  vector<fastjet::PseudoJet> particles = theClusterSequence.constituents(ja);
370  leftovers.insert(leftovers.end(),particles.begin(),particles.end());
371  }
372  }
373 
374  if ( verbose_ )cout<<"did not decluster."<<endl; // did not decluster into hard subjets
375 
376  ja.reset(0,0,0,0);
377  jb.reset(0,0,0,0);
378  leftovers.clear();
379  return false;
380 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
bool greaterByEtPseudoJet(fastjet::PseudoJet const &j1, fastjet::PseudoJet const &j2)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
bool adjacentCells(const fastjet::PseudoJet &jet1, const fastjet::PseudoJet &jet2, const std::vector< fastjet::PseudoJet > &cell_particles, const fastjet::ClusterSequence &theClusterSequence, double nCellMin) const
bool decomposeJet(const fastjet::PseudoJet &theJet, const fastjet::ClusterSequence &theClusterSequence, const std::vector< fastjet::PseudoJet > &cell_particles, double ptHard, double nCellMin, double deltarcut, fastjet::PseudoJet &ja, fastjet::PseudoJet &jb, std::vector< fastjet::PseudoJet > &leftovers) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
static std::string join(char **cmd)
Definition: RemoteFile.cc:18
T perp() const
Magnitude of transverse component.
fixed size matrix
HLT enums.
void run(const std::vector< fastjet::PseudoJet > &cell_particles, std::vector< fastjet::PseudoJet > &hardjetsOutput, boost::shared_ptr< fastjet::ClusterSequence > &fjClusterSeq)
Find the ProtoJets from the collection of input Candidates.