]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/tests/TestMakeELossFits.C
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / tests / TestMakeELossFits.C
1 #ifndef __CINT__
2 #include "AliForwardCorrectionManager.h"
3 #include "AliFMDCorrELossFit.h"
4 #include <TAxis.h>
5 #include <TFile.h>
6 #include <TList.h>
7 #include <TError.h>
8 #include <TH1.h>
9 #include <TF1.h>
10 #include <TClonesArray.h>
11 #else
12 class TAxis;
13 class AliFMDCirrELossFit;
14 class TH1;
15
16 #endif
17
18 /**
19  * Class to make correction object and save to file 
20  *
21  * Run like                                     
22  *
23  * @verbatim 
24  * gROOT->Macro("$ALICE_ROOT/PWGLF/FORWARD/analysis2/scripts/LoadLibs.C");
25  * gROOT->LoadMacro("$ALICE_ROOT/PWGLF/FORWARD/analysis2/scripts/Compile.C");
26  * Compile("$ALICE_ROOT/PWGLF/FORWARD/analysis2/MakeELossFit.C"); 
27  * MakeELossFit mef(sys, cms, field, mc, "AnalysisResults.root"); 
28  * mef.Run();
29  * @endverbatim 
30  * where @e sys, the collision system, is one of 
31  * - 1: pp
32  * - 2: PbPb
33  * 
34  * @e cms is the center of mass energy per nuclean in GeV (e.g., 2760
35  * for first PbPb data), @e field is (signed) L3 magnetic in kG, and 
36  * @e mc is wether this correction applies to MC data or not. 
37  * 
38  * The class generates a file like 
39  * @verbatim
40  * elossfits<sys>_<cms>GeV_<field>kG_<realmc>.root
41  * @endverbatim 
42  * in the working directory. This file can be moved to 
43  * @verbatim
44  * $(ALICE_ROOT)/PWGLF/FORWARD/corrections/ELossFits
45  * @endverbatim 
46  * and stored in SVN for later use. 
47  *
48  * @depcrecated 
49  * The class AliFMDELossFitter automatically generates the
50  * AliFMDCorrELossFit object.
51  *
52  * @ingroup pwglf_forward_scripts_tests
53  */
54
55 class MakeELossFit 
56 {
57 protected: 
58 public:
59   TList* fFitter;
60   TAxis* fAxis;
61   TClonesArray fFits;
62   UShort_t     fSys;
63   UShort_t     fCMS;
64   Short_t      fField;
65   Bool_t       fMC;
66
67   //__________________________________________________________________
68   MakeELossFit(UShort_t    sys, 
69                UShort_t    cms, 
70                Short_t     field, 
71                Bool_t      mc, 
72                const char* filename="forward_eloss.root") 
73     : fFitter(0), 
74       fAxis(0),
75       fFits("AliFMDCorrELossFit::ELossFit"),
76       fSys(sys), 
77       fCMS(cms), 
78       fField(field), 
79       fMC(mc)
80   {
81     TFile* file = TFile::Open(filename, "READ");
82     if (!file) { 
83       Error("MakeELossFit", "Failed to open file %s", filename);
84       return;
85     }
86     TList* forward = static_cast<TList*>(file->Get("Forward"));
87     // static_cast<TList*>(file->Get("PWGLFforwardDnDeta/Forward"));
88     if (!forward) { 
89       Error("MakeELossFit", "Couldn't get forward list from %s", filename);
90       return;
91     }
92
93     fFitter = static_cast<TList*>(forward->FindObject("fmdEnergyFitter"));
94     if (!fFitter) { 
95       Error("MakeELossFit", "Couldn't get fitter folder");
96       return;
97     }
98
99     fAxis = static_cast<TAxis*>(fFitter->FindObject("etaAxis"));
100     if (!fAxis) { 
101       Error("MakeELossFit", "Couldn't get eta axis");
102       return;
103     }
104     file->Close();
105
106 #if 0
107     AliForwardCorrectionManager& mgr = AliForwardCorrectionManager::Instance();
108     mgr.Init(sys, cms, field, mc, 0, true);
109 #endif
110   }
111   //__________________________________________________________________
112   TList* GetRing(UShort_t d, Char_t r) const
113   {
114     TList* rL = static_cast<TList*>(fFitter->FindObject(Form("FMD%d%c",d,r)));
115     if (!rL) { 
116       Warning("DrawFits", "List FMD%d%c not found", d, r);
117       return 0;
118     }
119     return rL;
120   }
121   //__________________________________________________________________
122   TList* GetEDists(UShort_t d, Char_t r) const
123   {
124     TList* rList = GetRing(d, r);
125     if (!rList) 
126       return 0;
127     // rList->ls();
128
129     TList* edists = static_cast<TList*>(rList->FindObject("EDists"));
130     if (!edists) { 
131       Error("PrintFits", "Couldn't get EDists list for FMD%d%c", d,r);
132       return 0; 
133     }
134     return edists;
135   }
136   //__________________________________________________________________
137   TH1* GetEDist(UShort_t d, Char_t r, UShort_t etabin)
138   {
139     TList* eList = GetEDists(d, r);
140     if (!eList) { 
141       Warning("GetEDist", "No list for FMD%d%c", d, r);
142       return 0;
143     }
144     // eList->ls();
145
146     TString cmp(Form("FMD%d%c_etabin%03d", d, r, etabin));
147
148     TIter next(eList);
149     TObject* o = 0;
150     while ((o = next())) { 
151       if (!cmp.CompareTo(o->GetName())) {
152         return static_cast<TH1*>(o);
153       }
154     }
155     return 0;
156   }
157 #ifndef __CINT__
158   //__________________________________________________________________
159   AliFMDCorrELossFit::ELossFit* FindBestFit(TH1* dist) 
160   {
161     TList* funcs = dist->GetListOfFunctions();
162     TIter  next(funcs);
163     TF1*   func = 0;
164     fFits.Clear();
165     Int_t  i = 0;
166     Info("FindBestFit", "%s", dist->GetName());
167     while ((func = static_cast<TF1*>(next()))) { 
168       AliFMDCorrELossFit::ELossFit* fit = 
169         new(fFits[i++]) AliFMDCorrELossFit::ELossFit(0,*func);
170       fit->CalculateQuality(10, .20, 1e-7);
171     }
172     fFits.Sort(false);
173     // fFits.Print();
174     return static_cast<AliFMDCorrELossFit::ELossFit*>(fFits.At(0));
175   }
176 #endif
177   //__________________________________________________________________
178   Bool_t Run()
179   {
180     if (!fFitter || !fAxis) { 
181       Error("Run", "Missing objects");
182       return kFALSE;
183     }
184     AliFMDCorrELossFit* obj = new AliFMDCorrELossFit;
185     obj->SetEtaAxis(*fAxis);
186
187     for (UShort_t d = 1; d <= 3; d++) { 
188       Info("Run", "detector is FMD%d", d);
189       UShort_t nQ = (d == 1 ? 1 : 2);
190       for (UShort_t q = 0; q < nQ; q++) { 
191         Char_t r = (q == 0 ? 'I' : 'O');
192         Int_t nBin = fAxis->GetNbins();
193         Info("Run", " ring is FMD%d%c - %d bins", d, r, nBin);
194         
195         Bool_t oneSeen = kFALSE;
196         for (UShort_t b = 1; b <= nBin; b++) { 
197           TH1* h = GetEDist(d, r, b);
198           if (oneSeen && !h) break;
199           if (!h)            continue;
200           if (!oneSeen)      oneSeen = true;
201
202           AliFMDCorrELossFit::ELossFit* best = FindBestFit(h);
203           best->Print("");
204           best->fDet  = d; 
205           best->fRing = r;
206           best->fBin  = b;
207           // Double_t eta = fAxis->GetBinCenter(b);
208           Info("Run", "    bin=%3d->%6.4f", b, eta);
209           obj->SetFit(d, r, b, new AliFMDCorrELossFit::ELossFit(*best));
210         }
211       }
212     }
213     AliForwardCorrectionManager& mgr = AliForwardCorrectionManager::Instance();
214     TString fname(mgr.GetFileName(AliForwardCorrectionManager::kELossFits,
215                                   fSys, fCMS, fField, fMC));
216     TFile* output = TFile::Open(fname.Data(), "RECREATE");
217     if (!output) { 
218       Warning("Run", "Failed to open output file %s", fname.Data());
219       return kFALSE;
220     }
221     obj->Write(mgr.GetObjectName(AliForwardCorrectionManager::kELossFits));
222     output->Write();
223     output->Close();
224     Info("Run", "File %s created.  It should be copied to %s and stored in SVN",
225          fname.Data(),mgr.GetFileDir(AliForwardCorrectionManager::kELossFits));
226     
227     return kTRUE;
228   }
229 };
230 //
231 // EOF
232 //