File indexing completed on 2024-04-06 12:25:22
0001 #ifndef RecoJets_FFTJetAlgorithms_ScaleCalculators_h
0002 #define RecoJets_FFTJetAlgorithms_ScaleCalculators_h
0003
0004 #include <vector>
0005 #include <algorithm>
0006
0007 #include "fftjet/SimpleFunctors.hh"
0008 #include "fftjet/RecombinedJet.hh"
0009 #include "fftjet/LinearInterpolator2d.hh"
0010 #include "fftjet/JetMagnitudeMapper2d.hh"
0011
0012 #include "RecoJets/FFTJetAlgorithms/interface/fftjetTypedefs.h"
0013
0014 namespace fftjetcms {
0015
0016 template <typename Arg1>
0017 class ConstDouble : public fftjet::Functor1<double, Arg1> {
0018 public:
0019 inline ConstDouble(const double value) : c_(value) {}
0020 ConstDouble() = delete;
0021 inline double operator()(const Arg1&) const override { return c_; }
0022
0023 private:
0024 double c_;
0025 };
0026
0027
0028 template <class T>
0029 class ProportionalToScale : public fftjet::Functor1<double, T> {
0030 public:
0031 inline ProportionalToScale(const double value) : c_(value) {}
0032 ProportionalToScale() = delete;
0033 inline double operator()(const T& r) const override { return r.scale() * c_; }
0034
0035 private:
0036 double c_;
0037 };
0038
0039
0040 template <class T>
0041 class MultiplyByConst : public fftjet::Functor1<double, T> {
0042 public:
0043 inline MultiplyByConst(const double factor, const fftjet::Functor1<double, T>* f, const bool takeOwnership = false)
0044 : c_(factor), func_(f), ownsPointer_(takeOwnership) {}
0045 MultiplyByConst() = delete;
0046
0047 inline ~MultiplyByConst() override {
0048 if (ownsPointer_)
0049 delete func_;
0050 }
0051
0052 inline double operator()(const T& r) const override { return (*func_)(r)*c_; }
0053
0054 private:
0055 double c_;
0056 const fftjet::Functor1<double, T>* func_;
0057 const bool ownsPointer_;
0058 };
0059
0060
0061 template <class T>
0062 class CompositeFunctor : public fftjet::Functor1<double, T> {
0063 public:
0064 inline CompositeFunctor(const fftjet::Functor1<double, double>* f1,
0065 const fftjet::Functor1<double, T>* f2,
0066 const bool takeOwnership = false)
0067 : f1_(f1), f2_(f2), ownsPointers_(takeOwnership) {}
0068 CompositeFunctor() = delete;
0069
0070 inline ~CompositeFunctor() override {
0071 if (ownsPointers_) {
0072 delete f1_;
0073 delete f2_;
0074 }
0075 }
0076
0077 inline double operator()(const T& r) const override { return (*f1_)((*f2_)(r)); }
0078
0079 private:
0080 const fftjet::Functor1<double, double>* f1_;
0081 const fftjet::Functor1<double, T>* f2_;
0082 const bool ownsPointers_;
0083 };
0084
0085
0086 template <class T>
0087 class ProductFunctor : public fftjet::Functor1<double, T> {
0088 public:
0089 inline ProductFunctor(const fftjet::Functor1<double, T>* f1,
0090 const fftjet::Functor1<double, T>* f2,
0091 const bool takeOwnership = false)
0092 : f1_(f1), f2_(f2), ownsPointers_(takeOwnership) {}
0093 ProductFunctor() = delete;
0094
0095 inline ~ProductFunctor() override {
0096 if (ownsPointers_) {
0097 delete f1_;
0098 delete f2_;
0099 }
0100 }
0101
0102 inline double operator()(const T& r) const override { return (*f1_)(r) * (*f2_)(r); }
0103
0104 private:
0105 const fftjet::Functor1<double, T>* f1_;
0106 const fftjet::Functor1<double, T>* f2_;
0107 const bool ownsPointers_;
0108 };
0109
0110
0111 template <class T>
0112 class MagnitudeDependent : public fftjet::Functor1<double, T> {
0113 public:
0114 inline MagnitudeDependent(const fftjet::Functor1<double, double>* f1, const bool takeOwnership = false)
0115 : f1_(f1), ownsPointer_(takeOwnership) {}
0116 MagnitudeDependent() = delete;
0117
0118 inline ~MagnitudeDependent() override {
0119 if (ownsPointer_)
0120 delete f1_;
0121 }
0122
0123 inline double operator()(const T& r) const override { return (*f1_)(r.magnitude()); }
0124
0125 private:
0126 const fftjet::Functor1<double, double>* f1_;
0127 const bool ownsPointer_;
0128 };
0129
0130
0131 class PeakEtaDependent : public fftjet::Functor1<double, fftjet::Peak> {
0132 public:
0133 inline PeakEtaDependent(const fftjet::Functor1<double, double>* f1, const bool takeOwnership = false)
0134 : f1_(f1), ownsPointer_(takeOwnership) {}
0135 PeakEtaDependent() = delete;
0136
0137 inline ~PeakEtaDependent() override {
0138 if (ownsPointer_)
0139 delete f1_;
0140 }
0141
0142 inline double operator()(const fftjet::Peak& r) const override { return (*f1_)(r.eta()); }
0143
0144 private:
0145 const fftjet::Functor1<double, double>* f1_;
0146 const bool ownsPointer_;
0147 };
0148
0149 class PeakEtaMagSsqDependent : public fftjet::Functor1<double, fftjet::Peak> {
0150 public:
0151 inline PeakEtaMagSsqDependent(const fftjet::LinearInterpolator2d* f1,
0152 const bool takeOwnership,
0153 const fftjet::JetMagnitudeMapper2d<fftjet::Peak>* jmmp,
0154 const bool ownjmp,
0155 const double fc)
0156 : f1_(f1), ownsPointer_(takeOwnership), jmmp_(jmmp), ownsjmmpPointer_(ownjmp), factor_(fc) {}
0157 PeakEtaMagSsqDependent() = delete;
0158
0159 inline ~PeakEtaMagSsqDependent() override {
0160 if (ownsPointer_)
0161 delete f1_;
0162 if (ownsjmmpPointer_)
0163 delete jmmp_;
0164 }
0165
0166 inline double operator()(const fftjet::Peak& r) const override {
0167 const double scale = r.scale();
0168 const double magnitude = r.magnitude();
0169 const double pt = scale * scale * factor_ * magnitude;
0170 const double partonpt = (*jmmp_)(pt, r);
0171 return (*f1_)(std::abs(r.eta()), partonpt);
0172 }
0173
0174 private:
0175 const fftjet::LinearInterpolator2d* f1_;
0176 const bool ownsPointer_;
0177 const fftjet::JetMagnitudeMapper2d<fftjet::Peak>* jmmp_;
0178 const bool ownsjmmpPointer_;
0179 const double factor_;
0180 };
0181
0182
0183 class JetEtaDependent : public fftjet::Functor1<double, fftjet::RecombinedJet<VectorLike> > {
0184 public:
0185 inline JetEtaDependent(const fftjet::Functor1<double, double>* f1, const bool takeOwnership = false)
0186 : f1_(f1), ownsPointer_(takeOwnership) {}
0187 JetEtaDependent() = delete;
0188
0189 inline ~JetEtaDependent() override {
0190 if (ownsPointer_)
0191 delete f1_;
0192 }
0193
0194 inline double operator()(const fftjet::RecombinedJet<VectorLike>& r) const override {
0195 return (*f1_)(r.vec().eta());
0196 }
0197
0198 private:
0199 const fftjet::Functor1<double, double>* f1_;
0200 const bool ownsPointer_;
0201 };
0202
0203
0204 class Polynomial : public fftjet::Functor1<double, double> {
0205 public:
0206 inline Polynomial(const std::vector<double>& coeffs) : coeffs_(nullptr), nCoeffs(coeffs.size()) {
0207 if (nCoeffs) {
0208 coeffs_ = new double[nCoeffs];
0209 std::copy(coeffs.begin(), coeffs.end(), coeffs_);
0210 }
0211 }
0212 Polynomial() = delete;
0213 inline ~Polynomial() override { delete[] coeffs_; }
0214
0215 inline double operator()(const double& x) const override {
0216 double sum = 0.0;
0217 const double* p = coeffs_ + nCoeffs - 1;
0218 for (unsigned i = 0; i < nCoeffs; ++i) {
0219 sum *= x;
0220 sum += *p--;
0221 }
0222 return sum;
0223 }
0224
0225 private:
0226 double* coeffs_;
0227 const unsigned nCoeffs;
0228 };
0229 }
0230
0231 #endif