[45e58be] | 1 | // $Id: ValenciaPlugin.cc 1209 2018-12-05 16:18:01Z vos $
|
---|
[f319c1d] | 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 |
|
---|
| 29 | FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
|
---|
| 30 |
|
---|
| 31 | namespace contrib{
|
---|
| 32 |
|
---|
| 33 | //----------------------------------------------------------------------
|
---|
| 34 | /// class that contains the algorithm parameters R and beta
|
---|
| 35 | // needed to run the Valencia algorithm
|
---|
| 36 | class ValenciaInfo {
|
---|
| 37 |
|
---|
| 38 | public:
|
---|
| 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 |
|
---|
| 46 | private:
|
---|
| 47 |
|
---|
| 48 | double R_, beta_, gamma_;
|
---|
| 49 | };
|
---|
| 50 |
|
---|
| 51 | class ValenciaBriefJet {
|
---|
| 52 | public:
|
---|
| 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.;
|
---|
[45e58be] | 68 | // modified diB in release 2.0.1
|
---|
| 69 | diB = pow(E,2*beta) * pow(jet.perp()/sqrt(jet.perp2()+jet.pz()*jet.pz()),2*info->gamma());
|
---|
[f319c1d] | 70 | }
|
---|
| 71 |
|
---|
| 72 | double distance(const ValenciaBriefJet * jet) const {
|
---|
| 73 | // inter-particle distance is similar to Durham and e+e- kt
|
---|
| 74 | // two differences:
|
---|
| 75 | // R is redefined for greater flexibility
|
---|
| 76 | // the energy is elevated to 2*beta
|
---|
| 77 | double dij = 1 - nx*jet->nx
|
---|
| 78 | - ny*jet->ny
|
---|
| 79 | - nz*jet->nz;
|
---|
| 80 |
|
---|
[45e58be] | 81 | if (pow(jet->E,2*beta) < pow(E,2*beta))
|
---|
[f319c1d] | 82 | dij *= 2 * pow(jet->E,2*beta);
|
---|
| 83 | else
|
---|
| 84 | dij *= 2 * pow(E,2*beta);
|
---|
| 85 |
|
---|
| 86 | dij/=pow(R,2);
|
---|
| 87 |
|
---|
| 88 | return dij;
|
---|
| 89 | }
|
---|
| 90 |
|
---|
| 91 | double beam_distance() const {
|
---|
| 92 | return diB;
|
---|
| 93 | }
|
---|
| 94 |
|
---|
| 95 |
|
---|
| 96 | double E, nx, ny, nz;
|
---|
| 97 | double diB;
|
---|
| 98 | double R, beta;
|
---|
| 99 | };
|
---|
| 100 |
|
---|
| 101 |
|
---|
| 102 | std::string ValenciaPlugin::description () const {
|
---|
| 103 | std::ostringstream desc;
|
---|
| 104 | desc << "Valencia plugin with R = " << R() << ", beta = " << beta() << " and gamma = " << gamma();
|
---|
| 105 | return desc.str();
|
---|
| 106 | }
|
---|
| 107 |
|
---|
| 108 | void ValenciaPlugin::run_clustering(fastjet::ClusterSequence & cs) const {
|
---|
| 109 | int njets = cs.jets().size();
|
---|
| 110 | ValenciaInfo vinfo(R(),beta(),gamma());
|
---|
| 111 |
|
---|
| 112 | NNH<ValenciaBriefJet,ValenciaInfo> nnh(cs.jets(),&vinfo);
|
---|
| 113 |
|
---|
| 114 | while (njets > 0) {
|
---|
| 115 | int i, j, k;
|
---|
| 116 | double dij = nnh.dij_min(i, j);
|
---|
| 117 |
|
---|
| 118 | if (j >= 0) {
|
---|
| 119 | cs.plugin_record_ij_recombination(i, j, dij, k);
|
---|
| 120 | nnh.merge_jets(i, j, cs.jets()[k], k);
|
---|
| 121 | } else {
|
---|
| 122 |
|
---|
| 123 | cs.plugin_record_iB_recombination(i, dij);
|
---|
| 124 | nnh.remove_jet(i);
|
---|
| 125 | }
|
---|
| 126 |
|
---|
| 127 | njets--;
|
---|
| 128 | }
|
---|
| 129 | }
|
---|
| 130 |
|
---|
| 131 |
|
---|
| 132 | } // namespace contrib
|
---|
| 133 |
|
---|
| 134 | FASTJET_END_NAMESPACE
|
---|