]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/EMCAL/AliFJWrapper.h
changed HLT "fastjet" trigger: new set methods in AliFJWrapper using chars.... re...
[u/mrichter/AliRoot.git] / HLT / EMCAL / AliFJWrapper.h
1 #ifndef AliFJWrapper_HH
2 #define AliFJWrapper_HH
3
4 #include <vector>
5 #include <fastjet/PseudoJet.hh>
6 #include <fastjet/JetDefinition.hh>
7 #include <fastjet/ClusterSequence.hh>
8 #include <fastjet/ClusterSequenceArea.hh>
9 #include <fastjet/AreaDefinition.hh>
10 #include <fastjet/SISConePlugin.hh>
11
12 #include "TNamed.h"
13
14 class AliFastJetHeader;
15
16 class AliFJWrapper : public TNamed
17 {
18  public:
19
20   AliFJWrapper(const char *_name, const char *_title);
21
22   virtual ~AliFJWrapper();
23   
24   virtual void CopySettingsFrom (const AliFJWrapper&     _wrapper);
25   virtual void SetupFromHeader  (const AliFastJetHeader* _header);
26
27   virtual void AddInputVector (double _px, double _py, double _pz, double _E, int _index = -99999);
28
29   virtual void AddInputVector (const fastjet::PseudoJet& _vec,                int _index = -99999);
30   virtual void AddInputVectors(const std::vector<fastjet::PseudoJet>& _vecs,  int _offsetIndex = -99999);
31
32   virtual Int_t  Run();
33   virtual void   Clear(const Option_t* _opt = "");
34
35   void SetStrategy     (const fastjet::Strategy _strat)             { fStrategy_ = _strat;  }
36   void SetAlgorithm    (const fastjet::JetAlgorithm _algor)         { fAlgor_    = _algor;  }
37   void SetRecombScheme (const fastjet::RecombinationScheme _scheme) { fScheme_   = _scheme; }
38   void SetAreaType     (const fastjet::AreaType _atype)             { fAreaType_ = _atype;  }
39
40   void SetupAlgorithmfromOpt (const char *option);
41   void SetupStrategyfromOpt  (const char *option);
42   void SetupSchemefromOpt    (const char *option);
43   void SetupAreaTypefromOpt  (const char *option);
44
45   void SetNRepeats      (const int    _nrepeat) { fNGhostRepeats_ = _nrepeat; }
46   void SetGhostArea     (const double _gharea)  { fGhostArea_     = _gharea;  }
47   void SetMaxRap        (const double _maxrap)  { fMaxRap_        = _maxrap;  }
48   void SetR             (const double _r)       { fR_             = _r;       }
49   void SetGridScatter   (const double _gridSc)  { fGridScatter_   = _gridSc;  }
50   void SetKtScatter     (const double _ktSc)    { fKtScatter_     = _ktSc;    }
51   void SetMeanGhostKt   (const double _meankt)  { fMeanGhostKt_   = _meankt;  }
52   void SetPluginAlgor   (const int    _plugin)  { fPluginAlgor_   = _plugin;  }
53   void SetUseArea4Vector(const bool   _useA4v)  { fUseArea4Vector_ = _useA4v; }
54
55   virtual void   GetMedianAndSigma(double& _median, double& _sigma, int _remove = 0);
56   virtual double GetMedianUsedForBgSubtraction() {return fMedianUsedForBgSubtraction_;}
57
58   std::vector<fastjet::PseudoJet> GetInputVectors()    { return fInputVectors_;}
59   std::vector<fastjet::PseudoJet> GetInclusiveJets()   { return fInclusiveJets_; }
60
61   std::vector<fastjet::PseudoJet> GetJetConstituents(const unsigned int _idx);
62
63   fastjet::ClusterSequenceArea*   GetClusterSequence() { return fClustSeq_;      }
64
65   std::vector<double>             GetSubtractedJetsPts(const double _median_pt = -1, const bool _sorted = false);
66
67   double                          GetJetArea         (const unsigned int _idx);
68   double                          GetJetSubtractedPt (const unsigned int _idx);
69
70  protected:
71
72   AliFJWrapper();
73   AliFJWrapper(const AliFJWrapper& _wrapper);
74   AliFJWrapper& operator = (const AliFJWrapper& _wrapper);// {return *this;}
75
76   std::vector<fastjet::PseudoJet> fInputVectors_;
77   std::vector<fastjet::PseudoJet> fInclusiveJets_;
78
79   std::vector<double> fSubtractedJetsPt_;
80   
81   fastjet::AreaDefinition        *fAreaDef_;//!
82   fastjet::VoronoiAreaSpec       *fVorAreaSpec_;//!
83   fastjet::GhostedAreaSpec       *fGhostedAreaSpec_;//!
84   fastjet::JetDefinition         *fJetDef_;//!
85   fastjet::JetDefinition::Plugin *fPlugin_;//!
86   fastjet::RangeDefinition       *fRange_;//!
87   fastjet::ClusterSequenceArea   *fClustSeq_;//!
88
89   // same as in the alifjheader
90   fastjet::Strategy               fStrategy_;
91   fastjet::JetAlgorithm           fAlgor_;
92   fastjet::RecombinationScheme    fScheme_;
93   fastjet::AreaType               fAreaType_;
94
95   int                             fNGhostRepeats_;
96   double                          fGhostArea_;
97   double                          fMaxRap_;
98   double                          fR_;
99
100   // not included in the header (within AliFastJetHeader)
101   // no setters for the moment - used default values in the constructor
102   double                          fGridScatter_;
103   double                          fKtScatter_;
104   double                          fMeanGhostKt_;
105   int                             fPluginAlgor_;
106
107   // extra parameters
108   double                          fMedianUsedForBgSubtraction_;
109   bool                            fUseArea4Vector_;
110
111   virtual void   SubtractBackground(const double _median_pt = -1);
112
113   ClassDef(AliFJWrapper, 0)
114 };
115
116 #endif