File: SISConeSphericalPlugin.cc

package info (click to toggle)
fastjet 3.4.0%2Bdfsg-1
  • links: PTS, VCS
  • area: main
  • in suites: bookworm, forky, sid, trixie
  • size: 9,540 kB
  • sloc: cpp: 78,628; python: 6,112; sh: 1,038; fortran: 673; makefile: 635; ansic: 161
file content (302 lines) | stat: -rw-r--r-- 11,625 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302

// fastjet stuff
#include "fastjet/ClusterSequence.hh"
#include "fastjet/SISConeSphericalPlugin.hh"

//// siscone stuff
#include "siscone/spherical/momentum.h"
#include "siscone/spherical/siscone.h"

// other stuff
#include<sstream>

FASTJET_BEGIN_NAMESPACE      // defined in fastjet/internal/base.hh

using namespace std;
using namespace siscone_spherical;

/// shortcut for converting siscone CSphmomentum into PseudoJet
template<> PseudoJet::PseudoJet(const siscone_spherical::CSphmomentum & four_vector) {
  (*this) = PseudoJet(four_vector.px,four_vector.py,four_vector.pz,
                      four_vector.E);
}


//======================================================================
// wrap-up around siscone's user-defined scales
namespace siscone_plugin_internal{
  /// @ingroup internal
  /// \class SISConeUserScale
  /// class that makes the transition between the internal SISCone
  /// user-defined scale choice (using SISCone's Cjet) and
  /// user-defined scale choices in the plugn above (using FastJet's
  /// PseudoJets)
  class SISConeSphericalUserScale  : public siscone_spherical::CSphsplit_merge::Cuser_scale_base{
  public:
    /// ctor takes the "fastjet-style" user-defined scale as well as a
    /// reference to the current cluster sequence (to access the
    /// particles if needed)
    SISConeSphericalUserScale(const SISConeSphericalPlugin::UserScaleBase *user_scale,
		     const ClusterSequence &cs)
      : _user_scale(user_scale), _cs(cs){}

    /// returns the scale associated to a given jet
    virtual double operator()(const siscone_spherical::CSphjet &jet) const{
      return _user_scale->result(_build_PJ_from_Cjet(jet));
    }

    /// returns true id the scasle associated to jet a is larger than
    /// the scale associated to jet b
    virtual bool is_larger(const siscone_spherical::CSphjet &a, const siscone_spherical::CSphjet &b) const{
      return _user_scale->is_larger(_build_PJ_from_Cjet(a), _build_PJ_from_Cjet(b));
    }

  private:
    /// constructs a PseudoJet from a siscone::Cjet
    ///
    /// Note that it is tempting to overload the PseudoJet ctor. This
    /// would not work because down the line we need to access the
    /// original PseudoJet through the ClusterSequence and therefore
    /// the PseudoJet structure needs to be aware of the
    /// ClusterSequence.
    PseudoJet _build_PJ_from_Cjet(const siscone_spherical::CSphjet &jet) const{
      PseudoJet j(jet.v.px, jet.v.py, jet.v.pz, jet.v.E);
      j.set_structure_shared_ptr(SharedPtr<PseudoJetStructureBase>(
                                   new SISConeSphericalPlugin::UserScaleBaseStructureType<siscone_spherical::CSphjet>(jet,_cs)));      
      return j;
    }

    const SISConeSphericalPlugin::UserScaleBase *_user_scale;
    const ClusterSequence & _cs;
  };
}

// end of the internal material
//======================================================================


/////////////////////////////////////////////
// static members declaration              //
/////////////////////////////////////////////
SharedPtr<SISConeSphericalPlugin>  SISConeSphericalPlugin::stored_plugin;
SharedPtr<std::vector<PseudoJet> > SISConeSphericalPlugin::stored_particles;
SharedPtr<CSphsiscone>             SISConeSphericalPlugin::stored_siscone;


