source: HiSusy/trunk/Delphes/Delphes-3.0.9/external/fastjet/plugins/ATLASCone/ATLASConePlugin.cc @ 5

Last change on this file since 5 was 5, checked in by zerwas, 11 years ago

update to Delphes-3.0.9

File size: 6.8 KB
Line 
1//STARTHEADER
2// $Id: ATLASConePlugin.cc 859 2012-11-28 01:49:23Z pavel $
3//
4// Copyright (c) 2007-2011, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
5//
6//----------------------------------------------------------------------
7// This file is part of FastJet.
8//
9//  FastJet is free software; you can redistribute it and/or modify
10//  it under the terms of the GNU General Public License as published by
11//  the Free Software Foundation; either version 2 of the License, or
12//  (at your option) any later version.
13//
14//  The algorithms that underlie FastJet have required considerable
15//  development and are described in hep-ph/0512210. If you use
16//  FastJet as part of work towards a scientific publication, please
17//  include a citation to the FastJet paper.
18//
19//  FastJet is distributed in the hope that it will be useful,
20//  but WITHOUT ANY WARRANTY; without even the implied warranty of
21//  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
22//  GNU General Public License for more details.
23//
24//  You should have received a copy of the GNU General Public License
25//  along with FastJet. If not, see <http://www.gnu.org/licenses/>.
26//----------------------------------------------------------------------
27//ENDHEADER
28
29// fastjet stuff
30#include "fastjet/ClusterSequence.hh"
31#include "fastjet/ATLASConePlugin.hh"
32
33// SpartyJet stuff
34#include "CommonUtils.hh"
35#include "JetConeFinderTool.hh"
36#include "JetSplitMergeTool.hh"
37
38// other stuff
39#include <vector>
40#include <sstream>
41
42FASTJET_BEGIN_NAMESPACE      // defined in fastjet/internal/base.hh
43
44using namespace std;
45
46bool ATLASConePlugin::_first_time = true;
47
48string ATLASConePlugin::description () const {
49  ostringstream desc;
50  desc << "ATLASCone plugin with R = "<< _radius
51       << ", seed threshold = " << _seedPt
52       << ", overlap threshold f = " << _f;
53  return desc.str();
54}
55
56void ATLASConePlugin::run_clustering(ClusterSequence & clust_seq) const {
57  // print a banner if we run this for the first time
58  _print_banner(clust_seq.fastjet_banner_stream());
59
60
61  // transfer the list of PseudoJet into a atlas::Jet::jet_list_t
62  //  jet_list_t is a vector<Jet*>
63  // We set the index of the 4-vect to trace the constituents at the end
64  //------------------------------------------------------------------
65  // cout << "ATLASConePlugin: transferring vectors from ClusterSequence" << endl;
66  atlas::JetConeFinderTool::jetcollection_t jets_ptr;
67  vector<atlas::Jet*> particles_ptr;
68
69  for (unsigned int i=0 ; i<clust_seq.jets().size() ; i++) {
70    const PseudoJet & mom = clust_seq.jets()[i];
71   
72    // first create the particle
73    atlas::Jet *particle = new atlas::Jet(mom.px(), mom.py(), mom.pz(), mom.E(), i);
74    particles_ptr.push_back(particle);
75
76    // then add it to the list of particles we'll use for the clustering
77    atlas::Jet *jet = new atlas::Jet;
78    jet->set_index(particle->index());
79    jet->addConstituent(particle);
80
81    // and finally add that one to the list of jets
82    jets_ptr.push_back(jet);
83  }
84  // cout << "ATLASCone: " << jets_ptr.size() << " particles to cluster" << endl;
85
86  // search the stable cones
87  //------------------------------------------------------------------
88  // cout << "ATLASConePlugin: searching for stable cones" << endl;
89  atlas::JetConeFinderTool stable_cone_finder;
90
91  // set the parameters
92  stable_cone_finder.m_coneR  = _radius;
93  stable_cone_finder.m_seedPt = _seedPt;
94
95  // really do the search.
96  // Note that the list of protocones is returned
97  // through the argument
98  stable_cone_finder.execute(jets_ptr);
99  // cout << "ATLASCone: " << jets_ptr.size() << " stable cones found" << endl;
100
101  // perform the split-merge
102  //------------------------------------------------------------------
103  // cout << "ATLASConePlugin: running the split-merge" << endl;
104  atlas::JetSplitMergeTool split_merge;
105 
106  // set the parameters
107  split_merge.m_f = _f;
108
109  // do the work
110  // again, the list of jets is returned through the argument
111  split_merge.execute(&jets_ptr);
112  // cout << "ATLASCone: " << jets_ptr.size() << " jets after split--merge" << endl;
113
114  // build the FastJet jets (a la SISConePlugin)
115  //------------------------------------------------------------------
116  // cout << "ATLASConePlugin: backporting jets to the ClusterSequence" << endl;
117  for (atlas::Jet::jet_list_t::iterator jet_it = jets_ptr.begin() ;
118         jet_it != jets_ptr.end(); jet_it++){
119    // iterate over the constituents, starting from the first one
120    // that we just take as a reference
121    atlas::Jet::constit_vect_t::iterator constit_it = (*jet_it)->firstConstituent();
122    // cout << " atlas: jet has " << (*jet_it)->getConstituentNum() << " constituents" << endl;
123    int jet_k = (*constit_it)->index();
124    constit_it++;
125   
126    // loop over the remaining particles
127    while (constit_it != (*jet_it)->lastConstituent()){
128      // take the last result of the merge
129      int jet_i = jet_k;
130      // and the next element of the jet
131      int jet_j = (*constit_it)->index();
132      // and merge them (with a fake dij)
133      double dij = 0.0;
134
135      // create the new jet by hand so that we can adjust its user index
136      // Note again the use of the E-scheme recombination here!
137      PseudoJet newjet = clust_seq.jets()[jet_i] + clust_seq.jets()[jet_j];
138      clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, newjet, jet_k);
139
140      // jump to the next constituent
141      constit_it++;
142    }
143
144    // we have merged all the jet's particles into a single object, so now
145    // "declare" it to be a beam (inclusive) jet.
146    // [NB: put a sensible looking d_iB just to be nice...]
147    double d_iB = clust_seq.jets()[jet_k].perp2();
148    clust_seq.plugin_record_iB_recombination(jet_k, d_iB);
149  }
150   
151  // cout << "ATLASConePlugin: Bye" << endl;
152  clear_list(particles_ptr);
153  clear_list(jets_ptr);
154}
155
156// print a banner for reference to the 3rd-party code
157void ATLASConePlugin::_print_banner(ostream *ostr) const{
158  if (! _first_time) return;
159  _first_time=false;
160
161  // make sure the user has not set the banner stream to NULL
162  if (!ostr) return; 
163
164  (*ostr) << "#-------------------------------------------------------------------------" << endl;
165  (*ostr) << "# You are running the ATLAS Cone plugin for FastJet                       " << endl;
166  (*ostr) << "# Original code from SpartyJet; interface by the FastJet authors          " << endl;
167  (*ostr) << "# If you use this plugin, please cite                                     " << endl;
168  (*ostr) << "#   P.A. Delsart, K. Geerlings, J. Huston, B. Martin and C. Vermilion,    " << endl;
169  (*ostr) << "#   SpartyJet, http://projects.hepforge.org/spartyjet                     " << endl;
170  (*ostr) << "# in addition to the usual FastJet reference.                             " << endl;
171  (*ostr) << "#-------------------------------------------------------------------------" << endl;
172
173  // make sure we really have the output done.
174  ostr->flush();
175}
176
177
178FASTJET_END_NAMESPACE      // defined in fastjet/internal/base.hh
Note: See TracBrowser for help on using the repository browser.