]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliTagAnalysis.cxx
AliAnalysisTaskMuonTrackingEff: removal of warnings, passing the AliCheckDetElementRe...
[u/mrichter/AliRoot.git] / ANALYSIS / AliTagAnalysis.cxx
1 /**************************************************************************
2  * Author: Panos Christakoglou.                                           *
3  * Contributors are mentioned in the code where appropriate.              *
4  *                                                                        *
5  * Permission to use, copy, modify and distribute this software and its   *
6  * documentation strictly for non-commercial purposes is hereby granted   *
7  * without fee, provided that the above copyright notice appears in all   *
8  * copies and that both the copyright notice and this permission notice   *
9  * appear in the supporting documentation. The authors make no claims     *
10  * about the suitability of this software for any purpose. It is          *
11  * provided "as is" without express or implied warranty.                  *
12  **************************************************************************/
13
14 /* $Id$ */
15
16 //-----------------------------------------------------------------
17 //           AliTagAnalysis class
18 //   This is the class to deal with the tag analysis
19 //   Origin: Panos Christakoglou, UOA-CERN, Panos.Christakoglou@cern.ch
20 //-----------------------------------------------------------------
21
22 //ROOT
23 #include <Riostream.h>
24 #include <TSystem.h>
25 #include <TChain.h>
26 #include <TFile.h>
27 #include <TEventList.h>
28 #include <TEntryList.h>
29 #include <TTreeFormula.h>
30 #include <TMap.h>
31
32 //ROOT-AliEn
33 #include <TGridResult.h>
34
35 #include "AliLog.h"
36
37 #include "AliRunTag.h"
38 #include "AliEventTag.h"
39 #include "AliTagAnalysis.h"
40 #include "AliEventTagCuts.h"
41 #include "AliDetectorTagCuts.h"
42 #include "AliLHCTagCuts.h"
43 #include "AliRunTagCuts.h"
44 #include "AliXMLCollection.h"
45
46 class TTree;
47
48 ClassImp(AliTagAnalysis)
49
50 //___________________________________________________________________________
51 AliTagAnalysis::AliTagAnalysis(): 
52   TObject(),
53   ftagresult(0x0),
54   fTagDirName(),
55   fChain(0x0),
56   fAnalysisType(),
57   fGlobalList(0) {
58   //Default constructor for a AliTagAnalysis
59 }
60
61 //___________________________________________________________________________
62 AliTagAnalysis::AliTagAnalysis(const char* type): 
63   TObject(),
64   ftagresult(0x0),
65   fTagDirName(),
66   fChain(0x0),
67   fAnalysisType(type),
68   fGlobalList(0) {
69   //constructor for a AliTagAnalysis
70 }
71
72 //___________________________________________________________________________
73 AliTagAnalysis::~AliTagAnalysis() {
74   //Default destructor for a AliTagAnalysis
75   if(ftagresult) delete ftagresult;
76   if(fChain) delete fChain;
77   if(fGlobalList) delete fGlobalList;
78 }
79
80 //___________________________________________________________________________
81 Bool_t  AliTagAnalysis::AddTagsFile(const char *alienUrl) {
82   // Add a single tags file to the chain
83
84   Bool_t rv = kTRUE ;
85
86   if (! fChain) fChain = new TChain("T");
87
88   TFile *f = TFile::Open(alienUrl,"READ");
89   fChain->Add(alienUrl);
90   AliInfo(Form("Chained tag files: %d ",fChain->GetEntries()));
91   delete f;
92
93   if (fChain->GetEntries() == 0 )
94     rv = kFALSE ;
95
96   return rv ;
97 }
98
99 //___________________________________________________________________________
100 void AliTagAnalysis::ChainLocalTags(const char *dirname) {
101   //Searches the entries of the provided direcory
102   //Chains the tags that are stored locally
103   fTagDirName = dirname;
104   TString fTagFilename;
105   
106   if (! fChain)  fChain = new TChain("T");
107   const char * tagPattern = 0x0;
108   if(fAnalysisType == "ESD") tagPattern = "ESD.tag.root";
109   else if(fAnalysisType == "AOD") tagPattern = "AOD.tag.root";
110   else AliFatal("Only ESD and AOD type is implemented!!!");
111
112   // Open the working directory
113   void * dirp = gSystem->OpenDirectory(fTagDirName);
114   const char * name = 0x0;
115   // Add all files matching *pattern* to the chain
116   while((name = gSystem->GetDirEntry(dirp))) {
117     if (strstr(name,tagPattern)) { 
118       fTagFilename = fTagDirName;
119       fTagFilename += "/";
120       fTagFilename += name;
121                 
122       fChain->Add(fTagFilename);  
123       printf("Tag file %s\n", fTagFilename.Data());
124       
125     }//pattern check
126   }//directory loop
127   AliInfo(Form("Chained tag files: %d ",fChain->GetEntries()));
128   fChain->ls();
129   
130 }
131
132
133 //___________________________________________________________________________
134 TChain * AliTagAnalysis::ChainGridTags(TGridResult *res) {
135   //Loops overs the entries of the TGridResult
136   //Chains the tags that are stored in the GRID
137   ftagresult = res;
138   Int_t nEntries = ftagresult->GetEntries();
139  
140   if (! fChain)  fChain = new TChain("T");
141
142   TString gridname = "alien://";
143   TString alienUrl;
144  
145   for(Int_t i = 0; i < nEntries; i++) {
146     alienUrl = ftagresult->GetKey(i,"turl");
147     fChain->Add(alienUrl);
148   }//grid result loop
149   return fChain;
150 }
151
152
153 //___________________________________________________________________________
154 TChain *AliTagAnalysis::QueryTags(AliRunTagCuts *runTagCuts, 
155                                   AliLHCTagCuts *lhcTagCuts, 
156                                   AliDetectorTagCuts *detTagCuts, 
157                                   AliEventTagCuts *evTagCuts) {
158   //Queries the tag chain using the defined 
159   //event tag cuts from the AliEventTagCuts object
160   //and returns a TChain along with the associated TEventList
161   AliInfo(Form("Querying the tags........"));
162
163   Bool_t aod = kFALSE;
164   TString aliceFile;
165   if(fAnalysisType == "ESD") aliceFile = "esdTree";
166   else if(fAnalysisType == "AOD") {
167       aliceFile = "aodTree";
168       aod = kTRUE;
169   }
170   else AliFatal("Only ESD and AOD type is implemented!!!");
171
172   //ESD file chain
173   TChain *esdChain = new TChain(aliceFile.Data());
174   //global entry list
175   fGlobalList = new TEntryList();
176   
177   //Defining tag objects
178   AliRunTag   *tag     = new AliRunTag;
179   AliEventTag *evTag   = new AliEventTag;
180   fChain->SetBranchAddress("AliTAG",&tag);
181
182   TString guid;
183   TString turl;
184   TString path;
185
186   TTree*      cTree     = 0; 
187   TEntryList* localList = new TEntryList();
188
189   Int_t iAccepted = 0;
190   Int_t iev       = 0;
191   Int_t ientry    = 0;
192   Int_t cEntries  = 0;
193   
194   for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
195     fChain->GetEntry(iTagFiles);
196     TTree* tree = fChain->GetTree();
197     if (cTree != tree) {
198         // Fix for aod tags: for each tree in the chain, merge the entries
199         cTree    = tree;
200         cEntries = tree->GetEntries();
201         iev      = 0;
202         ientry   = 0;
203     }
204
205     if(runTagCuts->IsAccepted(tag)) {
206       if(lhcTagCuts->IsAccepted(tag->GetLHCTag())) {
207         if(detTagCuts->IsAccepted(tag->GetDetectorTags())) {
208           if ((iev == 0) || !aod) localList->Reset();
209           Int_t iEvents = tag->GetNEvents();
210           const TClonesArray *tagList = tag->GetEventTags();
211           for(Int_t i = 0; i < iEvents; i++) {
212             evTag = (AliEventTag *) tagList->At(i);
213             guid = evTag->GetGUID(); 
214             turl = evTag->GetTURL(); 
215             path = evTag->GetPath();
216             localList->SetTreeName(aliceFile.Data());
217             if(turl!="") localList->SetFileName(turl.Data());
218             else localList->SetFileName(path.Data());
219
220             if(evTagCuts->IsAccepted(evTag)) {
221                 if(aod) localList->Enter(iev);
222                 else localList->Enter(i);
223             }
224             iev++;
225           }//event loop
226           if ((ientry == cEntries-1) || !aod) {
227               iAccepted += localList->GetN();
228               if(path != "") esdChain->AddFile(path);
229               else if(turl != "") esdChain->AddFile(turl);
230               fGlobalList->Add(localList);
231           }
232         }//detector tag cuts
233       }//lhc tag cuts
234     }//run tags cut
235     tag->Clear();
236     ientry++;
237   }//tag file loop
238   AliInfo(Form("Accepted events: %d",iAccepted));
239   esdChain->ls();
240   esdChain->SetEntryList(fGlobalList,"ne");
241   delete tag;
242   
243   return esdChain;
244 }
245
246 //___________________________________________________________________________
247 TChain *AliTagAnalysis::QueryTags(const char *fRunCut, 
248                                   const char *fLHCCut, 
249                                   const char *fDetectorCut, 
250                                   const char *fEventCut) {       
251   //Queries the tag chain using the defined      
252   //event tag cuts from the AliEventTagCuts object       
253   //and returns a TChain along with the associated TEventList    
254   AliInfo(Form("Querying the tags........"));    
255   
256   Bool_t aod = kFALSE;
257   TString aliceFile;
258   if(fAnalysisType == "ESD") aliceFile = "esdTree";
259   else if(fAnalysisType == "AOD") {
260       aliceFile = "aodTree";
261       aod = kTRUE;
262   }
263   else AliFatal("Only ESD and AOD type is implemented!!!");
264
265   //ESD file chain
266   TChain *esdChain = new TChain(aliceFile.Data());
267   //global entry list
268   fGlobalList = new TEntryList();
269   
270   //Defining tag objects         
271   AliRunTag   *tag   = new AliRunTag;    
272   AliEventTag *evTag = new AliEventTag;          
273   fChain->SetBranchAddress("AliTAG",&tag);       
274   
275   TString guid;          
276   TString turl;          
277   TString path;          
278   
279   TTreeFormula *fRunFormula = new TTreeFormula("fRun",fRunCut,fChain);   
280   TTreeFormula *fLHCFormula = new TTreeFormula("fLHC",fLHCCut,fChain);   
281   TTreeFormula *fDetectorFormula = new TTreeFormula("fDetector",fDetectorCut,fChain);
282   TTreeFormula *fEventFormula = new TTreeFormula("fEvent",fEventCut,fChain);
283   
284   TEntryList* localList = new TEntryList();
285
286   Int_t iev       = 0;
287   Int_t ientry    = 0;
288   Int_t cEntries  = 0;
289   Int_t current   = -1;          
290   Int_t iAccepted = 0;   
291
292   for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
293     fChain->GetEntry(iTagFiles);         
294     if (current != fChain->GetTreeNumber()) {    
295       fRunFormula->UpdateFormulaLeaves();        
296       fLHCFormula->UpdateFormulaLeaves();        
297       fDetectorFormula->UpdateFormulaLeaves();   
298       fEventFormula->UpdateFormulaLeaves();      
299       // Fix for aod tags: for each tree in the chain, merge the entries
300       cEntries = (fChain->GetTree())->GetEntries();
301       iev      = 0;
302       ientry   = 0;
303       //
304       current = fChain->GetTreeNumber();         
305     }    
306     if(fRunFormula->EvalInstance(iTagFiles) == 1) {      
307       if(fLHCFormula->EvalInstance(iTagFiles) == 1) {    
308         if(fDetectorFormula->EvalInstance(iTagFiles) == 1) {
309           if ((iev == 0) || !aod) localList->Reset();    
310           Int_t iEvents = fEventFormula->GetNdata();     
311           const TClonesArray *tagList = tag->GetEventTags();     
312           for(Int_t i = 0; i < iEvents; i++) {   
313             evTag = (AliEventTag *) tagList->At(i);      
314             guid = evTag->GetGUID();     
315             turl = evTag->GetTURL();     
316             path = evTag->GetPath();     
317             localList->SetTreeName(aliceFile.Data());
318             localList->SetFileName(turl.Data());
319             if(fEventFormula->EvalInstance(i) == 1) {
320                 if(aod) localList->Enter(iev);
321                 else localList->Enter(i);
322             }
323             iev++;
324           }//event loop          
325
326           if ((ientry == cEntries-1) || !aod) {  
327               if(path != "") esdChain->AddFile(path);    
328               else if(turl != "") esdChain->AddFile(turl);       
329               fGlobalList->Add(localList);
330               iAccepted += localList->GetN();
331           }
332         }//detector tag cuts
333       }//lhc tag cuts
334     }//run tag cut       
335     tag->Clear();
336     ientry++;
337   }//tag file loop       
338   AliInfo(Form("Accepted events: %d",iAccepted));        
339   esdChain->SetEntryList(fGlobalList,"ne");      
340
341   delete tag;
342   return esdChain;       
343 }
344
345 //___________________________________________________________________________
346 Bool_t AliTagAnalysis::CreateXMLCollection(const char* name, 
347                                            AliRunTagCuts *runTagCuts, 
348                                            AliLHCTagCuts *lhcTagCuts, 
349                                            AliDetectorTagCuts *detTagCuts, 
350                                            AliEventTagCuts *evTagCuts) {
351   //Queries the tag chain using the defined 
352   //event tag cuts from the AliEventTagCuts object
353   //and returns a XML collection
354   AliInfo(Form("Creating the collection........"));
355
356   Bool_t aod = kFALSE;
357   if(fAnalysisType == "AOD") aod = kTRUE;
358
359
360   AliXMLCollection *collection = new AliXMLCollection();
361   collection->SetCollectionName(name);
362   collection->WriteHeader();
363
364   TString guid;
365   TString turl;
366   TString lfn;
367   
368   TTree*      cTree = 0; 
369   TEntryList* localList = new TEntryList();
370   Int_t iAccepted = 0;
371   Int_t iev       = 0;
372   Int_t ientry    = 0;
373   Int_t cEntries  = 0;
374
375   Int_t iRejectedRun = 0;
376   Int_t iRejectedLHC = 0;
377   Int_t iRejectedDet = 0;
378   Int_t iRejectedEvt = 0;
379
380   Int_t iTotalEvents = 0;
381
382   Int_t iAcceptedEvtInFile = 0;
383   Int_t iRejectedEvtInFile = 0;
384
385   //Defining tag objects
386   AliRunTag   *tag   = new AliRunTag;
387   AliEventTag *evTag = new AliEventTag;
388   fChain->SetBranchAddress("AliTAG",&tag);
389
390   for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
391
392     fChain->GetEntry(iTagFiles);
393     TTree* tree = fChain->GetTree();
394     if (cTree != tree) {
395         // Fix for aod tags: for each tree in the chain, merge the entries
396         cTree    = tree;
397         cEntries = tree->GetEntries();
398         iev      = 0;
399         ientry   = 0;
400     }
401     //Event list
402     iTotalEvents += tag->GetNEvents();
403     if ((iev == 0) || !aod) localList->Reset();
404     if(runTagCuts->IsAccepted(tag)) {
405       if(lhcTagCuts->IsAccepted(tag->GetLHCTag())) {
406         if(detTagCuts->IsAccepted(tag->GetDetectorTags())) {
407           Int_t iEvents = tag->GetNEvents();
408           const TClonesArray *tagList = tag->GetEventTags();
409           iRejectedEvtInFile = 0;
410           iAcceptedEvtInFile = 0;
411           for(Int_t i = 0; i < iEvents; i++) {
412             evTag = (AliEventTag *) tagList->At(i);
413             guid = evTag->GetGUID(); 
414             turl = evTag->GetTURL(); 
415             lfn = turl(8,turl.Length());
416             if(evTagCuts->IsAccepted(evTag)) {
417                 if(aod) localList->Enter(iev);
418                 else localList->Enter(i);
419                 iAcceptedEvtInFile++;
420             }
421             else {
422               iRejectedEvt++;
423               iRejectedEvtInFile++;
424             }
425             iev++;
426           }//event loop
427           if ((ientry == cEntries-1) || !aod) {
428             iAccepted += localList->GetN();
429             collection->WriteBody(iTagFiles+1,guid,lfn,turl,localList,iAcceptedEvtInFile,iRejectedEvtInFile);
430           }
431         }//detector tag cuts
432         else {
433           iRejectedDet += tag->GetNEvents();
434         }
435       }//lhc tag cuts 
436       else {
437         iRejectedLHC += tag->GetNEvents();
438       }
439     }//run tag cuts
440     else {
441       iRejectedRun += tag->GetNEvents();
442     }
443     tag->Clear();
444     ientry++;
445   }//tag file loop
446   collection->WriteSummary(iTotalEvents, iAccepted, iRejectedRun, iRejectedLHC, iRejectedDet, iRejectedEvt);
447   collection->Export();
448   
449   delete tag;
450   return kTRUE;
451 }
452
453 //___________________________________________________________________________
454 Bool_t AliTagAnalysis::CreateXMLCollection(const char* name, 
455                                            const char *fRunCut, 
456                                            const char *fLHCCut, 
457                                            const char *fDetectorCut, 
458                                            const char *fEventCut) {
459   //Queries the tag chain using the defined 
460   //event tag cuts from the AliEventTagCuts object
461   //and returns a XML collection
462   AliInfo(Form("Creating the collection........"));
463
464   Bool_t aod = kFALSE;
465   if(fAnalysisType == "AOD") aod = kTRUE;
466
467   AliXMLCollection *collection = new AliXMLCollection();
468   collection->SetCollectionName(name);
469   collection->WriteHeader();
470
471   TString guid;
472   TString turl;
473   TString lfn;
474   TEntryList* localList = new TEntryList();
475   
476   Int_t iAccepted = 0;
477   Int_t iev       = 0;
478   Int_t ientry    = 0;
479   Int_t cEntries  = 0;
480
481   Int_t iRejectedRun = 0;
482   Int_t iRejectedLHC = 0;
483   Int_t iRejectedDet = 0;
484   Int_t iRejectedEvt = 0;
485
486   Int_t iTotalEvents = 0;
487
488   Int_t iAcceptedEvtInFile = 0;
489   Int_t iRejectedEvtInFile = 0;
490
491   //Defining tag objects
492   AliRunTag *tag     = new AliRunTag;
493   AliEventTag *evTag = new AliEventTag;
494   fChain->SetBranchAddress("AliTAG",&tag);
495
496   TTreeFormula *fRunFormula = new TTreeFormula("fRun",fRunCut,fChain);
497   TTreeFormula *fLHCFormula = new TTreeFormula("fLHC",fLHCCut,fChain);   
498   TTreeFormula *fDetectorFormula = new TTreeFormula("fDetector",fDetectorCut,fChain);
499   TTreeFormula *fEventFormula = new TTreeFormula("fEvent",fEventCut,fChain);
500
501   Int_t current = -1;
502   for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
503
504     fChain->GetEntry(iTagFiles);
505     if (current != fChain->GetTreeNumber()) {
506       fRunFormula->UpdateFormulaLeaves();
507       fLHCFormula->UpdateFormulaLeaves();
508       fDetectorFormula->UpdateFormulaLeaves();
509       fEventFormula->UpdateFormulaLeaves();
510       // Fix for aod tags: for each tree in the chain, merge the entries
511       cEntries = (fChain->GetTree())->GetEntries();
512       iev      = 0;
513       ientry   = 0;
514       //
515       current = fChain->GetTreeNumber();
516     }
517     //Event list
518     iTotalEvents += tag->GetNEvents();
519     if ((iev == 0) || !aod) localList->Reset();
520     if(fRunFormula->EvalInstance(iTagFiles) == 1) {
521       if(fLHCFormula->EvalInstance(iTagFiles) == 1) {    
522         if(fDetectorFormula->EvalInstance(iTagFiles) == 1) {     
523           Int_t iEvents = fEventFormula->GetNdata();
524           const TClonesArray *tagList = tag->GetEventTags();
525           iRejectedEvtInFile = 0;
526           iAcceptedEvtInFile = 0;
527           for(Int_t i = 0; i < iEvents; i++) {
528             evTag = (AliEventTag *) tagList->At(i);
529             guid = evTag->GetGUID(); 
530             turl = evTag->GetTURL(); 
531             lfn = turl(8,turl.Length());
532             if(fEventFormula->EvalInstance(i) == 1) {
533                 if(aod) localList->Enter(iev);
534                 else localList->Enter(i);
535                 iAcceptedEvtInFile++;
536             }
537             else {
538               iRejectedEvt++;
539               iRejectedEvtInFile++;
540             }
541             iev++;
542           }//event loop
543           if ((ientry == cEntries-1) || !aod) {
544             collection->WriteBody(iTagFiles+1,guid,lfn,turl,localList,iAcceptedEvtInFile, iRejectedEvtInFile);
545             iAccepted += localList->GetN();
546           }
547         }//detector tag cuts
548         else {
549           iRejectedDet += tag->GetNEvents();
550         }
551       }//lhc tag cuts 
552       else {
553         iRejectedLHC += tag->GetNEvents();
554       }
555     }//run tag cuts
556     else {
557       iRejectedRun += tag->GetNEvents();
558     }
559     ientry++;
560   }//tag file loop
561   collection->WriteSummary(iTotalEvents, iAccepted, iRejectedRun, iRejectedLHC, iRejectedDet, iRejectedEvt);
562   collection->Export();
563
564   delete tag;
565   return kTRUE;
566 }
567
568 //___________________________________________________________________________
569 TChain *AliTagAnalysis::GetInputChain(const char* system, const char *wn) {
570   //returns the chain+event list - used in batch sessions
571   // this function will be removed once the new root 
572   // improvements are committed
573   TString fsystem = system;
574   Int_t iAccepted = 0;
575
576   TChain *fAnalysisChain = 0;
577   if(fAnalysisType == "ESD") fAnalysisChain = new TChain("esdTree");
578   else if(fAnalysisType == "AOD") fAnalysisChain = new TChain("aodTree");
579   else AliFatal("Only ESD and AOD type is implemented!!!");
580   
581   //Event list
582   TEventList *fEventList = new TEventList();
583   AliXMLCollection *collection = AliXMLCollection::Open(wn);
584
585   collection->Reset();
586   while (collection->Next()) {
587     AliInfo(Form("Adding: %s",collection->GetTURL("")));
588     fAnalysisChain->Add(collection->GetTURL(""));
589     TEntryList *list = (TEntryList *)collection->GetEventList("");
590     for(Int_t i = 0; i < list->GetN(); i++) fEventList->Enter(iAccepted+list->GetEntry(i));
591
592     if(fsystem == "pp") iAccepted += 100;
593     else if(fsystem == "PbPb") iAccepted += 1;
594   }
595
596   fAnalysisChain->SetEventList(fEventList);
597   
598   AliInfo(Form("Number of selected events: %d",fEventList->GetN()));
599
600   return fAnalysisChain;
601 }
602
603 //___________________________________________________________________________
604 TChain *AliTagAnalysis::GetChainFromCollection(const char* collectionname, 
605                                                const char* treename) {
606   //returns the TChain+TEntryList object- used in batch sessions
607   TString aliceFile = treename;
608   Int_t iAccepted = 0;
609   TChain *fAnalysisChain = 0;
610   if(aliceFile == "esdTree") fAnalysisChain = new TChain("esdTree");
611   else if(aliceFile == "aodTree") fAnalysisChain = new TChain("aodTree");
612   else AliFatal("Inconsistent tree name - use esdTree or aodTree!");
613
614   //Event list
615   fGlobalList = new TEntryList();
616   AliXMLCollection *collection = AliXMLCollection::Open(collectionname);
617
618   // Tag selection summary per file
619   TMap *tagCutSummary = new TMap();
620   tagCutSummary->SetName("TagCutSumm");
621
622   collection->Reset();
623   while (collection->Next()) {
624     AliInfo(Form("Adding: %s",collection->GetTURL("")));
625     fAnalysisChain->Add(collection->GetTURL(""));
626     TEntryList *list = (TEntryList *)collection->GetEventList("");
627     list->SetTreeName(aliceFile.Data());
628     list->SetFileName(collection->GetTURL(""));
629     fGlobalList->Add(list);
630     iAccepted += list->GetN();
631     if (collection->GetCutSumm())
632       tagCutSummary->Add(new TObjString(collection->GetTURL("")), new TObjString(collection->GetCutSumm()));
633   }
634
635   fAnalysisChain->SetEntryList(fGlobalList,"ne");
636   
637   AliInfo(Form("Number of selected events: %d",iAccepted));
638
639   TList *aUserInfo = fAnalysisChain->GetUserInfo();
640   aUserInfo->Add(tagCutSummary);
641
642   Int_t iAccEv;
643   Int_t iTotalEvents;
644   Int_t iRejRun;
645   Int_t iRejLHC;
646   Int_t iRejDet;
647   Int_t iRejEvt;
648
649   collection->GetCollectionSummary(&iTotalEvents, &iAccEv, &iRejRun, &iRejLHC, &iRejDet, &iRejEvt);
650  
651   char nstr[2000];
652
653   sprintf(nstr, "TotalEvents=%i", iTotalEvents);
654   TObjString *iTotStr = new TObjString(nstr);
655   aUserInfo->Add(iTotStr);
656
657   sprintf(nstr, "AcceptedEvents=%i", iAccepted);
658   TObjString *iAccStr = new TObjString(nstr);
659   aUserInfo->Add(iAccStr);
660
661   sprintf(nstr, "RejectedRun=%i", iRejRun);
662   TObjString *iRejRunStr = new TObjString(nstr);
663   aUserInfo->Add(iRejRunStr);
664
665   sprintf(nstr, "RejectedLHC=%i", iRejLHC);
666   TObjString *iRejLHCStr = new TObjString(nstr);
667   aUserInfo->Add(iRejLHCStr);
668
669   sprintf(nstr, "RejectedDet=%i", iRejDet);
670   TObjString *iRejDetStr = new TObjString(nstr);
671   aUserInfo->Add(iRejDetStr);
672
673   sprintf(nstr, "RejectedEvt=%i", iRejEvt);
674   TObjString *iRejEvtStr = new TObjString(nstr);
675   aUserInfo->Add(iRejEvtStr);
676
677   return fAnalysisChain;
678 }