]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliTagAnalysis.cxx
Leak corrected (Rossella Romita)
[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  // AliDebug(Form("Chained tag files: %d ",fChain->GetEntries()));
129   fChain->ls();
130   
131 }
132
133
134 //___________________________________________________________________________
135 TChain * AliTagAnalysis::ChainGridTags(TGridResult *res) {
136   //Loops overs the entries of the TGridResult
137   //Chains the tags that are stored in the GRID
138   ftagresult = res;
139   Int_t nEntries = ftagresult->GetEntries();
140  
141   if (! fChain)  fChain = new TChain("T");
142
143   TString gridname = "alien://";
144   TString alienUrl;
145  
146   for(Int_t i = 0; i < nEntries; i++) {
147     alienUrl = ftagresult->GetKey(i,"turl");
148     fChain->Add(alienUrl);
149   }//grid result loop
150   return fChain;
151 }
152
153
154 //___________________________________________________________________________
155 TChain *AliTagAnalysis::QueryTags(AliRunTagCuts *runTagCuts, 
156                                   AliLHCTagCuts *lhcTagCuts, 
157                                   AliDetectorTagCuts *detTagCuts, 
158                                   AliEventTagCuts *evTagCuts) {
159   //Queries the tag chain using the defined 
160   //event tag cuts from the AliEventTagCuts object
161   //and returns a TChain along with the associated TEventList
162   AliInfo(Form("Querying the tags........"));
163
164   Bool_t aod = kFALSE;
165   TString aliceFile;
166   if(fAnalysisType == "ESD") aliceFile = "esdTree";
167   else if(fAnalysisType == "AOD") {
168       aliceFile = "aodTree";
169       aod = kTRUE;
170   }
171   else AliFatal("Only ESD and AOD type is implemented!!!");
172
173   //ESD file chain
174   TChain *esdChain = new TChain(aliceFile.Data());
175   //global entry list
176   fGlobalList = new TEntryList();
177   
178   //Defining tag objects
179   AliRunTag   *tag     = new AliRunTag;
180   AliEventTag *evTag   = new AliEventTag;
181   fChain->SetBranchAddress("AliTAG",&tag);
182
183   TString guid;
184   TString turl;
185   TString path;
186
187   TTree*      cTree     = 0; 
188   TEntryList* localList = new TEntryList();
189
190   Int_t iAccepted = 0;
191   Int_t iev       = 0;
192   Int_t ientry    = 0;
193   Int_t cEntries  = 0;
194   
195   for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
196     fChain->GetEntry(iTagFiles);
197     TTree* tree = fChain->GetTree();
198     if (cTree != tree) {
199         // Fix for aod tags: for each tree in the chain, merge the entries
200         cTree    = tree;
201         cEntries = tree->GetEntries();
202         iev      = 0;
203         ientry   = 0;
204     }
205
206     if(runTagCuts->IsAccepted(tag)) {
207       if(lhcTagCuts->IsAccepted(tag->GetLHCTag())) {
208         if(detTagCuts->IsAccepted(tag->GetDetectorTags())) {
209           if ((iev == 0) || !aod) localList->Reset();
210           Int_t iEvents = tag->GetNEvents();
211           const TClonesArray *tagList = tag->GetEventTags();
212           for(Int_t i = 0; i < iEvents; i++) {
213             evTag = (AliEventTag *) tagList->At(i);
214             guid = evTag->GetGUID(); 
215             turl = evTag->GetTURL(); 
216             path = evTag->GetPath();
217             localList->SetTreeName(aliceFile.Data());
218             if(turl!="") localList->SetFileName(turl.Data());
219             else localList->SetFileName(path.Data());
220
221             if(evTagCuts->IsAccepted(evTag)) {
222                 if(aod) localList->Enter(iev);
223                 else localList->Enter(i);
224             }
225             iev++;
226           }//event loop
227           if ((ientry == cEntries-1) || !aod) {
228               iAccepted += localList->GetN();
229               if(path != "") esdChain->AddFile(path);
230               else if(turl != "") esdChain->AddFile(turl);
231               fGlobalList->Add(localList);
232           }
233         }//detector tag cuts
234       }//lhc tag cuts
235     }//run tags cut
236     tag->Clear();
237     ientry++;
238   }//tag file loop
239   AliInfo(Form("Accepted events: %d",iAccepted));
240   esdChain->ls();
241   esdChain->SetEntryList(fGlobalList,"ne");
242   delete tag;
243   delete localList;
244   
245   return esdChain;
246 }
247
248 //___________________________________________________________________________
249 TChain *AliTagAnalysis::QueryTags(const char *fRunCut, 
250                                   const char *fLHCCut, 
251                                   const char *fDetectorCut, 
252                                   const char *fEventCut) {       
253   //Queries the tag chain using the defined      
254   //event tag cuts from the AliEventTagCuts object       
255   //and returns a TChain along with the associated TEventList    
256   AliInfo(Form("Querying the tags........"));    
257   
258   Bool_t aod = kFALSE;
259   TString aliceFile;
260   if(fAnalysisType == "ESD") aliceFile = "esdTree";
261   else if(fAnalysisType == "AOD") {
262       aliceFile = "aodTree";
263       aod = kTRUE;
264   }
265   else AliFatal("Only ESD and AOD type is implemented!!!");
266
267   //ESD file chain
268   TChain *esdChain = new TChain(aliceFile.Data());
269   //global entry list
270   fGlobalList = new TEntryList();
271   
272   //Defining tag objects         
273   AliRunTag   *tag   = new AliRunTag;    
274   AliEventTag *evTag = new AliEventTag;          
275   fChain->SetBranchAddress("AliTAG",&tag);       
276   
277   TString guid;          
278   TString turl;          
279   TString path;          
280   
281   TTreeFormula *fRunFormula = new TTreeFormula("fRun",fRunCut,fChain);   
282   TTreeFormula *fLHCFormula = new TTreeFormula("fLHC",fLHCCut,fChain);   
283   TTreeFormula *fDetectorFormula = new TTreeFormula("fDetector",fDetectorCut,fChain);
284   TTreeFormula *fEventFormula = new TTreeFormula("fEvent",fEventCut,fChain);
285   
286   TEntryList* localList = new TEntryList();
287
288   Int_t iev       = 0;
289   Int_t ientry    = 0;
290   Int_t cEntries  = 0;
291   Int_t current   = -1;          
292   Int_t iAccepted = 0;   
293
294   for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
295     fChain->GetEntry(iTagFiles);         
296     if (current != fChain->GetTreeNumber()) {    
297       fRunFormula->UpdateFormulaLeaves();        
298       fLHCFormula->UpdateFormulaLeaves();        
299       fDetectorFormula->UpdateFormulaLeaves();   
300       fEventFormula->UpdateFormulaLeaves();      
301       // Fix for aod tags: for each tree in the chain, merge the entries
302       cEntries = (fChain->GetTree())->GetEntries();
303       iev      = 0;
304       ientry   = 0;
305       //
306       current = fChain->GetTreeNumber();         
307     }    
308     if(fRunFormula->EvalInstance(iTagFiles) == 1) {      
309       if(fLHCFormula->EvalInstance(iTagFiles) == 1) {    
310         if(fDetectorFormula->EvalInstance(iTagFiles) == 1) {
311           if ((iev == 0) || !aod) localList->Reset();    
312           Int_t iEvents = fEventFormula->GetNdata();     
313           const TClonesArray *tagList = tag->GetEventTags();     
314           for(Int_t i = 0; i < iEvents; i++) {   
315             evTag = (AliEventTag *) tagList->At(i);      
316             guid = evTag->GetGUID();     
317             turl = evTag->GetTURL();     
318             path = evTag->GetPath();     
319             localList->SetTreeName(aliceFile.Data());
320             localList->SetFileName(turl.Data());
321             if(fEventFormula->EvalInstance(i) == 1) {
322                 if(aod) localList->Enter(iev);
323                 else localList->Enter(i);
324             }
325             iev++;
326           }//event loop          
327
328           if ((ientry == cEntries-1) || !aod) {  
329               if(path != "") esdChain->AddFile(path);    
330               else if(turl != "") esdChain->AddFile(turl);       
331               fGlobalList->Add(localList);
332               iAccepted += localList->GetN();
333           }
334         }//detector tag cuts
335       }//lhc tag cuts
336     }//run tag cut       
337     tag->Clear();
338     ientry++;
339   }//tag file loop       
340   AliInfo(Form("Accepted events: %d",iAccepted));        
341   esdChain->SetEntryList(fGlobalList,"ne");      
342
343   delete tag;
344   delete localList;
345   return esdChain;       
346 }
347
348 //___________________________________________________________________________
349 Bool_t AliTagAnalysis::CreateXMLCollection(const char* name, 
350                                            AliRunTagCuts *runTagCuts, 
351                                            AliLHCTagCuts *lhcTagCuts, 
352                                            AliDetectorTagCuts *detTagCuts, 
353                                            AliEventTagCuts *evTagCuts) {
354   //Queries the tag chain using the defined 
355   //event tag cuts from the AliEventTagCuts object
356   //and returns a XML collection
357   AliInfo(Form("Creating the collection........"));
358
359   Bool_t aod = kFALSE;
360   if(fAnalysisType == "AOD") aod = kTRUE;
361
362
363   AliXMLCollection *collection = new AliXMLCollection();
364   collection->SetCollectionName(name);
365   collection->WriteHeader();
366
367   TString guid;
368   TString turl;
369   TString lfn;
370   
371   TTree*      cTree = 0; 
372   TEntryList* localList = new TEntryList();
373   Int_t iAccepted = 0;
374   Int_t iev       = 0;
375   Int_t ientry    = 0;
376   Int_t cEntries  = 0;
377
378   Int_t iRejectedRun = 0;
379   Int_t iRejectedLHC = 0;
380   Int_t iRejectedDet = 0;
381   Int_t iRejectedEvt = 0;
382
383   Int_t iTotalEvents = 0;
384
385   Int_t iAcceptedEvtInFile = 0;
386   Int_t iRejectedEvtInFile = 0;
387
388   //Defining tag objects
389   AliRunTag   *tag   = new AliRunTag;
390   AliEventTag *evTag = new AliEventTag;
391   fChain->SetBranchAddress("AliTAG",&tag);
392
393   for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
394
395     fChain->GetEntry(iTagFiles);
396     TTree* tree = fChain->GetTree();
397     if (cTree != tree) {
398         // Fix for aod tags: for each tree in the chain, merge the entries
399         cTree    = tree;
400         cEntries = tree->GetEntries();
401         iev      = 0;
402         ientry   = 0;
403     }
404     //Event list
405     iTotalEvents += tag->GetNEvents();
406     if ((iev == 0) || !aod) localList->Reset();
407     if(runTagCuts->IsAccepted(tag)) {
408       if(lhcTagCuts->IsAccepted(tag->GetLHCTag())) {
409         if(detTagCuts->IsAccepted(tag->GetDetectorTags())) {
410           Int_t iEvents = tag->GetNEvents();
411           const TClonesArray *tagList = tag->GetEventTags();
412           iRejectedEvtInFile = 0;
413           iAcceptedEvtInFile = 0;
414           for(Int_t i = 0; i < iEvents; i++) {
415             evTag = (AliEventTag *) tagList->At(i);
416             guid = evTag->GetGUID(); 
417             turl = evTag->GetTURL(); 
418             lfn = turl(8,turl.Length());
419             if(evTagCuts->IsAccepted(evTag)) {
420                 if(aod) localList->Enter(iev);
421                 else localList->Enter(i);
422                 iAcceptedEvtInFile++;
423             }
424             else {
425               iRejectedEvt++;
426               iRejectedEvtInFile++;
427             }
428             iev++;
429           }//event loop
430           if ((ientry == cEntries-1) || !aod) {
431             iAccepted += localList->GetN();
432             collection->WriteBody(iTagFiles+1,guid,lfn,turl,localList,iAcceptedEvtInFile,iRejectedEvtInFile);
433           }
434         }//detector tag cuts
435         else {
436           iRejectedDet += tag->GetNEvents();
437         }
438       }//lhc tag cuts 
439       else {
440         iRejectedLHC += tag->GetNEvents();
441       }
442     }//run tag cuts
443     else {
444       iRejectedRun += tag->GetNEvents();
445     }
446     tag->Clear();
447     ientry++;
448   }//tag file loop
449   collection->WriteSummary(iTotalEvents, iAccepted, iRejectedRun, iRejectedLHC, iRejectedDet, iRejectedEvt);
450   collection->Export();
451   
452   delete tag;
453   delete localList;
454   return kTRUE;
455 }
456
457 //___________________________________________________________________________
458 Bool_t AliTagAnalysis::CreateXMLCollection(const char* name, 
459                                            const char *fRunCut, 
460                                            const char *fLHCCut, 
461                                            const char *fDetectorCut, 
462                                            const char *fEventCut) {
463   //Queries the tag chain using the defined 
464   //event tag cuts from the AliEventTagCuts object
465   //and returns a XML collection
466   AliInfo(Form("Creating the collection........"));
467
468   Bool_t aod = kFALSE;
469   if(fAnalysisType == "AOD") aod = kTRUE;
470
471   AliXMLCollection *collection = new AliXMLCollection();
472   collection->SetCollectionName(name);
473   collection->WriteHeader();
474
475   TString guid;
476   TString turl;
477   TString lfn;
478   TEntryList* localList = new TEntryList();
479   
480   Int_t iAccepted = 0;
481   Int_t iev       = 0;
482   Int_t ientry    = 0;
483   Int_t cEntries  = 0;
484
485   Int_t iRejectedRun = 0;
486   Int_t iRejectedLHC = 0;
487   Int_t iRejectedDet = 0;
488   Int_t iRejectedEvt = 0;
489
490   Int_t iTotalEvents = 0;
491
492   Int_t iAcceptedEvtInFile = 0;
493   Int_t iRejectedEvtInFile = 0;
494
495   //Defining tag objects
496   AliRunTag *tag     = new AliRunTag;
497   AliEventTag *evTag = new AliEventTag;
498   fChain->SetBranchAddress("AliTAG",&tag);
499
500   TTreeFormula *fRunFormula = new TTreeFormula("fRun",fRunCut,fChain);
501   TTreeFormula *fLHCFormula = new TTreeFormula("fLHC",fLHCCut,fChain);   
502   TTreeFormula *fDetectorFormula = new TTreeFormula("fDetector",fDetectorCut,fChain);
503   TTreeFormula *fEventFormula = new TTreeFormula("fEvent",fEventCut,fChain);
504
505   Int_t current = -1;
506   for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
507
508     fChain->GetEntry(iTagFiles);
509     if (current != fChain->GetTreeNumber()) {
510       fRunFormula->UpdateFormulaLeaves();
511       fLHCFormula->UpdateFormulaLeaves();
512       fDetectorFormula->UpdateFormulaLeaves();
513       fEventFormula->UpdateFormulaLeaves();
514       // Fix for aod tags: for each tree in the chain, merge the entries
515       cEntries = (fChain->GetTree())->GetEntries();
516       iev      = 0;
517       ientry   = 0;
518       //
519       current = fChain->GetTreeNumber();
520     }
521     //Event list
522     iTotalEvents += tag->GetNEvents();
523     if ((iev == 0) || !aod) localList->Reset();
524     if(fRunFormula->EvalInstance(iTagFiles) == 1) {
525       if(fLHCFormula->EvalInstance(iTagFiles) == 1) {    
526         if(fDetectorFormula->EvalInstance(iTagFiles) == 1) {     
527           Int_t iEvents = fEventFormula->GetNdata();
528           const TClonesArray *tagList = tag->GetEventTags();
529           iRejectedEvtInFile = 0;
530           iAcceptedEvtInFile = 0;
531           for(Int_t i = 0; i < iEvents; i++) {
532             evTag = (AliEventTag *) tagList->At(i);
533             guid = evTag->GetGUID(); 
534             turl = evTag->GetTURL(); 
535             lfn = turl(8,turl.Length());
536             if(fEventFormula->EvalInstance(i) == 1) {
537                 if(aod) localList->Enter(iev);
538                 else localList->Enter(i);
539                 iAcceptedEvtInFile++;
540             }
541             else {
542               iRejectedEvt++;
543               iRejectedEvtInFile++;
544             }
545             iev++;
546           }//event loop
547           if ((ientry == cEntries-1) || !aod) {
548             collection->WriteBody(iTagFiles+1,guid,lfn,turl,localList,iAcceptedEvtInFile, iRejectedEvtInFile);
549             iAccepted += localList->GetN();
550           }
551         }//detector tag cuts
552         else {
553           iRejectedDet += tag->GetNEvents();
554         }
555       }//lhc tag cuts 
556       else {
557         iRejectedLHC += tag->GetNEvents();
558       }
559     }//run tag cuts
560     else {
561       iRejectedRun += tag->GetNEvents();
562     }
563     ientry++;
564   }//tag file loop
565   collection->WriteSummary(iTotalEvents, iAccepted, iRejectedRun, iRejectedLHC, iRejectedDet, iRejectedEvt);
566   collection->Export();
567
568   delete tag;
569   return kTRUE;
570 }
571
572 //___________________________________________________________________________
573 TChain *AliTagAnalysis::GetInputChain(const char* system, const char *wn) {
574   //returns the chain+event list - used in batch sessions
575   // this function will be removed once the new root 
576   // improvements are committed
577   TString fsystem = system;
578   Int_t iAccepted = 0;
579
580   TChain *fAnalysisChain = 0;
581   if(fAnalysisType == "ESD") fAnalysisChain = new TChain("esdTree");
582   else if(fAnalysisType == "AOD") fAnalysisChain = new TChain("aodTree");
583   else AliFatal("Only ESD and AOD type is implemented!!!");
584   
585   //Event list
586   TEventList *fEventList = new TEventList();
587   AliXMLCollection *collection = AliXMLCollection::Open(wn);
588
589   collection->Reset();
590   while (collection->Next()) {
591     AliInfo(Form("Adding: %s",collection->GetTURL("")));
592     fAnalysisChain->Add(collection->GetTURL(""));
593     TEntryList *list = (TEntryList *)collection->GetEventList("");
594     for(Int_t i = 0; i < list->GetN(); i++) fEventList->Enter(iAccepted+list->GetEntry(i));
595
596     if(fsystem == "pp") iAccepted += 100;
597     else if(fsystem == "PbPb") iAccepted += 1;
598   }
599
600   fAnalysisChain->SetEventList(fEventList);
601   
602   AliInfo(Form("Number of selected events: %d",fEventList->GetN()));
603
604   return fAnalysisChain;
605 }
606
607 //___________________________________________________________________________
608 TChain *AliTagAnalysis::GetChainFromCollection(const char* collectionname, 
609                                                const char* treename) {
610   //returns the TChain+TEntryList object- used in batch sessions
611   TString aliceFile = treename;
612   Int_t iAccepted = 0;
613   TChain *fAnalysisChain = 0;
614   if(aliceFile == "esdTree") fAnalysisChain = new TChain("esdTree");
615   else if(aliceFile == "aodTree") fAnalysisChain = new TChain("aodTree");
616   else AliFatal("Inconsistent tree name - use esdTree or aodTree!");
617
618   //Event list
619   fGlobalList = new TEntryList();
620   AliXMLCollection *collection = AliXMLCollection::Open(collectionname);
621
622   // Tag selection summary per file
623   TMap *tagCutSummary = new TMap();
624   tagCutSummary->SetName("TagCutSumm");
625
626   collection->Reset();
627   while (collection->Next()) {
628     AliInfo(Form("Adding: %s",collection->GetTURL("")));
629     fAnalysisChain->Add(collection->GetTURL(""));
630     TEntryList *list = (TEntryList *)collection->GetEventList("");
631     list->SetTreeName(aliceFile.Data());
632     list->SetFileName(collection->GetTURL(""));
633     fGlobalList->Add(list);
634     iAccepted += list->GetN();
635     if (collection->GetCutSumm())
636       tagCutSummary->Add(new TObjString(collection->GetTURL("")), new TObjString(collection->GetCutSumm()));
637   }
638
639   fAnalysisChain->SetEntryList(fGlobalList,"ne");
640   
641   AliInfo(Form("Number of selected events: %d",iAccepted));
642
643   TList *aUserInfo = fAnalysisChain->GetUserInfo();
644   aUserInfo->Add(tagCutSummary);
645
646   Int_t iAccEv;
647   Int_t iTotalEvents;
648   Int_t iRejRun;
649   Int_t iRejLHC;
650   Int_t iRejDet;
651   Int_t iRejEvt;
652
653   collection->GetCollectionSummary(&iTotalEvents, &iAccEv, &iRejRun, &iRejLHC, &iRejDet, &iRejEvt);
654  
655   char nstr[2000];
656
657   sprintf(nstr, "TotalEvents=%i", iTotalEvents);
658   TObjString *iTotStr = new TObjString(nstr);
659   aUserInfo->Add(iTotStr);
660
661   sprintf(nstr, "AcceptedEvents=%i", iAccepted);
662   TObjString *iAccStr = new TObjString(nstr);
663   aUserInfo->Add(iAccStr);
664
665   sprintf(nstr, "RejectedRun=%i", iRejRun);
666   TObjString *iRejRunStr = new TObjString(nstr);
667   aUserInfo->Add(iRejRunStr);
668
669   sprintf(nstr, "RejectedLHC=%i", iRejLHC);
670   TObjString *iRejLHCStr = new TObjString(nstr);
671   aUserInfo->Add(iRejLHCStr);
672
673   sprintf(nstr, "RejectedDet=%i", iRejDet);
674   TObjString *iRejDetStr = new TObjString(nstr);
675   aUserInfo->Add(iRejDetStr);
676
677   sprintf(nstr, "RejectedEvt=%i", iRejEvt);
678   TObjString *iRejEvtStr = new TObjString(nstr);
679   aUserInfo->Add(iRejEvtStr);
680
681   return fAnalysisChain;
682 }