File indexing completed on 2025-01-18 09:57:14
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023 #ifndef __FASTJET_CONTRIB_RECURSIVELUNDEEGENERATOR_HH__
0024 #define __FASTJET_CONTRIB_RECURSIVELUNDEEGENERATOR_HH__
0025
0026 #include "LundEEHelpers.hh"
0027
0028 #include <fastjet/internal/base.hh>
0029 #include "fastjet/tools/Recluster.hh"
0030 #include "fastjet/JetDefinition.hh"
0031 #include "fastjet/PseudoJet.hh"
0032 #include <string>
0033 #include <vector>
0034 #include <utility>
0035 #include <queue>
0036
0037 using namespace std;
0038
0039 FASTJET_BEGIN_NAMESPACE
0040
0041 namespace contrib{
0042
0043
0044
0045
0046
0047 class LundEEDeclustering {
0048 public:
0049
0050
0051 const PseudoJet & pair() const {return pair_;}
0052
0053 const PseudoJet & harder() const {return harder_;}
0054
0055 const PseudoJet & softer() const {return softer_;}
0056
0057
0058
0059 double m() const {return m_;}
0060
0061
0062 double eta() const {return eta_;}
0063
0064
0065 double sin_theta() const {return sin_theta_;}
0066
0067
0068 double z() const {return z_;}
0069
0070
0071 double kt() const {return kt_;}
0072
0073
0074 double lnkt() const {return lnkt_;}
0075
0076
0077 double kappa() const {return kappa_;}
0078
0079
0080 int iplane() const {return iplane_;}
0081
0082
0083
0084 int depth() const {return depth_;}
0085
0086
0087
0088
0089 int leaf_iplane() const {return leaf_iplane_;}
0090
0091
0092 int sign_s() const {return sign_s_;}
0093
0094
0095
0096
0097
0098 double psi() const {return psi_;}
0099
0100
0101 void set_psi(double psi) {psi_ = psi;}
0102
0103
0104 double psibar() const {return psibar_;}
0105
0106
0107
0108 std::pair<double,double> const lund_coordinates() const {
0109 return std::pair<double,double>(eta_,lnkt_);
0110 }
0111
0112 virtual ~LundEEDeclustering() {}
0113
0114 private:
0115 int iplane_;
0116 double psi_, psibar_, lnkt_, eta_;
0117 double m_, z_, kt_, kappa_, sin_theta_;
0118 PseudoJet pair_, harder_, softer_;
0119 int depth_ = -1, leaf_iplane_ = -1;
0120 int sign_s_;
0121
0122 protected:
0123
0124
0125
0126 LundEEDeclustering() {}
0127
0128
0129
0130 LundEEDeclustering(const PseudoJet& pair,
0131 const PseudoJet& j1, const PseudoJet& j2,
0132 int iplane = -1, double psi = 0.0, double psibar = 0.0, int depth = -1, int leaf_iplane = -1, int sign_s = 1);
0133
0134 friend class RecursiveLundEEGenerator;
0135
0136 };
0137
0138
0139
0140
0141 inline bool operator<(const LundEEDeclustering& d1, const LundEEDeclustering& d2) {
0142 return d1.kt() < d2.kt();
0143 }
0144
0145
0146
0147
0148
0149 class RecursiveLundEEGenerator {
0150 public:
0151
0152
0153
0154
0155
0156 RecursiveLundEEGenerator(int max_depth = 0, bool dynamical_psi_ref = false) :
0157 max_depth_(max_depth), nx_(1,0,0,0), ny_(0,1,0,0), dynamical_psi_ref_(dynamical_psi_ref)
0158 {}
0159
0160
0161 virtual ~RecursiveLundEEGenerator() {}
0162
0163
0164
0165
0166
0167
0168 virtual std::vector<LundEEDeclustering> result(const ClusterSequence & cs) const {
0169 std::vector<PseudoJet> exclusive_jets = cs.exclusive_jets(2);
0170 assert(exclusive_jets.size() == 2);
0171
0172
0173 if (exclusive_jets[0].pz() < exclusive_jets[1].pz()) {
0174 std::swap(exclusive_jets[0],exclusive_jets[1]);
0175 }
0176
0177 PseudoJet d_ev = exclusive_jets[0] - exclusive_jets[1];
0178 lund_plane::Matrix3 rotmat = lund_plane::Matrix3::from_direction(d_ev);
0179
0180 std::vector<LundEEDeclustering> declusterings;
0181 int depth = 0;
0182 int max_iplane_sofar = 1;
0183 for (unsigned ijet = 0; ijet < exclusive_jets.size(); ijet++) {
0184
0185
0186 PseudoJet axis = d_ev/sqrt(d_ev.modp2());
0187 PseudoJet ref_plane = axis;
0188
0189 int sign_s = ijet == 0? +1 : -1;
0190
0191 append_to_vector(declusterings, exclusive_jets[ijet], depth, ijet, max_iplane_sofar,
0192 rotmat, sign_s, exclusive_jets[0], exclusive_jets[1], ref_plane, 0., true);
0193 }
0194
0195
0196 typedef LundEEDeclustering LD;
0197
0198
0199
0200 sort(declusterings.begin(), declusterings.end(),
0201 [](const LD & d1, LD & d2){return d1.kt() > d2.kt();});
0202
0203 return declusterings;
0204 }
0205
0206 private:
0207
0208
0209
0210
0211
0212 void append_to_vector(std::vector<LundEEDeclustering> & declusterings,
0213 const PseudoJet & jet, int depth,
0214 int iplane, int & max_iplane_sofar,
0215 const lund_plane::Matrix3 & rotmat, int sign_s,
0216 const PseudoJet & harder,
0217 const PseudoJet & softer,
0218 const PseudoJet & psibar_ref_plane,
0219 const double & last_psibar, bool first_time) const {
0220 PseudoJet j1, j2;
0221 if (!jet.has_parents(j1, j2)) return;
0222 if (j1.modp2() < j2.modp2()) std::swap(j1,j2);
0223
0224
0225 lund_plane::Matrix3 new_rotmat;
0226 if (dynamical_psi_ref_) {
0227 new_rotmat = lund_plane::Matrix3::from_direction(rotmat.transpose()*(sign_s*jet)) * rotmat;
0228 } else {
0229 new_rotmat = rotmat;
0230 }
0231 PseudoJet rx = new_rotmat * nx_;
0232 PseudoJet ry = new_rotmat * ny_;
0233 PseudoJet u1 = j1/j1.modp(), u2 = j2/j2.modp();
0234 PseudoJet du = u2 - u1;
0235 double x = du.px() * rx.px() + du.py() * rx.py() + du.pz() * rx.pz();
0236 double y = du.px() * ry.px() + du.py() * ry.py() + du.pz() * ry.pz();
0237 double psi = atan2(y,x);
0238
0239
0240 double psibar = 0.;
0241 PseudoJet n1, n2;
0242
0243
0244 if (first_time) {
0245
0246
0247 n1 = lund_plane::cross_product(psibar_ref_plane,j1);
0248 n2 = lund_plane::cross_product(j1,j2);
0249
0250 double signed_angle = 0.;
0251 n2 /= n2.modp();
0252 if (n1.modp() != 0) {
0253 n1 /= n1.modp();
0254 signed_angle = lund_plane::signed_angle_between_planes(n1,n2,j1);
0255 }
0256
0257 psibar = lund_plane::map_to_pi(j1.phi() + signed_angle);
0258 }
0259
0260 else {
0261 n2 = lund_plane::cross_product(j1,j2);
0262 n2 /= n2.modp();
0263 psibar = lund_plane::map_to_pi(last_psibar + lund_plane::signed_angle_between_planes(psibar_ref_plane, n2, j1));
0264 }
0265
0266 int leaf_iplane = -1;
0267
0268
0269 bool recurse_into_softer = (depth < max_depth_ || max_depth_ < 0);
0270 if (recurse_into_softer) {
0271 max_iplane_sofar += 1;
0272 leaf_iplane = max_iplane_sofar;
0273 }
0274
0275 LundEEDeclustering declust(jet, j1, j2, iplane, psi, psibar, depth, leaf_iplane, sign_s);
0276 declusterings.push_back(declust);
0277
0278
0279
0280
0281 append_to_vector(declusterings, j1, depth, iplane, max_iplane_sofar, new_rotmat, sign_s, u1, u2, n2, psibar, false);
0282 if (recurse_into_softer) {
0283 append_to_vector(declusterings, j2, depth+1, leaf_iplane, max_iplane_sofar, new_rotmat, sign_s, u1, u2, n2, psibar, false);
0284 }
0285 }
0286
0287 int max_depth_ = 0;
0288
0289 PseudoJet nx_;
0290 PseudoJet ny_;
0291 bool dynamical_psi_ref_;
0292 };
0293
0294 }
0295
0296 FASTJET_END_NAMESPACE
0297
0298
0299 std::ostream & operator<<(std::ostream & ostr, const fastjet::contrib::LundEEDeclustering & d);
0300
0301 #endif