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