]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/analysis/AliFMDAnaParameters.h
Added extra indices to the correction files
[u/mrichter/AliRoot.git] / FMD / analysis / AliFMDAnaParameters.h
1 #ifndef ALIFMDANAPARAMETERS_H
2 #define ALIFMDANAPARAMETERS_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights
4  * reserved. 
5  *
6  * Latest changes by Hans Hjersing Dalsgaard, NBI, hans.dalsgaard@cern.ch
7  *
8  * See cxx source for full Copyright notice                               
9  */
10 //
11 //The design of this class is based on the AliFMDParameters class. Its purpose
12 //is to hold parameters for the analysis such as background correction and 
13 //fit functions.
14 //
15 //Author: Hans Hjersing Dalsgaard, NBI, hans.dalsgaard@cern.ch
16 //
17 //____________________________________________________________________
18
19 #ifndef ROOT_TNamed
20 # include <TNamed.h>
21 #endif
22 #ifndef ROOT_TArrayI
23 # include <TArrayI.h>
24 #endif
25 #ifndef ALIFMDUSHORTMAP_H
26 # include <AliFMDUShortMap.h>
27 #endif
28 #ifndef ALIFMDBOOLMAP_H
29 # include <AliFMDBoolMap.h>
30 #endif
31 #include "AliCDBEntry.h"
32
33 #include "TFile.h"
34 #include "TObjArray.h"
35 #include "TH2F.h"
36 #include "TAxis.h"
37 #include "TH1F.h"
38 #include "AliFMDAnaCalibBackgroundCorrection.h"
39 #include "AliFMDAnaCalibEnergyDistribution.h"
40 #include "AliFMDAnaCalibEventSelectionEfficiency.h"
41 #include <TVector2.h>
42 #include <TString.h>
43 //#include "AliPWG0Helper.h"
44 class AliESDEvent;
45
46 /**
47  * @ingroup FMD_ana
48  */
49 class AliFMDAnaParameters : public TNamed
50 {
51 public:
52   /** Enumeration of things to initialize */ 
53   enum What { 
54     /** Pulser gain */ 
55     kBackgroundCorrection         = 0x1, // Background Correction 
56     kEnergyDistributions          = 0x2, // Energy Distributions
57     kEventSelectionEfficiency     = 0x4  // Event Selection Efficiency
58   };
59   
60   enum Trigger { kMB1 = 0, kMB2, kSPDFASTOR };
61   
62   enum Energy { k900 , k10000, k14000 };
63   
64   enum MagField {k0G, k5G};
65   
66   enum Species {kPP, kPbPb};
67   
68   /** Singleton access
69       @return  single to */
70   static AliFMDAnaParameters* Instance();
71   
72   void Init(Bool_t forceReInit=kTRUE, UInt_t what=kBackgroundCorrection|kEnergyDistributions|kEventSelectionEfficiency);
73   Float_t GetVtxCutZ();
74   Int_t GetNvtxBins();
75   Int_t GetNetaBins();
76   Float_t GetEtaMin();  
77   Float_t GetEtaMax();
78   Float_t GetMPV(Int_t det, Char_t ring, Float_t eta);
79   Float_t GetSigma(Int_t det, Char_t ring, Float_t eta);
80   Float_t Get2MIPWeight(Int_t det, Char_t ring, Float_t eta);
81   Float_t Get3MIPWeight(Int_t det, Char_t ring, Float_t eta);
82   //static const char* GetBackgroundPath() { return fgkBackgroundCorrection;}
83   // static const char* GetEdistPath()      { return fgkEnergyDists;}
84   static const char* GetBackgroundID() { return fgkBackgroundID;}
85   static const char* GetEdistID()      { return fgkEnergyDistributionID;}
86   static const char* GetEventSelectionEffID()      { return fgkEventSelectionEffID;}
87   TH2F* GetBackgroundCorrection(Int_t det, Char_t ring, Int_t vtxbin);
88   TH1F* GetDoubleHitCorrection(Int_t det, Char_t ring);
89   
90   Float_t GetEventSelectionEfficiency(Int_t vtxbin);
91   Float_t GetPhiFromSector(UShort_t det, Char_t ring, UShort_t sec) const;
92   Float_t GetEtaFromStrip(UShort_t det, Char_t ring, UShort_t sec, UShort_t strip, Float_t zvtx) const;
93   Float_t  GetStripLength(Char_t ring, UShort_t strip)  ;
94   Float_t  GetBaseStripLength(Char_t ring, UShort_t strip)  ;
95   Float_t  GetMaxR(Char_t ring) const;
96   Float_t  GetMinR(Char_t ring) const;
97   void     SetBackgroundPath(const Char_t* bgpath) {fBackgroundPath.Form(bgpath);}
98   void     SetEnergyPath(const Char_t* epath) {fEnergyPath.Form(epath);}
99   void     SetEventSelectionPath(const Char_t* evpath) {fEventSelectionEffPath.Form(evpath);}
100   void     SetProcessPrimary(Bool_t prim=kTRUE) {fProcessPrimary = prim;}
101   void     SetProcessHits(Bool_t hits=kTRUE) {fProcessHits = hits;}
102   Bool_t   GetProcessPrimary() {return fProcessPrimary;}
103   Bool_t   GetProcessHits() {return fProcessHits;}
104   void     GetVertex(AliESDEvent* esd, Double_t* vertexXYZ);
105   void     SetTriggerDefinition(Trigger trigger) {fTrigger = trigger;}
106   Trigger  GetTriggerDefinition() {return fTrigger;}
107   Bool_t   IsEventTriggered(AliESDEvent* esd);
108   void     SetEnergy(Energy energy) {fEnergy = energy;}
109   void     SetMagField(MagField magfield) {fMagField = magfield;}
110   char*    GetPath(const char* species);
111 protected:
112   
113   AliFMDAnaParameters();
114   
115   AliFMDAnaParameters(const AliFMDAnaParameters& o) 
116     : TNamed(o),
117       fIsInit(o.fIsInit),
118       fBackground(o.fBackground),
119       fEnergyDistribution(o.fEnergyDistribution),
120       fEventSelectionEfficiency(o.fEventSelectionEfficiency),
121       fCorner1(o.fCorner1),
122       fCorner2(o.fCorner2),
123       fEnergyPath(o.fEnergyPath),
124       fBackgroundPath(o.fBackgroundPath),
125       fEventSelectionEffPath(o.fEventSelectionEffPath),
126       fProcessPrimary(o.fProcessPrimary),
127       fProcessHits(o.fProcessHits),
128       fTrigger(o.fTrigger),
129       fEnergy(o.fEnergy),
130       fMagField(o.fMagField),
131       fSpecies(o.fSpecies)
132   {}
133   AliFMDAnaParameters& operator=(const AliFMDAnaParameters&) { return *this; }
134   virtual ~AliFMDAnaParameters() {}
135   
136   static AliFMDAnaParameters* fgInstance;   // Static singleton instance
137   
138   //  AliCDBEntry* GetEntry(const char* path, Bool_t fatal=kTRUE) const ;
139   void InitBackground();
140   void InitEnergyDists();
141   void InitEventSelectionEff();
142   
143   TH1F* GetEnergyDistribution(Int_t det, Char_t ring, Float_t eta);
144   TObjArray* GetBackgroundArray();
145   
146   TAxis* GetRefAxis();
147   void SetCorners(Char_t ring) ;
148   
149   Bool_t fIsInit;
150   //TObjArray*  fBackgroundArray;
151   // TObjArray*  fEdistArray;
152   AliFMDAnaCalibBackgroundCorrection*         fBackground;
153   AliFMDAnaCalibEnergyDistribution*           fEnergyDistribution;
154   AliFMDAnaCalibEventSelectionEfficiency*     fEventSelectionEfficiency;
155   
156   //static const char* fgkBackgroundCorrection;
157   //static const char* fgkEnergyDists;
158   static const char* fgkBackgroundID;
159   static const char* fgkEnergyDistributionID ;
160   static const char* fgkEventSelectionEffID ;
161   
162   TVector2 fCorner1;
163   TVector2 fCorner2;
164   TString  fEnergyPath;
165   TString  fBackgroundPath;
166   TString  fEventSelectionEffPath;
167   Bool_t   fProcessPrimary;
168   Bool_t   fProcessHits; 
169   Trigger  fTrigger;
170   Energy   fEnergy;
171   MagField fMagField;
172   Species  fSpecies;
173   
174   ClassDef(AliFMDAnaParameters,0) // Manager of parameters
175 };
176
177 #endif
178 //____________________________________________________________________
179 //
180 // Local Variables:
181 //   mode: C++
182 // End:
183 //
184 // EOF
185 //
186