]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/CalibMacros/CPass0/makeOCDB.C
Merge branch 'feature-movesplit'
[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://", Int_t detectorBitsQualityFlag = -1)
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   // Quality flags
82   Bool_t TPC_qf = kTRUE;
83   Bool_t TOF_qf = kTRUE;
84   Bool_t TRD_qf = kTRUE;
85   Bool_t T0_qf  = kTRUE;
86   Bool_t SDD_qf = kTRUE;
87   Bool_t SPD_qf = kTRUE;
88
89   if (detectorBitsQualityFlag != -1){
90     TPC_qf = ((detectorBitsQualityFlag & AliDAQ::kTPC_QF) == AliDAQ::kTPC_QF)? kTRUE : kFALSE;
91     TOF_qf = ((detectorBitsQualityFlag & AliDAQ::kTOF_QF) == AliDAQ::kTOF_QF)? kTRUE : kFALSE;
92     TRD_qf = ((detectorBitsQualityFlag & AliDAQ::kTRD_QF) == AliDAQ::kTRD_QF)? kTRUE : kFALSE;
93     T0_qf  = ((detectorBitsQualityFlag & AliDAQ::kT0_QF)  == AliDAQ::kT0_QF)?  kTRUE : kFALSE;
94     SDD_qf = ((detectorBitsQualityFlag & AliDAQ::kSDD_QF) == AliDAQ::kSDD_QF)? kTRUE : kFALSE;
95     SPD_qf = ((detectorBitsQualityFlag & AliDAQ::kSPD_QF) == AliDAQ::kSPD_QF)? kTRUE : kFALSE;
96   }    
97
98   Printf("Quality flags: detectorBitsQualityFlag = %d, TPC = %d, TOF = %d, TRD = %d, T0 = %d, SDD = %d, SPD = %d", detectorBitsQualityFlag, (Int_t)TPC_qf, (Int_t)TOF_qf, (Int_t)TRD_qf, (Int_t)T0_qf, (Int_t)SDD_qf, (Int_t)SPD_qf);
99
100   // TPC part
101   AliTPCPreprocessorOffline *procesTPC = 0;
102   if (detStr.Contains("TPC") && TPC_qf){
103     Printf("\n******* Calibrating TPC *******");
104     procesTPC = new AliTPCPreprocessorOffline;
105     // switch on parameter validation
106     procesTPC->SetTimeGainRange(0.5,5.0);
107     procesTPC->SetMaxVDriftCorr(0.2); 
108     //procesTPC->SetMinTracksVdrift(100000);
109     procesTPC->SwitchOnValidation();
110
111     // Make timegain calibration
112     //proces.CalibTimeGain("CalibObjects.root", runNumber,AliCDBRunRange::Infinity(),targetOCDBstorage);
113     if (isMagFieldON) procesTPC->CalibTimeGain("CalibObjects.root", runNumber,runNumber,targetStorage);
114     
115     // Make vdrift calibration
116     //proces.CalibTimeVdrift("CalibObjects.root",runNumber,AliCDBRunRange::Infinity(),targetOCDBstorage);
117     procesTPC->CalibTimeVdrift("CalibObjects.root",runNumber,runNumber,targetStorage);
118   }
119   else {
120     Printf("\n******* NOT Calibrating TPC: detStr = %s, TPC_qf = %d *******", detStr.Data(), (Int_t)TPC_qf);
121   }
122
123   // TOF part
124   AliTOFAnalysisTaskCalibPass0 *procesTOF=0;
125   if (detStr.Contains("TOF") && detStr.Contains("TPC") && TOF_qf){
126     procesTOF = new AliTOFAnalysisTaskCalibPass0;
127     Printf("\n******* Calibrating TOF *******");
128     if (isMagFieldON) procesTOF->ProcessOutput("CalibObjects.root", targetStorage);
129     else {
130       printf("Not calibrating TOF in case of mag field OFF\n");
131     }
132   }
133   else {
134     Printf("\n******* NOT Calibrating TOF: detStr = %s, TOF_qf = %d *******", detStr.Data(), (Int_t)TOF_qf);
135   }
136
137   // T0 part
138   AliT0PreprocessorOffline *procesT0= 0;
139   if (detStr.Contains("T0") && T0_qf) {
140     Printf("\n******* Calibrating T0 *******");
141     // Make  calibration of channels offset
142     procesT0 = new AliT0PreprocessorOffline;
143     if(isLHC10)
144       procesT0->CalibOffsetChannels("CalibObjects.root",runNumber, runNumber, targetStorage);
145     else 
146       procesT0->Process("CalibObjects.root",runNumber, runNumber, targetStorage);
147   }
148   else {
149     Printf("\n******* NOT Calibrating T0: detStr = %s, T0_qf = %d *******", detStr.Data(), (Int_t)T0_qf);
150   }
151
152   //TRD part
153   AliTRDPreprocessorOffline *procesTRD = 0;
154   if (detStr.Contains("TRD") && detStr.Contains("TPC") && TRD_qf){
155     Printf("\n******* Calibrating TRD *******");
156     procesTRD = new  AliTRDPreprocessorOffline;
157     if(isLHC10) procesTRD->SetSwitchOnChamberStatus(kFALSE);
158     procesTRD->SetLinearFitForVdrift(kTRUE);
159     procesTRD->SetMinStatsVdriftT0PH(600*10);
160     procesTRD->SetMinStatsVdriftLinear(50);
161     procesTRD->SetMinStatsGain(600);
162     procesTRD->SetLimitValidateNoData(60);
163     procesTRD->SetLimitValidateBadCalib(60);
164     procesTRD->SetAlternativeDriftVelocityFit(kTRUE);
165     procesTRD->Init("CalibObjects.root");
166     Int_t versionVdriftUsed = procesTRD->GetVersionVdriftUsed();
167     Int_t subversionVdriftUsed = procesTRD->GetSubVersionVdriftUsed();
168     Int_t versionGainUsed = procesTRD->GetVersionGainUsed();
169     Int_t subversionGainUsed = procesTRD->GetSubVersionGainUsed();
170     Int_t versionExBUsed = procesTRD->GetVersionExBUsed();
171     Int_t subversionExBUsed = procesTRD->GetSubVersionExBUsed();
172     printf("version and subversion vdrift %d and %d\n",versionVdriftUsed,subversionVdriftUsed);
173     printf("version and subversion gain %d and %d\n",versionGainUsed,subversionGainUsed);
174     printf("version and subversion exb %d and %d\n",versionExBUsed,subversionExBUsed);
175     procesTRD->Process("CalibObjects.root",runNumber,runNumber,targetStorage);
176   }
177   else {
178     Printf("\n******* NOT Calibrating TRD: detStr = %s, TRD_qf = %d *******", detStr.Data(), (Int_t)TRD_qf);
179   }
180   
181   TF1 *gsf = (TF1 *)gROOT->GetFunction("gaus");
182   if (gsf) for (int i=gsf->GetNpar();i--;) gsf->SetParError(i,0); // reset errors from previous fits
183   
184   //Mean Vertex
185   AliMeanVertexPreprocessorOffline * procesMeanVtx=0;
186   if (detStr.Contains("ITSSPD") && SPD_qf) {
187     Printf("\n******* Calibrating MeanVertex *******");
188     procesMeanVtx = new AliMeanVertexPreprocessorOffline;
189     procesMeanVtx->ProcessOutput("CalibObjects.root", targetStorage, runNumber);
190   }
191   else {
192     Printf("\n******* NOT Calibrating MeanVertex: detStr = %s, SPD_qf = %d *******", detStr.Data(), (Int_t)SPD_qf);
193   }
194
195   //
196   // Print calibration status into the stdout
197   //
198   Int_t trdStatus = (procesTRD) ?  procesTRD->GetStatus():0;
199   Int_t tofStatus = (procesTOF) ?  procesTOF->GetStatus():0;
200   Int_t t0Status  = (procesT0)  ?  procesT0->GetStatus():0;
201   Int_t tpcStatus = (procesTPC) ?  procesTPC->GetStatus():0;
202   Int_t meanVtxStatus = (procesMeanVtx) ? procesMeanVtx->GetStatus():0;
203   //
204   printf("\n");
205   printf("******* CPass0 calibration status *******\n");
206   printf("TRD calibration status=%d\n",trdStatus);
207   printf("TOF calibration status=%d\n",tofStatus);
208   printf("TPC calibration status=%d\n",tpcStatus);
209   printf("T0  calibration status=%d\n",t0Status);
210   printf("MeanVertex  calibration status=%d\n",meanVtxStatus);
211   //
212   TTreeSRedirector *pcstream = new TTreeSRedirector("cpassStat.root","recreate");
213   printCalibStat(runNumber, "CalibObjects.root",pcstream);
214   delete pcstream;
215   return;
216 }
217
218
219 // function to print statistics used to calibrate the various detectors
220
221 void printCalibStat(Int_t run, const char * fname,  TTreeSRedirector * pcstream){
222
223   //
224   // Dump the statistical information about all histograms in the calibration files 
225   //    into the statistical tree, print on the screen (log files) as well 
226   //
227   //
228   // 1. Default dump for all histograms
229   //    Information to dump:
230   //    stat =Entries, Mean, MeanError,  RMS, MaxBin
231   //    Branch naming convention:
232   //    <detName>_<hisName><statName>
233   //
234   // 2. Detector statistical information  - to be implemented by expert
235   //                                      - First version implemented by MI 
236   //  
237   // 
238
239   TFile *fin = TFile::Open(fname);
240   if (!fin) return;
241   const Double_t kMaxHis=10000;
242   
243   TList * keyList = fin->GetListOfKeys();
244   Int_t nkeys=keyList->GetEntries();
245   Double_t *hisEntries = new Double_t[kMaxHis];
246   Double_t *hisMean = new Double_t[kMaxHis];
247   Double_t *hisMeanError = new Double_t[kMaxHis];
248   Double_t *hisRMS = new Double_t[kMaxHis];
249   Double_t *hisMaxBin = new Double_t[kMaxHis];
250   Int_t counter=0;
251   
252   if (pcstream) (*pcstream)<<"calibStatAll"<<"run="<<run;
253   for (Int_t ikey=0; ikey<nkeys; ikey++){
254     TObject * object = fin->Get(keyList->At(ikey)->GetName());
255     if (!object) continue;
256     if (object->InheritsFrom("TCollection")==0) continue;
257     TSeqCollection *collection  = (TSeqCollection*)object; 
258     Int_t nentries= collection->GetEntries();
259     for (Int_t ihis=0; ihis<nentries; ihis++){
260       TObject * ohis = collection->At(ihis);
261       if (!ohis) continue;
262       if (ohis->InheritsFrom("TH1")==0) continue;
263       TH1* phis = (TH1*)ohis;
264       hisEntries[counter]=phis->GetEntries();   
265       Int_t idim=1;
266       if (ohis->InheritsFrom("TH2")) idim=2;
267       if (ohis->InheritsFrom("TH3")) idim=3;        
268       hisMean[counter]=phis->GetMean(idim);     
269       hisMeanError[counter]=phis->GetMeanError(idim);   
270       hisRMS[counter]=phis->GetRMS(idim);       
271       hisMaxBin[counter]=phis->GetBinCenter(phis->GetMaximumBin());     
272       if (pcstream) (*pcstream)<<"calibStatAll"<<
273                       Form("%s_%sEntries=",keyList->At(ikey)->GetName(), phis->GetName())<<hisEntries[counter]<<        
274                       Form("%s_%sMean=",keyList->At(ikey)->GetName(), phis->GetName())<<hisMean[counter]<<      
275                       Form("%s_%sMeanError=",keyList->At(ikey)->GetName(), phis->GetName())<<hisMeanError[counter]<<    
276                       Form("%s_%sRMS=",keyList->At(ikey)->GetName(), phis->GetName())<<hisRMS[counter]<<        
277                       Form("%s_%sMaxBin=",keyList->At(ikey)->GetName(), phis->GetName())<<hisMaxBin[counter];   
278       //printf("Histo:\t%s_%s\t%f\t%d\n",keyList->At(ikey)->GetName(), phis->GetName(), hisEntries[counter],idim);
279       counter++;
280     }
281     delete object;
282   }    
283   
284   //
285   // Expert dump example (MI first iteration):
286   //
287   // 0.)  TOF dump
288   //
289
290   Int_t tofEvents=0;
291   Int_t tofTracks=0;
292   TList * TOFCalib = (TList*)fin->Get("TOFHistos");      
293   if (TOFCalib) {
294     TH1 *histoEvents = (TH1*)TOFCalib->FindObject("hHistoVertexTimestamp");
295     TH1 *histoTracks = (TH1*)TOFCalib->FindObject("hHistoDeltatTimestamp");
296     if (histoEvents && histoTracks){
297       tofEvents = TMath::Nint(histoEvents->GetEntries());
298       tofTracks = TMath::Nint(histoTracks->GetEntries());
299     }
300     delete TOFCalib;
301   }
302   printf("Monalisa TOFevents\t%d\n",tofEvents);
303   if (pcstream) (*pcstream)<<"calibStatAll"<<"TOFevents="<<tofEvents;
304   printf("Monalisa TOFtracks\t%d\n",tofTracks);
305   if (pcstream) (*pcstream)<<"calibStatAll"<<"TOFtracks="<<tofTracks;
306
307   //
308   // 1.)  TPC  dump - usefull events/tracks  for the calibration
309   //
310   Int_t tpcEvents=0;
311   Int_t tpcTracks=0;
312   TObject* obj = dynamic_cast<TObject*>(fin->Get("TPCCalib"));
313   TObjArray* array = dynamic_cast<TObjArray*>(obj);
314   TDirectory* dir = dynamic_cast<TDirectory*>(obj);
315   AliTPCcalibTime  * calibTime = NULL;
316   if (dir) {
317     calibTime = dynamic_cast<AliTPCcalibTime*>(dir->Get("calibTime"));
318   }
319   else if (array){
320     calibTime = (AliTPCcalibTime *)array->FindObject("calibTime");
321   }
322   if (calibTime) {
323       tpcEvents = TMath::Nint(calibTime->GetTPCVertexHisto(0)->GetEntries());
324       tpcTracks = TMath::Nint(calibTime->GetResHistoTPCITS(0)->GetEntries());
325   }
326   printf("Monalisa TPCevents\t%d\n",tpcEvents);
327   if (pcstream) (*pcstream)<<"calibStatAll"<<"TPCevents="<<tpcEvents;
328   printf("Monalisa TPCtracks\t%d\n",tpcTracks);
329   if (pcstream) (*pcstream)<<"calibStatAll"<<"TPCtracks="<<tpcTracks;
330
331   //
332   // 2. TRD dump 
333   //
334   Int_t trdEvents=0;
335   Int_t trdTracks=0;
336   TList * TRDCalib = (TList*)fin->Get("TRDCalib");      
337   if (TRDCalib) {
338     TH1  *histoEvents = (TH1*)TRDCalib->FindObject("NEventsInput_AliTRDCalibTask");
339     TH1  *histoTracks = (TH1*)TRDCalib->FindObject("AbsoluteGain_AliTRDCalibTask");
340     if (histoEvents && histoTracks){
341       trdEvents= TMath::Nint(histoEvents->GetEntries());
342       trdTracks= TMath::Nint(histoTracks->GetEntries());
343     }
344     delete TRDCalib;
345   }
346   printf("Monalisa TRDevents\t%d\n",trdEvents);
347   if (pcstream) (*pcstream)<<"calibStatAll"<<"TRDevents="<<trdEvents;
348   printf("Monalisa TRDtracks\t%d\n",trdTracks);
349   if (pcstream) (*pcstream)<<"calibStatAll"<<"TRDtracks="<<trdTracks;
350
351   //
352   // 3. T0 dump 
353   //
354   Int_t T0Events=0;
355   TList * T0Calib = (TList*)fin->Get("T0Calib");      
356   if (T0Calib) {
357     TH1  *histoEvents = (TH1*) T0Calib->FindObject("fTzeroORAplusORC");
358     if (histoEvents){
359       T0Events= TMath::Nint(histoEvents->GetEntries());
360     }
361     delete T0Calib;
362   }
363   printf("Monalisa T0events\t%d\n",T0Events);
364   if (pcstream) (*pcstream)<<"calibStatAll"<<"T0events="<<T0Events;
365
366   //
367   // 4. Mean vertex -   dump 
368     Int_t meanVertexEvents=0;
369   TList * meanVertexCalib = (TList*)fin->Get("MeanVertex");      
370   if (meanVertexCalib) {
371     TH1  *histoEvents = (TH1*) meanVertexCalib->FindObject("hTRKVertexX");
372     if (histoEvents){
373       meanVertexEvents = TMath::Nint(histoEvents->GetEntries());
374     }
375     delete meanVertexCalib;
376   }
377   printf("Monalisa MeanVertexevents\t%d\n",meanVertexEvents);
378   if (pcstream) (*pcstream)<<"calibStatAll"<<"MeanVertexevents="<<meanVertexEvents;
379
380   //
381   // 5. SDD dump 
382   //
383   Int_t sddEvents=0;
384   Int_t sddTracks=0;
385   TList * SDDCalib = (TList*)fin->Get("clistSDDCalib");      
386   if (SDDCalib) {
387     TH1  *histoEvents = (TH1*) SDDCalib->FindObject("hNEvents");
388     if (histoEvents ){
389       sddEvents = TMath::Nint(histoEvents->GetBinContent(4));
390       sddTracks = TMath::Nint(histoEvents->GetBinContent(5));
391     }
392     delete SDDCalib;
393   }
394   printf("Monalisa SDDevents\t%d\n",sddEvents);
395   if (pcstream) (*pcstream)<<"calibStatAll"<<"SDDevents="<<sddEvents;
396   printf("Monalisa SDDtracks\t%d\n",sddTracks);
397   if (pcstream) (*pcstream)<<"calibStatAll"<<"SDDtracks="<<sddTracks;
398
399   //
400   if (pcstream) (*pcstream)<<"calibStatAll"<<"\n";
401   delete fin;
402
403 }
404