]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/CalibMacros/CalibTimeGain.C
Macro to create TimeGain OCDB entry from the calibration data
[u/mrichter/AliRoot.git] / TPC / CalibMacros / CalibTimeGain.C
1
2 // small macro to generate and update OCDB entries for a given run:
3 //
4 // this is a TObjArray which has at 0 the MIP-Spline and at 1 the Fermi-Plateau-Spline ...
5
6 /* How to use it...
7
8 gSystem->Load("libANALYSIS");
9 gSystem->Load("libSTAT");
10 gSystem->Load("libTPCcalib");
11 gSystem->AddIncludePath("-I$ALICE_ROOT/TPC");
12
13 .L $ALICE_ROOT/TPC/CalibMacros/CalibTimeGain.C++
14
15 TFile fcalib("CalibObjectsTrain1.root");
16 //TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib");
17 //AliTPCcalibTimeGain * gain = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGain");
18 AliTPCcalibTimeGain * gain = ( AliTPCcalibTimeGain *)fcalib.Get("calibTimeGain")
19
20 ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
21 CalibTimeGain(gain, ocdbStorage.Data(),70000,100000)
22
23 */
24
25 #include "TObjArray.h"
26 #include "TGraphErrors.h"
27 #include "THnSparse.h"
28 #include "TCanvas.h"
29
30 #include "AliTPCcalibTimeGain.h"
31 #include "AliSplineFit.h"
32 #include "AliCDBMetaData.h"
33 #include "AliCDBId.h"
34 #include "AliCDBManager.h"
35 #include "AliCDBStorage.h"
36
37
38 Bool_t CalibTimeGain(AliTPCcalibTimeGain * gain, Char_t * storagePath, Int_t startRunNumber, Int_t endRunNumber, Bool_t updateOCDB = kFALSE, Int_t minEntriesGaussFit = 500, Bool_t makeQAplot=kTRUE){
39
40   TObjArray * splineArray = new TObjArray(4);
41   gain->GetHistGainTime()->GetAxis(5)->SetRangeUser(startRunNumber, endRunNumber);
42
43   // 1.) try to create MIP spline
44   gain->GetHistGainTime()->GetAxis(2)->SetRangeUser(1.51,2.49); // only beam data
45   gain->GetHistGainTime()->GetAxis(4)->SetRangeUser(0.39,0.51); // only MIP pions
46   //
47   TGraphErrors * graphMIP = AliTPCcalibBase::FitSlices(gain->GetHistGainTime(),0,1,minEntriesGaussFit,10);
48   if (graphMIP->GetN()==0) graphMIP = 0x0;
49   AliSplineFit * fitMIP = 0;
50   if (graphMIP) fitMIP = AliTPCcalibTimeGain::MakeSplineFit(graphMIP);
51   if (graphMIP) graphMIP->SetName("TGRAPHERRORS_MEAN_GAIN_BEAM_ALL");// set proper names according to naming convention
52   (*splineArray)[0] = fitMIP;
53   
54   // 2.) try to create Cosmic spline
55   gain->GetHistGainTime()->GetAxis(2)->SetRangeUser(0.51,1.49); // only cosmics
56   gain->GetHistGainTime()->GetAxis(4)->SetRangeUser(20,100);    // only Fermi-Plateau muons
57   //
58   TGraphErrors * graphCosmic = AliTPCcalibBase::FitSlices(gain->GetHistGainTime(),0,1,minEntriesGaussFit,10);
59   if (graphCosmic->GetN()==0) graphCosmic = 0x0;
60   AliSplineFit * fitCosmic = 0;
61   if (graphCosmic) fitCosmic = AliTPCcalibTimeGain::MakeSplineFit(graphCosmic);
62   if (graphCosmic) graphCosmic->SetName("TGRAPHERRORS_MEAN_GAIN_COSMIC_ALL"); // set proper names according to naming convention
63   (*splineArray)[1] = fitCosmic;
64
65   // with naming convention and backward compatibility
66   (*splineArray)[2] = graphMIP;
67   (*splineArray)[3] = graphCosmic;
68   cout << "graphCosmic: " << graphCosmic << " graphMIP " << graphMIP << endl;
69   //
70   // store OCDB entry
71   //
72   if (!fitCosmic && !fitMIP) return kFALSE;
73   if (!graphCosmic && !graphMIP) return kFALSE;
74   //
75   AliCDBMetaData *metaData= new AliCDBMetaData();
76   metaData->SetObjectClassName("TObjArray");
77   metaData->SetResponsible("Alexander Kalweit");
78   metaData->SetBeamPeriod(1);
79   metaData->SetAliRootVersion("05-24-00"); //root version
80   metaData->SetComment("Calibration of the time dependence of the gain due to pressure and temperature changes.");
81   AliCDBId id1("TPC/Calib/TimeGain", startRunNumber, endRunNumber);
82   AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
83   if (updateOCDB) gStorage->Put(splineArray, id1, metaData);
84
85   if (makeQAplot) {
86     TCanvas * canvQA = new TCanvas("canvQA", "time dependent gain QA histogram");
87     canvQA->cd();
88     TGraphErrors * gr = gain->GetGraphGainVsTime(0,minEntriesGaussFit);
89     TH2D * GainTime = gain->GetHistGainTime()->Projection(0,1);
90     GainTime->GetXaxis()->SetTimeDisplay(kTRUE);
91     GainTime->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
92     GainTime->Draw("colz");
93     //
94     if (graphCosmic) {
95       graphCosmic->SetMarkerStyle(25);
96       graphCosmic->Draw("lp");
97       TGraph * grfitCosmic = fitCosmic->MakeGraph(graphCosmic->GetX()[0],graphCosmic->GetX()[graphCosmic->GetN()-1],50000,0);
98       grfitCosmic->SetLineColor(2);
99       grfitCosmic->Draw("lu");
100     }
101     if (fitMIP) {
102       graphMIP->SetMarkerStyle(25);
103       graphMIP->Draw("lp");
104       TGraph * grfitMIP = fitMIP->MakeGraph(graphMIP->GetX()[0],graphMIP->GetX()[graphMIP->GetN()-1],50000,0);
105       grfitMIP->SetLineColor(2);
106       grfitMIP->Draw("lu");
107
108     }
109  
110   }
111
112   return kTRUE;
113
114 }
115
116
117 /* Make dummy entry - pseudo code 
118
119  TObjArray * splineArray = new TObjArray(2);
120
121  AliSplineFit * fitMIP = 0;
122  AliSplineFit * fitCosmic = 0;
123  (*splineArray)[0] = fitMIP;
124  (*splineArray)[1] = fitCosmic;
125  
126   AliCDBMetaData *metaData= new AliCDBMetaData();
127   metaData->SetObjectClassName("TObjArray");
128   metaData->SetResponsible("Alexander Kalweit");
129   metaData->SetBeamPeriod(1);
130   metaData->SetAliRootVersion("05-23-02"); //root version
131   metaData->SetComment("Calibration of the time dependence of the gain due to pressure and temperature changes");
132   AliCDBId id1("TPC/Calib/TimeGain", 0, AliCDBRunRange::Infinity());
133   AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage("local:///d/alice05/akalweit/projects/TimeGainCalibration");
134   gStorage->Put(splineArray, id1, metaData);
135
136
137 */
138
139
140
141