/////////////////////////////////////////////
// now comes the implementation itself     //
/////////////////////////////////////////////
string SISConeSphericalPlugin::description () const {
  ostringstream desc;
  
  const string on = "on";
  const string off = "off";

  string sm_scale_string = "split-merge uses " + 
    split_merge_scale_name(Esplit_merge_scale(split_merge_scale()));

  desc << "Spherical SISCone jet algorithm with " ;
  desc << "cone_radius = "       << cone_radius        () << ", ";
  if (_progressive_removal)
    desc << "progressive-removal mode, ";
  else 
    desc << "overlap_threshold = " << overlap_threshold  () << ", ";
  desc << "n_pass_max = "        << n_pass_max         () << ", ";
  desc << "protojet_Emin = "     << protojet_Emin()      << ", ";
  if (_progressive_removal && _user_scale) {
    desc << "using a user-defined scale for ordering of stable cones";
    string user_scale_desc = _user_scale->description();
    if (user_scale_desc != "") desc << " (" << user_scale_desc << ")";
  } else {
    desc <<  sm_scale_string;
  }
  if (!_progressive_removal){
    desc << "caching turned "      << (caching() ? on : off);
    desc << ", SM stop scale = "     << _split_merge_stopping_scale;
  }

  // add a note to the description if we use the pt-weighted splitting
  if (_use_E_weighted_splitting){
    desc << ", using E-weighted splitting";
  }

  if (_use_jet_def_recombiner){
    desc << ", using jet-definition's own recombiner";
  }

  // create a fake siscone object so that we can find out more about it
  CSphsiscone siscone;
  if (siscone.merge_identical_protocones) {
    desc << ", and (IR unsafe) merge_indentical_protocones=true" ;
  }

  desc << ", SISCone code v" << siscone_version();

  return desc.str();
}


