]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/macros/simpleTrending.C
5272f84103a9aa645d70c273ace35decfceb4f1a
[u/mrichter/AliRoot.git] / PWGPP / macros / simpleTrending.C
1 void processContainer(TObject* object, TTreeSRedirector* debugStreamer, TString name);
2 void loadLibraries();
3 void simpleTrending(TString inputFileName, Int_t run, TString filterExpr=".*", TString trendingFileName="trending.root", TString treeName="trending", TString fileOpenMode="update" );
4
5 TString treeName;
6
7 void simpleTrending(TString inputFileName, Int_t run, TString filterExpr, TString trendingFileName, TString debugTreeName, TString fileOpenMode )
8 {
9
10   // Dump the statistical information about all histograms in the file
11   //    into a tree
12   //
13   //
14   // 1. Default dump for all histograms
15   //    Information to dump:
16   //    stat =Entries, Mean, MeanError,  RMS, MaxBin
17   //    Branch naming convention:
18   //    <detName>_<hisName><statName>
19   //
20   // 2. Detector statistical information  - to be implemented by expert
21   //                                      - First version implemented by MI
22   //                                      - updated for QA by MK
23   //  
24   // 
25
26   loadLibraries();
27
28   treeName=debugTreeName;
29
30   TFile *inputFile = TFile::Open(inputFileName.Data());
31   if (!inputFile) return;
32   
33   TList * keyList = inputFile->GetListOfKeys();
34   Int_t nkeys=keyList->GetEntries();
35  
36   TRegexp filterRegexp=filterExpr.Data();
37
38   //check if we have a matching container, only then create the output file
39   TList* keyList=inputFile->GetListOfKeys();
40   Int_t nkeys=keyList->GetEntries();
41   Bool_t containerExists=kFALSE;
42   for (Int_t i=0; i<nkeys; i++)
43   {
44     TObject* object=keyList->At(i);
45     if (!object) continue;
46     TString name=object->GetName();
47     if (name.Contains(filterRegexp))
48     {
49       containerExists=kTRUE;
50       break;
51     }
52   }
53   if (!containerExists) 
54   {
55     printf("container %s does not exist in %s\n",filterExpr.Data(),inputFileName.Data());
56     return;
57   }
58
59   TTreeSRedirector *pcstream = new TTreeSRedirector(trendingFileName,fileOpenMode);
60   (*pcstream)<<treeName.Data()<<"run="<<run;
61   
62   //main loop over the top level objects, filtering is done here
63   for (Int_t i=0; i<nkeys; i++)
64   {
65     TObject * object = inputFile->Get(keyList->At(i)->GetName());
66     if (!object) continue;
67     TString name=object->GetName();
68     if (!name.Contains(filterRegexp)) continue;
69     processContainer(object,pcstream,name);
70   }
71   
72   //
73   // Expert dump example (MI first iteration):
74   //
75   // 0.)  TOF dump
76   //
77
78   Int_t tofEvents=0;
79   Int_t tofTracks=0;
80   TList * TOFCalib = (TList*)inputFile->Get("TOFHistos");      
81   if (TOFCalib) {
82     TH1 *histoEvents = (TH1*)TOFCalib->FindObject("hHistoVertexTimestamp");
83     TH1 *histoTracks = (TH1*)TOFCalib->FindObject("hHistoDeltatTimestamp");
84     if (histoEvents && histoTracks){
85       tofEvents = TMath::Nint(histoEvents->GetEntries());
86       tofTracks = TMath::Nint(histoTracks->GetEntries());
87     }
88     delete TOFCalib;
89   }
90   printf("Monalisa TOFevents\t%d\n",tofEvents);
91   if (pcstream) (*pcstream)<<treeName.Data()<<"TOFevents="<<tofEvents;
92   printf("Monalisa TOFtracks\t%d\n",tofTracks);
93   if (pcstream) (*pcstream)<<treeName.Data()<<"TOFtracks="<<tofTracks;
94
95   //
96   // 1.)  TPC  dump - usefull events/tracks  for the calibration
97   //
98   Int_t tpcEvents=0;
99   Int_t tpcTracks=0;
100   TObject* obj = dynamic_cast<TObject*>(inputFile->Get("TPCCalib"));
101   TObjArray* array = dynamic_cast<TObjArray*>(obj);
102   TDirectory* dir = dynamic_cast<TDirectory*>(obj);
103   AliTPCcalibTime  * calibTime = NULL;
104   if (dir) {
105     calibTime = dynamic_cast<AliTPCcalibTime*>(dir->Get("calibTime"));
106   }
107   else if (array){
108     calibTime = (AliTPCcalibTime *)array->FindObject("calibTime");
109   }
110   if (calibTime) {
111       tpcEvents = TMath::Nint(calibTime->GetTPCVertexHisto(0)->GetEntries());
112       tpcTracks = TMath::Nint(calibTime->GetResHistoTPCITS(0)->GetEntries());
113   }
114   printf("Monalisa TPCevents\t%d\n",tpcEvents);
115   if (pcstream) (*pcstream)<<treeName.Data()<<"TPCevents="<<tpcEvents;
116   printf("Monalisa TPCtracks\t%d\n",tpcTracks);
117   if (pcstream) (*pcstream)<<treeName.Data()<<"TPCtracks="<<tpcTracks;
118
119   //
120   // 2. TRD dump 
121   //
122   Int_t trdEvents=0;
123   Int_t trdTracks=0;
124   TList * TRDCalib = (TList*)inputFile->Get("TRDCalib");      
125   if (TRDCalib) {
126     TH1  *histoEvents = (TH1*)TRDCalib->FindObject("NEventsInput_AliTRDCalibTask");
127     TH1  *histoTracks = (TH1*)TRDCalib->FindObject("AbsoluteGain_AliTRDCalibTask");
128     if (histoEvents && histoTracks){
129       trdEvents= TMath::Nint(histoEvents->GetEntries());
130       trdTracks= TMath::Nint(histoTracks->GetEntries());
131     }
132     delete TRDCalib;
133   }
134   printf("Monalisa TRDevents\t%d\n",trdEvents);
135   if (pcstream) (*pcstream)<<treeName.Data()<<"TRDevents="<<trdEvents;
136   printf("Monalisa TRDtracks\t%d\n",trdTracks);
137   if (pcstream) (*pcstream)<<treeName.Data()<<"TRDtracks="<<trdTracks;
138
139   //
140   // 3. T0 dump 
141   //
142   Int_t T0Events=0;
143   TList * T0Calib = (TList*)inputFile->Get("T0Calib");      
144   if (T0Calib) {
145     TH1  *histoEvents = (TH1*) T0Calib->FindObject("fTzeroORAplusORC");
146     if (histoEvents){
147       T0Events= TMath::Nint(histoEvents->GetEntries());
148     }
149     delete T0Calib;
150   }
151   printf("Monalisa T0events\t%d\n",T0Events);
152   if (pcstream) (*pcstream)<<treeName.Data()<<"T0events="<<T0Events;
153
154   //
155   // 4. Mean vertex -   dump 
156   // Not present in CPass1
157   /*
158     Int_t meanVertexEvents=0;
159   TList * meanVertexCalib = (TList*)inputFile->Get("MeanVertex");      
160   if (meanVertexCalib) {
161     TH1  *histoEvents = (TH1*) meanVertexCalib->FindObject("hTRKVertexX");
162     if (histoEvents){
163       meanVertexEvents = TMath::Nint(histoEvents->GetEntries());
164     }
165     delete meanVertexCalib;
166   }
167   printf("Monalisa MeanVertexevents\t%d\n",meanVertexEvents);
168   if (pcstream) (*pcstream)<<treeName.Data()<<"MeanVertexevents="<<meanVertexEvents;
169   */
170
171   //
172   // 5. SDD dump 
173   //
174   Int_t sddEvents=0;
175   Int_t sddTracks=0;
176   TList * SDDCalib = (TList*)inputFile->Get("clistSDDCalib");      
177   if (SDDCalib) {
178     TH1  *histoEvents = (TH1*) SDDCalib->FindObject("hNEvents");
179     if (histoEvents ){
180       sddEvents = TMath::Nint(histoEvents->GetBinContent(4));
181       sddTracks = TMath::Nint(histoEvents->GetBinContent(5));
182     }
183     delete SDDCalib;
184   }
185   printf("Monalisa SDDevents\t%d\n",sddEvents);
186   if (pcstream) (*pcstream)<<treeName.Data()<<"SDDevents="<<sddEvents;
187   printf("Monalisa SDDtracks\t%d\n",sddTracks);
188   if (pcstream) (*pcstream)<<treeName.Data()<<"SDDtracks="<<sddTracks;
189
190   //
191   if (pcstream) (*pcstream)<<treeName.Data()<<"\n";
192   delete pcstream;
193   delete inputFile;
194
195 }
196
197 void processContainer(TObject* inputObject, TTreeSRedirector* pcstream, TString parentname)
198 {
199   //recursively process the contents of an object:
200   //might be a TDirectory
201   //or a TCollection
202   //pus information about the contained histograms in the debugStreamer
203
204   TDirectory* inputDir=NULL;
205   TSeqCollection* inputCollection=NULL;
206   TH1* inputHistogram=NULL;
207   
208   TString inputObjectName=inputObject->GetName();
209
210   //TDirectory* inputDir=dynamic_cast<TDirectory*>(inputObject);
211   //TSeqCollection* inputCollection=dynamic_cast<TSeqCollection*>(inputObject);
212   //TH1* inputHistogram=dynamic_cast<TH1*>(inputObject);
213   if (inputObject->InheritsFrom("TDirectory")) inputDir=dynamic_cast<TDirectory*>(inputObject);
214   if (inputObject->InheritsFrom("TSeqCollection")) inputCollection=dynamic_cast<TSeqCollection*>(inputObject);
215   if (inputObject->InheritsFrom("TH1")) inputHistogram=dynamic_cast<TH1*>(inputObject);
216
217   if (inputCollection)
218   {
219     printf("processing collection: %s\n",inputCollection->GetName());
220     Int_t nentries= inputCollection->GetEntries();
221     for (Int_t i=0; i<nentries; i++)
222     {
223       TObject * object = inputCollection->At(i);
224       if (!object) continue;
225       TString name=parentname+"/"+object->GetName();
226       processContainer(object,pcstream,name);
227     }
228   } 
229   else if (inputDir)
230   {
231     printf("processing directory: %s\n",inputDir->GetName());
232     TList* keyList=inputDir->GetListOfKeys();
233     Int_t nkeys=keyList->GetEntries();
234     for (Int_t i=0; i<nkeys; i++)
235     {
236       TObject * object = inputDir->Get(keyList->At(i)->GetName());
237       if (!object) continue;
238       TString name=parentname+"/"+object->GetName();
239       processContainer(object,pcstream,name);
240     }
241
242   }
243   else if (inputHistogram)
244   {
245     Double_t hisEntries;
246     Double_t hisMean;
247     Double_t hisMeanError;
248     Double_t hisRMS;
249     Double_t hisMaxBin;
250     if (inputHistogram->InheritsFrom("TH1")==0) continue;
251     TH1* phis = (TH1*)inputHistogram;
252     hisEntries=phis->GetEntries();      
253     Int_t idim=1;
254     if (inputHistogram->InheritsFrom("TH2")) idim=2;
255     if (inputHistogram->InheritsFrom("TH3")) idim=3;        
256     hisMean=phis->GetMean(idim);        
257     hisMeanError=phis->GetMeanError(idim);      
258     hisRMS=phis->GetRMS(idim);  
259     hisMaxBin=phis->GetBinCenter(phis->GetMaximumBin());
260     TString name=parentname;
261     printf("histogram: %s\n",name.Data());
262     if (pcstream) (*pcstream)<<treeName.Data()<<
263         Form("%s_Entries=",name.Data())<<hisEntries<<   
264         Form("%s_Mean=",name.Data())<<hisMean<< 
265         Form("%s_MeanError=",name.Data())<<hisMeanError<<       
266         Form("%s_RMS=",name.Data())<<hisRMS<<   
267         Form("%s_MaxBin=",name)<<hisMaxBin;     
268   }
269
270 }
271
272 void loadLibraries()
273 {
274   gSystem->SetIncludePath("-I. -I$ROOTSYS/include -I$ALICE_ROOT/include -I$ALICE_ROOT -I$ALICE_ROOT/ITS -I$ALICE_ROOT/TRD -I$ALICE_ROOT/PWGPP -     I$ALICE_ROOT/PWGPP/TRD");
275   gSystem->Load("libANALYSIS");
276   gSystem->Load("libANALYSISalice");
277   gSystem->Load("libANALYSIScalib");
278   gSystem->Load("libESDfilter");
279   gSystem->Load("libCORRFW");
280   gSystem->Load("libTender");
281   gSystem->Load("libPWGPP.so");
282   gSystem->Load("libAliHLTTrigger.so");
283   gSystem->Load("libPWGTools");
284   gSystem->Load("libEMCALUtils");
285   gSystem->Load("libPHOSUtils");
286   gSystem->Load("libPWGCaloTrackCorrBase");
287   gSystem->Load("libPWGGACaloTrackCorrelations");
288   gSystem->Load("libPWGGACaloTasks");
289   gSystem->Load("libPWGGAPHOSTasks");
290   gSystem->Load("libPWGTools");
291   gSystem->Load("libPWGEMCAL");
292   gSystem->Load("libPWGGAEMCALTasks");
293   gSystem->Load("libPWGmuon");
294   gSystem->Load("libPWGPPMUONlite");
295   gSystem->Load("libPWGmuondep");
296   gSystem->Load("libPWGLFforward2");
297 }
298