]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/macros/CalibrationDB/AliEMCALSetCDB.C
0a6358050f0f92212116d06dcb5698d829465024
[u/mrichter/AliRoot.git] / EMCAL / macros / CalibrationDB / AliEMCALSetCDB.C
1
2 // Script to create calibration parameters and store them into CDB
3 // Two sets of calibration parameters can be created:
4 // 1) equal parameters
5 // 2) randomly distributed parameters for decalibrated detector silumations
6 // Modified from PHOS script for EMCAL by Gustavo Conesa
7
8 #if !defined(__CINT__)
9 #include "TControlBar.h"
10 #include "TString.h"
11 #include "TRandom.h"
12 #include "TH2F.h"
13 #include "TCanvas.h"
14
15 #include "AliRun.h"
16 #include "AliEMCALCalibData.h"
17 #include "AliCDBMetaData.h"
18 #include "AliCDBId.h"
19 #include "AliCDBEntry.h"
20 #include "AliCDBManager.h"
21 #include "AliCDBStorage.h"
22 #endif
23
24
25 void AliEMCALSetCDB()
26 {
27   TControlBar *menu = new TControlBar("vertical","EMCAL CDB");
28   menu->AddButton("Help to run EMCAL CDB","Help()",
29                   "Explains how to use EMCAL CDS menus");
30   menu->AddButton("Equal CC","SetCC(0)",
31                   "Set equal calibration coefficients");
32   menu->AddButton("Decalibrate","SetCC(1)",
33                   "Set random decalibration calibration coefficients");
34   menu->AddButton("Read equal CC","GetCC(0)",
35                   "Read initial equal calibration coefficients");
36   menu->AddButton("Read random CC","GetCC(1)",
37                   "Read random decalibration calibration coefficients");
38   menu->Show();
39 }
40
41 //------------------------------------------------------------------------
42 void Help()
43 {
44   char *string =
45     "\nSet calibration parameters and write them into ALICE CDB. Press button \"Equal CC\" to create equal pedestals and gain factors. Press button \"Decalibrate\" to create random pedestals and gain factors to imitate decalibrated detector\n";
46   printf(string);
47 }
48
49 //------------------------------------------------------------------------
50 void SetCC(Int_t flag=0)
51 {
52   // Writing calibration coefficients into the Calibration DB
53   // Arguments:
54   //   flag=0: all calibration coefficients are equal
55   //   flag=1: all calibration coefficients random (decalibration)
56   // Author: Boris Polishchuk (Boris.Polichtchouk at cern.ch)
57
58   TString DBFolder;
59   Int_t firstRun   =  0;
60   Int_t lastRun    = 10;
61   Int_t beamPeriod =  1;
62   char* objFormat  = "";
63
64   if      (flag == 0) {
65     DBFolder  ="local://InitCalibDB";
66     firstRun  =  0;
67     lastRun   =  0;
68     objFormat = "EMCAL initial gain factors and pedestals";
69   }
70   else if (flag == 1) {
71     DBFolder  ="local://DeCalibDB";
72     firstRun  =  0;
73     lastRun   = 12;
74     objFormat = "EMCAL random pedestals and ADC gain factors (12x48x24)";
75   }
76   
77   AliEMCALCalibData *calibda=new AliEMCALCalibData("EMCAL");
78   
79   Float_t fADCpedestal = 0.0;
80   Float_t fADCchannel  = 0.0153;  // 250 GeV / (16*1024)
81
82   TRandom rn;
83   Int_t nSMod  = 12;
84   Int_t nCol   = 48;
85   Int_t nRow   = 24;
86   Int_t nRow2  = 12; //Modules 11 and 12 are half modules
87
88   for(Int_t supermodule=0; supermodule < nSMod; supermodule++) {
89     for(Int_t column=0; column< nCol; column++) {
90       if(supermodule >= 10)
91         nRow = nRow2;
92       for(Int_t row=0; row< nRow; row++) {
93         if (flag == 1) {
94           // Decalibration:
95           // Spread calibration coefficients uniformly around central value 
96           fADCchannel  = rn.Uniform(0.00140,0.00160);
97           fADCpedestal = 0;
98         }
99         calibda->SetADCchannel (supermodule,column,row,fADCchannel);
100         calibda->SetADCpedestal(supermodule,column,row,fADCpedestal);
101         cout<<"Set SM: "<<supermodule<<" col "<<column<<" row "<<row
102             <<" chan "<<fADCchannel<<" ped fADCpedestal"<<endl;
103       }
104     }
105   }
106
107   //Store calibration data into database
108   
109   AliCDBMetaData md;
110   md.SetComment(objFormat);
111   md.SetBeamPeriod(beamPeriod);
112   md.SetResponsible("Boris Polichtchouk");
113   
114   AliCDBId id("EMCAL/Calib/Data",firstRun,lastRun);
115
116   AliCDBManager* man = AliCDBManager::Instance();  
117   AliCDBStorage* loc = man->GetStorage(DBFolder.Data());
118   loc->Put(calibda, id, &md);
119
120 }
121
122 //------------------------------------------------------------------------
123 void GetCC(Int_t flag=0)
124 {
125   // Read calibration coefficients into the Calibration DB
126   // Arguments:
127   //   flag=0: all calibration coefficients are equal
128   //   flag=1: all calibration coefficients random (decalibration)
129   // Author: Yuri.Kharlov at cern.ch
130
131   TString DBFolder;
132
133   if      (flag == 0) {
134     DBFolder  ="local://InitCalibDB";
135   }
136   else if (flag == 1) {
137     DBFolder  ="local://DeCalibDB";
138   }
139
140   AliEMCALCalibData* clb  = (AliEMCALCalibData*)
141     (AliCDBManager::Instance()
142      ->GetStorage(DBFolder.Data())
143      ->Get("EMCAL/Calib/Data",
144            gAlice->GetRunNumber())->GetObject());
145
146   static const Int_t nSMod = 12;
147   static const Int_t nCol  = 48;
148   Int_t nRow  = 24;
149   Int_t nRow2 = 12; //Modules 11 and 12 are half modules
150
151   TH2F *hPed[nSMod], *hGain[nSMod];
152   TCanvas *cPed   = new TCanvas("cPed" ,"Pedestals Mod 1-6"   , 10,10,400,800);
153   TCanvas *cGain  = new TCanvas("cGain","Gain factors Mod 1-6", 410,10,400,800);
154   TCanvas *cPed2  = new TCanvas("cPed2","Pedestals SMod 7-12", 10,10,400,800);
155   TCanvas *cGain2 = new TCanvas("cGain2","Gain factors SMod 7-12", 410,10,400,800);
156   cPed   ->Divide(2,3);
157   cGain  ->Divide(2,3);
158   cPed2  ->Divide(2,3);
159   cGain2 ->Divide(2,3);
160   for (Int_t supermodule=0; supermodule<nSMod; supermodule++) {
161
162     if(supermodule >= 10)
163       nRow = nRow2;
164
165     TString namePed="hPed";
166     namePed+=supermodule;
167     TString titlePed="Pedestals in supermodule ";
168     titlePed+=supermodule;
169     hPed[supermodule] = new TH2F(namePed.Data(),titlePed.Data(),
170                             nCol,1.,1.*nCol,nRow,1.,1.*nRow);
171
172     TString nameGain="hGain";
173     nameGain+=supermodule;
174     TString titleGain="Gain factors in supermodule ";
175     titleGain+=supermodule;
176     hGain[supermodule] = new TH2F(nameGain.Data(),titleGain.Data(),
177                                     nCol,1.,1.*nCol,nRow,1.,1.*nRow);
178     for (Int_t column=0; column<nCol; column++) {
179       for (Int_t row=0; row<nRow; row++) {
180         Float_t ped  = clb->GetADCpedestal(supermodule,column,row);
181         Float_t gain = clb->GetADCchannel (supermodule,column,row);
182         //cout<<"Get SM: "<<supermodule<<" col "<<column<<" row "<<row
183         //  <<" chan "<<gain<<endl;
184         hPed[supermodule] ->SetBinContent(column,row,ped);
185         hGain[supermodule]->SetBinContent(column,row,gain);
186       }
187     }
188     if(supermodule < 7){
189       cPed ->cd(supermodule);
190       hPed[supermodule] ->Draw("lego2");
191       cGain->cd(supermodule);
192       hGain[supermodule]->Draw("lego2");
193     }
194     else{
195       cPed2 ->cd(supermodule-6);
196       hPed[supermodule] ->Draw("lego2");
197       cGain2->cd(supermodule-6);
198       hGain[supermodule]->Draw("lego2");
199     }
200
201   }
202   cPed   ->Print("pedestals_SM_1_6.eps");
203   cGain  ->Print("gains_SM_1-6.eps");
204   cPed2  ->Print("pedestals_SM_7-12.eps");
205   cGain2 ->Print("gains_SM_7-12.eps");
206 }