]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/CalibMacros/CPass1/makeOCDB.C
Update, according to macros on alien from latest running (LHC13g + update from T0...
[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://")
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
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     Printf("TPC won't be calibrated at CPass1 for the time being... Doing nothing here");
71     //procesTPC = new AliTPCPreprocessorOffline;
72     // switch on parameter validation
73     //procesTPC->SetTimeGainRange(0.5,4.0);
74     //procesTPC->SetMinTracksVdrift(100000);
75     //procesTPC->SwitchOnValidation();
76     // Make timegain calibration
77     //proces.CalibTimeGain("CalibObjects.root", runNumber,AliCDBRunRange::Infinity(),targetStorage);
78     // Make vdrift calibration
79     //proces.CalibTimeVdrift("CalibObjects.root",runNumber,AliCDBRunRange::Infinity(),targetStorage);
80   }
81
82   // TOF part
83   AliTOFAnalysisTaskCalibPass0 *procesTOF=0;
84   if (detStr.Contains("TOF") && detStr.Contains("TPC")){
85     procesTOF = new AliTOFAnalysisTaskCalibPass0;
86     Printf("\n******* Calibrating TOF *******");
87     procesTOF->ProcessOutput("CalibObjects.root", targetStorage);
88   }
89
90   // T0 part
91   AliT0PreprocessorOffline *procesT0= 0;
92   if (detStr.Contains("T0")) {
93     Printf("\n******* Calibrating T0 *******");
94     procesT0 = new AliT0PreprocessorOffline;
95     // Make  calibration of channels offset
96     procesT0->setDArun(100000);
97     procesT0->Process("CalibObjects.root",runNumber, runNumber, targetStorage);
98   }
99
100   //TRD part
101   AliTRDPreprocessorOffline *procesTRD = 0;
102   if (detStr.Contains("TRD") && detStr.Contains("TPC")){
103     Printf("\n******* Calibrating TRD *******");
104     procesTRD = new  AliTRDPreprocessorOffline;
105     procesTRD->SetLinearFitForVdrift(kTRUE);
106     procesTRD->SetMinStatsVdriftT0PH(600*10);
107     procesTRD->SetMinStatsVdriftLinear(50);
108     procesTRD->SetMinStatsGain(600);
109     procesTRD->SetLimitValidateNoData(60);
110     procesTRD->SetLimitValidateBadCalib(60);
111     procesTRD->SetAlternativeDriftVelocityFit(kTRUE);
112     procesTRD->Init("CalibObjects.root");
113     Int_t versionVdriftUsed = procesTRD->GetVersionVdriftUsed();
114     Int_t subversionVdriftUsed = procesTRD->GetSubVersionVdriftUsed();
115     Int_t versionGainUsed = procesTRD->GetVersionGainUsed();
116     Int_t subversionGainUsed = procesTRD->GetSubVersionGainUsed();
117     Int_t versionExBUsed = procesTRD->GetVersionExBUsed();
118     Int_t subversionExBUsed = procesTRD->GetSubVersionExBUsed();
119     printf("version and subversion vdrift %d and %d\n",versionVdriftUsed,subversionVdriftUsed);
120     printf("version and subversion gain %d and %d\n",versionGainUsed,subversionGainUsed);
121     printf("version and subversion exb %d and %d\n",versionExBUsed,subversionExBUsed);
122     procesTRD->Process("CalibObjects.root",runNumber,runNumber,targetStorage);
123   }
124   
125   // switched OFF at CPass1 in any case
126   /*
127   //Mean Vertex
128   AliMeanVertexPreprocessorOffline * procesMeanVtx=0;
129   if (detStr.Contains("ITSSPD")) {
130     Printf("\n******* Calibrating MeanVertex *******");
131     procesMeanVtx = new AliMeanVertexPreprocessorOffline;
132     procesMeanVtx->ProcessOutput("CalibObjects.root", targetStorage, runNumber);
133   }
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("******* CPass1 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   TObject* obj = dynamic_cast<TObject*>(fin->Get("TPCCalib"));
252   TObjArray* array = dynamic_cast<TObjArray*>(obj);
253   TDirectory* dir = dynamic_cast<TDirectory*>(obj);
254   AliTPCcalibTime  * calibTime = NULL;
255   if (dir) {
256     calibTime = dynamic_cast<AliTPCcalibTime*>(dir->Get("calibTime"));
257   }
258   else if (array){
259     calibTime = (AliTPCcalibTime *)array->FindObject("calibTime");
260   }
261   if (calibTime) {
262       tpcEvents = TMath::Nint(calibTime->GetTPCVertexHisto(0)->GetEntries());
263       tpcTracks = TMath::Nint(calibTime->GetResHistoTPCITS(0)->GetEntries());
264   }
265   printf("Monalisa TPCevents\t%d\n",tpcEvents);
266   if (pcstream) (*pcstream)<<"calibStatAll"<<"TPCevents="<<tpcEvents;
267   printf("Monalisa TPCtracks\t%d\n",tpcTracks);
268   if (pcstream) (*pcstream)<<"calibStatAll"<<"TPCtracks="<<tpcTracks;
269
270   //
271   // 2. TRD dump 
272   //
273   Int_t trdEvents=0;
274   Int_t trdTracks=0;
275   TList * TRDCalib = (TList*)fin->Get("TRDCalib");      
276   if (TRDCalib) {
277     TH1  *histoEvents = (TH1*)TRDCalib->FindObject("NEventsInput_AliTRDCalibTask");
278     TH1  *histoTracks = (TH1*)TRDCalib->FindObject("AbsoluteGain_AliTRDCalibTask");
279     if (histoEvents && histoTracks){
280       trdEvents= TMath::Nint(histoEvents->GetEntries());
281       trdTracks= TMath::Nint(histoTracks->GetEntries());
282     }
283     delete TRDCalib;
284   }
285   printf("Monalisa TRDevents\t%d\n",trdEvents);
286   if (pcstream) (*pcstream)<<"calibStatAll"<<"TRDevents="<<trdEvents;
287   printf("Monalisa TRDtracks\t%d\n",trdTracks);
288   if (pcstream) (*pcstream)<<"calibStatAll"<<"TRDtracks="<<trdTracks;
289
290   //
291   // 3. T0 dump 
292   //
293   Int_t T0Events=0;
294   TList * T0Calib = (TList*)fin->Get("T0Calib");      
295   if (T0Calib) {
296     TH1  *histoEvents = (TH1*) T0Calib->FindObject("fTzeroORAplusORC");
297     if (histoEvents){
298       T0Events= TMath::Nint(histoEvents->GetEntries());
299     }
300     delete T0Calib;
301   }
302   printf("Monalisa T0events\t%d\n",T0Events);
303   if (pcstream) (*pcstream)<<"calibStatAll"<<"T0events="<<T0Events;
304
305   //
306   // 4. Mean vertex -   dump 
307   // Not present in CPass1
308   /*
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   //
323   // 5. SDD dump 
324   //
325   Int_t sddEvents=0;
326   Int_t sddTracks=0;
327   TList * SDDCalib = (TList*)fin->Get("clistSDDCalib");      
328   if (SDDCalib) {
329     TH1  *histoEvents = (TH1*) SDDCalib->FindObject("hNEvents");
330     if (histoEvents ){
331       sddEvents = TMath::Nint(histoEvents->GetBinContent(4));
332       sddTracks = TMath::Nint(histoEvents->GetBinContent(5));
333     }
334     delete SDDCalib;
335   }
336   printf("Monalisa SDDevents\t%d\n",sddEvents);
337   if (pcstream) (*pcstream)<<"calibStatAll"<<"SDDevents="<<sddEvents;
338   printf("Monalisa SDDtracks\t%d\n",sddTracks);
339   if (pcstream) (*pcstream)<<"calibStatAll"<<"SDDtracks="<<sddTracks;
340
341   //
342   if (pcstream) (*pcstream)<<"calibStatAll"<<"\n";
343   delete fin;
344
345 }
346
347