]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/AliCorrectionManagerBase.cxx
merging trunk to TPCdev
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / AliCorrectionManagerBase.cxx
1 #include "AliCorrectionManagerBase.h"
2 #include "AliOADBForward.h"
3 #include "AliForwardUtil.h"
4 #include <AliLog.h>
5 #include <TMath.h>
6 #include <iostream>
7 #include <TROOT.h>
8 #include <TSystem.h>
9 #include <TBrowser.h>
10 #include <TParameter.h>
11 #include <TFileMerger.h>
12
13 //____________________________________________________________________
14 AliCorrectionManagerBase::AliCorrectionManagerBase()
15   : fCorrections(),
16     fIsInit(false),
17     fRun(kIgnoreValue), 
18     fSys(kIgnoreValue), 
19     fSNN(kIgnoreValue), 
20     fField(kIgnoreField), 
21     fMC(false), 
22     fSatellite(false), 
23     fDB(0),
24     fDebug(false),
25     fFallBack(false)
26 {
27 }
28
29 //____________________________________________________________________
30 AliCorrectionManagerBase::AliCorrectionManagerBase(Bool_t)
31   : fCorrections(16),
32     fIsInit(false),
33     fRun(kIgnoreValue), 
34     fSys(kIgnoreValue), 
35     fSNN(kIgnoreValue), 
36     fField(kIgnoreField), 
37     fMC(false), 
38     fSatellite(false), 
39     fDB(0),
40     fDebug(false),
41     fFallBack(false)
42 {
43   fCorrections.SetOwner(false);
44   fCorrections.SetName("corrections");
45 }
46 //____________________________________________________________________
47 AliCorrectionManagerBase::AliCorrectionManagerBase(const 
48                                                    AliCorrectionManagerBase& o)
49   : TObject(o),
50     fCorrections(),
51     fIsInit(o.fIsInit),
52     fRun(o.fRun), 
53     fSys(o.fSys), 
54     fSNN(o.fSNN), 
55     fField(o.fField), 
56     fMC(o.fMC), 
57     fSatellite(o.fSatellite), 
58     fDB(o.fDB),
59     fDebug(o.fDebug),
60     fFallBack(o.fFallBack)
61 {
62   fCorrections.SetOwner(false);
63   Int_t n = o.fCorrections.GetEntriesFast();
64   for (Int_t i = 0; i < n; i++) { 
65     fCorrections.AddAt(o.fCorrections.At(i), i);
66   }
67 }
68 //____________________________________________________________________
69 AliCorrectionManagerBase&
70 AliCorrectionManagerBase::operator=(const AliCorrectionManagerBase& o)
71 {
72   if (&o == this) return *this;
73
74   fIsInit       = o.fIsInit;
75   fRun          = o.fRun; 
76   fSys          = o.fSys; 
77   fSNN          = o.fSNN; 
78   fField        = o.fField; 
79   fMC           = o.fMC; 
80   fSatellite    = o.fSatellite;
81   fDB           = o.fDB;
82   fDebug        = o.fDebug;
83   fFallBack     = o.fFallBack;
84
85   fCorrections.Clear();
86   Int_t n = o.fCorrections.GetEntriesFast();
87   for (Int_t i = 0; i < n; i++) { 
88     fCorrections.AddAt(o.fCorrections.At(i), i);
89   }
90   return *this;
91 }
92
93 //____________________________________________________________________
94 AliCorrectionManagerBase::~AliCorrectionManagerBase()
95 {
96   // fCorrections.Delete();
97 }
98
99 //____________________________________________________________________
100 void
101 AliCorrectionManagerBase::Print(Option_t* option) const
102 {
103   char ind[gROOT->GetDirLevel()+1];
104   for (Int_t i = 0; i < gROOT->GetDirLevel(); i++) ind[i] = ' ';
105   ind[gROOT->GetDirLevel()] = '\0';
106
107   std::cout << ind << GetName() << ":\n"
108             << ind << "  Initialised:      " 
109             << (fIsInit ? "yes" : "no") << std::endl;
110   if (fIsInit) 
111     std::cout << ind << "  Run number:       " << fRun << "\n"
112               << ind << "  Collision system: " 
113               << AliForwardUtil::CollisionSystemString(fSys) << "\n"
114               << ind << "  Sqrt(s_NN):       "
115               << AliForwardUtil::CenterOfMassEnergyString(fSNN) << "\n"
116               << ind << "  Magnetic field:   " 
117               << AliForwardUtil::MagneticFieldString(fField) << "\n"
118               << ind << "  For simulations:  " << (fMC ? "yes" : "no") << "\n"
119               << ind << "  For satellites:   " 
120               << (fSatellite ? "yes" : "no") << std::endl;
121
122   TString opt(option);
123   opt.ToUpper();
124   if (!opt.Contains("R")) return;
125   
126   gROOT->IncreaseDirLevel();
127   Int_t n = fCorrections.GetEntriesFast();
128   for (Int_t id = 0; id < n; id++) { 
129     const Correction* c = GetCorrection(id);
130     c->Print(option);
131   }
132   gROOT->DecreaseDirLevel();  
133   
134 }
135
136 //____________________________________________________________________
137 void
138 AliCorrectionManagerBase::Browse(TBrowser* b)
139 {
140   b->Add(&fCorrections);
141 }
142
143 //____________________________________________________________________
144 void
145 AliCorrectionManagerBase::SetPrefix(const TString& prefix)
146 {
147   Int_t n = fCorrections.GetEntriesFast();
148   for (Int_t id = 0; id < n; id++) { 
149     Correction* c = GetCorrection(id);
150     const char* old = c->GetTitle();
151     TString     oldf(gSystem->BaseName(old));
152     c->SetFile(gSystem->ConcatFileName(prefix, oldf));
153   }
154 }
155
156 //____________________________________________________________________
157 Bool_t
158 AliCorrectionManagerBase::Store(TObject*     o,
159                                 ULong_t     runNo,
160                                 UShort_t    sys, 
161                                 UShort_t    sNN, 
162                                 Short_t     field, 
163                                 Bool_t      mc,
164                                 Bool_t      sat, 
165                                 const char* file,
166                                 const char* meth) const
167 {
168   Bool_t ret = false;
169   Int_t n = fCorrections.GetEntriesFast();
170   for (Int_t id = 0; id < n; id++) { 
171     const Correction* c = GetCorrection(id);
172     if (!o->IsA()->InheritsFrom(c->fCls)) continue;
173
174     ret = c->StoreIt(fDB, o, runNo, sys, sNN, field, mc, sat, file, meth);
175     break;
176   }
177   return ret;
178 }
179     
180 //____________________________________________________________________
181 Bool_t
182 AliCorrectionManagerBase::Append(const TString& addition, 
183                                  const TString& destination) const
184 {
185   if (addition.IsNull()) {
186     AliWarning("No addition specified");
187     return false;
188   }
189   if (destination.IsNull()) { 
190     AliWarning("No destination storage specified");
191     return false;
192   }
193   TFileMerger merger;
194   merger.SetPrintLevel(1);
195   merger.OutputFile(destination, "UPDATE");
196   merger.AddFile(addition);
197   if (!merger.PartialMerge()) {
198     AliInfoF("Failed to merge %s with %s", 
199              addition.Data(), destination.Data());
200     return false;
201   }
202   if (destination.BeginsWith("$OADB_PATH") ||
203       destination.BeginsWith("$ALICE_ROOT"))
204     AliInfoF("Now commit %s to subversion", destination.Data());
205   return true;
206 }
207   
208 //____________________________________________________________________
209 void
210 AliCorrectionManagerBase::RegisterCorrection(Int_t id, Correction* corr)
211 {
212   fCorrections.AddAtAndExpand(corr, id);
213 }
214
215 //____________________________________________________________________
216 void
217 AliCorrectionManagerBase::RegisterCorrection(Int_t id, 
218                                              const TString& tableName, 
219                                              const TString& fileName, 
220                                              TClass*        cls, 
221                                              UShort_t       fields,
222                                              Bool_t         enabled)
223 {
224   RegisterCorrection(id,new Correction(tableName,fileName,cls,fields,enabled));
225 }
226
227 //____________________________________________________________________
228 AliCorrectionManagerBase::Correction*
229 AliCorrectionManagerBase::GetCorrection(Int_t id)
230 {
231   if (id < 0 || id > fCorrections.GetEntriesFast()) return 0;
232   return static_cast<Correction*>(fCorrections.At(id));
233 }
234
235 //____________________________________________________________________
236 const AliCorrectionManagerBase::Correction*
237 AliCorrectionManagerBase::GetCorrection(Int_t id) const
238 {
239   if (id < 0 || id > fCorrections.GetEntriesFast()) return 0;
240   return static_cast<Correction*>(fCorrections.At(id));
241 }
242
243 //____________________________________________________________________
244 void 
245 AliCorrectionManagerBase::SetCorrectionFile(Int_t id, const TString& fileName)
246 {
247   Correction* c = GetCorrection(id);
248   if (!c) return;
249   c->SetFile(fileName);
250 }
251
252 //____________________________________________________________________
253 Int_t
254 AliCorrectionManagerBase::GetId(const TString& what) const
255 {
256   Int_t n = fCorrections.GetEntriesFast();
257   for (Int_t id = 0; id < n; id++) { 
258     const Correction* c = GetCorrection(id);
259     if (what.EqualTo(c->GetName(), TString::kIgnoreCase)) return id;
260   }
261   return -1;
262 }
263
264 //____________________________________________________________________
265 void
266 AliCorrectionManagerBase::EnableCorrection(Int_t id, Bool_t enable)
267 {
268   Correction* c = GetCorrection(id);
269   if (!c) { 
270     AliWarningF("Cannot enable non-existing correction at %d", id);
271     return;
272   }
273   c->fEnabled = enable;
274 }
275
276 //____________________________________________________________________
277 Int_t
278 AliCorrectionManagerBase::GetId(const TObject* obj) const
279 {
280   Int_t   n   = fCorrections.GetEntriesFast();
281   TClass* ocl = obj->IsA();
282   for (Int_t id = 0; id < n; id++) { 
283     const Correction* c = GetCorrection(id);
284     if (ocl->InheritsFrom(c->fCls)) return id;
285   }
286   return -1;
287 }
288 //____________________________________________________________________
289 TObject*
290 AliCorrectionManagerBase::Get(Int_t id)
291 {
292   Correction* c = GetCorrection(id);
293   if (!c) {
294     AliWarningF("Cannot find correction with id %d", id);
295     return 0;
296   }
297   return c->Get();
298 }
299 //____________________________________________________________________
300 const TObject*
301 AliCorrectionManagerBase::Get(Int_t id) const
302 {
303   const Correction* c = GetCorrection(id);
304   if (!c) {
305     AliWarningF("Cannot find correction with id %d", id);
306     return 0;
307   }
308   return c->Get();
309 }
310
311 //____________________________________________________________________
312 Bool_t
313 AliCorrectionManagerBase::InitCorrections(ULong_t    run, 
314                                           UShort_t   sys, 
315                                           UShort_t   sNN, 
316                                           Short_t    fld, 
317                                           Bool_t     mc, 
318                                           Bool_t     sat,
319                                           Bool_t     force)
320 {
321   if (force) fIsInit = false;
322   if (!CheckConditions(run, sys, sNN, fld, mc, sat)) return false;
323   if (!ReadCorrections(run, sys, sNN, fld, mc, sat)) return false;
324   fIsInit = true;
325
326   if (fDB) {
327     delete fDB;
328     fDB = 0;
329   }
330
331   return true;
332 }
333
334 //____________________________________________________________________
335 Bool_t
336 AliCorrectionManagerBase::CheckConditions(ULong_t    run, 
337                                           UShort_t   sys, 
338                                           UShort_t   sNN, 
339                                           Short_t    fld, 
340                                           Bool_t     mc, 
341                                           Bool_t     sat)
342 {
343   if (!fIsInit) return true;
344
345   AliInfo("We are already initialised - checking settings...");
346   Bool_t same = true;
347   if (fRun != run) {
348     same = false;
349   }
350   if (fSys != sys) { 
351     AliWarningF("Initialised collision system %s (%d) and "
352                 "passed same %s (%d) does not match", 
353                 AliForwardUtil::CollisionSystemString(fSys), fSys,
354                 AliForwardUtil::CollisionSystemString(sys), sys);
355     same = false;
356   }
357   if (TMath::Abs(fSNN - sNN) >= 10) {
358     AliWarningF("Initialised center of mass energy per nuclean "
359                 "%s (%d) and passed same %s (%d) does not match",
360                 AliForwardUtil::CenterOfMassEnergyString(fSNN), fSNN,
361                 AliForwardUtil::CenterOfMassEnergyString(sNN), sNN);
362     same = false;
363   }
364   if (fField != fld) {
365       AliWarningF("Initialied L3 magnetic field %s (%d) and passed "
366                   "same %s (%d) does not match", 
367                   AliForwardUtil::MagneticFieldString(fField), fField,
368                   AliForwardUtil::MagneticFieldString(fld), fld);
369       same = false;
370   }
371   if (fMC != mc) {
372     AliWarningF("Initialied data type (%s) and passed "
373                 "same (%s) does not match", 
374                 (fMC ? "MC" : "real"), (mc ? "MC" : "real"));
375     same = false;
376   }
377   if (fSatellite != sat) {
378     AliWarningF("Initialied collision ip type (%s) and passed "
379                 "same (%s) does not match", 
380                 (fSatellite ? "satellite" : "nominal"), 
381                 (sat ? "satellite" : "nominal"));
382     same = false;
383   }
384   if (!same) {
385     AliWarning("Intialised parameters and these are not the same " 
386                "- PROCEED WITH CAUTION!");
387   }
388   else
389     AliInfo("Initialized values consistent with data");
390   
391   return true;
392
393 }
394
395 //____________________________________________________________________
396 Bool_t
397 AliCorrectionManagerBase::ReadCorrection(Int_t      id,
398                                          ULong_t    run, 
399                                          UShort_t   sys, 
400                                          UShort_t   sNN, 
401                                          Short_t    fld, 
402                                          Bool_t     mc, 
403                                          Bool_t     sat)
404 {
405   if (!fDB) {
406     // We should always open the database, since we're not
407     // streamingthat object to disk.
408     fDB = new AliOADBForward;
409   }
410
411   Correction* c = GetCorrection(id);
412   if (!c->fEnabled) return true;
413   return c->ReadIt(fDB, run, sys, sNN, fld, mc, sat, fDebug, fFallBack);
414 }
415
416 //____________________________________________________________________
417 Bool_t
418 AliCorrectionManagerBase::ReadCorrections(ULong_t    run, 
419                                           UShort_t   sys, 
420                                           UShort_t   sNN, 
421                                           Short_t    fld, 
422                                           Bool_t     mc, 
423                                           Bool_t     sat)
424 {
425   if (fIsInit) return true;
426   if (fRun       == run && 
427       fSys       == sys && 
428       fField     == fld && 
429       fMC        == mc  && 
430       fSatellite == sat &&
431       TMath::Abs(fSNN - sNN) < 11) { 
432     // Already initialized for this - return
433     fIsInit = true;
434     return true;
435   }
436   if (!fDB) {
437     // We should always open the database, since we're not
438     // streamingthat object to disk.
439     fDB = new AliOADBForward;
440   }
441
442   fRun       = run;
443   fSys       = sys; 
444   fSNN       = sNN;
445   fField     = fld;
446   fMC        = mc;
447   fSatellite = sat;
448   Int_t  n   = fCorrections.GetEntriesFast();
449   Bool_t ret = true;
450   for (Int_t id = 0; id < n; id++) 
451     if (!ReadCorrection(id, run, sys, sNN, fld, mc, sat)) ret = false;
452   return ret;
453 }
454
455 //====================================================================
456 AliCorrectionManagerBase::Correction::Correction() 
457   : TNamed(), 
458     fCls(0), 
459     fClientCls(""),
460     fQueryFields(0), 
461     fEnabled(false), 
462     fLastEntry(),
463     fObject(0)
464 {}
465
466 //____________________________________________________________________
467 AliCorrectionManagerBase::Correction::Correction(const TString& tableName, 
468                                                  const TString& fileName, 
469                                                  TClass*        cls,
470                                                  UShort_t       fields,
471                                                  Bool_t         enabled) 
472   : TNamed(tableName, fileName), 
473     fCls(cls), 
474     fClientCls(cls->GetName()),
475     fQueryFields(fields), 
476     fEnabled(enabled), 
477     fLastEntry(""),
478     fObject(0)
479 {}
480
481 //____________________________________________________________________
482 AliCorrectionManagerBase::Correction::Correction(const Correction& o)
483   : TNamed(o), 
484     fCls(o.fCls), 
485     fClientCls(o.fClientCls),
486     fQueryFields(o.fQueryFields),
487     fEnabled(o.fEnabled), 
488     fLastEntry(o.fLastEntry),
489     fObject(o.fObject)
490 {}
491
492 //____________________________________________________________________
493 AliCorrectionManagerBase::Correction&
494 AliCorrectionManagerBase::Correction::operator=(const Correction& o)
495 {
496   if (&o == this) return *this;
497   SetName(o.GetName());
498   SetTitle(o.GetTitle());
499   fCls             = o.fCls;
500   //fClientCls       = o.fClientCls;
501   fQueryFields     = o.fQueryFields;
502   fEnabled         = o.fEnabled;
503   fLastEntry       = o.fLastEntry;
504   fObject          = o.fObject;
505   return *this;
506 }
507
508 //____________________________________________________________________
509 Bool_t
510 AliCorrectionManagerBase::Correction::ReadIt(AliOADBForward* db, 
511                                              ULong_t         run, 
512                                              UShort_t        sys, 
513                                              UShort_t        sNN, 
514                                              Short_t         fld, 
515                                              Bool_t          mc, 
516                                              Bool_t          sat,
517                                              Bool_t          vrb,
518                                              Bool_t          fallback)
519 {
520   if (!fEnabled) {
521     AliWarningF("Correction %s not enabled", GetName());
522     return 0;
523   }
524
525   // Assume failure 
526   fObject = 0;
527
528   // Massage fields according to settings 
529   if (!(fQueryFields & kRun))       run = kIgnoreValue;
530   if (!(fQueryFields & kSys))       sys = kIgnoreValue;
531   if (!(fQueryFields & kSNN))       sNN = kIgnoreValue;
532   if (!(fQueryFields & kField))     fld = AliOADBForward::kInvalidField; // kIgnoreField;
533   if (!(fQueryFields & kMC))        mc  = false;
534   if (!(fQueryFields & kSatellite)) sat = false;
535
536   // Check if table is open, and if not try to open it 
537   if (!db->FindTable(fName, true)) {
538     if (!db->Open(fTitle, fName, false, vrb, fallback)) {
539       AliWarningF("Failed to open table %s from %s", GetName(), GetTitle());
540       AliWarningF("content of %s for %s:", 
541                   gSystem->WorkingDirectory(), GetName());
542       gSystem->Exec("pwd; ls -l");
543       return false;
544     }
545   }
546   
547   // Query database 
548   AliOADBForward::Entry* e = db->Get(fName, run, AliOADBForward::kDefault, 
549                                      sys, sNN, fld, mc, sat);
550   // Check return value 
551   if (!e || !e->fData) {
552     AliWarningF("Failed to get %s from database in %s with "
553                 "run=%lu sys=%hu sNN=%hu fld=%hd %s %s", 
554                 GetName(), GetTitle(), run, sys, sNN, fld, 
555                 (mc ? "MC" : "real"), (sat ? "satellite" : "nominal"));
556     return false;
557   }
558
559   // Ge the returned data
560   TObject* o = e->fData;
561
562   const TClass* cl = TheClass();
563   // Check return class 
564   if (!o->IsA()->InheritsFrom(cl)) { 
565     AliWarningF("%p is not pointer to a %s object but a %s", 
566                 o, fCls->GetName(), o->ClassName());
567     return false;
568   }
569
570   // Success 
571   fObject    = o;
572   fLastEntry = e->GetTitle();
573
574   return true;
575 }
576
577 //____________________________________________________________________
578 Bool_t
579 AliCorrectionManagerBase::Correction::StoreIt(AliOADBForward* db, 
580                                               TObject*        obj,
581                                               ULong_t         run, 
582                                               UShort_t        sys, 
583                                               UShort_t        sNN, 
584                                               Short_t         fld, 
585                                               Bool_t          mc, 
586                                               Bool_t          sat,
587                                               const char*     file, 
588                                               const char*     meth) const
589 {
590   // Info("StoreIt", "Storing run=%lu sys=%hy sNN=%d fld=%d mc=%d sat=%d", 
591   //       run, sys, sNN, fld, mc, sat);
592   const TClass* cl = TheClass();
593
594   // Check value class 
595   if (!obj->IsA()->InheritsFrom(cl)) { 
596     AliWarningF("%p is not pointer to a %s object but a %s", 
597                 obj, cl->GetName(), obj->ClassName());
598     return false;
599   }
600
601   Bool_t          local    = file || !db;
602   TString         fileName = (local ? file : fTitle.Data());
603   AliOADBForward* tdb      = (local ? new AliOADBForward : db);
604   
605   // Try to open the table read/write 
606   if (!tdb->Open(fileName, Form("%s/%s", GetName(), meth), true, true)) {
607     AliWarningF("Failed to open table %s in %s", GetName(), fileName.Data());
608     return false;
609   }
610
611   // Massage fields according to settings 
612   if (!(fQueryFields & kRun))       run = kIgnoreValue;
613   if (!(fQueryFields & kSys))       sys = kIgnoreValue;
614   if (!(fQueryFields & kSNN))       sNN = kIgnoreValue;
615   if (!(fQueryFields & kField))     fld = AliOADBForward::kInvalidField; // kIgnoreField;
616   if (!(fQueryFields & kMC))        mc  = false;
617   if (!(fQueryFields & kSatellite)) sat = false;
618   
619   // Try to insert the object 
620   if (!tdb->Insert(fName, obj, run, sys, sNN, fld, mc, sat)) { 
621     AliWarningF("Failed to insert into %s off database in %s with "
622                 "run=%lu sys=%hu sNN=%hu fld=%hd %s %s", 
623                 GetName(), GetTitle(), run, sys, sNN, fld, 
624                 (mc ? "MC" : "real"), (sat ? "satellite" : "nominal"));
625     return false;
626   }
627
628   if (local) { 
629     tdb->Close();
630     delete tdb;
631
632     AliInfoF("Correction object %s written to DB in %s - merge this with "
633              "%s to store for good", obj->GetName(), fileName.Data(), 
634              GetTitle());
635   }
636
637   // Success 
638   return true;
639 }
640 //____________________________________________________________________
641 TObject*
642 AliCorrectionManagerBase::Correction::Get()   
643 {
644   if (!fEnabled) {
645     AliWarningF("Correction %s not enabled", GetName());
646     return 0;
647   }
648   return fObject;
649 }
650 //____________________________________________________________________
651 const TObject*
652 AliCorrectionManagerBase::Correction::Get() const
653 {
654   if (!fEnabled) {
655     AliWarningF("Correction %s not enabled", GetName());
656     return 0;
657   }
658   return fObject;
659 }
660
661 //____________________________________________________________________
662 const TClass*
663 AliCorrectionManagerBase::Correction::TheClass() const
664 {
665   if (fCls) return fCls;
666   if (fClientCls.IsNull()) { 
667     AliErrorF("No class name set for correction %s", GetName());
668     return 0;
669   }
670   fCls = gROOT->GetClass(fClientCls);
671   if (!fCls) { 
672     AliErrorF("Couldn't get class %s for correction %s", 
673               fClientCls.Data(), GetName());
674     return 0;
675   }
676   return fCls;
677 }
678
679 //____________________________________________________________________
680 void
681 AliCorrectionManagerBase::Correction::Print(Option_t* option) const
682 {
683   char ind[gROOT->GetDirLevel()+1];
684   for (Int_t i = 0; i < gROOT->GetDirLevel(); i++) ind[i] = ' ';
685   ind[gROOT->GetDirLevel()] = '\0';
686
687   std::cout << ind << GetName() << ":  " << (fEnabled ? "en" : "dis") 
688             << "abled" << std::endl;
689   if (!fEnabled) return;
690
691   TString flds;
692   if (fQueryFields & kRun)       flds.Append("run");
693   if (fQueryFields & kSys)       flds.Append("|sys");
694   if (fQueryFields & kSNN)       flds.Append("|sNN");
695   if (fQueryFields & kField)     flds.Append("|field");
696   if (fQueryFields & kMC)        flds.Append("|MC");
697   if (fQueryFields & kSatellite) flds.Append("|Satellite");
698   if (flds.BeginsWith("|")) flds.Remove(0,1);
699
700   const TClass* cl = TheClass();
701
702   std::cout << "   Path:            " << GetTitle() << "\n"
703             << "   Data class:      " << cl->GetName() << "\n"
704             << "   Query fields:    " << flds << std::endl;
705   
706   if (fObject && !fLastEntry.IsNull()) 
707     std::cout << "   Entry:           " << fLastEntry << std::endl;
708   
709   TString opt(option);
710   opt.ToUpper();
711   if (!opt.Contains("D") || !fObject) return;
712
713   gROOT->IncreaseDirLevel();
714   fObject->Print();
715   gROOT->DecreaseDirLevel();
716 }
717
718 //____________________________________________________________________
719 void
720 AliCorrectionManagerBase::Correction::Browse(TBrowser* b)
721 {
722   b->Add(const_cast<TClass*>(fCls), "Class");
723   TString flds;
724   if (fQueryFields & kRun)       flds.Append("run");
725   if (fQueryFields & kSys)       flds.Append("|sys");
726   if (fQueryFields & kSNN)       flds.Append("|sNN");
727   if (fQueryFields & kField)     flds.Append("|field");
728   if (fQueryFields & kMC)        flds.Append("|MC");
729   if (fQueryFields & kSatellite) flds.Append("|Satellite");
730   if (flds.BeginsWith("|")) flds.Remove(0,1);
731
732   b->Add(new TObjString(flds), "Query fields");
733   b->Add(new TParameter<bool>("Enabled", fEnabled));
734   b->Add(new TObjString(fLastEntry), "Entry");
735   if (fObject) b->Add(fObject);
736 }
737 //
738 // EOF
739 //