File indexing completed on 2026-06-21 08:23:36
0001 #ifndef __FJCONTRIB_IFNPLUGIN_HH__
0002 #define __FJCONTRIB_IFNPLUGIN_HH__
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013 #include "fastjet/contrib/FlavNeutraliser.hh"
0014
0015 FASTJET_BEGIN_NAMESPACE
0016 namespace contrib{
0017
0018
0019
0020 class IFNPlugin : public JetDefinition::Plugin {
0021 public:
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032 IFNPlugin(
0033 const JetDefinition & jet_def,
0034 double alpha,
0035 double omega = -1,
0036 FlavRecombiner::FlavSummation flav_summation = FlavRecombiner::net,
0037 bool use_mass_flav = false
0038 ) : _jet_def(jet_def),
0039 _p(0.5*alpha),
0040 _q(0.5*(2-alpha)),
0041 _a(omega > 0 ? omega : 3 - alpha),
0042 _modulo_2(flav_summation == FlavRecombiner::modulo_2),
0043 _measure_in(FlavNeutraliser::general),
0044 _use_mass_flav(use_mass_flav),
0045 _spherical_algo(_jet_def.is_spherical()),
0046 _pp(1),
0047 _recursive(true) {
0048 if (flav_summation == FlavRecombiner::any_abs) {
0049 throw Error("IFNPlugin: FlavRecombiner::any_abs is not supported");
0050 }
0051 check_mod2_consistency();
0052 }
0053
0054
0055
0056
0057
0058
0059
0060
0061
0062 IFNPlugin(
0063 JetDefinition jet_def,
0064 bool modulo_2,
0065 FlavNeutraliser::measure measure_in,
0066 bool use_mass_flav = false,
0067 double pp = 1,
0068 bool recursive_neutralisation = true)
0069 : _jet_def(jet_def),
0070 _modulo_2(modulo_2),
0071 _measure_in(measure_in),
0072 _use_mass_flav(use_mass_flav),
0073 _spherical_algo(_jet_def.is_spherical()),
0074 _pp(pp),
0075 _recursive(recursive_neutralisation) {check_mod2_consistency();}
0076
0077
0078
0079
0080
0081
0082
0083
0084
0085
0086
0087
0088
0089
0090
0091
0092
0093 IFNPlugin(
0094 JetDefinition jet_def,
0095 double p, double q, double a,
0096 bool modulo_2 = false,
0097 FlavNeutraliser::measure measure_in = FlavNeutraliser::general,
0098 bool use_mass_flav = false,
0099 double pp = 1,
0100 bool recursive_neutralisation = true)
0101 : _jet_def(jet_def),
0102 _p(p), _q(q), _a(a),
0103 _modulo_2(modulo_2),
0104 _measure_in(measure_in),
0105 _use_mass_flav(use_mass_flav),
0106 _spherical_algo(_jet_def.is_spherical()),
0107 _pp(pp),
0108 _recursive(recursive_neutralisation)
0109 {check_mod2_consistency();}
0110
0111 typedef FlavNeutraliser::measure Measure;
0112
0113 void set_recursive(bool val) {_recursive = val;}
0114 bool recursive() const {return _recursive;}
0115
0116 void set_measure(FlavNeutraliser::measure measure) {_measure_in = measure;}
0117 FlavNeutraliser::measure measure() const {return _measure_in;}
0118
0119
0120 double R() const FASTJET_OVERRIDE { return _jet_def.R(); }
0121
0122
0123 std::string description() const FASTJET_OVERRIDE;
0124
0125 void run_clustering(ClusterSequence &) const FASTJET_OVERRIDE;
0126
0127 bool is_spherical() const FASTJET_OVERRIDE {return _spherical_algo;}
0128
0129 private:
0130
0131
0132
0133 void check_mod2_consistency() const;
0134
0135
0136 JetDefinition _jet_def;
0137
0138
0139 double _p;
0140 double _q;
0141
0142
0143 double _a;
0144
0145
0146 bool _modulo_2;
0147 FlavNeutraliser::measure _measure_in;
0148 bool _use_mass_flav;
0149 bool _spherical_algo;
0150 double _pp;
0151 bool _recursive;
0152 };
0153
0154
0155 typedef IFNPlugin FlavNeutraliserPlugin;
0156
0157 }
0158
0159
0160 FASTJET_END_NAMESPACE
0161
0162 #endif