CMS 3D CMS Logo

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