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