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