CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch2/src/RecoJets/FFTJetAlgorithms/interface/ScaleCalculators.h

Go to the documentation of this file.
00001 #ifndef RecoJets_FFTJetAlgorithms_ScaleCalculators_h
00002 #define RecoJets_FFTJetAlgorithms_ScaleCalculators_h
00003 
00004 #include <vector>
00005 #include <algorithm>
00006 
00007 #include "fftjet/SimpleFunctors.hh"
00008 #include "fftjet/RecombinedJet.hh"
00009 
00010 #include "RecoJets/FFTJetAlgorithms/interface/fftjetTypedefs.h"
00011 
00012 namespace fftjetcms {
00013     // Return a predefined constant
00014     template <typename Arg1>
00015     class ConstDouble : public fftjet::Functor1<double,Arg1>
00016     {
00017     public:
00018         inline ConstDouble(const double value) : c_(value) {}
00019         inline double operator()(const Arg1&) const {return c_;}
00020 
00021     private:
00022         ConstDouble();
00023         double c_;
00024     };
00025 
00026 
00027     // Multiplication of a constant by the "scale()" member
00028     template <class T>
00029     class ProportionalToScale : public fftjet::Functor1<double,T>
00030     {
00031     public:
00032         inline ProportionalToScale(const double value) : c_(value) {}
00033         inline double operator()(const T& r) const {return r.scale()*c_;}
00034 
00035     private:
00036         ProportionalToScale();
00037         double c_;
00038     };
00039 
00040 
00041     // Multiplication by a constant
00042     template <class T>
00043     class MultiplyByConst : public fftjet::Functor1<double,T>
00044     {
00045     public:
00046         inline MultiplyByConst(const double factor,
00047                                const fftjet::Functor1<double,T>* f,
00048                                const bool takeOwnership=false)
00049             : c_(factor), func_(f), ownsPointer_(takeOwnership) {}
00050 
00051         inline ~MultiplyByConst() {if (ownsPointer_) delete func_;}
00052 
00053         inline double operator()(const T& r) const {return (*func_)(r)*c_;}
00054 
00055     private:
00056         MultiplyByConst();
00057         double c_;
00058         const fftjet::Functor1<double,T>* func_;
00059         const bool ownsPointer_;
00060     };
00061 
00062 
00063     // Function composition
00064     template <class T>
00065     class CompositeFunctor : public fftjet::Functor1<double,T>
00066     {
00067     public:
00068         inline CompositeFunctor(const fftjet::Functor1<double,double>* f1,
00069                                 const fftjet::Functor1<double,T>* f2,
00070                                 const bool takeOwnership=false)
00071             : f1_(f1), f2_(f2), ownsPointers_(takeOwnership) {}
00072 
00073         inline ~CompositeFunctor()
00074             {if (ownsPointers_) {delete f1_; delete f2_;}}
00075 
00076         inline double operator()(const T& r) const {return (*f1_)((*f2_)(r));}
00077 
00078     private:
00079         CompositeFunctor();
00080         const fftjet::Functor1<double,double>* f1_;
00081         const fftjet::Functor1<double,T>* f2_;
00082         const bool ownsPointers_;
00083     };
00084 
00085 
00086     // Product of two functors
00087     template <class T>
00088     class ProductFunctor : public fftjet::Functor1<double,T>
00089     {
00090     public:
00091         inline ProductFunctor(const fftjet::Functor1<double,T>* f1,
00092                               const fftjet::Functor1<double,T>* f2,
00093                               const bool takeOwnership=false)
00094             : f1_(f1), f2_(f2), ownsPointers_(takeOwnership) {}
00095 
00096         inline ~ProductFunctor()
00097             {if (ownsPointers_) {delete f1_; delete f2_;}}
00098 
00099         inline double operator()(const T& r) const
00100             {return (*f1_)(r) * (*f2_)(r);}
00101 
00102     private:
00103         ProductFunctor();
00104         const fftjet::Functor1<double,T>* f1_;
00105         const fftjet::Functor1<double,T>* f2_;
00106         const bool ownsPointers_;
00107     };
00108 
00109 
00110     // Function dependent on magnitude
00111     template <class T>
00112     class MagnitudeDependent : public fftjet::Functor1<double,T>
00113     {
00114     public:
00115         inline MagnitudeDependent(const fftjet::Functor1<double,double>* f1,
00116                                   const bool takeOwnership=false)
00117             : f1_(f1), ownsPointer_(takeOwnership) {}
00118 
00119         inline ~MagnitudeDependent() {if (ownsPointer_) delete f1_;}
00120 
00121         inline double operator()(const T& r) const
00122             {return (*f1_)(r.magnitude());}
00123 
00124     private:
00125         MagnitudeDependent();
00126         const fftjet::Functor1<double,double>* f1_;
00127         const bool ownsPointer_;
00128     };
00129 
00130 
00131     // Functions dependent on peak eta
00132     class PeakEtaDependent : public fftjet::Functor1<double,fftjet::Peak>
00133     {
00134     public:
00135         inline PeakEtaDependent(const fftjet::Functor1<double,double>* f1,
00136                                 const bool takeOwnership=false)
00137             : f1_(f1), ownsPointer_(takeOwnership) {}
00138 
00139         inline ~PeakEtaDependent() {if (ownsPointer_) delete f1_;}
00140 
00141         inline double operator()(const fftjet::Peak& r) const
00142             {return (*f1_)(r.eta());}
00143 
00144     private:
00145         PeakEtaDependent();
00146         const fftjet::Functor1<double,double>* f1_;
00147         const bool ownsPointer_;
00148     };
00149 
00150 
00151     class PeakEtaMagSsqDependent : public fftjet::Functor1<double,fftjet::Peak>
00152     {
00153     public:
00154         inline PeakEtaMagSsqDependent(
00155             const fftjet::LinearInterpolator2d* f1,
00156             const bool takeOwnership,
00157             const fftjet::JetMagnitudeMapper2d<fftjet::Peak>* jmmp,
00158             const bool ownjmp, const double fc)
00159             : f1_(f1), ownsPointer_(takeOwnership), jmmp_(jmmp),
00160               ownsjmmpPointer_(ownjmp), factor_(fc) {}
00161 
00162         inline ~PeakEtaMagSsqDependent()
00163         {
00164             if (ownsPointer_) delete f1_;
00165             if (ownsjmmpPointer_) delete jmmp_;
00166         }
00167 
00168         inline double operator()(const fftjet::Peak& r) const
00169         {
00170             const double scale = r.scale();
00171             const double magnitude = r.magnitude();
00172             const double pt = scale*scale*factor_*magnitude;
00173             const double partonpt = (*jmmp_)(pt,r);
00174             return (*f1_)(std::abs(r.eta()),partonpt)  ;
00175         }
00176 
00177     private:
00178         PeakEtaMagSsqDependent();
00179         const fftjet::LinearInterpolator2d* f1_;
00180         const bool ownsPointer_;
00181         const fftjet::JetMagnitudeMapper2d <fftjet::Peak>* jmmp_;
00182         const bool ownsjmmpPointer_;
00183         const double factor_;
00184     };
00185 
00186 
00187     // Functions dependent on jet eta
00188     class JetEtaDependent : 
00189         public fftjet::Functor1<double,fftjet::RecombinedJet<VectorLike> >
00190     {
00191     public:
00192         inline JetEtaDependent(const fftjet::Functor1<double,double>* f1,
00193                                const bool takeOwnership=false)
00194             : f1_(f1), ownsPointer_(takeOwnership) {}
00195 
00196         inline ~JetEtaDependent() {if (ownsPointer_) delete f1_;}
00197 
00198         inline double operator()(
00199             const fftjet::RecombinedJet<VectorLike>& r) const
00200             {return (*f1_)(r.vec().eta());}
00201 
00202     private:
00203         JetEtaDependent();
00204         const fftjet::Functor1<double,double>* f1_;
00205         const bool ownsPointer_;
00206     };
00207 
00208 
00209     // A simple polynomial. Coefficients are in the order c0, c1, c2, ...
00210     class Polynomial : public fftjet::Functor1<double,double>
00211     {
00212     public:
00213         inline Polynomial(const std::vector<double>& coeffs)
00214             : coeffs_(0), nCoeffs(coeffs.size())
00215         {
00216             if (nCoeffs)
00217             {
00218                 coeffs_ = new double[nCoeffs];
00219                 std::copy(coeffs.begin(), coeffs.end(), coeffs_);
00220             }
00221         }
00222         inline ~Polynomial() {delete [] coeffs_;}
00223 
00224         inline double operator()(const double& x) const
00225         {
00226             double sum = 0.0;
00227             const double* p = coeffs_ + nCoeffs - 1;
00228             for (unsigned i=0; i<nCoeffs; ++i)
00229             {
00230                 sum *= x;
00231                 sum += *p--;
00232             }
00233             return sum;
00234         }
00235 
00236     private:
00237         Polynomial();
00238         double* coeffs_;
00239         const unsigned nCoeffs;
00240     };
00241 }
00242 
00243 #endif // RecoJets_FFTJetAlgorithms_ScaleCalculators_h