CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
ScaleCalculators.h
Go to the documentation of this file.
1 #ifndef RecoJets_FFTJetAlgorithms_ScaleCalculators_h
2 #define RecoJets_FFTJetAlgorithms_ScaleCalculators_h
3 
4 #include <vector>
5 #include <algorithm>
6 
7 #include "fftjet/SimpleFunctors.hh"
8 #include "fftjet/RecombinedJet.hh"
9 #include "fftjet/LinearInterpolator2d.hh"
10 #include "fftjet/JetMagnitudeMapper2d.hh"
11 
13 
14 namespace fftjetcms {
15  // Return a predefined constant
16  template <typename Arg1>
17  class ConstDouble : public fftjet::Functor1<double, Arg1> {
18  public:
19  inline ConstDouble(const double value) : c_(value) {}
20  ConstDouble() = delete;
21  inline double operator()(const Arg1&) const override { return c_; }
22 
23  private:
24  double c_;
25  };
26 
27  // Multiplication of a constant by the "scale()" member
28  template <class T>
29  class ProportionalToScale : public fftjet::Functor1<double, T> {
30  public:
31  inline ProportionalToScale(const double value) : c_(value) {}
32  ProportionalToScale() = delete;
33  inline double operator()(const T& r) const override { return r.scale() * c_; }
34 
35  private:
36  double c_;
37  };
38 
39  // Multiplication by a constant
40  template <class T>
41  class MultiplyByConst : public fftjet::Functor1<double, T> {
42  public:
43  inline MultiplyByConst(const double factor, const fftjet::Functor1<double, T>* f, const bool takeOwnership = false)
44  : c_(factor), func_(f), ownsPointer_(takeOwnership) {}
45  MultiplyByConst() = delete;
46 
47  inline ~MultiplyByConst() override {
48  if (ownsPointer_)
49  delete func_;
50  }
51 
52  inline double operator()(const T& r) const override { return (*func_)(r)*c_; }
53 
54  private:
55  double c_;
56  const fftjet::Functor1<double, T>* func_;
57  const bool ownsPointer_;
58  };
59 
60  // Function composition
61  template <class T>
62  class CompositeFunctor : public fftjet::Functor1<double, T> {
63  public:
64  inline CompositeFunctor(const fftjet::Functor1<double, double>* f1,
65  const fftjet::Functor1<double, T>* f2,
66  const bool takeOwnership = false)
67  : f1_(f1), f2_(f2), ownsPointers_(takeOwnership) {}
68  CompositeFunctor() = delete;
69 
70  inline ~CompositeFunctor() override {
71  if (ownsPointers_) {
72  delete f1_;
73  delete f2_;
74  }
75  }
76 
77  inline double operator()(const T& r) const override { return (*f1_)((*f2_)(r)); }
78 
79  private:
80  const fftjet::Functor1<double, double>* f1_;
81  const fftjet::Functor1<double, T>* f2_;
82  const bool ownsPointers_;
83  };
84 
85  // Product of two functors
86  template <class T>
87  class ProductFunctor : public fftjet::Functor1<double, T> {
88  public:
89  inline ProductFunctor(const fftjet::Functor1<double, T>* f1,
90  const fftjet::Functor1<double, T>* f2,
91  const bool takeOwnership = false)
92  : f1_(f1), f2_(f2), ownsPointers_(takeOwnership) {}
93  ProductFunctor() = delete;
94 
95  inline ~ProductFunctor() override {
96  if (ownsPointers_) {
97  delete f1_;
98  delete f2_;
99  }
100  }
101 
102  inline double operator()(const T& r) const override { return (*f1_)(r) * (*f2_)(r); }
103 
104  private:
105  const fftjet::Functor1<double, T>* f1_;
106  const fftjet::Functor1<double, T>* f2_;
107  const bool ownsPointers_;
108  };
109 
110  // Function dependent on magnitude
111  template <class T>
112  class MagnitudeDependent : public fftjet::Functor1<double, T> {
113  public:
114  inline MagnitudeDependent(const fftjet::Functor1<double, double>* f1, const bool takeOwnership = false)
115  : f1_(f1), ownsPointer_(takeOwnership) {}
116  MagnitudeDependent() = delete;
117 
118  inline ~MagnitudeDependent() override {
119  if (ownsPointer_)
120  delete f1_;
121  }
122 
123  inline double operator()(const T& r) const override { return (*f1_)(r.magnitude()); }
124 
125  private:
126  const fftjet::Functor1<double, double>* f1_;
127  const bool ownsPointer_;
128  };
129 
130  // Functions dependent on peak eta
131  class PeakEtaDependent : public fftjet::Functor1<double, fftjet::Peak> {
132  public:
133  inline PeakEtaDependent(const fftjet::Functor1<double, double>* f1, const bool takeOwnership = false)
134  : f1_(f1), ownsPointer_(takeOwnership) {}
135  PeakEtaDependent() = delete;
136 
137  inline ~PeakEtaDependent() override {
138  if (ownsPointer_)
139  delete f1_;
140  }
141 
142  inline double operator()(const fftjet::Peak& r) const override { return (*f1_)(r.eta()); }
143 
144  private:
145  const fftjet::Functor1<double, double>* f1_;
146  const bool ownsPointer_;
147  };
148 
149  class PeakEtaMagSsqDependent : public fftjet::Functor1<double, fftjet::Peak> {
150  public:
151  inline PeakEtaMagSsqDependent(const fftjet::LinearInterpolator2d* f1,
152  const bool takeOwnership,
153  const fftjet::JetMagnitudeMapper2d<fftjet::Peak>* jmmp,
154  const bool ownjmp,
155  const double fc)
156  : f1_(f1), ownsPointer_(takeOwnership), jmmp_(jmmp), ownsjmmpPointer_(ownjmp), factor_(fc) {}
157  PeakEtaMagSsqDependent() = delete;
158 
159  inline ~PeakEtaMagSsqDependent() override {
160  if (ownsPointer_)
161  delete f1_;
162  if (ownsjmmpPointer_)
163  delete jmmp_;
164  }
165 
166  inline double operator()(const fftjet::Peak& r) const override {
167  const double scale = r.scale();
168  const double magnitude = r.magnitude();
169  const double pt = scale * scale * factor_ * magnitude;
170  const double partonpt = (*jmmp_)(pt, r);
171  return (*f1_)(std::abs(r.eta()), partonpt);
172  }
173 
174  private:
175  const fftjet::LinearInterpolator2d* f1_;
176  const bool ownsPointer_;
177  const fftjet::JetMagnitudeMapper2d<fftjet::Peak>* jmmp_;
178  const bool ownsjmmpPointer_;
179  const double factor_;
180  };
181 
182  // Functions dependent on jet eta
183  class JetEtaDependent : public fftjet::Functor1<double, fftjet::RecombinedJet<VectorLike> > {
184  public:
185  inline JetEtaDependent(const fftjet::Functor1<double, double>* f1, const bool takeOwnership = false)
186  : f1_(f1), ownsPointer_(takeOwnership) {}
187  JetEtaDependent() = delete;
188 
189  inline ~JetEtaDependent() override {
190  if (ownsPointer_)
191  delete f1_;
192  }
193 
194  inline double operator()(const fftjet::RecombinedJet<VectorLike>& r) const override {
195  return (*f1_)(r.vec().eta());
196  }
197 
198  private:
199  const fftjet::Functor1<double, double>* f1_;
200  const bool ownsPointer_;
201  };
202 
203  // A simple polynomial. Coefficients are in the order c0, c1, c2, ...
204  class Polynomial : public fftjet::Functor1<double, double> {
205  public:
206  inline Polynomial(const std::vector<double>& coeffs) : coeffs_(nullptr), nCoeffs(coeffs.size()) {
207  if (nCoeffs) {
208  coeffs_ = new double[nCoeffs];
209  std::copy(coeffs.begin(), coeffs.end(), coeffs_);
210  }
211  }
212  Polynomial() = delete;
213  inline ~Polynomial() override { delete[] coeffs_; }
214 
215  inline double operator()(const double& x) const override {
216  double sum = 0.0;
217  const double* p = coeffs_ + nCoeffs - 1;
218  for (unsigned i = 0; i < nCoeffs; ++i) {
219  sum *= x;
220  sum += *p--;
221  }
222  return sum;
223  }
224 
225  private:
226  double* coeffs_;
227  const unsigned nCoeffs;
228  };
229 } // namespace fftjetcms
230 
231 #endif // RecoJets_FFTJetAlgorithms_ScaleCalculators_h
const fftjet::Functor1< double, double > * f1_
const fftjet::Functor1< double, T > * func_
double operator()(const fftjet::Peak &r) const override
MultiplyByConst(const double factor, const fftjet::Functor1< double, T > *f, const bool takeOwnership=false)
double operator()(const double &x) const override
JetEtaDependent(const fftjet::Functor1< double, double > *f1, const bool takeOwnership=false)
const fftjet::Functor1< double, T > * f1_
double operator()(const T &r) const override
const fftjet::Functor1< double, double > * f1_
ProportionalToScale(const double value)
double operator()(const T &r) const override
MagnitudeDependent(const fftjet::Functor1< double, double > *f1, const bool takeOwnership=false)
double operator()(const Arg1 &) const override
const fftjet::Functor1< double, T > * f2_
CompositeFunctor(const fftjet::Functor1< double, double > *f1, const fftjet::Functor1< double, T > *f2, const bool takeOwnership=false)
double operator()(const fftjet::RecombinedJet< VectorLike > &r) const override
ProductFunctor(const fftjet::Functor1< double, T > *f1, const fftjet::Functor1< double, T > *f2, const bool takeOwnership=false)
Polynomial(const std::vector< double > &coeffs)
PeakEtaMagSsqDependent(const fftjet::LinearInterpolator2d *f1, const bool takeOwnership, const fftjet::JetMagnitudeMapper2d< fftjet::Peak > *jmmp, const bool ownjmp, const double fc)
double operator()(const T &r) const override
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const fftjet::LinearInterpolator2d * f1_
PeakEtaDependent(const fftjet::Functor1< double, double > *f1, const bool takeOwnership=false)
const fftjet::Functor1< double, double > * f1_
const fftjet::Functor1< double, T > * f2_
double operator()(const T &r) const override
double operator()(const T &r) const override
const fftjet::JetMagnitudeMapper2d< fftjet::Peak > * jmmp_
double operator()(const fftjet::Peak &r) const override
long double T
const fftjet::Functor1< double, double > * f1_
tuple size
Write out results.
ConstDouble(const double value)