CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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  GreaterByEtPseudoJet compEt;
93  sort( centralJets.begin(), centralJets.end(), compEt );
94 
95  // These will store the 4-vectors of each hard jet
96  vector<math::XYZTLorentzVector> p4_hardJets;
97 
98  // These will store the indices of each subjet that
99  // are present in each jet
100  vector<vector<int> > indices( centralJets.size() );
101 
102  // Loop over central jets, attempt to find substructure
103  vector<fastjet::PseudoJet>::iterator jetIt = centralJets.begin(),
104  centralJetsEnd = centralJets.end();
105  if ( verbose_ )cout<<"Loop over jets"<<endl;
106  int i=0;
107  for ( ; jetIt != centralJetsEnd; ++jetIt ) {
108  if ( verbose_ )cout<<"\nJet "<<i<<endl;
109  i++;
110  fastjet::PseudoJet localJet = *jetIt;
111 
112  // Get the 4-vector for this jet
113  p4_hardJets.push_back( math::XYZTLorentzVector(localJet.px(), localJet.py(), localJet.pz(), localJet.e() ));
114 
115  // jet decomposition. try to find 3 or 4 hard, well-localized subjets, characteristic of a boosted top.
116  if ( verbose_ )cout<<"local jet pt = "<<localJet.perp()<<endl;
117  if ( verbose_ )cout<<"deltap = "<<ptFracBins_[sumEtBinId]<<endl;
118 
119  double ptHard = ptFracBins_[sumEtBinId]*localJet.perp();
120  vector<fastjet::PseudoJet> leftoversAll;
121 
122  // stage 1: primary decomposition. look for when the jet declusters into two hard subjets
123  if ( verbose_ ) cout << "Doing decomposition 1" << endl;
124  fastjet::PseudoJet ja, jb;
125  vector<fastjet::PseudoJet> leftovers1;
126  bool hardBreak1 = decomposeJet(localJet,*fjClusterSeq,cell_particles,ptHard,nCellMin,deltarcut,ja,jb,leftovers1);
127  leftoversAll.insert(leftoversAll.end(),leftovers1.begin(),leftovers1.end());
128 
129  // stage 2: secondary decomposition. look for when the hard subjets found above further decluster into two hard sub-subjets
130  //
131  // ja -> jaa+jab ?
132  if ( verbose_ ) cout << "Doing decomposition 2. ja->jaa+jab?" << endl;
133  fastjet::PseudoJet jaa, jab;
134  vector<fastjet::PseudoJet> leftovers2a;
135  bool hardBreak2a = false;
136  if (hardBreak1) hardBreak2a = decomposeJet(ja,*fjClusterSeq,cell_particles,ptHard,nCellMin,deltarcut,jaa,jab,leftovers2a);
137  leftoversAll.insert(leftoversAll.end(),leftovers2a.begin(),leftovers2a.end());
138  // jb -> jba+jbb ?
139  if ( verbose_ ) cout << "Doing decomposition 2. ja->jba+jbb?" << endl;
140  fastjet::PseudoJet jba, jbb;
141  vector<fastjet::PseudoJet> leftovers2b;
142  bool hardBreak2b = false;
143  if (hardBreak1) hardBreak2b = decomposeJet(jb,*fjClusterSeq,cell_particles,ptHard,nCellMin,deltarcut,jba,jbb,leftovers2b);
144  leftoversAll.insert(leftoversAll.end(),leftovers2b.begin(),leftovers2b.end());
145 
146  // NOTE: it might be good to consider some checks for whether these subjets can be further decomposed. e.g., the above procedure leaves
147  // open the possibility of "subjets" that actually consist of two or more distinct hard clusters. however, this kind of thing
148  // is a rarity for the simulations so far considered.
149 
150  // proceed if one or both of the above hard subjets successfully decomposed
151  if ( verbose_ ) cout << "Done with decomposition" << endl;
152 
153  if ( verbose_ ) cout<<"hardBreak1 = "<<hardBreak1<<endl;
154  if ( verbose_ ) cout<<"hardBreak2a = "<<hardBreak2a<<endl;
155  if ( verbose_ ) cout<<"hardBreak2b = "<<hardBreak2b<<endl;
156 
157  fastjet::PseudoJet hardA = blankJet, hardB = blankJet, hardC = blankJet, hardD = blankJet;
158  if (!hardBreak1) {
159  hardA = localJet;
160  hardB = blankJet;
161  hardC = blankJet;
162  hardD = blankJet;
163  if(verbose_)cout<<"Hardbreak failed. Save subjet1=localJet"<<endl;
164  }
165  if (hardBreak1 && !hardBreak2a && !hardBreak2b) {
166  hardA = ja;
167  hardB = jb;
168  hardC = blankJet;
169  hardD = blankJet;
170  if(verbose_)cout<<"First decomposition succeeded, both second decompositions failed. Save subjet1=ja subjet2=jb"<<endl;
171  }
172  if (hardBreak1 && hardBreak2a && !hardBreak2b) {
173  hardA = jaa;
174  hardB = jab;
175  hardC = jb;
176  hardD = blankJet;
177  if(verbose_)cout<<"First decomposition succeeded, ja split succesfully, jb did not split. Save subjet1=jaa subjet2=jab subjet3=jb"<<endl;
178  }
179  if (hardBreak1 && !hardBreak2a && hardBreak2b) {
180  hardA = jba;
181  hardB = jbb;
182  hardC = ja;
183  hardD = blankJet;
184  if(verbose_)cout<<"First decomposition succeeded, jb split succesfully, ja did not split. Save subjet1=jba subjet2=jbb subjet3=ja"<<endl;
185  }
186  if (hardBreak1 && hardBreak2a && hardBreak2b) {
187  hardA = jaa;
188  hardB = jab;
189  hardC = jba;
190  hardD = jbb;
191  if(verbose_)cout<<"First decomposition and both secondary decompositions succeeded. Save subjet1=jaa subjet2=jab subjet3=jba subjet4=jbb"<<endl;
192  }
193 
194  // check if we are left with >= 3 hard subjets
195  fastjet::PseudoJet subjet1 = blankJet;
196  fastjet::PseudoJet subjet2 = blankJet;
197  fastjet::PseudoJet subjet3 = blankJet;
198  fastjet::PseudoJet subjet4 = blankJet;
199  subjet1 = hardA; subjet2 = hardB; subjet3 = hardC; subjet4 = hardD;
200 
201  // record the hard subjets
202  vector<fastjet::PseudoJet> hardSubjets;
203 
204  if ( verbose_ ) {
205  std::cout << "HardA : user_index = " << hardA.user_index() << ", (Pt,Y,Phi,M) = ("
206  << hardA.pt() << ", " << hardA.rapidity() << ", "
207  << hardA.phi() << ", " << hardA.m() << ")" << std::endl;
208 
209  std::cout << "HardB : user_index = " << hardB.user_index() << ", (Pt,Y,Phi,M) = ("
210  << hardB.pt() << ", " << hardB.rapidity() << ", "
211  << hardB.phi() << ", " << hardB.m() << ")" << std::endl;
212 
213  std::cout << "HardC : user_index = " << hardC.user_index() << ", (Pt,Y,Phi,M) = ("
214  << hardC.pt() << ", " << hardC.rapidity() << ", "
215  << hardC.phi() << ", " << hardC.m() << ")" << std::endl;
216 
217  std::cout << "HardD : user_index = " << hardD.user_index() << ", (Pt,Y,Phi,M) = ("
218  << hardD.pt() << ", " << hardD.rapidity() << ", "
219  << hardD.phi() << ", " << hardD.m() << ")" << std::endl;
220  }
221 
222  // Check to see if any subjects are counted amongst the "hard" subjets from previous
223  // lines. NOTE: In Fastjet 3.0, the default "user_index" changed from 0 to -1, so
224  // this can no longer be used as a designator for the veto of "blankJet" subjets,
225  // and now switch to pt > some small value.
226  if ( subjet1.pt() > 0.0001 )
227  hardSubjets.push_back(subjet1);
228  if ( subjet2.pt() > 0.0001 )
229  hardSubjets.push_back(subjet2);
230  if ( subjet3.pt() > 0.0001 )
231  hardSubjets.push_back(subjet3);
232  if ( subjet4.pt() > 0.0001 )
233  hardSubjets.push_back(subjet4);
234  sort(hardSubjets.begin(), hardSubjets.end(), compEt );
235 
236  // Use new fastjet functionality to create a Pseudojet from constituents
237  fastjet::PseudoJet candidate = join(hardSubjets);
238  // Reset the jet's 4-vector to the "ungroomed" value
239  candidate.reset_momentum( jetIt->px(), jetIt->py(), jetIt->pz(), jetIt->e() );
240 
241  if ( verbose_ ) {
242  std::cout << "Final top-jet candidate: (Pt,Y,Phi,M) = ("
243  << candidate.pt() << ", " << candidate.rapidity() << ", "
244  << candidate.phi() << ", " << candidate.m() << ")" << std::endl;
245  std::vector<fastjet::PseudoJet> pieces = candidate.pieces();
246  std::cout << "Number of pieces = " << pieces.size() << std::endl;
247  for ( std::vector<fastjet::PseudoJet>::const_iterator ibegin = pieces.begin(), iend = pieces.end(), i = ibegin;
248  i != iend; ++i ) {
249  std::cout << " Piece : " << i - ibegin << ", (Pt,Y,Phi,M) = ("
250  << i->pt() << ", " << i->rapidity() << ", "
251  << i->phi() << ", " << i->m() << ")" << std::endl;
252  }
253  }
254  // Add to the list
255  hardjetsOutput.push_back( candidate );
256  }
257 }
258 
259 
260 
261 
262 //-----------------------------------------------------------------------
263 // determine whether two clusters (made of calorimeter towers) are living on "adjacent" cells. if they are, then
264 // we probably shouldn't consider them to be independent objects!
265 //
266 // From Sal: Ignoring genjet case
267 //
268 bool CATopJetAlgorithm::adjacentCells(const fastjet::PseudoJet & jet1, const fastjet::PseudoJet & jet2,
269  const vector<fastjet::PseudoJet> & cell_particles,
270  const fastjet::ClusterSequence & theClusterSequence,
271  double nCellMin ) const {
272 
273 
274  double eta1 = jet1.rapidity();
275  double phi1 = jet1.phi();
276  double eta2 = jet2.rapidity();
277  double phi2 = jet2.phi();
278 
279  double deta = abs(eta2 - eta1) / 0.087;
280  double dphi = fabs( reco::deltaPhi(phi2,phi1) ) / 0.087;
281 
282  return ( ( deta + dphi ) <= nCellMin );
283 }
284 
285 
286 
287 //-------------------------------------------------------------------------
288 // attempt to decompose a jet into "hard" subjets, where hardness is set by ptHard
289 //
290 bool CATopJetAlgorithm::decomposeJet(const fastjet::PseudoJet & theJet,
291  const fastjet::ClusterSequence & theClusterSequence,
292  const vector<fastjet::PseudoJet> & cell_particles,
293  double ptHard, double nCellMin, double deltarcut,
294  fastjet::PseudoJet & ja, fastjet::PseudoJet & jb,
295  vector<fastjet::PseudoJet> & leftovers) const {
296 
297  bool goodBreak;
298  fastjet::PseudoJet j = theJet;
299  double InputObjectPt = j.perp();
300  if ( verbose_ )cout<<"Input Object Pt = "<<InputObjectPt<<endl;
301  if ( verbose_ )cout<<"ptHard = "<<ptHard<<endl;
302  leftovers.clear();
303  if ( verbose_ )cout<<"start while loop"<<endl;
304 
305  while (1) { // watch out for infinite loop!
306  goodBreak = theClusterSequence.has_parents(j,ja,jb);
307  if (!goodBreak){
308  if ( verbose_ )cout<<"bad break. this is one cell. can't decluster anymore."<<endl;
309  break; // this is one cell, can't decluster anymore
310  }
311 
312  if ( verbose_ )cout<<"good break. ja Pt = "<<ja.perp()<<" jb Pt = "<<jb.perp()<<endl;
313 
315 
316  // check if clusters are adjacent using a constant deltar adjacency.
317  double clusters_deltar=fabs(ja.eta()-jb.eta())+fabs(deltaPhi(ja.phi(),jb.phi()));
318 
319  if ( verbose_ && useAdjacency_ ==1)cout<<"clusters_deltar = "<<clusters_deltar<<endl;
320  if ( verbose_ && useAdjacency_ ==1)cout<<"deltar cut = "<<deltarcut<<endl;
321 
322  if ( useAdjacency_==1 && clusters_deltar < deltarcut){
323  if ( verbose_ )cout<<"clusters too close. consant adj. break."<<endl;
324  break;
325  }
326 
327  // Check if clusters are adjacent using a DeltaR adjacency which is a function of pT.
328  double clusters_deltaR=deltaR( ja.rapidity(), ja.phi(), jb.rapidity(), jb.phi() );
329 
330  if ( verbose_ && useAdjacency_ ==2)cout<<"clusters_deltaR = "<<clusters_deltaR<<endl;
331  if ( verbose_ && useAdjacency_ ==2)cout<<"0.4-0.0004*InputObjectPt = "<<0.4-0.0004*InputObjectPt<<endl;
332 
333  if ( useAdjacency_==2 && clusters_deltaR < 0.4-0.0004*InputObjectPt)
334  {
335  if ( verbose_ )cout<<"clusters too close. modified adj. break."<<endl;
336  break;
337  }
338 
339  // Check if clusters are adjacent in the calorimeter.
340  if ( useAdjacency_==3 && adjacentCells(ja,jb,cell_particles,theClusterSequence,nCellMin) ){
341  if ( verbose_ )cout<<"clusters too close in the calorimeter. calorimeter adj. break."<<endl;
342  break; // the clusters are "adjacent" in the calorimeter => shouldn't have decomposed
343  }
344 
345  if ( verbose_ )cout<<"clusters pass distance cut"<<endl;
346 
348 
349  if ( verbose_ )cout<<"ptHard = "<<ptHard<<endl;
350 
351  if (ja.perp() < ptHard && jb.perp() < ptHard){
352  if ( verbose_ )cout<<"two soft clusters. dead end"<<endl;
353  break; // broke into two soft clusters, dead end
354  }
355 
356  if (ja.perp() > ptHard && jb.perp() > ptHard){
357  if ( verbose_ )cout<<"two hard clusters. done"<<endl;
358  return true; // broke into two hard clusters, we're done!
359  }
360 
361  else if (ja.perp() > jb.perp()) { // broke into one hard and one soft, ditch the soft one and try again
362  if ( verbose_ )cout<<"ja hard jb soft. try to split hard. j = ja"<<endl;
363  j = ja;
364  vector<fastjet::PseudoJet> particles = theClusterSequence.constituents(jb);
365  leftovers.insert(leftovers.end(),particles.begin(),particles.end());
366  }
367  else {
368  if ( verbose_ )cout<<"ja hard jb soft. try to split hard. j = jb"<<endl;
369  j = jb;
370  vector<fastjet::PseudoJet> particles = theClusterSequence.constituents(ja);
371  leftovers.insert(leftovers.end(),particles.begin(),particles.end());
372  }
373  }
374 
375  if ( verbose_ )cout<<"did not decluster."<<endl; // did not decluster into hard subjets
376 
377  ja.reset(0,0,0,0);
378  jb.reset(0,0,0,0);
379  leftovers.clear();
380  return false;
381 }
int i
Definition: DBlmapReader.cc:9
tuple pieces
Definition: csv2json.py:31
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
int j
Definition: DBlmapReader.cc:9
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
static std::string join(char **cmd)
Definition: RemoteFile.cc:18
T perp() const
Magnitude of transverse component.
tuple cout
Definition: gather_cfg.py:145
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.