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