]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/CalibMacros/CPass0/makeOCDB.C
#97528: Commit to the trunk changes needed for mirroring OCDB entries during CPass
[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
16 void makeOCDB(Int_t runNumber, TString  targetOCDBstorage="", TString sourceOCDBstorage="raw://")
17 {
18
19   //
20   // extract OCDB entries for detectors participating in the calibration for the current run
21   //
22
23   gROOT->Macro("$ALICE_ROOT/PWGPP/CalibMacros/CPass0/LoadLibraries.C");
24   gROOT->LoadMacro("$ALICE_ROOT/PWGPP/CalibMacros/CPass0/ConfigCalibTrain.C");
25
26   // switch off log info
27   AliLog::SetClassDebugLevel("AliESDEvent",0);
28
29   // config GRP
30   printf("runNumber from runCalibTrain = %d\n",runNumber);
31   ConfigCalibTrain(runNumber, sourceOCDBstorage.Data());
32
33   // check the presence of the detectors
34   AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
35   AliGRPObject* grpData = dynamic_cast<AliGRPObject*>(entry->GetObject());
36   if (!grpData) {printf("Failed to get GRP data for run",runNumber); return;}
37   Int_t activeDetectors = grpData->GetDetectorMask();
38   TString detStr = AliDAQ::ListOfTriggeredDetectors(activeDetectors);
39   printf("Detectors in the data:\n%s\n",detStr.Data());
40
41   // Steering Tasks - set output storage
42   // DefaultStorage set already before - in ConfigCalibTrain.C
43
44   // Setting the mirror SEs for the default storage
45   TString mirrorsStr("ALICE::CERN::OCDB,ALICE::FZK::SE,ALICE::LLNL::SE");
46   AliCDBManager::Instance()->SetMirrorSEs(mirrorsStr.Data());
47   printf("List of mirror SEs set to: \"%s\"\n",mirrorsStr.Data());
48
49   // activate target OCDB storage
50   AliCDBStorage* targetStorage = 0x0;
51   if (targetOCDBstorage.Length()==0) targetOCDBstorage+="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
52   else if (targetOCDBstorage.CompareTo("same",TString::kIgnoreCase) == 0 ){
53     targetStorage = AliCDBManager::Instance()->GetDefaultStorage();
54   }
55   else {
56     targetStorage = AliCDBManager::Instance()->GetStorage(targetOCDBstorage.Data());
57     targetStorage->SetMirrorSEs(mirrorsStr.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   //
144   printf("\n");
145   printf("******* CPass0 calibration status *******\n");
146   printf("TRD calibration status=%d\n",trdStatus);
147   printf("TOF calibration status=%d\n",tofStatus);
148   printf("TPC calibration status=%d\n",tpcStatus);
149   printf("T0  calibration status=%d\n",t0Status);
150   //
151   TTreeSRedirector *pcstream = new TTreeSRedirector("cpassStat.root","recreate");
152   printCalibStat(runNumber, "CalibObjects.root",pcstream);
153   delete pcstream;
154   return;
155 }
156
157
158 // function to print statistics used to calibrate the various detectors
159
160 void printCalibStat(Int_t run, const char * fname,  TTreeSRedirector * pcstream){
161
162   //
163   // Dump the statistical information about all histograms in the calibration files 
164   //    into the statistical tree, print on the screen (log files) as well 
165   //
166   //
167   // 1. Default dump for all histograms
168   //    Information to dump:
169   //    stat =Entries, Mean, MeanError,  RMS, MaxBin
170   //    Branch naming convention:
171   //    <detName>_<hisName><statName>
172   //
173   // 2. Detector statistical information  - to be implemented by expert
174   //                                      - First version implemented by MI 
175   //  
176   // 
177
178   TFile *fin = TFile::Open(fname);
179   if (!fin) return;
180   const Double_t kMaxHis=10000;
181   
182   TList * keyList = fin->GetListOfKeys();
183   Int_t nkeys=keyList->GetEntries();
184   Double_t *hisEntries = new Double_t[kMaxHis];
185   Double_t *hisMean = new Double_t[kMaxHis];
186   Double_t *hisMeanError = new Double_t[kMaxHis];
187   Double_t *hisRMS = new Double_t[kMaxHis];
188   Double_t *hisMaxBin = new Double_t[kMaxHis];
189   Int_t counter=0;
190   
191   if (pcstream) (*pcstream)<<"calibStatAll"<<"run="<<run;
192   for (Int_t ikey=0; ikey<nkeys; ikey++){
193     TObject * object = fin->Get(keyList->At(ikey)->GetName());
194     if (!object) continue;
195     if (object->InheritsFrom("TCollection")==0) continue;
196     TSeqCollection *collection  = (TSeqCollection*)object; 
197     Int_t nentries= collection->GetEntries();
198     for (Int_t ihis=0; ihis<nentries; ihis++){
199       TObject * ohis = collection->At(ihis);
200       if (!ohis) continue;
201       if (ohis->InheritsFrom("TH1")==0) continue;
202       TH1* phis = (TH1*)ohis;
203       hisEntries[counter]=phis->GetEntries();   
204       Int_t idim=1;
205       if (ohis->InheritsFrom("TH2")) idim=2;
206       if (ohis->InheritsFrom("TH3")) idim=3;        
207       hisMean[counter]=phis->GetMean(idim);     
208       hisMeanError[counter]=phis->GetMeanError(idim);   
209       hisRMS[counter]=phis->GetRMS(idim);       
210       hisMaxBin[counter]=phis->GetBinCenter(phis->GetMaximumBin());     
211       if (pcstream) (*pcstream)<<"calibStatAll"<<
212                       Form("%s_%sEntries=",keyList->At(ikey)->GetName(), phis->GetName())<<hisEntries[counter]<<        
213                       Form("%s_%sMean=",keyList->At(ikey)->GetName(), phis->GetName())<<hisMean[counter]<<      
214                       Form("%s_%sMeanError=",keyList->At(ikey)->GetName(), phis->GetName())<<hisMeanError[counter]<<    
215                       Form("%s_%sRMS=",keyList->At(ikey)->GetName(), phis->GetName())<<hisRMS[counter]<<        
216                       Form("%s_%sMaxBin=",keyList->At(ikey)->GetName(), phis->GetName())<<hisMaxBin[counter];   
217       //printf("Histo:\t%s_%s\t%f\t%d\n",keyList->At(ikey)->GetName(), phis->GetName(), hisEntries[counter],idim);
218       counter++;
219     }
220     delete object;
221   }    
222   
223   //
224   // Expert dump example (MI first iteration):
225   //
226   // 0.)  TOF dump
227   //
228
229   Int_t tofEvents=0;
230   Int_t tofTracks=0;
231   TList * TOFCalib = (TList*)fin->Get("TOFHistos");      
232   if (TOFCalib) {
233     TH1 *histoEvents = (TH1*)TOFCalib->FindObject("hHistoVertexTimestamp");
234     TH1 *histoTracks = (TH1*)TOFCalib->FindObject("hHistoDeltatTimestamp");
235     if (histoEvents && histoTracks){
236       tofEvents = TMath::Nint(histoEvents->GetEntries());
237       tofTracks = TMath::Nint(histoTracks->GetEntries());
238     }
239     delete TOFCalib;
240   }
241   printf("Monalisa TOFevents\t%d\n",tofEvents);
242   if (pcstream) (*pcstream)<<"calibStatAll"<<"TOFevents="<<tofEvents;
243   printf("Monalisa TOFtracks\t%d\n",tofTracks);
244   if (pcstream) (*pcstream)<<"calibStatAll"<<"TOFtracks="<<tofTracks;
245
246   //
247   // 1.)  TPC  dump - usefull events/tracks  for the calibration
248   //
249   Int_t tpcEvents=0;
250   Int_t tpcTracks=0;
251   TList * TPCCalib = (TList*)fin->Get("TPCCalib");      
252   if (TPCCalib) {
253     AliTPCcalibTime  * calibTime = (AliTPCcalibTime  *)TPCCalib->FindObject("calibTime");
254     if (calibTime){
255       tpcEvents = TMath::Nint(calibTime->GetTPCVertexHisto(0)->GetEntries());
256       tpcTracks = TMath::Nint(calibTime->GetResHistoTPCITS(0)->GetEntries());
257     }
258   }
259   printf("Monalisa TPCevents\t%d\n",tpcEvents);
260   if (pcstream) (*pcstream)<<"calibStatAll"<<"TPCevents="<<tpcEvents;
261   printf("Monalisa TPCtracks\t%d\n",tpcTracks);
262   if (pcstream) (*pcstream)<<"calibStatAll"<<"TPCtracks="<<tpcTracks;
263
264   //
265   // 2. TRD dump 
266   //
267   Int_t trdEvents=0;
268   Int_t trdTracks=0;
269   TList * TRDCalib = (TList*)fin->Get("TRDCalib");      
270   if (TRDCalib) {
271     TH1  *histoEvents = (TH1*)TRDCalib->FindObject("NEventsInput_AliTRDCalibTask");
272     TH1  *histoTracks = (TH1*)TRDCalib->FindObject("AbsoluteGain_AliTRDCalibTask");
273     if (histoEvents && histoTracks){
274       trdEvents= TMath::Nint(histoEvents->GetEntries());
275       trdTracks= TMath::Nint(histoTracks->GetEntries());
276     }
277     delete TRDCalib;
278   }
279   printf("Monalisa TRDevents\t%d\n",trdEvents);
280   if (pcstream) (*pcstream)<<"calibStatAll"<<"TRDevents="<<trdEvents;
281   printf("Monalisa TRDtracks\t%d\n",trdTracks);
282   if (pcstream) (*pcstream)<<"calibStatAll"<<"TRDtracks="<<trdTracks;
283
284   //
285   // 3. T0 dump 
286   //
287   Int_t T0Events=0;
288   TList * T0Calib = (TList*)fin->Get("T0Calib");      
289   if (T0Calib) {
290     TH1  *histoEvents = (TH1*) T0Calib->FindObject("fTzeroORAplusORC");
291     if (histoEvents){
292       T0Events= TMath::Nint(histoEvents->GetEntries());
293     }
294     delete T0Calib;
295   }
296   printf("Monalisa T0events\t%d\n",T0Events);
297   if (pcstream) (*pcstream)<<"calibStatAll"<<"T0events="<<T0Events;
298
299   //
300   // 4. Mean vertex -   dump 
301     Int_t meanVertexEvents=0;
302   TList * meanVertexCalib = (TList*)fin->Get("MeanVertex");      
303   if (meanVertexCalib) {
304     TH1  *histoEvents = (TH1*) meanVertexCalib->FindObject("hTRKVertexX");
305     if (histoEvents){
306       meanVertexEvents = TMath::Nint(histoEvents->GetEntries());
307     }
308     delete meanVertexCalib;
309   }
310   printf("Monalisa MeanVertexevents\t%d\n",meanVertexEvents);
311   if (pcstream) (*pcstream)<<"calibStatAll"<<"MeanVertexevents="<<meanVertexEvents;
312
313   //
314   // 5. SDD dump 
315   //
316   Int_t sddEvents=0;
317   Int_t sddTracks=0;
318   TList * SDDCalib = (TList*)fin->Get("clistSDDCalib");      
319   if (SDDCalib) {
320     TH1  *histoEvents = (TH1*) SDDCalib->FindObject("hNEvents");
321     if (histoEvents ){
322       sddEvents = TMath::Nint(histoEvents->GetBinContent(4));
323       sddTracks = TMath::Nint(histoEvents->GetBinContent(5));
324     }
325     delete SDDCalib;
326   }
327   printf("Monalisa SDDevents\t%d\n",sddEvents);
328   if (pcstream) (*pcstream)<<"calibStatAll"<<"SDDevents="<<sddEvents;
329   printf("Monalisa SDDtracks\t%d\n",sddTracks);
330   if (pcstream) (*pcstream)<<"calibStatAll"<<"SDDtracks="<<sddTracks;
331
332   //
333   if (pcstream) (*pcstream)<<"calibStatAll"<<"\n";
334   delete fin;
335
336 }
337