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