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