]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/CalibMacros/CalibTimeGain.C
Update variable description file
[u/mrichter/AliRoot.git] / TPC / CalibMacros / CalibTimeGain.C
1 // Macro to generate and update OCDB entries for a given run:
2 // this is a TObjArray which has at 0 the MIP-Spline and at 1 the Fermi-Plateau-Spline ...
3 // Responsible: marian.ivanov@cern.ch
4 // Responsible: A.Kalweit@gsi.de
5
6 /* How to use it locally:
7 //
8 // Load libraries
9 gSystem->Load("libANALYSIS");
10 gSystem->Load("libSTAT");
11 gSystem->Load("libTPCcalib");
12 gSystem->AddIncludePath("-I$ALICE_ROOT/TPC");
13 .L $ALICE_ROOT/TPC/CalibMacros/CalibTimeGain.C+
14
15 //Make calibration
16 CalibTimeGain("CalibObjectsTrain1.root",0,120000);
17 TBrowser b;
18 b.Add(gainArray);
19
20 */
21 #if !defined(__CINT__) || defined(__MAKECINT__)
22 #include "TObjArray.h"
23 #include "TGraphErrors.h"
24 #include "THnSparse.h"
25 #include "TCanvas.h"
26 #include "TFile.h"
27
28 #include "AliTPCcalibTimeGain.h"
29 #include "AliSplineFit.h"
30 #include "AliCDBMetaData.h"
31 #include "AliCDBId.h"
32 #include "AliCDBManager.h"
33 #include "AliCDBStorage.h"
34 #endif
35
36
37
38
39 TGraphErrors * graphMIP    = 0;       // graph time dependence of MIP
40 TGraphErrors * graphCosmic = 0;       // graph time dependence at Plateu
41 AliSplineFit * fitMIP = 0;            // fit of dependence - MIP
42 AliSplineFit * fitCosmic = 0;         // fit of dependence - Plateu
43 TObjArray    * gainArray = new TObjArray(4); // array to be stored in the OCDB
44 AliTPCcalibTimeGain * gainMIP =0;     // calibration component for MIP
45 AliTPCcalibTimeGain * gainCosmic =0;  // calibration component for cosmic
46
47
48 void UpdateOCDBGain(Int_t  startRunNumber, Int_t endRunNumber, const char* storagePath);
49 void ReadGainGlobal(Char_t* fileName="CalibObjectsTrain1.root");
50 void MakeQAPlot(Float_t  FPtoMIPratio);
51 Bool_t AnalyzeGain(Int_t startRunNumber, Int_t endRunNumber, Int_t minEntriesGaussFit = 500, Float_t FPtoMIPratio = 1.43); 
52
53
54
55 void CalibTimeGain(Char_t* fileName="CalibObjectsTrain1.root", Int_t startRunNumber=0, Int_t endRunNumber=AliCDBRunRange::Infinity(),  TString  ocdbStorage=""){
56   //
57   // Update OCDB gain
58   //
59   ReadGainGlobal(fileName);
60   AnalyzeGain(startRunNumber,endRunNumber, 1000,1.43);
61   MakeQAPlot(1.43);  
62   if (ocdbStorage.Length()==0) ocdbStorage+="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
63   UpdateOCDBGain( startRunNumber, endRunNumber, ocdbStorage.Data());
64 }
65
66
67
68
69 void ReadGainGlobal(Char_t* fileName){
70   //
71   // read calibration entries from file
72   // 
73   TFile fcalib(fileName);
74   TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib");
75   if (array){
76     gainMIP    = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGain");
77     gainCosmic = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGainCosmic");
78   }else{
79     gainMIP    = ( AliTPCcalibTimeGain *)fcalib.Get("calibTimeGain");
80     gainCosmic = ( AliTPCcalibTimeGain *)fcalib.Get("calibTimeGainCosmic");
81   }
82   TH1 * hisT=0;
83   Int_t firstBinA =0, lastBinA=0;
84
85   if (gainCosmic){ 
86     hisT= gainCosmic->GetHistGainTime()->Projection(1);
87     firstBinA = hisT->FindFirstBinAbove(2);
88     lastBinA  = hisT->FindLastBinAbove(2);    
89     gainCosmic->GetHistGainTime()->GetAxis(1)->SetRange(firstBinA,lastBinA);
90     delete hisT;
91   }
92
93   if (gainMIP){ 
94     hisT= gainCosmic->GetHistGainTime()->Projection(1);
95     firstBinA = hisT->FindFirstBinAbove(2);
96     lastBinA  = hisT->FindLastBinAbove(2);    
97     gainMIP->GetHistGainTime()->GetAxis(1)->SetRange(firstBinA,lastBinA);
98     delete hisT;
99   }
100
101 }
102
103
104
105 Bool_t AnalyzeGain(Int_t startRunNumber, Int_t endRunNumber, Int_t minEntriesGaussFit,  Float_t FPtoMIPratio){
106   //
107   //
108   //
109   gainMIP->GetHistGainTime()->GetAxis(5)->SetRangeUser(startRunNumber, endRunNumber);
110   // 1.) try to create MIP spline
111   gainMIP->GetHistGainTime()->GetAxis(2)->SetRangeUser(1.51,2.49); // only beam data
112   gainMIP->GetHistGainTime()->GetAxis(4)->SetRangeUser(0.39,0.51); // only MIP pions
113   //
114   graphMIP = AliTPCcalibBase::FitSlices(gainMIP->GetHistGainTime(),0,1,minEntriesGaussFit,10,0.1,0.7);
115   if (graphMIP->GetN()==0) graphMIP = 0x0;
116   if (graphMIP) fitMIP = AliTPCcalibTimeGain::MakeSplineFit(graphMIP);
117   if (graphMIP) graphMIP->SetName("TGRAPHERRORS_MEAN_GAIN_BEAM_ALL");// set proper names according to naming convention
118   gainArray->AddAt(fitMIP,0);
119   
120
121   // 2.) try to create Cosmic spline
122   gainCosmic->GetHistGainTime()->GetAxis(2)->SetRangeUser(0.51,1.49); // only cosmics
123   gainCosmic->GetHistGainTime()->GetAxis(4)->SetRangeUser(20,100);    // only Fermi-Plateau muons
124   //
125   graphCosmic = AliTPCcalibBase::FitSlices(gainCosmic->GetHistGainTime(),0,1,minEntriesGaussFit,10);
126   if (graphCosmic->GetN()==0) graphCosmic = 0x0;
127   //
128   if (graphCosmic) {
129     for(Int_t i=0; i < graphCosmic->GetN(); i++) {
130       graphCosmic->GetY()[i] /= FPtoMIPratio;
131       graphCosmic->GetEY()[i] /= FPtoMIPratio;
132     }
133   }
134   //
135   if (graphCosmic) fitCosmic = AliTPCcalibTimeGain::MakeSplineFit(graphCosmic);
136   if (graphCosmic) graphCosmic->SetName("TGRAPHERRORS_MEAN_GAIN_COSMIC_ALL"); // set proper names according to naming convention
137   gainArray->AddAt(fitCosmic,1);
138   // with naming convention and backward compatibility
139   gainArray->AddAt(graphMIP,2);
140   gainArray->AddAt(graphCosmic,3);
141   cout << "graphCosmic: " << graphCosmic << " graphMIP " << graphMIP << endl;
142   return kTRUE;
143
144 }
145
146
147
148 void UpdateOCDBGain(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
149   //
150   // Update OCDB entry
151   //
152   AliCDBMetaData *metaData= new AliCDBMetaData();
153   metaData->SetObjectClassName("TObjArray");
154   metaData->SetResponsible("Alexander Kalweit");
155   metaData->SetBeamPeriod(1);
156   metaData->SetAliRootVersion("05-24-00"); //root version
157   metaData->SetComment("Calibration of the time dependence of the gain due to pressure and temperature changes.");
158   AliCDBId id1("TPC/Calib/TimeGain", startRunNumber, endRunNumber);
159   AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
160   gStorage->Put(gainArray, id1, metaData);    
161 }
162
163 void MakeQAPlot(Float_t  FPtoMIPratio) {
164   //
165   // Make QA plot to visualize results
166   //
167   //
168   //
169   if (graphCosmic) {
170     TCanvas * canvasCosmic = new TCanvas("gain Cosmic", "time dependent gain QA histogram cosmic");
171     canvasCosmic->cd();
172     TH2D * gainHistoCosmic = gainCosmic->GetHistGainTime()->Projection(0,1);
173     gainHistoCosmic->SetDirectory(0);
174     gainHistoCosmic->SetName("GainHistoCosmic");
175     gainHistoCosmic->GetXaxis()->SetTimeDisplay(kTRUE);
176     gainHistoCosmic->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
177     gainHistoCosmic->Draw("colz");
178     graphCosmic->SetMarkerStyle(25);
179     graphCosmic->Draw("lp");
180     graphCosmic->SetMarkerStyle(25);
181     TGraph * grfitCosmic = fitCosmic->MakeGraph(graphCosmic->GetX()[0],graphCosmic->GetX()[graphCosmic->GetN()-1],50000,0);
182     if (grfitCosmic) {
183       for(Int_t i=0; i < grfitCosmic->GetN(); i++) {
184         grfitCosmic->GetY()[i] *= FPtoMIPratio; 
185       }
186       for(Int_t i=0; i < graphCosmic->GetN(); i++) {
187         graphCosmic->GetY()[i] *= FPtoMIPratio; 
188       }
189     }
190     graphCosmic->Draw("lp");
191     grfitCosmic->SetLineColor(2);
192     grfitCosmic->Draw("lu");
193     gainArray->AddLast(gainHistoCosmic);
194     gainArray->AddLast(canvasCosmic->Clone());
195     delete canvasCosmic;    
196   }
197   if (fitMIP) {
198     TCanvas * canvasMIP = new TCanvas("gain MIP", "time dependent gain QA histogram MIP");
199     canvasMIP->cd();
200     TH2D * gainHistoMIP    = gainMIP->GetHistGainTime()->Projection(0,1);
201     gainHistoMIP->SetName("GainHistoCosmic");
202     gainHistoMIP->SetDirectory(0);
203     gainHistoMIP->GetXaxis()->SetTimeDisplay(kTRUE);
204     gainHistoMIP->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
205     gainHistoMIP->Draw("colz");
206     graphMIP->SetMarkerStyle(25);
207     graphMIP->Draw("lp");
208     TGraph * grfitMIP = fitMIP->MakeGraph(graphMIP->GetX()[0],graphMIP->GetX()[graphMIP->GetN()-1],50000,0);
209     grfitMIP->SetLineColor(2);
210     grfitMIP->Draw("lu");    
211     gainArray->AddLast(gainHistoMIP);
212     gainArray->AddLast(canvasMIP->Clone());
213     delete canvasMIP;    
214   }  
215 }
216
217