]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/muon/AliMergeableCollection.cxx
Updates to run with deltas (L. Cunqueiro)
[u/mrichter/AliRoot.git] / PWG3 / muon / AliMergeableCollection.cxx
1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 *                                                                        *
4 * Author: The ALICE Off-line Project.                                    *
5 * Contributors are mentioned in the code where appropriate.              *
6 *                                                                        *
7 * Permission to use, copy, modify and distribute this software and its   *
8 * documentation strictly for non-commercial purposes is hereby granted   *
9 * without fee, provided that the above copyright notice appears in all   *
10 * copies and that both the copyright notice and this permission notice   *
11 * appear in the supporting documentation. The authors make no claims     *
12 * about the suitability of this software for any purpose. It is          *
13 * provided "as is" without express or implied warranty.                  *
14 **************************************************************************/
15
16 // $Id: AliMergeableCollection.cxx 50593 2011-07-14 17:42:28Z martinez $
17
18 ///
19 /// A mergeable object container. 
20 ///
21 /// For each tuple (key1,key2,..,keyN) a (hash)list of mergeable objects is associated.
22 /// Note that key1, key2 (optional), ..., keyN (optional) are strings. 
23 /// Those strings should not contain "/" themselves.
24 ///
25 /// More helper functions might be added in the future (e.g. Project, etc...)
26
27 #include "AliMergeableCollection.h"
28
29 ClassImp(AliMergeableCollection)
30
31 #include "AliLog.h"
32 #include "Riostream.h"
33 #include "TError.h"
34 #include "THashList.h"
35 #include "TKey.h"
36 #include "TMap.h"
37 #include "TObjArray.h"
38 #include "TObjString.h"
39 #include "TRegexp.h"
40 #include "TROOT.h"
41 #include "TSystem.h"
42 #include "TH1.h"
43 //#include "TH2.h"
44
45 //_____________________________________________________________________________
46 AliMergeableCollection::AliMergeableCollection(const char* name, const char* title) 
47 : TNamed(name,title), fMap(0x0), fMustShowEmptyObject(0), fMapVersion(0), fMessages()
48 {
49   /// Ctor
50 }
51
52 //_____________________________________________________________________________
53 AliMergeableCollection::~AliMergeableCollection()
54 {
55   /// dtor. Note that the map is owner
56   if ( fMap ) fMap->DeleteAll();
57   delete fMap;
58 }
59
60 //_____________________________________________________________________________
61 Bool_t 
62 AliMergeableCollection::Adopt(TObject* obj)
63 {
64   /// Adopt a given object at top level (i.e. no key)
65   return InternalAdopt("",obj);
66 }
67
68 //_____________________________________________________________________________
69 Bool_t 
70 AliMergeableCollection::Adopt(const char* identifier, TObject* obj)
71 {
72   /// Adopt a given object, and associate it with pair key
73   TString sidentifier(identifier);
74   if ( ! sidentifier.IsNull() ){
75     if ( ! sidentifier.EndsWith("/") ) sidentifier.Append("/");
76     if ( ! sidentifier.BeginsWith("/") ) sidentifier.Prepend("/");
77   }
78   return InternalAdopt(sidentifier.Data(),obj);
79 }
80
81 //_____________________________________________________________________________
82 void AliMergeableCollection::ClearMessages()
83 {
84   /// clear pending messages
85   fMessages.clear();
86 }
87
88 //_____________________________________________________________________________
89 TIterator*
90 AliMergeableCollection::CreateIterator(Bool_t direction) const
91 {
92   /// Create an iterator (must be deleted by the client)
93   return fMap ? new AliMergeableCollectionIterator(this,direction) : 0x0;
94 }
95
96 //_____________________________________________________________________________
97 AliMergeableCollection*
98 AliMergeableCollection::Clone(const char* name) const
99 {
100   /// Clone this collection.
101   /// We loose the messages.
102   
103   AliMergeableCollection* newone = new AliMergeableCollection(name,GetTitle());
104   
105   newone->fMap = static_cast<TMap*>(fMap->Clone());
106   newone->fMustShowEmptyObject = fMustShowEmptyObject;
107   newone->fMapVersion = fMapVersion;  
108   
109   return newone;
110 }
111
112 //_____________________________________________________________________________
113 void 
114 AliMergeableCollection::Delete(Option_t*)
115 {
116   /// Delete all the objects
117   fMap->DeleteAll();
118   delete fMap;
119   fMap=0x0;
120 }
121
122 //_____________________________________________________________________________
123 TObject* 
124 AliMergeableCollection::FindObject(const char* fullIdentifier) const
125 {
126   /// Find an object by its full identifier.
127   
128   return GetObject(fullIdentifier);
129 }
130
131 //_____________________________________________________________________________
132 TObject* 
133 AliMergeableCollection::FindObject(const TObject *object) const
134 {
135   /// Find an object 
136   AliWarning("This method is awfully inefficient. Please improve it or use FindObject(const char*)");
137   TIter next(CreateIterator());
138   TObject* obj;
139   while ( ( obj=next() ) )
140   {
141     if ( obj->IsEqual(object) ) return obj;
142   }
143   return 0x0;
144 }
145
146
147 //_____________________________________________________________________________
148 TList*
149 AliMergeableCollection::CreateListOfKeys(Int_t index) const
150 {
151   /// Create the list of keys at level index
152   
153   TList* list = new TList;
154   list->SetOwner(kTRUE);
155   
156   TObjArray* ids = SortAllIdentifiers();
157   TIter next(ids);
158   TObjString* str;
159   
160   while ( ( str = static_cast<TObjString*>(next()) ) )
161   {
162     TString oneid = GetKey(str->String().Data(),index,kFALSE);
163     if (oneid.Length()>0 && !list->Contains(oneid))
164     {
165       list->Add(new TObjString(oneid));
166     }
167   }
168   
169   delete ids;
170   return list;
171 }
172
173 //_____________________________________________________________________________
174 TList* 
175 AliMergeableCollection::CreateListOfObjectNames(const char* identifier) const
176 {
177   /// Create list of object names for /key1/key2/key...
178   /// Returned list must be deleted by client
179   
180   TList* listOfNames = new TList;
181   listOfNames->SetOwner(kTRUE);
182   
183   TIter next(Map());
184   TObjString* str;
185   
186   while ( ( str = static_cast<TObjString*>(next()) ) )
187   {
188     TString currIdentifier = str->String();
189     if ( currIdentifier.CompareTo(identifier) ) continue;
190     
191     THashList* list = static_cast<THashList*>(Map()->GetValue(identifier));
192     
193     TIter nextObject(list);
194     TObject* obj;
195     
196     while ( ( obj = nextObject() ) )
197     {
198       listOfNames->Add(new TObjString(obj->GetName()));
199     }    
200   }
201   
202   return listOfNames;
203 }
204
205
206 //_____________________________________________________________________________
207 TString
208 AliMergeableCollection::GetIdentifier(const char* fullIdentifier) const
209 {
210   /// Extract the identifier from the fullIdentifier
211   TString sfullIdentifier(fullIdentifier);
212   TObjArray* arr = sfullIdentifier.Tokenize("/");
213   TString identifier = "";
214   for ( Int_t istr=0; istr<arr->GetLast(); istr++ ) {
215     identifier += "/" + GetKey(fullIdentifier, istr, kTRUE);
216   }
217   delete arr;
218   identifier.Append("/");
219   return identifier;
220 }
221
222 //_____________________________________________________________________________
223 TString
224 AliMergeableCollection::GetKey(const char* identifier, Int_t index, Bool_t idContainsObjName) const
225 {
226   /// Extract the index element of the key pair from the fullIdentifier
227   TString sidentifier(identifier);
228   if ( ! idContainsObjName ) sidentifier.Append("/dummy");
229   return InternalDecode(sidentifier.Data(),index);
230 }
231
232 //_____________________________________________________________________________
233 TString
234 AliMergeableCollection::GetObjectName(const char* fullIdentifier) const
235 {
236   /// Extract the object name from an identifier
237   
238   return InternalDecode(fullIdentifier,-1);  
239 }
240
241
242 //_____________________________________________________________________________
243 TObject* 
244 AliMergeableCollection::GetObject(const char* fullIdentifier) const
245 {
246   /// Get object key1/key2/.../objectName:action
247   
248   TObjArray* arr = TString(fullIdentifier).Tokenize(":");
249   
250   TString fullIdWithoutAction(static_cast<TObjString*>(arr->At(0))->String());
251   TString action;
252   
253   if ( arr->GetLast() > 0 ) 
254   {
255     action = static_cast<TObjString*>(arr->At(1))->String();
256     action.ToUpper();
257   }
258   
259   delete arr;
260   
261   return GetObject(GetIdentifier(fullIdWithoutAction).Data(), GetObjectName(fullIdWithoutAction));
262   
263   
264   //  if (obj)
265   //  {
266   //    TH2* h2(0x0);
267   //    
268   //    if ( action == "PX" && ( (h2 = dynamic_cast<TH2*>(obj)) ) ) 
269   //    {
270   //      return h2->ProjectionX(NormalizeName(identifier.Data(),action.Data()).Data());
271   //    }
272   //    else if ( action == "PY" && ( (h2 = dynamic_cast<TH2*>(obj)) ) ) 
273   //    {
274   //      return h2->ProjectionY(NormalizeName(identifier.Data(),action.Data()).Data());
275   //    }
276   //    else if ( action == "PFX" && ( (h2 = dynamic_cast<TH2*>(obj)) ) ) 
277   //    {
278   //      return h2->ProfileX(NormalizeName(identifier.Data(),action.Data()).Data());
279   //    }
280   //    else if ( action == "PFY" && ( (h2 = dynamic_cast<TH2*>(obj)) ) ) 
281   //    {
282   //      return h2->ProfileY(NormalizeName(identifier.Data(),action.Data()).Data());
283   //    }
284   //    
285   //  }
286   //  else
287   //  {
288   //    AliDebug(1,Form("Object %s not found",sidentifier));
289   //  }
290   //  return obj;
291 }
292
293
294 //_____________________________________________________________________________
295 TObject* 
296 AliMergeableCollection::GetObject(const char* identifier, 
297                                   const char* objectName) const
298 {
299   /// Get object for (identifier,objectName) triplet
300   
301   TString sidentifier(identifier);
302   if ( ! sidentifier.IsNull() ) {
303     if ( ! sidentifier.BeginsWith("/") ) sidentifier.Prepend("/");
304     if ( ! sidentifier.EndsWith("/") ) sidentifier.Append("/");
305   }
306   return InternalObject(sidentifier.Data(),objectName);
307 }
308
309 //_____________________________________________________________________________
310 TObject* AliMergeableCollection::GetSum(const char* idPattern)
311 {
312   /// Sum objects
313   /// The pattern must be in the form:
314   /// /key1_pattern1 key1_pattern2,key1_pattern,.../key2_pattern1,key2_pattern,.../.../objectName_pattern1,objectName_pattern...
315   /// The logical or between patterns separated by commas is taken
316   
317   TObject* sumObject = 0x0;
318   
319   // Build array of lists of pattern
320   TString idPatternString(idPattern);
321   TObjArray* keyList = idPatternString.Tokenize("/");
322   TObjArray keyMatrix(keyList->GetEntries());
323   keyMatrix.SetOwner();
324   for ( Int_t ikey=0; ikey<keyList->GetEntries(); ikey++ ) {
325     TObjArray* subKeyList = ((TObjString*)keyList->At(ikey))->GetString().Tokenize(",");
326     keyMatrix.AddAt(subKeyList, ikey);
327   }
328   delete keyList;
329   
330   TString debugMsg = "Adding objects:";
331   
332   TIter next(Map());
333   TObjString* str;
334   while ( ( str = static_cast<TObjString*>(next()) ) )
335   {
336     TString identifier = str->String();
337     
338     Bool_t listMatchPattern = kTRUE;
339     for ( Int_t ikey=0; ikey<keyMatrix.GetEntries()-1; ikey++ ) {
340       TString currKey = GetKey(identifier, ikey, kFALSE);
341       Bool_t matchKey = kFALSE;
342       TObjArray* subKeyList = static_cast<TObjArray*> ( keyMatrix.At(ikey) );
343       for ( Int_t isub=0; isub<subKeyList->GetEntries(); isub++ ) {
344         TString subKeyString = static_cast<TObjString*> (subKeyList->At(isub))->GetString();
345         if ( currKey.Contains(subKeyString.Data()) ) {
346           matchKey = kTRUE;
347           break;
348         }
349       } // loop on the list of patterns of each key
350       if ( ! matchKey ) {
351         listMatchPattern = kFALSE;
352         break;
353       }
354     } // loop on keys in the idPattern
355     if ( ! listMatchPattern ) continue;
356     
357     THashList* list = static_cast<THashList*>(Map()->GetValue(identifier.Data()));
358     
359     TIter nextObj(list);
360     TObject* obj;
361     
362     while ( ( obj = nextObj()) )
363     {
364       TString currKey = obj->GetName();
365       Bool_t matchKey = kFALSE;
366       TObjArray* subKeyList = static_cast<TObjArray*> ( keyMatrix.Last() );
367       for ( Int_t isub=0; isub<subKeyList->GetEntries(); isub++ ) {
368         TString subKeyString = static_cast<TObjString*> (subKeyList->At(isub))->GetString();
369         if ( currKey.Contains(subKeyString.Data()) ) {
370           matchKey = kTRUE;
371           break;
372         }
373       }
374       if ( ! matchKey ) continue;
375       if ( ! sumObject ) sumObject = obj->Clone();
376       else MergeObject(sumObject, obj);
377       debugMsg += Form(" %s%s",identifier.Data(),obj->GetName());
378     } // loop on objects in list
379   } // loop on identifiers in map
380   
381   AliDebug(1,debugMsg.Data());
382   
383   return sumObject;
384 }
385
386 //_____________________________________________________________________________
387 Bool_t AliMergeableCollection::InternalAdopt(const char* identifier, TObject* obj)
388 {
389   /// Adopt an obj
390   
391   if (!obj)
392   {
393     Error("Adopt","Cannot adopt a null object");
394     return kFALSE;
395   }
396   
397   if ( ! obj->IsA()->InheritsFrom(TObject::Class()) ||
398         ! obj->IsA()->GetMethodWithPrototype("Merge", "TCollection*") ) {
399     Error("Adopt","Cannot adopt an object which is not mergeable!"); 
400   }
401   
402   THashList* hlist = 0x0;
403   
404   hlist = static_cast<THashList*>(Map()->GetValue(identifier));
405   
406   if (!hlist)
407   {
408     hlist = new THashList;
409     hlist->SetOwner(kTRUE);
410     Map()->Add(new TObjString(identifier),hlist);
411     hlist->SetName(identifier);
412   }
413   
414   TObject* existingObj = hlist->FindObject(obj->GetName());
415   
416   if (existingObj)
417   {
418     AliError(Form("Cannot adopt an already existing object : %s -> %s",identifier,existingObj->GetName()));
419     return kFALSE;
420   }
421   
422   if ( obj->IsA()->InheritsFrom(TH1::Class()) ) (static_cast<TH1*> ( obj ))->SetDirectory(0);  
423   
424   hlist->AddLast(obj);
425   
426   return kTRUE;
427   
428 }
429
430 //_____________________________________________________________________________
431 TString
432 AliMergeableCollection::InternalDecode(const char* identifier, Int_t index) const
433 {
434   /// Extract the index-th element of the identifier (/key1/key2/.../keyN/objectName)
435   /// object is index=-1 (i.e. last)
436   
437   if ( strlen(identifier) > 0 && identifier[0] != '/' ) 
438   {    
439     AliError(Form("identifier %s is malformed (should start with /)",identifier));
440     return "";
441   }
442   
443   TObjArray* array = TString(identifier).Tokenize("/");
444
445   if ( index >= array->GetLast() ) 
446   {
447     AliError(Form("Requiring index %i of identifier %s which only have %i",index, identifier, array->GetLast()));
448     delete array;
449     return "";
450   }
451
452   TString value("");
453   
454   if ( index < 0 ) 
455   {
456     value = static_cast<TObjString*>(array->Last())->String();    
457   }
458   else if ( index <= array->GetLast() ) 
459   {
460     value = static_cast<TObjString*>(array->At(index))->String();
461   }
462   
463   delete array;
464   
465   return value;
466 }
467
468 //_____________________________________________________________________________
469 TObject* 
470 AliMergeableCollection::InternalObject(const char* identifier,
471                                        const char* objectName) const
472 {
473   /// Get object for (identifier,objectName) 
474   
475   if (!fMap) 
476   {
477     return 0x0;
478   }
479   
480   THashList* hlist = static_cast<THashList*>(Map()->GetValue(identifier));
481   if (!hlist) 
482   {
483     TString msg(Form("Did not find hashlist for identifier=%s dir=%s",identifier,gDirectory ? gDirectory->GetName() : "" ));
484     fMessages[msg.Data()]++;
485     return 0x0;
486   }
487   
488   TObject* obj = hlist->FindObject(objectName);
489   if (!obj)
490   {
491     TString msg(Form("Did not find objectName=%s in %s",objectName,identifier));
492     fMessages[msg.Data()]++;
493   }
494   return obj;
495 }
496
497
498 //_____________________________________________________________________________
499 Bool_t AliMergeableCollection::IsEmptyObject(TObject* obj) const
500 {
501   /// Check if object is empty
502   /// (done only for TH1, so far)
503     
504   if ( obj->IsA()->InheritsFrom(TH1::Class()) ) {
505     TH1* histo = static_cast<TH1*> (obj);
506     if ( histo->GetEntries() == 0 ) return kTRUE;
507   }
508
509   return kFALSE;
510
511 }
512
513
514 //_____________________________________________________________________________
515 TMap* AliMergeableCollection::Map() const
516 {
517   /// Wrapper to insure proper key formats (i.e. new vs old)
518   
519   if (!fMap)
520   {
521     fMap = new TMap;
522     fMap->SetOwnerKeyValue(kTRUE,kTRUE);
523     fMapVersion = 1;
524   }
525   else
526   {
527     if ( fMapVersion < 1 ) 
528     {
529       AliInfo("Remapping");
530       // change the keys
531       TIter next(fMap);
532       TObjString* str;
533       
534       while ( ( str = static_cast<TObjString*>(next()) ) )
535       {
536         if ( str->String().Contains("./") )
537         {
538           TString newkey(str->String());
539           
540           newkey.ReplaceAll("./","");
541           
542           TObject* o = fMap->GetValue(str);
543           
544           TPair* p = fMap->RemoveEntry(str);
545           if (!p)
546           {
547             AliError("oups oups oups");
548             return 0x0;
549           }
550           
551           fMap->Add(new TObjString(newkey.Data()),o);
552           
553           delete p;
554         }
555       }
556       
557       fMapVersion = 1;
558     }
559   }
560   
561   return fMap;
562 }
563
564 //_____________________________________________________________________________
565 Long64_t
566 AliMergeableCollection::Merge(TCollection* list)
567 {
568   // Merge a list of AliMergeableCollection objects with this
569   // Returns the number of merged objects (including this).
570   
571   if (!list) return 0;
572   
573   if (list->IsEmpty()) return 1;
574   
575   TIter next(list);
576   TObject* currObj;
577   TList mapList;
578   Int_t count(0);
579   
580   while ( ( currObj = next() ) )
581   {
582     AliMergeableCollection* mergeCol = dynamic_cast<AliMergeableCollection*>(currObj);
583     if (!mergeCol) {
584       AliFatal(Form("object named \"%s\" is a %s instead of an AliMergeableCollection!", currObj->GetName(), currObj->ClassName()));
585       continue;
586     }
587     
588     ++count;
589     
590     if ( mergeCol->fMap ) mergeCol->Map(); // to insure keys in the new format
591     
592     TIter nextIdentifier(mergeCol->fMap);
593     TObjString* identifier;
594
595     while ( ( identifier = static_cast<TObjString*>(nextIdentifier()) ) )
596     {
597       THashList* otherList = static_cast<THashList*>(mergeCol->fMap->GetValue(identifier->String().Data()));
598
599       TIter nextObject(otherList);
600       TObject* obj;
601       
602       while ( ( obj = nextObject() ) )
603       {
604         TString newid(Form("%s%s",identifier->String().Data(),obj->GetName()));
605         
606         TObject* thisObject = GetObject(newid.Data());
607         
608         if (!thisObject)
609         {
610           AliDebug(1,Form("Adopting a new object = %s%s",identifier->String().Data(),obj->GetName()));
611           
612           Bool_t ok = Adopt(identifier->String(), obj->Clone());
613           
614           if (!ok)
615           {
616             AliError(Form("Adoption of object %s failed",obj->GetName()));
617           }
618         }
619         else
620         {
621           // add it...
622           AliDebug(1,Form("Merging object = %s%s",
623                           identifier->String().Data(),
624                           obj->GetName()));
625           
626           MergeObject(thisObject, obj);
627         }
628       } // loop on objects in map
629     } // loop on identifiers
630   } // loop on collections in list
631          
632   AliDebug(1,Form("count=%d",count));
633   
634   return count+1;
635 }
636
637 //_____________________________________________________________________________
638 Bool_t AliMergeableCollection::MergeObject(TObject* baseObject, TObject* objToAdd)
639 {
640   /// Add objToAdd to baseObject
641   
642   if ( baseObject->IsA()->Class() != objToAdd->IsA()->Class() ) {
643     printf("MergeObject: Cannot add %s to %s", objToAdd->ClassName(), baseObject->ClassName());
644     return kFALSE;
645   }
646   if ( ! baseObject->IsA()->InheritsFrom(TObject::Class()) ||
647       ! baseObject->IsA()->GetMethodWithPrototype("Merge", "TCollection*") ) {
648     printf("MergeObject: Objects are not mergeable!");
649     return kFALSE;
650   }  
651   
652   TList list;
653   list.Add(objToAdd);
654   
655   TString listArgs = Form("((TCollection*)0x%lx)", (ULong_t)&list);
656   Int_t error = 0;
657   baseObject->Execute("Merge", listArgs.Data(), &error);
658   return kTRUE;
659 }
660
661 //_____________________________________________________________________________
662 TString AliMergeableCollection::NormalizeName(const char* identifier,const char* action) const
663 {
664   /// Replace / by _ to build a root-compliant histo name
665   TString name(GetName());
666   
667   name += "_";
668   name += identifier;
669   name += "_";
670   name += action;
671   name.ReplaceAll("/","_");
672   name.ReplaceAll("-","_");
673   return name;
674 }
675
676 //_____________________________________________________________________________
677 Int_t 
678 AliMergeableCollection::NumberOfObjects() const
679 {
680   /// Get the number of objects we hold
681   TIter next(CreateIterator(this));
682   Int_t count(0);
683   while ( next() ) ++count;
684   return count;
685 }
686
687 //_____________________________________________________________________________
688 Int_t 
689 AliMergeableCollection::NumberOfKeys() const
690 {
691   /// Get the number of keys we have
692   return fMap ? fMap->GetSize() : 0;
693 }
694
695 //_____________________________________________________________________________
696 void 
697 AliMergeableCollection::Print(Option_t* option) const
698 {
699   /// Print all the objects we hold, in a hopefully visually pleasing
700   /// way.
701   ///
702   /// Option can be used to select given part only, using the schema :
703   /// /*/*/*/*/*
704   /// Where the stars are wilcards for /key1/key2/.../objectName
705   ///
706   /// if * is used it is assumed to be a wildcard for objectName
707   ///
708   /// For other selections the full syntax /*/*/*/*/* must be used.
709   /// 
710   /// Use "-" as objectName to disable object's name output
711   
712   cout << Form("AliMergeableCollection(%s,%s) : %d keys and %d objects",
713                GetName(),GetTitle(),
714                NumberOfKeys(), NumberOfObjects()) << endl;
715   
716   if (!strlen(option)) return;
717     
718   TObjArray* select = TString(option).Tokenize("/");
719   
720   TString sreObjectName(select->Last()->GetName());
721   TRegexp reObjectName(sreObjectName.Data(),kTRUE);
722   
723   TObjArray* identifiers = SortAllIdentifiers();
724   
725   printf("identifiers entries %i\n", identifiers->GetEntries());
726     
727   TIter nextIdentifier(identifiers);
728   
729   TObjString* sid(0x0);
730   
731   while ( ( sid = static_cast<TObjString*>(nextIdentifier()) ) )
732   {
733     Bool_t identifierPrinted(kFALSE);
734
735     TString identifier(sid->String());
736     
737     Bool_t matchPattern = kTRUE;
738     for ( Int_t isel=0; isel<select->GetLast(); isel++ ) {
739       if ( ! GetKey(identifier.Data(), isel, kFALSE).Contains(TRegexp(select->At(isel)->GetName(),kTRUE)) ) {
740         matchPattern = kFALSE;
741         break;
742       }
743     }
744     if ( ! matchPattern ) continue;
745     
746     if ( sreObjectName == "*" ) {
747       identifierPrinted = kTRUE;
748       cout << identifier.Data() << endl;
749     }
750       
751     THashList * list = static_cast<THashList*>(Map()->GetValue(sid->String().Data()));      
752     TObjArray names;
753     names.SetOwner(kTRUE);
754     TIter nextUnsortedObj(list);
755     TObject* obj;
756     while ( ( obj = nextUnsortedObj() ) )
757     {
758       names.Add(new TObjString(obj->GetName()));
759     }
760     names.Sort();
761     TIter nextObjName(&names);
762     TObjString* oname;
763     while ( ( oname = static_cast<TObjString*>(nextObjName()) ) )
764     {
765       TString objName(oname->String());
766       if (objName.Contains(reObjectName) )
767       {
768         obj = list->FindObject(objName.Data());
769         if ( IsEmptyObject(obj) && ! fMustShowEmptyObject ) continue;
770         if (!identifierPrinted)
771         {
772           cout << identifier.Data() << endl;
773           identifierPrinted = kTRUE;
774         }
775         cout << Form("    (%s) %s", obj->ClassName(), obj->GetName());
776         if ( obj->IsA()->InheritsFrom(TH1::Class()) ) {
777           TH1* histo = static_cast<TH1*> (obj);
778           cout << Form(" %s Entries=%d Sum=%g",histo->GetTitle(),Int_t(histo->GetEntries()),histo->GetSumOfWeights());
779         }
780         cout << endl;
781       }
782     }
783     if (!identifierPrinted && sreObjectName=="-" )
784     { 
785       // to handle the case where we used objectName="-" to disable showing the objectNames,
786       // but we still want to see the matching keys maybe...
787       cout << identifier.Data() << endl;
788     }
789   }
790   
791   delete select;
792   
793   delete identifiers;
794 }
795
796 //_____________________________________________________________________________
797 void 
798 AliMergeableCollection::PrintMessages(const char* prefix) const
799 {
800   /// Print pending messages
801   
802   std::map<std::string,int>::const_iterator it;
803   
804   for ( it = fMessages.begin(); it != fMessages.end(); ++it ) 
805   {
806     cout << Form("%s : message %s appeared %5d times",prefix,it->first.c_str(),it->second) << endl;
807   }
808 }
809
810
811 //_____________________________________________________________________________
812 UInt_t 
813 AliMergeableCollection::EstimateSize(Bool_t show) const
814 {
815   /// Estimate the memory (in kilobytes) used by some objects
816
817 //  For TH1:
818 //  sizeof(TH1) + (nbins+2)*(nbytes_per_bin) +name+title_sizes 
819 //  if you have errors add (nbins+2)*8 
820     
821   TIter next(CreateIterator());
822   
823   TObject* obj;
824   UInt_t size(0);
825   
826   while ( ( obj = next() ) )
827   {
828     UInt_t thissize=0;
829     if ( obj->IsA()->InheritsFrom(TH1::Class()) ) {
830       TH1* histo = static_cast<TH1*> (obj);
831       Int_t nbins = (histo->GetNbinsX()+2);
832     
833       if (histo->GetNbinsY()>1)
834       {
835         nbins *= (histo->GetNbinsY()+2);
836       }
837     
838       if (histo->GetNbinsZ()>1)
839       {
840         nbins *= (histo->GetNbinsZ()+2);
841       }
842       
843       Bool_t hasErrors = ( histo->GetSumw2N() > 0 );
844     
845       TString cname(histo->ClassName());
846     
847       Int_t nbytesPerBin(0);
848     
849       if (cname.Contains(TRegexp("C$")) ) nbytesPerBin = sizeof(Char_t);
850       if (cname.Contains(TRegexp("S$")) ) nbytesPerBin = sizeof(Short_t);
851       if (cname.Contains(TRegexp("I$")) ) nbytesPerBin = sizeof(Int_t);
852       if (cname.Contains(TRegexp("F$")) ) nbytesPerBin = sizeof(Float_t);
853       if (cname.Contains(TRegexp("D$")) ) nbytesPerBin = sizeof(Double_t);
854         
855       if (!nbytesPerBin)
856       {
857         AliError(Form("Could not get the number of bytes per bin for histo %s of class %s. Thus the size estimate will be wrong !",
858                       histo->GetName(),histo->ClassName()));
859         continue;
860       }
861     
862       thissize = sizeof(histo) + nbins*(nbytesPerBin) + strlen(histo->GetName())
863       + strlen(histo->GetTitle());
864       
865       if ( hasErrors) thissize += nbins*8;
866     }
867     else {
868       AliWarning(Form("Cannot estimate size of %s\n", obj->ClassName()));
869       continue;
870     }
871
872     size += thissize;
873     
874     if ( show ) 
875     {
876       AliInfo(Form("Size of %30s is %20d bytes",obj->GetName(),thissize));
877     }
878   } // loop on objects
879
880   return size;
881 }
882
883 //_____________________________________________________________________________
884 void AliMergeableCollection::PruneEmptyObjects()
885 {
886   /// Delete the empty objects
887   /// (Implemented for TH1 only)
888   TIter next(Map());
889   TObjString* key;
890   
891   TList toBeRemoved;
892   toBeRemoved.SetOwner(kTRUE);
893   
894   while ( ( key = static_cast<TObjString*>(next()) ) )
895   {
896     TString identifier(key->String());
897     THashList* hlist = static_cast<THashList*>(Map()->GetValue(identifier.Data()));
898     TIter nextObject(hlist);
899     TObject* obj;
900     while ( ( obj = nextObject() ) )
901     {
902       if ( IsEmptyObject(obj) ) toBeRemoved.Add(new TObjString(Form("%s%s",identifier.Data(),obj->GetName())));
903     }
904   }
905   
906   TIter nextTBR(&toBeRemoved);
907   while ( ( key = static_cast<TObjString*>(nextTBR()) ) )
908   {
909     Remove(key->GetString().Data());
910     AliDebug(2,Form("Removing %s", key->GetString().Data()));
911   }
912 }
913
914 //_____________________________________________________________________________
915 AliMergeableCollection* 
916 AliMergeableCollection::Project(const char* identifier) const
917 {
918   /// To be implemented : would create a new collection starting at /key1/key2/...
919   
920   if (!fMap) return 0x0;
921   
922   AliMergeableCollection* mergCol = new AliMergeableCollection(Form("%s %s",GetName(),identifier),
923                                                                GetTitle());
924   
925   TIter next(Map());
926   TObjString* str;
927   
928   while ( ( str = static_cast<TObjString*>(next()) ) )
929   {
930     TString currIdentifier = str->String();
931     if ( ! currIdentifier.Contains(identifier) ) continue;
932     
933     THashList* list = static_cast<THashList*>(Map()->GetValue(identifier));
934     
935     TIter nextObj(list);
936     TObject* obj;
937     
938     while ( ( obj = nextObj()) )
939     {
940       TObject* clone = obj->Clone();
941
942       TString newkey(currIdentifier.Data());
943       newkey.ReplaceAll(identifier,"");
944
945       if (newkey=="/") newkey="";
946       
947       mergCol->InternalAdopt(newkey.Data(),clone);
948     }    
949   }
950
951   return mergCol;
952 }
953
954 //_____________________________________________________________________________
955 TObject* 
956 AliMergeableCollection::Remove(const char* fullIdentifier)
957 {
958   ///
959   /// Remove a given object (given its fullIdentifier=/key1/key2/.../objectName)
960   ///
961   /// Note that we do *not* remove the /key1/key2/... entry even if there's no
962   /// more object for this triplet.
963   ///
964   /// Not very efficient. Could be improved ?
965   ///
966   
967   TString identifier = GetIdentifier(fullIdentifier);
968   
969   THashList* hlist = dynamic_cast<THashList*>(Map()->GetValue(identifier.Data()));
970   
971   if (!hlist)
972   {
973     AliWarning(Form("Could not get hlist for key=%s",identifier.Data()));
974     return 0x0;
975   }
976     
977   TObject* obj = GetObject(fullIdentifier);
978   if (!obj)
979   {
980     AliError(Form("Could not find object %s",fullIdentifier));
981     return 0x0;
982   }
983   
984   TObject* rmObj = hlist->Remove(obj);
985   if (!rmObj)
986   {
987     AliError("Remove failed");
988     return 0x0;
989   }
990   
991   return rmObj;
992 }
993
994 //_____________________________________________________________________________
995 TObjArray*
996 AliMergeableCollection::SortAllIdentifiers() const
997 {
998   /// Sort our internal identifiers. Returned array must be deleted.
999   TObjArray* identifiers = new TObjArray;
1000   identifiers->SetOwner(kFALSE); 
1001   TIter next(Map());
1002   TObjString* sid;
1003   
1004   while ( ( sid = static_cast<TObjString*>(next()) ) )
1005   {
1006     if ( !identifiers->FindObject(sid->String().Data()) )
1007     {
1008       identifiers->Add(sid);      
1009     }
1010   }
1011   identifiers->Sort();
1012   return identifiers;
1013 }
1014
1015
1016 ///////////////////////////////////////////////////////////////////////////////
1017 //
1018 // AliMergeableCollectionIterator
1019 //
1020 ///////////////////////////////////////////////////////////////////////////////
1021
1022 class AliMergeableCollectionIterator;
1023
1024 //_____________________________________________________________________________
1025 AliMergeableCollectionIterator::AliMergeableCollectionIterator(const AliMergeableCollection* mcol, Bool_t dir)
1026 : fkMergeableCollection(mcol), fMapIterator(0x0), fHashListIterator(0x0), fDirection(dir)
1027 {
1028   /// Default ctor
1029 }
1030
1031 //_____________________________________________________________________________
1032 AliMergeableCollectionIterator&
1033 AliMergeableCollectionIterator::operator=(const TIterator&)
1034 {
1035   /// Overriden operator= (imposed by Root's declaration of TIterator ?)
1036   Fatal("TIterator::operator=","Not implementeable"); // because there's no clone in TIterator :-(
1037   return *this;
1038 }
1039
1040 //_____________________________________________________________________________
1041 AliMergeableCollectionIterator::~AliMergeableCollectionIterator()
1042 {
1043   /// dtor
1044   Reset();
1045 }
1046
1047 //_____________________________________________________________________________
1048 TObject* AliMergeableCollectionIterator::Next()
1049 {
1050   /// Advance to next object in the collection
1051   
1052   if (!fHashListIterator)
1053   {
1054     if ( !fMapIterator ) 
1055     {
1056       fMapIterator = fkMergeableCollection->fMap->MakeIterator(fDirection);
1057     }
1058     TObjString* key = static_cast<TObjString*>(fMapIterator->Next());
1059     if (!key)
1060     {
1061       // we are done
1062       return 0x0;
1063     }      
1064     THashList* list = static_cast<THashList*>(fkMergeableCollection->Map()->GetValue(key->String().Data()));
1065     if (!list) return 0x0;
1066     fHashListIterator = list->MakeIterator(fDirection);
1067   }
1068
1069   TObject* obj = fHashListIterator->Next();
1070   
1071   if (!obj) 
1072   {
1073     delete fHashListIterator;
1074     fHashListIterator = 0x0;
1075     return Next();
1076   }
1077   
1078   return obj;
1079 }
1080
1081 //_____________________________________________________________________________
1082 void AliMergeableCollectionIterator::Reset()
1083 {
1084   /// Reset the iterator
1085   delete fHashListIterator;
1086   delete fMapIterator;
1087 }