// overloading the base class implementation
void SISConeSphericalPlugin::run_clustering(ClusterSequence & clust_seq) const {

  CSphsiscone::set_banner_stream(clust_seq.fastjet_banner_stream());

  CSphsiscone   local_siscone;
  CSphsiscone * siscone;

  unsigned n = clust_seq.jets().size();

  bool new_siscone = true; // by default we'll be running it

  if (caching() && !_progressive_removal) {

    // Establish if we have a cached run with the same R, npass and
    // particles. If not then do any tidying up / reallocation that's
    // necessary for the next round of caching, otherwise just set
    // relevant pointers so that we can reuse and old run.
    if (stored_siscone.get() != 0) {
      new_siscone = !(stored_plugin->cone_radius()   == cone_radius()
                      && stored_plugin->n_pass_max() == n_pass_max()
                      && stored_particles->size()    == n);
      if (!new_siscone) {
        for(unsigned i = 0; i < n; i++) {
          // only check momentum because indices will be correctly dealt
          // with anyway when extracting the clustering order.
          new_siscone |= !have_same_momentum(clust_seq.jets()[i],
                                             (*stored_particles)[i]);
        }
      }
    }

    // allocate the new siscone, etc., if need be
    if (new_siscone) {
      stored_siscone  .reset( new CSphsiscone );
      stored_particles.reset( new std::vector<PseudoJet>(clust_seq.jets()));
      reset_stored_plugin();
    }

    siscone = stored_siscone.get();
  } else {
    siscone = &local_siscone;
  }

  // make sure stopping scale is set in siscone
  siscone->SM_var2_hardest_cut_off = _split_merge_stopping_scale*_split_merge_stopping_scale;

  // set the specific parameters
  // when running with ghosts for passive areas, do not put the
  // ghosts into the stable-cone search (not relevant)
  siscone->stable_cone_soft_E2_cutoff = ghost_separation_scale()
                                      * ghost_separation_scale();
  // set the type of splitting we want (default=std one, true->pt-weighted split)
  siscone->set_E_weighted_splitting(_use_E_weighted_splitting);


  if (new_siscone) {
    // transfer fastjet initial particles into the siscone type
    std::vector<CSphmomentum> siscone_momenta(n);
    for(unsigned i = 0; i < n; i++) {
      const PseudoJet & p = clust_seq.jets()[i]; // shorthand
      siscone_momenta[i] = CSphmomentum(p.px(), p.py(), p.pz(), p.E());
    }

    // run the jet finding
    //cout << "plg sms: " << split_merge_scale() << endl;
    if (_progressive_removal){
      // handle the optional user-defined scale choice
      SharedPtr<siscone_plugin_internal::SISConeSphericalUserScale> internal_scale;
      if (_user_scale){
	internal_scale.reset(new siscone_plugin_internal::SISConeSphericalUserScale(_user_scale, clust_seq));
	siscone->set_user_scale(internal_scale.get());
      }
      siscone->compute_jets_progressive_removal(siscone_momenta, cone_radius(),
						n_pass_max(), protojet_or_ghost_Emin(), 
						Esplit_merge_scale(split_merge_scale()));
    } else {
      siscone->compute_jets(siscone_momenta, cone_radius(), overlap_threshold(),
			    n_pass_max(), protojet_or_ghost_Emin(), 
			    Esplit_merge_scale(split_merge_scale()));
    }
  } else {
    // rerun the jet finding
    // just run the overlap part of the jets.
    //cout << "plg rcmp sms: " << split_merge_scale() << endl;
    siscone->recompute_jets(overlap_threshold(), protojet_or_ghost_Emin(), 
			    Esplit_merge_scale(split_merge_scale()));
  }

  // extract the jets [in reverse order -- to get nice ordering in pt at end]
  int njet = siscone->jets.size();

  // allocate space for the extras object
  SISConeSphericalExtras * extras = new SISConeSphericalExtras(n);

  // the ordering in which the inclusive jets are transfered here is
  // deliberate and ensures that when a user asks for
  // inclusive_jets(), they are provided in the order in which SISCone
  // created them.
  for (int ijet = njet-1; ijet >= 0; ijet--) {
    const CSphjet & jet = siscone->jets[ijet]; // shorthand

    // Successively merge the particles that make up the cone jet
    // until we have all particles in it.  Start off with the zeroth
    // particle.
    int jet_k = jet.contents[0];
    for (unsigned ipart = 1; ipart < jet.contents.size(); ipart++) {
      // take the last result of the merge
      int jet_i = jet_k;
      // and the next element of the jet
      int jet_j = jet.contents[ipart];
      // and merge them (with a fake dij)
      double dij = 0.0;

      if (_use_jet_def_recombiner) {
       clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, jet_k);
      } else {
       // create the new jet by hand so that we can adjust its user index
       PseudoJet newjet = clust_seq.jets()[jet_i] + clust_seq.jets()[jet_j];
       // set the user index to be the pass in which the jet was discovered
       // *** The following line was commented for 3.0.1 ***
       //newjet.set_user_index(jet.pass);
       clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, newjet, jet_k);
      }

    }

    // we have merged all the jet's particles into a single object, so now
    // "declare" it to be a beam (inclusive) jet.
    // [NB: put a sensible looking d_iB just to be nice...]
    double d_iB = clust_seq.jets()[jet_k].perp2();
    clust_seq.plugin_record_iB_recombination(jet_k, d_iB);

    // now record the pass of the jet in the extras object
    extras->_pass[clust_seq.jets()[jet_k].cluster_hist_index()] = jet.pass;
  }

  // now copy the list of protocones into an "extras" objects
  for (unsigned ipass = 0; ipass < siscone->protocones_list.size(); ipass++) {
    for (unsigned ipc = 0; ipc < siscone->protocones_list[ipass].size(); ipc++) {
      PseudoJet protocone(siscone->protocones_list[ipass][ipc]);
      protocone.set_user_index(ipass);
      extras->_protocones.push_back(protocone);
    }
  }
  extras->_most_ambiguous_split = siscone->most_ambiguous_split;

  // tell it what the jet definition was
  extras->_jet_def_plugin = this;

  // give the extras object to the cluster sequence.
  // 
  // As of v3.1 of FastJet, extras are automatically owned (as
  // SharedPtr) by the ClusterSequence and auto_ptr is deprecated. So
  // we can use a simple pointer here
  clust_seq.plugin_associate_extras(extras);
}

void SISConeSphericalPlugin::reset_stored_plugin() const{
  stored_plugin.reset( new SISConeSphericalPlugin(*this));
}


FASTJET_END_NAMESPACE      // defined in fastjet/internal/base.hh