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