1 | //STARTHEADER |
---|
2 | // $Id: EECambridgePlugin.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/EECambridgePlugin.hh" |
---|
32 | #include "fastjet/NNH.hh" |
---|
33 | |
---|
34 | // other stuff |
---|
35 | #include <sstream> |
---|
36 | #include <limits> |
---|
37 | |
---|
38 | FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh |
---|
39 | |
---|
40 | using namespace std; |
---|
41 | |
---|
42 | //---------------------------------------------------------------------- |
---|
43 | /// class to help run an e+e- Cambridge algorithm |
---|
44 | class EECamBriefJet { |
---|
45 | public: |
---|
46 | void init(const PseudoJet & jet) { |
---|
47 | double norm = 1.0/sqrt(jet.modp2()); |
---|
48 | nx = jet.px() * norm; |
---|
49 | ny = jet.py() * norm; |
---|
50 | nz = jet.pz() * norm; |
---|
51 | } |
---|
52 | |
---|
53 | double distance(const EECamBriefJet * jet) const { |
---|
54 | double dij = 1 - nx*jet->nx |
---|
55 | - ny*jet->ny |
---|
56 | - nz*jet->nz; |
---|
57 | return dij; |
---|
58 | } |
---|
59 | |
---|
60 | double beam_distance() const { |
---|
61 | return numeric_limits<double>::max(); |
---|
62 | } |
---|
63 | |
---|
64 | private: |
---|
65 | double nx, ny, nz; |
---|
66 | }; |
---|
67 | |
---|
68 | |
---|
69 | string EECambridgePlugin::description () const { |
---|
70 | ostringstream desc; |
---|
71 | desc << "EECambridge plugin with ycut = " << ycut() ; |
---|
72 | return desc.str(); |
---|
73 | } |
---|
74 | |
---|
75 | void EECambridgePlugin::run_clustering(ClusterSequence & cs) const { |
---|
76 | int njets = cs.jets().size(); |
---|
77 | NNH<EECamBriefJet> nnh(cs.jets()); |
---|
78 | |
---|
79 | double Q2 = cs.Q2(); |
---|
80 | |
---|
81 | while (njets > 0) { |
---|
82 | int i, j, k; |
---|
83 | // here we get a minimum based on the purely angular variable from the NNH class |
---|
84 | // (called dij there, but vij in the Cambridge article (which uses dij for |
---|
85 | // a kt distance...) |
---|
86 | double vij = nnh.dij_min(i, j); // i,j are return values... |
---|
87 | |
---|
88 | // next we work out the dij (ee kt distance), and based on its |
---|
89 | // value decide whether we have soft-freezing (represented here by |
---|
90 | // a "Beam" clustering) or not |
---|
91 | double dij; |
---|
92 | if (j >= 0) { |
---|
93 | double scale = min(cs.jets()[i].E(), cs.jets()[j].E()); |
---|
94 | dij = 2 * vij * scale * scale; |
---|
95 | if (dij > Q2 * ycut()) { |
---|
96 | // we'll call the softer partner a "beam" jet |
---|
97 | if (cs.jets()[i].E() > cs.jets()[j].E()) std::swap(i,j); |
---|
98 | j = -1; |
---|
99 | } |
---|
100 | } else { |
---|
101 | // for the last particle left, just use yij = 1 |
---|
102 | dij = Q2; |
---|
103 | } |
---|
104 | |
---|
105 | if (j >= 0) { |
---|
106 | cs.plugin_record_ij_recombination(i, j, dij, k); |
---|
107 | nnh.merge_jets(i, j, cs.jets()[k], k); |
---|
108 | } else { |
---|
109 | cs.plugin_record_iB_recombination(i, dij); |
---|
110 | nnh.remove_jet(i); |
---|
111 | } |
---|
112 | njets--; |
---|
113 | } |
---|
114 | |
---|
115 | } |
---|
116 | |
---|
117 | FASTJET_END_NAMESPACE // defined in fastjet/internal/base.hh |
---|