]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/DEV/AliSISConeJetHeader.h
file with geometry definitions, added temporarily until the one in EVE/alice-data...
[u/mrichter/AliRoot.git] / JETAN / DEV / AliSISConeJetHeader.h
1 #ifndef ALISISCONEJETHEADER_H
2 #define ALISISCONEJETHEADER_H
3  
4 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice                               */
6
7 /* $Id$ */
8  
9 //---------------------------------------------------------------------
10 // SISCone (FastJet v2.3.4) finder algorithm interface
11 // Finder Header Class 
12 // Author: swensy.jangal@ires.in2p3.fr
13 //---------------------------------------------------------------------
14  
15 #include "fastjet/JetDefinition.hh"
16 #include "AliJetHeader.h"
17  
18 class AliSISConeJetHeader : public AliJetHeader
19 {
20  public:
21  
22   AliSISConeJetHeader();
23   virtual ~AliSISConeJetHeader() { }
24
25   // Getters
26   Bool_t                       GetBGMode()                     const           {return fBGMode;}
27   Int_t                        GetActiveAreaRepeats()          const           {return fActiveAreaRepeats;}
28   Int_t                        GetAreaTypeNumber()             const           {return fAreaTypeNumber;}
29   Int_t                        GetBGAlgorithm()                const           {return fBGAlgo;}        
30   Int_t                        GetNPassMax()                   const           {return fNPassMax;}
31   Int_t                        GetSplitMergeScale()            const           {return fSplitMergeScaleNumber;}
32   Double_t                     GetGhostEtaMax()                const           {return fGhostEtaMax;}
33   Double_t                     GetGhostArea()                  const           {return fGhostArea;}
34   Double_t                     GetEffectiveRFact()             const           {return fEffectiveRFact;}
35   Double_t                     GetRapMax()                     const           {return fRapMax;}
36   Double_t                     GetRapMin()                     const           {return fRapMin;}
37   Double_t                     GetPhiMax()                     const           {return fPhiMax;}
38   Double_t                     GetPhiMin()                     const           {return fPhiMin;}
39   Double_t                     GetConeRadius()                 const           {return fConeRadius;}
40   Double_t                     GetOverlapThreshold()           const           {return fOverlapThreshold;}
41   Double_t                     GetPtProtojetMin()              const           {return fPtProtoJetMin;}
42   Double_t                     GetRForRho()                    const           {return fRRho;}
43   Double_t                     GetCaching()                    const           {return fCaching;}
44   Double_t                     GetSplitMergeStoppingScale()    const           {return fSplitMergeStoppingScale;}
45   Double_t                     GetMinJetPt()                   const           {return fMinJetPt;}  
46   Double_t                     GetGridScatter()                const           {return fGridScatter;}
47   Double_t                     GetKtScatter()                  const           {return fKtScatter;}
48   Double_t                     GetMeanGhostKt()                const           {return fMeanGhostKt;}
49
50   // Setters
51   void                         SetBGAlgorithm(Int_t value)                     {fBGAlgo = value;}
52   void                         SetBGMode(Bool_t value)                         {fBGMode = value;}
53   void                         SetCaching(Bool_t value)                        {fCaching = value;}
54   void                         SetComment(TString com)                         {fComment=com;}
55   void                         SetComment(const char* com)                     {AliJetHeader::SetComment(com);}
56   void                         SetGhostEtaMax(Double_t f)                      {fGhostEtaMax = f;}
57   void                         SetGhostArea(Double_t f)                        {fGhostArea = f;}
58   void                         SetActiveAreaRepeats(Int_t f)                   {fActiveAreaRepeats =f;}
59   void                         SetAreaTypeNumber(Int_t f)                      {fAreaTypeNumber = f;}
60   void                         SetEffectiveRFact(Double_t value)               {fEffectiveRFact = value;}       
61   void                         SetConeRadius(Double_t value)                   {fConeRadius = value;}
62   void                         SetMinJetPt(Double_t value)                     {fMinJetPt = value;}
63   void                         SetNPassMax(Int_t value)                        {fNPassMax = value;}
64   void                         SetOverlapThreshold(Double_t value)             {fOverlapThreshold = value;}
65   void                         SetPhiRange(Double_t fmin, Double_t fmax)       {fPhiMin = fmin; fPhiMax = fmax;}
66   void                         SetPtProtojetMin(Double_t value)                {fPtProtoJetMin = value;}
67   void                         SetRapRange(Double_t fmin, Double_t fmax)       {fRapMin = fmin; fRapMax = fmax;}
68   void                         SetRForRho(Double_t value)                      {fRRho = value;}
69   void                         SetSplitMergeScale(Int_t value)                 {fSplitMergeScaleNumber = value;}
70   void                         SetSplitMergeStoppingScale(Double_t value)      {fSplitMergeStoppingScale = value;}        
71   void                         SetGridScatter(Double_t value)                  {fGridScatter = value;}
72   void                         SetKtScatter(Double_t value)                    {fKtScatter = value;}
73   void                         SetMeanGhostKt(Double_t value)                  {fMeanGhostKt = value;}
74
75   // Added for background stuff
76   fastjet::Strategy            GetStrategy() const                             {return fStrategy;}
77   fastjet::RecombinationScheme GetRecombScheme() const                         {return fRecombScheme;}
78   Double_t                     GetRparamBkg() const                            {return fRparamBkg;}
79   Bool_t                       Use4VectorArea() const                          {return fkUse4VectorArea;}
80
81   void                         SetUse4VectorArea()                             {fkUse4VectorArea = kTRUE;}
82   void                         SetStrategy(fastjet::Strategy f)                {fStrategy = f;}
83   void                         SetRecombScheme(fastjet::RecombinationScheme f) {fRecombScheme = f;}
84   void                         SetRparamBkg(Double_t f)                        {fRparamBkg = f;}
85
86   // others
87   void                         PrintParameters() const;
88
89  protected:
90
91
92   Int_t    fActiveAreaRepeats;                // How many times do you want to caculate active areas?
93   Int_t    fAreaTypeNumber;                   // Kind of area
94   Int_t    fBGAlgo;                           // Algorithm for rho calculus
95   Bool_t   fCaching;                          // Do we record found cones for this set of data?
96   Double_t fConeRadius;                       // Cone radius
97   Double_t fEffectiveRFact;                   // Radius for Voronoi diagram
98   Double_t fGhostEtaMax;                      // Maximum eta in which a ghost can be generated
99   Double_t fGhostArea;                        // Area of one ghost
100   Double_t fGridScatter;                      // fractional random fluctuations of the position of the ghosts on the y-phi grid
101   Double_t fKtScatter;                        // fractional random fluctuations of the tranverse momentum of the ghosts on the y-phi grid
102   Double_t fMeanGhostKt;                      // average transverse momentum of the ghosts.
103   Double_t fMinJetPt;                         // Jet minimum energy
104   Int_t    fNPassMax;                         // maximum number of passes
105   Double_t fOverlapThreshold;                 // overlap parameter
106   Double_t fPhiMax, fPhiMin;                  // Phi range
107   Double_t fPtProtoJetMin;                    // pT min of protojets
108   Double_t fRapMax, fRapMin;                  // Eta range
109   Double_t fRRho;                             // Radius to determine rho
110   Int_t    fSplitMergeScaleNumber;            // Kind of recombination in split/merge procedure, there's only one
111   Double_t fSplitMergeStoppingScale;          // Stopping scale for split/merge procedure in case of area calculus
112
113   // Added for background
114   Double_t fRparamBkg;                        //R param for bkg calculation
115   fastjet::Strategy fStrategy;                // fastjet::Best;
116   fastjet::RecombinationScheme fRecombScheme; // fastjet::BIpt_scheme;
117   Bool_t   fkUse4VectorArea;                  // Toggle use of 4-vector area
118
119   ClassDef(AliSISConeJetHeader,5)             // SISCONE header class
120 };
121  
122 #endif