]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/CalibMacros/CPass0/makeOCDB.C
Calibrating TOF only if TPC is in the data taking (see https://savannah.cern.ch/bugs...
[u/mrichter/AliRoot.git] / PWGPP / CalibMacros / CPass0 / makeOCDB.C
1 /*
2   macro to extract the OCDB entries
3
4   input: CalibObjects.root
5   ouput: TimeGain and TimeVdrift calibration objects for TPC and TRD
6
7   Example:
8   .L $ALICE_ROOT/PWGPP/CalibMacros/CPass0/makeOCDB.C
9   makeOCDB("105160");
10
11 */
12
13 void printCalibStat(Int_t run, const char * fname,  TTreeSRedirector * pcstream);
14
15
16 void makeOCDB(Int_t runNumber, TString  targetOCDBstorage="", TString sourceOCDBstorage="raw://")
17 {
18
19   //
20   // extract OCDB entries for detectors participating in the calibration for the current run
21   //
22
23   gROOT->Macro("$ALICE_ROOT/PWGPP/CalibMacros/CPass0/LoadLibraries.C");
24   gROOT->LoadMacro("$ALICE_ROOT/PWGPP/CalibMacros/CPass0/ConfigCalibTrain.C");
25
26   // switch off log info
27   AliLog::SetClassDebugLevel("AliESDEvent",0);
28
29   // config GRP
30   printf("runNumber from runCalibTrain = %d\n",runNumber);
31   ConfigCalibTrain(runNumber, sourceOCDBstorage.Data());
32
33   // check the presence of the detectors
34   AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
35   AliGRPObject* grpData = dynamic_cast<AliGRPObject*>(entry->GetObject()); 
36   if (!grpData) {printf("Failed to get GRP data for run",runNumber); return;}
37   Int_t activeDetectors = grpData->GetDetectorMask(); 
38   TString detStr = AliDAQ::ListOfTriggeredDetectors(activeDetectors);
39   printf("Detectors in the data:\n%s\n",detStr.Data());
40   printf("Monalisa: Detectors in the data:\t%s\n",detStr.Data());
41
42   // Steering Tasks - set output storage
43   // DefaultStorage set already before - in ConfigCalibTrain.C
44   //targetOCDBstorage+="?se=ALICE::CERN::SE";
45
46   // for tests, the target OCDB has to be different from raw://, but since the calibration procedure may need to read the OCDB,
47   // one cannot set the target OCDB as new default storage, since it could it that it does not contain all the necessary entries. 
48   // This is true only if the target OCDB storage is different from the one used as default source and to configure the train
49   if (targetOCDBstorage != sourceOCDBstorage) AliCDBManager::Instance()->SetSpecificStorage("*/*/*", targetOCDBstorage.Data());
50
51   if (gSystem->AccessPathName("TPC", kFileExists)==0) {  
52     AliCDBManager::Instance()->SetSpecificStorage("TPC/Calib/Correction","local://");
53   }
54   // set OCDB storage
55   if (targetOCDBstorage.Length()==0) targetOCDBstorage+="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
56
57   // TPC part
58   AliTPCPreprocessorOffline *procesTPC = 0;
59   if  (detStr.Contains("TPC")){ 
60     Printf("\n******* Calibrating TPC *******");
61     procesTPC = new AliTPCPreprocessorOffline; 
62     // switch on parameter validation
63     procesTPC->SetTimeGainRange(0.5,4.0);
64     procesTPC->SwitchOnValidation();
65
66     // Make timegain calibration
67     //proces.CalibTimeGain("CalibObjects.root", runNumber,AliCDBRunRange::Infinity(),targetOCDBstorage);
68     procesTPC->CalibTimeGain("CalibObjects.root", runNumber,runNumber,targetOCDBstorage);
69     
70     // Make vdrift calibration
71     //proces.CalibTimeVdrift("CalibObjects.root",runNumber,AliCDBRunRange::Infinity(),targetOCDBstorage);
72     procesTPC->CalibTimeVdrift("CalibObjects.root",runNumber,runNumber,targetOCDBstorage);
73   }
74
75   // TOF part
76   AliTOFAnalysisTaskCalibPass0 *procesTOF=0;
77   if ( detStr.Contains("TOF") && detStr.Contains("TPC")){    
78     procesTOF = new AliTOFAnalysisTaskCalibPass0;
79     Printf("\n******* Calibrating TOF *******");
80     procesTOF->ProcessOutput("CalibObjects.root", targetOCDBstorage);
81   }
82
83   // T0 part
84   AliT0PreprocessorOffline *procesT0= 0;
85   if ( detStr.Contains("T0")) {
86     Printf("\n******* Calibrating T0 *******");
87     // Make  calibration of channels offset
88     procesT0= new  AliT0PreprocessorOffline;
89     procesT0->Process("CalibObjects.root",runNumber, runNumber, targetOCDBstorage);
90   }
91
92   //TRD part
93   AliTRDPreprocessorOffline *procesTRD = 0;
94   if ( detStr.Contains("TRD")){
95     Printf("\n******* Calibrating TRD *******");
96     procesTRD = new  AliTRDPreprocessorOffline;
97     procesTRD->SetLinearFitForVdrift(kTRUE);
98     procesTRD->SetMinStatsVdriftT0PH(600*10);
99     procesTRD->SetMinStatsVdriftLinear(50);
100     procesTRD->SetMinStatsGain(600);
101     procesTRD->SetLimitValidateNoData(60);
102     procesTRD->SetLimitValidateBadCalib(60);
103     procesTRD->SetAlternativeDriftVelocityFit(kTRUE);
104     procesTRD->Init("CalibObjects.root");
105     Int_t versionVdriftUsed = procesTRD->GetVersionVdriftUsed();
106     Int_t subversionVdriftUsed = procesTRD->GetSubVersionVdriftUsed();
107     Int_t versionGainUsed = procesTRD->GetVersionGainUsed();
108     Int_t subversionGainUsed = procesTRD->GetSubVersionGainUsed();
109     Int_t versionExBUsed = procesTRD->GetVersionExBUsed();
110     Int_t subversionExBUsed = procesTRD->GetSubVersionExBUsed();
111     printf("version and subversion vdrift %d and %d\n",versionVdriftUsed,subversionVdriftUsed);
112     printf("version and subversion gain %d and %d\n",versionGainUsed,subversionGainUsed);
113     printf("version and subversion exb %d and %d\n",versionExBUsed,subversionExBUsed);
114     procesTRD->Process("CalibObjects.root",runNumber,runNumber,targetOCDBstorage);
115   }
116   
117   //Mean Vertex
118   AliMeanVertexPreprocessorOffline * procesMeanVtx=0;
119   if ( detStr.Contains("ITSSPD")) {
120     Printf("\n******* Calibrating MeanVertex *******");
121     procesMeanVtx =  new AliMeanVertexPreprocessorOffline;
122     procesMeanVtx->ProcessOutput("CalibObjects.root", targetOCDBstorage, runNumber);
123   }
124
125   //
126   // Print calibration status into the stdout
127   //
128   Int_t trdStatus = (procesTRD) ?  procesTRD->GetStatus():0;
129   Int_t tofStatus = (procesTOF) ?  procesTOF->GetStatus():0;
130   Int_t t0Status  = (procesT0)  ?  procesT0->GetStatus():0;
131   Int_t tpcStatus = (procesTPC) ? ((procesTPC->ValidateTimeDrift() || procesTPC->ValidateTimeGain())==kFALSE):0;
132   //
133   printf("\n");
134   printf("******* CPass0 calibration status *******\n");
135   printf("TRD calibration status=%d\n",trdStatus);
136   printf("TOF calibration status=%d\n",tofStatus);
137   printf("TPC calibration status=%d\n",tpcStatus);
138   printf("T0  calibration status=%d\n",t0Status);
139   //
140   TTreeSRedirector *pcstream = new TTreeSRedirector("cpassStat.root","recreate");
141   printCalibStat(runNumber, "CalibObjects.root",pcstream);
142   delete pcstream;
143   return;
144 }
145
146
147 // function to print statistics used to calibrate the various detectors
148
149 void printCalibStat(Int_t run, const char * fname,  TTreeSRedirector * pcstream){
150
151   //
152   // Dump the statistical information about all histograms in the calibration files 
153   //    into the statistical tree, print on the screen (log files) as well 
154   //
155   //
156   // 1. Default dump for all histograms
157   //    Information to dump:
158   //    stat =Entries, Mean, MeanError,  RMS, MaxBin
159   //    Branch naming convention:
160   //    <detName>_<hisName><statName>
161   //
162   // 2. Detector statistical information  - to be implemented by expert
163   //                                      - First version implemented by MI 
164   //  
165   // 
166
167   TFile *fin = TFile::Open(fname);
168   if (!fin) return;
169   const Double_t kMaxHis=10000;
170   
171   TList * keyList = fin->GetListOfKeys();
172   Int_t nkeys=keyList->GetEntries();
173   Double_t *hisEntries = new Double_t[kMaxHis];
174   Double_t *hisMean = new Double_t[kMaxHis];
175   Double_t *hisMeanError = new Double_t[kMaxHis];
176   Double_t *hisRMS = new Double_t[kMaxHis];
177   Double_t *hisMaxBin = new Double_t[kMaxHis];
178   Int_t counter=0;
179   
180   if (pcstream) (*pcstream)<<"calibStatAll"<<"run="<<run;
181   for (Int_t ikey=0; ikey<nkeys; ikey++){
182     TObject * object = fin->Get(keyList->At(ikey)->GetName());
183     if (!object) continue;
184     if (object->InheritsFrom("TCollection")==0) continue;
185     TSeqCollection *collection  = (TSeqCollection*)object; 
186     Int_t nentries= collection->GetEntries();
187     for (Int_t ihis=0; ihis<nentries; ihis++){
188       TObject * ohis = collection->At(ihis);
189       if (!ohis) continue;
190       if (ohis->InheritsFrom("TH1")==0) continue;
191       TH1* phis = (TH1*)ohis;
192       hisEntries[counter]=phis->GetEntries();   
193       Int_t idim=1;
194       if (ohis->InheritsFrom("TH2")) idim=2;
195       if (ohis->InheritsFrom("TH3")) idim=3;        
196       hisMean[counter]=phis->GetMean(idim);     
197       hisMeanError[counter]=phis->GetMeanError(idim);   
198       hisRMS[counter]=phis->GetRMS(idim);       
199       hisMaxBin[counter]=phis->GetBinCenter(phis->GetMaximumBin());     
200       if (pcstream) (*pcstream)<<"calibStatAll"<<
201                       Form("%s_%sEntries=",keyList->At(ikey)->GetName(), phis->GetName())<<hisEntries[counter]<<        
202                       Form("%s_%sMean=",keyList->At(ikey)->GetName(), phis->GetName())<<hisMean[counter]<<      
203                       Form("%s_%sMeanError=",keyList->At(ikey)->GetName(), phis->GetName())<<hisMeanError[counter]<<    
204                       Form("%s_%sRMS=",keyList->At(ikey)->GetName(), phis->GetName())<<hisRMS[counter]<<        
205                       Form("%s_%sMaxBin=",keyList->At(ikey)->GetName(), phis->GetName())<<hisMaxBin[counter];   
206       //printf("Histo:\t%s_%s\t%f\t%d\n",keyList->At(ikey)->GetName(), phis->GetName(), hisEntries[counter],idim);
207       counter++;
208     }
209     delete object;
210   }    
211   
212   //
213   // Expert dump example (MI first iteration):
214   //
215   // 0.)  TOF dump
216   //
217
218   Int_t tofEvents=0;
219   Int_t tofTracks=0;
220   TList * TOFCalib = (TList*)fin->Get("TOFHistos");      
221   if (TOFCalib) {
222     TH1 *histoEvents = (TH1*)TOFCalib->FindObject("hHistoVertexTimestamp");
223     TH1 *histoTracks = (TH1*)TOFCalib->FindObject("hHistoDeltatTimestamp");
224     if (histoEvents && histoTracks){
225       tofEvents = TMath::Nint(histoEvents->GetEntries());
226       tofTracks = TMath::Nint(histoTracks->GetEntries());
227     }
228     delete TOFCalib;
229   }
230   printf("Monalisa TOFevents\t%d\n",tofEvents);
231   if (pcstream) (*pcstream)<<"calibStatAll"<<"TOFevents="<<tofEvents;
232   printf("Monalisa TOFtracks\t%d\n",tofTracks);
233   if (pcstream) (*pcstream)<<"calibStatAll"<<"TOFtracks="<<tofTracks;
234
235   //
236   // 1.)  TPC  dump - usefull events/tracks  for the calibration
237   //
238   Int_t tpcEvents=0;
239   Int_t tpcTracks=0;
240   TList * TPCCalib = (TList*)fin->Get("TPCCalib");      
241   if (TPCCalib) {
242     AliTPCcalibTime  * calibTime = (AliTPCcalibTime  *)TPCCalib->FindObject("calibTime");
243     if (calibTime){
244       tpcEvents = TMath::Nint(calibTime->GetTPCVertexHisto(0)->GetEntries());
245       tpcTracks = TMath::Nint(calibTime->GetResHistoTPCITS(0)->GetEntries());
246     }
247   }
248   printf("Monalisa TPCevents\t%d\n",tpcEvents);
249   if (pcstream) (*pcstream)<<"calibStatAll"<<"TPCevents="<<tpcEvents;
250   printf("Monalisa TPCtracks\t%d\n",tpcTracks);
251   if (pcstream) (*pcstream)<<"calibStatAll"<<"TPCtracks="<<tpcTracks;
252
253   //
254   // 2. TRD dump 
255   //
256   Int_t trdEvents=0;
257   Int_t trdTracks=0;
258   TList * TRDCalib = (TList*)fin->Get("TRDCalib");      
259   if (TRDCalib) {
260     TH1  *histoEvents = (TH1*)TRDCalib->FindObject("NEventsInput_AliTRDCalibTask");
261     TH1  *histoTracks = (TH1*)TRDCalib->FindObject("AbsoluteGain_AliTRDCalibTask");
262     if (histoEvents && histoTracks){
263       trdEvents= TMath::Nint(histoEvents->GetEntries());
264       trdTracks= TMath::Nint(histoTracks->GetEntries());
265     }
266     delete TRDCalib;
267   }
268   printf("Monalisa TRDevents\t%d\n",trdEvents);
269   if (pcstream) (*pcstream)<<"calibStatAll"<<"TRDevents="<<trdEvents;
270   printf("Monalisa TRDtracks\t%d\n",trdTracks);
271   if (pcstream) (*pcstream)<<"calibStatAll"<<"TRDtracks="<<trdTracks;
272
273   //
274   // 3. T0 dump 
275   //
276   Int_t T0Events=0;
277   TList * T0Calib = (TList*)fin->Get("T0Calib");      
278   if (T0Calib) {
279     TH1  *histoEvents = (TH1*) T0Calib->FindObject("fTzeroORAplusORC");
280     if (histoEvents){
281       T0Events= TMath::Nint(histoEvents->GetEntries());
282     }
283     delete T0Calib;
284   }
285   printf("Monalisa T0events\t%d\n",T0Events);
286   if (pcstream) (*pcstream)<<"calibStatAll"<<"T0events="<<T0Events;
287
288   //
289   // 4. Mean vertex -   dump 
290     Int_t meanVertexEvents=0;
291   TList * meanVertexCalib = (TList*)fin->Get("MeanVertex");      
292   if (meanVertexCalib) {
293     TH1  *histoEvents = (TH1*) meanVertexCalib->FindObject("hTRKVertexX");
294     if (histoEvents){
295       meanVertexEvents = TMath::Nint(histoEvents->GetEntries());
296     }
297     delete meanVertexCalib;
298   }
299   printf("Monalisa MeanVertexevents\t%d\n",meanVertexEvents);
300   if (pcstream) (*pcstream)<<"calibStatAll"<<"MeanVertexevents="<<meanVertexEvents;
301
302   //
303   // 5. SDD dump 
304   //
305   Int_t sddEvents=0;
306   Int_t sddTracks=0;
307   TList * SDDCalib = (TList*)fin->Get("clistSDDCalib");      
308   if (SDDCalib) {
309     TH1  *histoEvents = (TH1*) SDDCalib->FindObject("hNEvents");
310     if (histoEvents ){
311       sddEvents = TMath::Nint(histoEvents->GetBinContent(4));
312       sddTracks = TMath::Nint(histoEvents->GetBinContent(5));
313     }
314     delete SDDCalib;
315   }
316   printf("Monalisa SDDevents\t%d\n",sddEvents);
317   if (pcstream) (*pcstream)<<"calibStatAll"<<"SDDevents="<<sddEvents;
318   printf("Monalisa SDDtracks\t%d\n",sddTracks);
319   if (pcstream) (*pcstream)<<"calibStatAll"<<"SDDtracks="<<sddTracks;
320
321   //
322   if (pcstream) (*pcstream)<<"calibStatAll"<<"\n";
323   delete fin;
324
325 }
326