Fork me on GitHub

source: git/external/fastjet/contribs/ValenciaPlugin/ValenciaPlugin.cc@ dda357d

ImprovedOutputFile Timing dual_readout llp
Last change on this file since dda357d was f319c1d, checked in by Ulrike Schnoor <schnooru@…>, 7 years ago

added Valencia jet algorithm as fastjet/contribs and exclusive jet clustering

  • Property mode set to 100644
File size: 3.5 KB
RevLine 
[f319c1d]1// $Id: ValenciaPlugin.cc 776 2015-02-24 17:53:27Z vos $
2//
3// Copyright (c) 2014, Marcel Vos and Ignacio Garcia
4//
5//----------------------------------------------------------------------
6// This file is part of FastJet contrib.
7//
8// It is free software; you can redistribute it and/or modify it under
9// the terms of the GNU General Public License as published by the
10// Free Software Foundation; either version 2 of the License, or (at
11// your option) any later version.
12//
13// It is distributed in the hope that it will be useful, but WITHOUT
14// ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
15// or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public
16// License for more details.
17//
18// You should have received a copy of the GNU General Public License
19// along with this code. If not, see <http://www.gnu.org/licenses/>.
20//----------------------------------------------------------------------
21
22#include "ValenciaPlugin.hh"
23#include "fastjet/NNH.hh"
24
25// strings and streams
26#include <sstream>
27#include <limits>
28
29FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
30
31namespace contrib{
32
33//----------------------------------------------------------------------
34/// class that contains the algorithm parameters R and beta
35// needed to run the Valencia algorithm
36class ValenciaInfo {
37
38public:
39 ValenciaInfo(double Ri, double betai, double gammai)
40 { R_ = Ri; beta_ = betai; gamma_ = gammai;}
41
42 double beta() { return beta_; }
43 double gamma() { return gamma_; }
44 double R() { return R_; }
45
46private:
47
48 double R_, beta_, gamma_;
49};
50
51class ValenciaBriefJet {
52public:
53 void init(const PseudoJet & jet, ValenciaInfo * info) {
54 double norm = 1.0/sqrt(jet.modp2());
55 // pseudo-jet information needed to calculate distance
56 nx = jet.px() * norm;
57 ny = jet.py() * norm;
58 nz = jet.pz() * norm;
59 E = jet.E();
60
61 // algorithm parameters needed in distance calculation
62 R = info->R();
63 beta = info->beta();
64
65 // beam distance
66 // E-to-the-2*beta times sin(polar angle)-to-the-2*gamma
67 if (E==0. || jet.perp()==0.) diB=0.;
68 diB = pow(E,2*beta) * pow(jet.perp()/E,2*info->gamma());
69 }
70
71 double distance(const ValenciaBriefJet * jet) const {
72 // inter-particle distance is similar to Durham and e+e- kt
73 // two differences:
74 // R is redefined for greater flexibility
75 // the energy is elevated to 2*beta
76 double dij = 1 - nx*jet->nx
77 - ny*jet->ny
78 - nz*jet->nz;
79
80 if (jet->E < E)
81 dij *= 2 * pow(jet->E,2*beta);
82 else
83 dij *= 2 * pow(E,2*beta);
84
85 dij/=pow(R,2);
86
87 return dij;
88 }
89
90 double beam_distance() const {
91 return diB;
92 }
93
94
95 double E, nx, ny, nz;
96 double diB;
97 double R, beta;
98};
99
100
101std::string ValenciaPlugin::description () const {
102 std::ostringstream desc;
103 desc << "Valencia plugin with R = " << R() << ", beta = " << beta() << " and gamma = " << gamma();
104 return desc.str();
105}
106
107void ValenciaPlugin::run_clustering(fastjet::ClusterSequence & cs) const {
108 int njets = cs.jets().size();
109 ValenciaInfo vinfo(R(),beta(),gamma());
110
111 NNH<ValenciaBriefJet,ValenciaInfo> nnh(cs.jets(),&vinfo);
112
113 while (njets > 0) {
114 int i, j, k;
115 double dij = nnh.dij_min(i, j);
116
117 if (j >= 0) {
118 cs.plugin_record_ij_recombination(i, j, dij, k);
119 nnh.merge_jets(i, j, cs.jets()[k], k);
120 } else {
121
122 cs.plugin_record_iB_recombination(i, dij);
123 nnh.remove_jet(i);
124 }
125
126 njets--;
127 }
128}
129
130
131} // namespace contrib
132
133FASTJET_END_NAMESPACE
Note: See TracBrowser for help on using the repository browser.