Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-30 10:06:42

0001 #ifndef __DISTANCEMEASURE_HH__
0002 #define __DISTANCEMEASURE_HH__
0003 
0004 #include <float.h>
0005 #include <algorithm>
0006 #include "fastjet/PseudoJet.hh"
0007 
0008 FASTJET_BEGIN_NAMESPACE      // defined in fastjet/internal/base.hh
0009 
0010 namespace contrib {
0011     namespace QCDAwarePlugin {
0012 
0013         /**
0014          * the DistanceMeasure abstract base class calculates the
0015          * distance between two PseudoJets. the default implemented
0016          * DistanceMeasures include kt, C/A, and anti-kt
0017          */
0018 
0019         class DistanceMeasure {
0020             public:
0021                 virtual double dij(const fastjet::PseudoJet& pji, const fastjet::PseudoJet& pjj) const = 0;
0022                 virtual double diB(const fastjet::PseudoJet& pji) const = 0;
0023                 virtual double R() const = 0;
0024                 virtual std::string algname() const = 0;
0025 
0026                 virtual ~DistanceMeasure() {};
0027         };
0028 
0029 
0030         class KtMeasure : public DistanceMeasure {
0031             private:
0032                 double _R;
0033 
0034             public:
0035                 KtMeasure(double R) : _R(R) {}
0036 
0037                 inline double dij(const fastjet::PseudoJet& pji, const fastjet::PseudoJet& pjj) const {
0038                     double drbyR2 = pji.squared_distance(pjj) / (_R * _R);
0039                     return std::min(pji.perp2(), pjj.perp2()) * drbyR2;
0040                 }
0041 
0042                 inline double diB(const fastjet::PseudoJet& pji) const {
0043                     return pji.perp2();
0044                 }
0045 
0046                 virtual double R() const {
0047                     return _R;
0048                 }
0049 
0050                 std::string algname() const {
0051                     return "kt";
0052                 }
0053         };
0054 
0055 
0056         class AntiKtMeasure : public DistanceMeasure {
0057             private:
0058                 double _R;
0059 
0060             public:
0061                 AntiKtMeasure(double R) : _R(R) {}
0062 
0063                 inline double dij(const fastjet::PseudoJet& pji, const fastjet::PseudoJet& pjj) const {
0064                     double drbyR2 = pji.squared_distance(pjj) / (_R * _R);
0065                     return 1.0 / std::max(pji.perp2(), pjj.perp2()) * drbyR2;
0066                 }
0067 
0068                 inline double diB(const fastjet::PseudoJet& pji) const {
0069                     return 1.0 / pji.perp2();
0070                 }
0071 
0072                 virtual double R() const {
0073                     return _R;
0074                 }
0075 
0076                 std::string algname() const {
0077                     return "anti-kt";
0078                 }
0079         };
0080 
0081 
0082         class CAMeasure : public DistanceMeasure {
0083             private:
0084                 double _R;
0085 
0086             public:
0087                 CAMeasure(double R) : _R(R) {}
0088 
0089                 inline double dij(const fastjet::PseudoJet& pji, const fastjet::PseudoJet& pjj) const {
0090                     return pji.squared_distance(pjj) / (_R * _R);
0091                 }
0092 
0093 
0094                 inline double diB(const fastjet::PseudoJet& pji) const {
0095                     return 1.0;
0096                 }
0097 
0098                 virtual double R() const {
0099                     return _R;
0100                 }
0101 
0102                 std::string algname() const {
0103                     return "Cambridge-Aachen";
0104                 }
0105         };
0106 
0107     } // QCDAware
0108 } // contrib
0109 
0110 FASTJET_END_NAMESPACE
0111 
0112 #endif