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