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