]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/keeps/AliForwardCorrectionManagerOADB.cxx
Mega commit of many changes to PWGLFforward
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / keeps / AliForwardCorrectionManagerOADB.cxx
1 //
2 // Manager (singleton) of corrections 
3 // 
4 #include "AliForwardCorrectionManagerOADB.h"
5 #include "AliFMDCorrSecondaryMap.h"
6 #include "AliFMDCorrDoubleHit.h"
7 #include "AliFMDCorrELossFit.h"
8 #include "AliFMDCorrVertexBias.h"
9 #include "AliFMDCorrMergingEfficiency.h"
10 #include "AliFMDCorrAcceptance.h"
11 #include "AliForwardUtil.h"
12 #include "AliOADBForward.h"
13 #include <TString.h>
14 #include <AliLog.h>
15 #include <TFile.h>
16 #include <TSystem.h>
17 #include <TBrowser.h>
18 #include <TROOT.h>
19 #include <iostream>
20 #include <iomanip>
21     
22 //____________________________________________________________________
23 AliForwardCorrectionManagerOADB* AliForwardCorrectionManagerOADB::fgInstance= 0;
24 const char* AliForwardCorrectionManagerOADB::fgkSecondaryMapSkel = "secondary";
25 const char* AliForwardCorrectionManagerOADB::fgkDoubleHitSkel    = "doublehit";
26 const char* AliForwardCorrectionManagerOADB::fgkELossFitsSkel    = "elossfits";
27 const char* AliForwardCorrectionManagerOADB::fgkVertexBiasSkel   = "vertexbias";
28 const char* AliForwardCorrectionManagerOADB::fgkMergingEffSkel   = "merging";
29 const char* AliForwardCorrectionManagerOADB::fgkAcceptanceSkel   = "acceptance";
30
31 #define PREFIX  "$(ALICE_ROOT)/OADB/PWGLF/FORWARD/CORRECTIONS/data/"
32 #define DB_NAME "fmd_corrections.root"
33
34 //____________________________________________________________________
35 AliForwardCorrectionManagerOADB& AliForwardCorrectionManagerOADB::Instance()
36 {
37   // 
38   // Access to the singleton object 
39   // 
40   // Return:
41   //    Reference to the singleton object 
42   //
43   if (!fgInstance) fgInstance= new AliForwardCorrectionManagerOADB;
44   return *fgInstance;
45 }
46
47 //____________________________________________________________________
48 AliForwardCorrectionManagerOADB::AliForwardCorrectionManagerOADB()
49   : TObject(), 
50     fInit(kFALSE),
51     fRunNo(0),
52     fSys(0),
53     fSNN(0),
54     fField(999),
55     fMC(false), 
56     fSat(false),
57     fELossFitsPath(PREFIX DB_NAME),
58     fMergingEffPath(PREFIX DB_NAME), 
59     fSecondaryMapPath(PREFIX DB_NAME),
60     fDoubleHitPath(PREFIX DB_NAME),
61     fVertexBiasPath(PREFIX DB_NAME),
62     fAcceptancePath(PREFIX DB_NAME),
63     fELossFit(0),
64     fSecondaryMap(0),
65     fDoubleHit(0),
66     fVertexBias(0),
67     fMergingEfficiency(0),
68     fAcceptance(0),
69     fDB(0)
70 {
71   // 
72   // Default constructor 
73   //
74 }
75 //____________________________________________________________________
76 AliForwardCorrectionManagerOADB::AliForwardCorrectionManagerOADB(const AliForwardCorrectionManagerOADB& o)
77   : TObject(o),
78     fInit(o.fInit),
79     fRunNo(o.fRunNo),
80     fSys(o.fSys),
81     fSNN(o.fSNN),
82     fField(o.fField),
83     fMC(o.fMC),
84     fSat(o.fSat),
85     fELossFitsPath(o.fELossFitsPath),
86     fMergingEffPath(o.fMergingEffPath), 
87     fSecondaryMapPath(o.fSecondaryMapPath),
88     fDoubleHitPath(o.fDoubleHitPath),
89     fVertexBiasPath(o.fVertexBiasPath),
90     fAcceptancePath(o.fAcceptancePath),
91     fELossFit(o.fELossFit),
92     fSecondaryMap(o.fSecondaryMap),
93     fDoubleHit(o.fDoubleHit),
94     fVertexBias(o.fVertexBias),
95     fMergingEfficiency(o.fMergingEfficiency),
96   fAcceptance(o.fAcceptance),
97   fDB(0)
98 {
99   // 
100   // Copy constructor 
101   // 
102   // Parameters:
103   //    o Object to copy from 
104   //
105 }
106 //____________________________________________________________________
107 AliForwardCorrectionManagerOADB&
108 AliForwardCorrectionManagerOADB::operator=(const AliForwardCorrectionManagerOADB& o)
109 {
110   // 
111   // Assignment operator 
112   // 
113   // Parameters:
114   //    o Object to assign from 
115   // 
116   // Return:
117   //    Reference to this object 
118   //
119   fInit             = o.fInit;
120   fRunNo            = o.fRunNo;
121   fSys              = o.fSys;
122   fSNN              = o.fSNN;
123   fField            = o.fField;
124   fMC               = o.fMC;
125   fSat              = o.fSat;
126   fELossFitsPath    = o.fELossFitsPath;
127   fMergingEffPath   = o.fMergingEffPath;
128   fSecondaryMapPath = o.fSecondaryMapPath;
129   fDoubleHitPath    = o.fDoubleHitPath;
130   fVertexBiasPath   = o.fVertexBiasPath;
131   fAcceptancePath   = o.fAcceptancePath;
132   fELossFit         = o.fELossFit;
133   fSecondaryMap     = o.fSecondaryMap;
134   fDoubleHit        = o.fDoubleHit;
135   fVertexBias       = o.fVertexBias;
136   fMergingEfficiency= o.fMergingEfficiency;
137   fAcceptance       = o.fAcceptance;
138   fDB               = o.fDB;
139   return *this;
140 }
141
142 //____________________________________________________________________
143 void
144 AliForwardCorrectionManagerOADB::SetPrefix(const char* prefix)
145 {
146   /** 
147    *
148    * @param prefix Prefix to correction objects. 
149    */
150   fELossFitsPath    = Form("%s/%s", prefix, DB_NAME);
151   fMergingEffPath   = Form("%s/%s", prefix, DB_NAME); 
152   fSecondaryMapPath = Form("%s/%s", prefix, DB_NAME);
153   fDoubleHitPath    = Form("%s/%s", prefix, DB_NAME);
154   fVertexBiasPath   = Form("%s/%s", prefix, DB_NAME);
155   fAcceptancePath   = Form("%s/%s", prefix, DB_NAME);
156   
157 }
158 //____________________________________________________________________
159 void
160 AliForwardCorrectionManagerOADB::SetFile(ECorrection what, const char* filename)
161 {
162   /** 
163    * Set the file directory for a type 
164    * 
165    * @param what     Type 
166    * @param dirname  Directory name 
167    */
168   TString *path = 0;
169   if      (what & kSecondaryMap)        path = &fSecondaryMapPath;
170   else if (what & kDoubleHit)           path = &fDoubleHitPath;
171   else if (what & kELossFits)           path = &fELossFitsPath;
172   else if (what & kVertexBias)          path = &fVertexBiasPath;
173   else if (what & kMergingEfficiency)   path = &fMergingEffPath;
174   else if (what & kAcceptance)          path = &fAcceptancePath;
175   else { 
176     AliWarning(Form("No such path defined for 0x%02x", what));
177     return;
178   }
179   if (!path) {
180     AliWarning(Form("Couldn't find string for path 0x%02x", what));
181     return;
182   }
183   *path = filename;
184 }
185
186 //____________________________________________________________________
187 Bool_t
188 AliForwardCorrectionManagerOADB::Init(ULong_t     runNo, 
189                                       const char* sys, 
190                                       Float_t     sNN, 
191                                       Float_t     field,
192                                       Bool_t      mc,
193                                       Bool_t      sat,
194                                       UInt_t      what,
195                                       Bool_t      force)
196 {
197   // 
198   // Read in correction based on passed parameters
199   // 
200   // Parameters:
201   //    collisionSystem Collision system string 
202   //    cmsNN           Center of mass energy per nucleon pair [GeV]
203   //    field           Magnetic field [kG]
204   //    mc              Monte-carlo switch
205   //    what            What to read in 
206   //    force           Force (re-)reading of specified things
207   // 
208   // Return:
209   //    true on success
210   //
211   UShort_t col = AliForwardUtil::ParseCollisionSystem(sys);
212   // AliInfo(Form("Initialising with cms='%s', sNN=%fGeV field=%fkG", 
213   //           cms, sNN, field));
214   return Init(runNo, col, 
215               AliForwardUtil::ParseCenterOfMassEnergy(col, sNN),
216               AliForwardUtil::ParseMagneticField(field), 
217               mc, sat, what, force);
218 }
219
220 //____________________________________________________________________
221 Bool_t
222 AliForwardCorrectionManagerOADB::Init(ULong_t  runNo, 
223                                       UShort_t sys, 
224                                       UShort_t sNN, 
225                                       Short_t  field,
226                                       Bool_t   mc,
227                                       Bool_t   sat,
228                                       UInt_t   what,
229                                       Bool_t   force)
230 {
231   // 
232   // Read in corrections based on the parameters given 
233   // 
234   // Parameters:
235   //    collisionSystem Collision system
236   //    cmsNN           Center of mass energy per nuclean pair [GeV]
237   //    field           Magnetic field setting [kG]
238   //    mc              Monte-carlo switch
239   //    what            What to read in. 
240   //    force           Force (re-)reading of specified things
241   // 
242   // Return:
243   //    
244   //
245   if (force) fInit = kFALSE;
246   if (!fDB) { 
247     // We should always open the database, since we're not
248     // streamingthat object to disk.
249     fDB = new AliOADBForward;
250   }
251   if (fInit) {
252     // Check that the initialisation and the passed parameters 
253     // match - if not give an error but continue - this allows 
254     // users to override particular settings. 
255     
256     AliInfo("We are already initialised - checking settings...");
257     Bool_t same = true;
258     if (fRunNo != runNo) {
259       same = false;
260     }
261     if (fSys != sys) { 
262       AliWarningF("Initialised collision system %s (%d) and "
263                   "passed same %s (%d) does not match", 
264                   AliForwardUtil::CollisionSystemString(fSys), fSys,
265                   AliForwardUtil::CollisionSystemString(sys), sys);
266       same = false;
267     }
268     if (TMath::Abs(fSNN - sNN) >= 10) {
269       AliWarningF("Initialised center of mass energy per nuclean "
270                   "%s (%d) and passed same %s (%d) does not match",
271                   AliForwardUtil::CenterOfMassEnergyString(fSNN), fSNN,
272                   AliForwardUtil::CenterOfMassEnergyString(sNN), sNN);
273       same = false;
274     }
275     if (fField != field) {
276       AliWarningF("Initialied L3 magnetic field %s (%d) and passed "
277                   "same %s (%d) does not match", 
278                   AliForwardUtil::MagneticFieldString(fField), fField,
279                   AliForwardUtil::MagneticFieldString(field), field);
280       same = false;
281     }
282     if (fMC != mc) {
283       AliWarningF("Initialied data type (%s) and passed "
284                   "same (%s) does not match", 
285                   (fMC ? "MC" : "real"), (mc ? "MC" : "real"));
286       same = false;
287     }
288     if (fSat != sat) {
289       AliWarningF("Initialied collision ip type (%s) and passed "
290                   "same (%s) does not match", 
291                   (fSat ? "satelitte" : "nominal"), 
292                   (sat ? "satellite" : "nominal"));
293       same = false;
294     }
295     if (!same) {
296       AliWarning("Intialised parameters and these are not the same " 
297                  "- PROCEED WITH CAUTION!");
298     }
299     else
300       AliInfo("Initialized values consistent with data");
301     
302     return kTRUE;
303   }
304
305   Bool_t ret = kTRUE;
306   if (fRunNo == runNo && 
307       fSys   == sys   && 
308       fField == field &&
309       fMC    == mc    && 
310       fSat   == sat   &&
311       TMath::Abs(fSNN - sNN) < 10) {
312     // We're already initialised for these settings - do nothing and return
313     fInit = kTRUE;
314     return ret;
315   }
316   // Set cached parameters 
317   fRunNo = runNo;
318   fSys   = sys;
319   fSNN   = sNN;
320   fField = field;
321   fMC    = mc;
322   fSat   = sat;
323
324   // AliInfo(Form("Initialising with cms=%d, sNN=%dGeV field=%dkG", 
325   //           cms, sNN, field));
326   // Read secondary map if requested 
327   if (what & kSecondaryMap) {
328     if (!ReadSecondaryMap(runNo, sys, sNN, field, sat)) {
329       AliWarningF("Failed to read in secondary map for "
330                   "run=%lu, sys=%hu, sNN=%huGeV, field=%hdkG, %s", 
331                   runNo, sys, sNN, field, (sat ? "satellite" : "nominal"));
332       ret = kFALSE;
333     }
334   }
335   // Read double hit if requested 
336   if (what & kDoubleHit) {
337     if (!ReadDoubleHit(runNo, sys, sNN, field, sat)) {
338       AliWarningF("Failed to read in double hit correction for "
339                   "run=%lu, sys=%hu, sNN=%huGeV, field=%hdkG, %s", 
340                   runNo, sys, sNN, field, (sat ? "satellite" : "nominal"));
341       ret = kFALSE;
342     }
343   }
344   // Read energy loss fits if requested 
345   if (what & kELossFits) {
346     if (!ReadELossFits(runNo, sys, sNN, field, mc, sat)) {
347       AliWarningF("Failed to read in energy loss fits for "
348                   "run=%lu, sys=%hu, sNN=%huGeV, field=%hdkG, %s, %s", 
349                   runNo, sys, sNN, field, mc ? "MC" : "real", 
350                   (sat ? "satellite" : "nominal"));
351       ret = kFALSE;
352     }
353   }
354   // Read acceptance correction if requested 
355   if (what & kAcceptance) {
356     if (!ReadAcceptance(runNo, sys, sNN, sat)) {
357       AliWarningF("Failed to read in acceptance for "
358                   "run=%lu, sys=%hu, sNN=%huGeV, %s", 
359                   runNo, sys, sNN, (sat ? "satellite" : "nominal"));
360       ret = kFALSE;
361     }
362   }
363   // Read event selection efficiencies if requested 
364   if (what & kVertexBias) {
365     if (!ReadVertexBias(runNo, sys, sNN, field, sat)) {
366       AliWarningF("Failed to read in vertex bias correction for "
367                   "run=%lu, sys=%hu, sNN=%huGeV, field=%hdkG, %s", 
368                   runNo, sys, sNN, field, (sat ? "satellite" : "nominal"));
369       ret = kFALSE;
370     }
371   }
372   // Read merging efficiencies if requested 
373   if (what & kMergingEfficiency) {
374     if (!ReadMergingEfficiency(runNo, sys, sNN, field, sat)) {
375       AliWarningF("Failed to read in hit merging efficiency for "
376                   "run=%lu, sys=%hu, sNN=%huGeV, field=%hdkG, %s", 
377                   runNo, sys, sNN, field, (sat ? "satellite" : "nominal"));
378       ret = kFALSE;
379     }
380   }
381   fInit = kTRUE;
382   return ret;
383 }
384 //____________________________________________________________________
385 const TAxis* 
386 AliForwardCorrectionManagerOADB::GetEtaAxis() const
387 {
388   if (!fSecondaryMap) return 0;
389   return &(fSecondaryMap->GetEtaAxis());
390 }
391 //____________________________________________________________________
392 const TAxis* 
393 AliForwardCorrectionManagerOADB::GetVertexAxis() const
394 {
395   if (!fSecondaryMap) return 0;
396   return &(fSecondaryMap->GetVertexAxis());
397 }
398
399 //____________________________________________________________________
400 const Char_t* 
401 AliForwardCorrectionManagerOADB::GetFileName(ECorrection what) const
402 {
403   // 
404   // Get the path to the specified object 
405   // 
406   // Parameters:
407   //    what  Which stuff to get the path for 
408   //    sys   Collision system
409   //    sNN   Center of mass energy [GeV]
410   //    field Magnetic field in the L3 magnet [kG]
411   //    mc    Whether the correction objects should be valid for MC
412   // 
413   // Return:
414   //    The full path or null 
415   //
416   return gSystem->BaseName(GetFilePath(what));
417 }
418 //____________________________________________________________________
419 const Char_t*
420 AliForwardCorrectionManagerOADB::GetFileDir(ECorrection what) const
421 {
422   // 
423   // Get the path to the specified object 
424   // 
425   // Parameters:
426   //    what  Which stuff to get the path for 
427   // 
428   // Return:
429   //    The full path or null 
430   //
431   return gSystem->DirName(GetFilePath(what));
432 }
433
434 //____________________________________________________________________
435 const TString&
436 AliForwardCorrectionManagerOADB::GetFilePath(ECorrection what) const
437 {
438   // 
439   // Get the full path to the object.  Note, the manager must be
440   // initialised for this to work
441   // 
442   // Parameters:
443   //    what Which stuff to get the path for 
444   // 
445   // Return:
446   //    The full path or null
447   //
448   switch (what) {    
449   case kSecondaryMap:        return fSecondaryMapPath;
450   case kDoubleHit:           return fDoubleHitPath;
451   case kELossFits:           return fELossFitsPath;
452   case kVertexBias:          return fVertexBiasPath;
453   case kMergingEfficiency:   return fMergingEffPath;
454   case kAcceptance:          return fAcceptancePath;
455   default: break;
456   }
457   static TString null;
458   return null;
459 }
460
461
462 //____________________________________________________________________
463 const Char_t*
464 AliForwardCorrectionManagerOADB::GetObjectName(ECorrection what) const
465 {
466   // 
467   // Get the object name corresponding to correction type 
468   // 
469   // Parameters:
470   //    what Correction 
471   // 
472   // Return:
473   //    Object name or null
474   //
475   switch (what) {
476   case kSecondaryMap:       return fgkSecondaryMapSkel;
477   case kDoubleHit:          return fgkDoubleHitSkel;
478   case kELossFits:          return fgkELossFitsSkel;
479   case kVertexBias:         return fgkVertexBiasSkel;
480   case kMergingEfficiency:  return fgkMergingEffSkel;
481   case kAcceptance:         return fgkAcceptanceSkel;
482   default: break;
483   }
484   return 0;
485 }
486
487 //____________________________________________________________________
488 const TClass*
489 AliForwardCorrectionManagerOADB::GetObjectClass(ECorrection what) const
490 {
491   // 
492   // Get the object name corresponding to correction type 
493   // 
494   // Parameters:
495   //    what Correction 
496   // 
497   // Return:
498   //    Object name or null
499   //
500   switch (what) {
501   case kSecondaryMap:       return AliFMDCorrSecondaryMap::Class();
502   case kDoubleHit:          return AliFMDCorrDoubleHit::Class();
503   case kELossFits:          return AliFMDCorrELossFit::Class();
504   case kVertexBias:         return AliFMDCorrVertexBias::Class();
505   case kMergingEfficiency:  return AliFMDCorrMergingEfficiency::Class();
506   case kAcceptance:         return AliFMDCorrAcceptance::Class();
507   default: break;
508   }
509   return 0;
510 }
511 //____________________________________________________________________
512 AliForwardCorrectionManagerOADB::ECorrection
513 AliForwardCorrectionManagerOADB::GetObjectType(const TString& what) const
514 {
515   TString w(what);
516   w.ToLower();
517   
518   if      (w.EqualTo(GetObjectName(kSecondaryMap)))       
519     return kSecondaryMap;
520   else if (w.EqualTo(GetObjectName(kDoubleHit)))          
521     return kDoubleHit;
522   else if (w.EqualTo(GetObjectName(kELossFits)))          
523     return kELossFits;
524   else if (w.EqualTo(GetObjectName(kVertexBias)))         
525     return kVertexBias;
526   else if (w.EqualTo(GetObjectName(kMergingEfficiency)))  
527     return kMergingEfficiency;
528   else if (w.EqualTo(GetObjectName(kAcceptance)))         
529     return kAcceptance;
530   
531   return kAll;
532 }
533 //____________________________________________________________________
534 AliForwardCorrectionManagerOADB::ECorrection
535 AliForwardCorrectionManagerOADB::GetObjectType(const TObject* obj) const
536 {
537   TClass* cl = obj->IsA();
538   
539   if      (cl->InheritsFrom(GetObjectClass(kSecondaryMap)))       
540     return kSecondaryMap;
541   else if (cl->InheritsFrom(GetObjectClass(kDoubleHit)))          
542     return kDoubleHit;
543   else if (cl->InheritsFrom(GetObjectClass(kELossFits)))          
544     return kELossFits;
545   else if (cl->InheritsFrom(GetObjectClass(kVertexBias)))         
546     return kVertexBias;
547   else if (cl->InheritsFrom(GetObjectClass(kMergingEfficiency)))  
548     return kMergingEfficiency;
549   else if (cl->InheritsFrom(GetObjectClass(kAcceptance)))
550     return kAcceptance;
551   
552   return kAll;
553 }
554
555 //____________________________________________________________________
556 TObject*
557 AliForwardCorrectionManagerOADB::GetObject(ECorrection what, 
558                                            ULong_t     runNo,
559                                            UShort_t    sys, 
560                                            UShort_t    sNN, 
561                                            Short_t     field,
562                                            Bool_t      mc,
563                                            Bool_t      sat) const
564 {
565   // 
566   // Get the path to the specified object 
567   // 
568   // Parameters:
569   //    what  Which stuff to get the path for 
570   //    sys   Collision system
571   //    sNN   Center of mass energy [GeV]
572   //    field Magnetic field in the L3 magnet [kG]
573   //    mc    Whether the correction objects should be valid for MC
574   // 
575   // Return:
576   //    The full path or null 
577   //
578   if (!fDB) {
579     AliWarning("Database not opened");
580     return 0;
581   }
582   TString tableName = GetObjectName(what);
583   if (!fDB->FindTable(tableName, true)) {
584     if (!fDB->Open(GetFilePath(what), tableName, false, true)) { 
585       AliWarningF("Failed to open table %s from %s", 
586                   tableName.Data(), GetFilePath(what).Data());
587       return 0;
588     }
589   }
590
591   TObject* o = fDB->GetData(GetObjectName(what), runNo,
592                             AliOADBForward::kDefault, 
593                             sys, sNN, field, mc, sat);
594   
595   return o;
596 }
597 //____________________________________________________________________
598 TObject*
599 AliForwardCorrectionManagerOADB::GetObject(ECorrection what) const
600 {
601   // 
602   // Get the object that contaisn the specified correction
603   // 
604   // Parameters:
605   //    what Which object to get
606   // 
607   // Return:
608   //    The object or null
609   //
610   if (!fInit) { 
611     AliWarning("Corrections manager initialised, do a forced Init(...)");
612     return 0;
613   }
614   
615   return GetObject(what, fRunNo, fSys, fSNN, fField, fMC, fSat);
616 }
617
618 #define CHECK_TYPE(O,RET,CL) do {       \
619   RET = dynamic_cast<CL*>(O);           \
620   if (!O) AliWarningF("%p is not a pointer to a %s object, but a %s", \
621                       O, #CL, O->ClassName()); } while (false)
622
623 //____________________________________________________________________
624 Bool_t 
625 AliForwardCorrectionManagerOADB::ReadSecondaryMap(ULong_t  runNo, 
626                                                   UShort_t sys, 
627                                                   UShort_t sNN, 
628                                                   Short_t  field,
629                                                   Bool_t   sat)
630 {
631   // 
632   // Read in the secondary map 
633   // 
634   // Parameters:
635   //    sys   Collision system
636   //    sNN   Center of mass energy [GeV]
637   //    field Magnetic field in the L3 magnet [kG]
638   // 
639   // Return:
640   //    True on success, false otherwise 
641   //
642   if (fInit) { 
643     AliWarning("Corrections manager initialised, do a forced Init(...)");
644     return kFALSE;
645   }
646
647   TObject* o = GetObject(kSecondaryMap, runNo, sys, sNN, field, false, sat);
648   if (!o) return kFALSE;
649
650   CHECK_TYPE(o, fSecondaryMap, AliFMDCorrSecondaryMap);
651   if (!fSecondaryMap) return kFALSE;
652
653   return kTRUE;
654 }
655 //____________________________________________________________________
656 Bool_t 
657 AliForwardCorrectionManagerOADB::ReadDoubleHit(ULong_t  runNo,
658                                                UShort_t sys, 
659                                                UShort_t sNN, 
660                                                Short_t  field, 
661                                                Bool_t   sat)
662 {
663   // 
664   // Read in the double hit correction
665   // 
666   // Parameters:
667   //    sys   Collision system
668   //    sNN   Center of mass energy [GeV]
669   //    field Magnetic field in the L3 magnet [kG]
670   // 
671   // Return:
672   //    True on success, false otherwise 
673   //
674   if (fInit) { 
675     AliWarning("Corrections manager initialised, do a forced Init(...)");
676     return kFALSE;
677   }
678
679   TObject* o = GetObject(kDoubleHit, runNo, sys, sNN, field, false, sat);
680   if (!o) return kFALSE;
681
682   CHECK_TYPE(o, fDoubleHit, AliFMDCorrDoubleHit);
683   if (!fDoubleHit) return kFALSE;
684
685   return kTRUE;
686 }
687
688 //____________________________________________________________________
689 Bool_t 
690 AliForwardCorrectionManagerOADB::ReadELossFits(ULong_t  runNo,
691                                                UShort_t sys, 
692                                                UShort_t sNN, 
693                                                Short_t  field, 
694                                                Bool_t   mc,
695                                                Bool_t   sat)
696 {
697   // 
698   // Read in the energy loss fits 
699   // 
700   // Parameters:
701   //    sys   Collision system
702   //    sNN   Center of mass energy [GeV]
703   //    field Magnetic field in the L3 magnet [kG]
704   //    mc    Whether the correction objects should be valid for MC
705   // 
706   // Return:
707   //    True on success, false otherwise 
708   //
709   if (fInit) { 
710     AliWarning("Corrections manager initialised, do a forced Init(...)");
711     return kFALSE;
712   }
713
714   TObject* o = GetObject(kELossFits, runNo, sys, sNN, field, mc, sat);
715   if (!o) return kFALSE;
716
717   CHECK_TYPE(o, fELossFit, AliFMDCorrELossFit);
718   if (!fELossFit) return kFALSE;
719
720   return kTRUE;
721 }
722
723 //____________________________________________________________________
724 Bool_t 
725 AliForwardCorrectionManagerOADB::ReadVertexBias(ULong_t  runNo, 
726                                                 UShort_t sys, 
727                                                 UShort_t sNN, 
728                                                 Short_t  field, 
729                                                 Bool_t   sat)
730 {
731   // 
732   // Read in the event selection efficiency 
733   // 
734   // Parameters:
735   //    sys   Collision system
736   //    sNN   Center of mass energy [GeV]
737   //    field Magnetic field in the L3 magnet [kG]
738   // 
739   // Return:
740   //    True on success, false otherwise 
741   //
742   if (fInit) { 
743     AliWarning("Corrections manager initialised, do a forced Init(...)");
744     return kFALSE;
745   }
746
747   TObject* o = GetObject(kVertexBias, runNo, sys, sNN, field, false, sat);
748   if (!o) return kFALSE;
749
750   CHECK_TYPE(o, fVertexBias, AliFMDCorrVertexBias);
751   if (!fVertexBias) return kFALSE;
752
753   return kTRUE;
754 }
755
756 //____________________________________________________________________
757 Bool_t 
758 AliForwardCorrectionManagerOADB::ReadMergingEfficiency(ULong_t runNo, 
759                                                        UShort_t sys, 
760                                                        UShort_t sNN, 
761                                                        Short_t  field, 
762                                                        Bool_t   sat)
763 {
764   // 
765   // Read in the merging efficiency 
766   // 
767   // Parameters:
768   //    sys   Collision system
769   //    sNN   Center of mass energy [GeV]
770   //    field Magnetic field in the L3 magnet [kG]
771   // 
772   // Return:
773   //    True on success, false otherwise 
774   //
775   if (fInit) { 
776     AliWarning("Corrections manager initialised, do a forced Init(...)");
777     return kFALSE;
778   }
779
780   TObject* o = GetObject(kMergingEfficiency, runNo, sys, sNN, field, false,sat);
781   if (!o) return kFALSE;
782
783   CHECK_TYPE(o, fMergingEfficiency, AliFMDCorrMergingEfficiency);
784   if (!fMergingEfficiency) return kFALSE;
785
786   return kTRUE;
787 }
788
789 //____________________________________________________________________
790 Bool_t 
791 AliForwardCorrectionManagerOADB::ReadAcceptance(ULong_t  runNo, 
792                                                 UShort_t sys, 
793                                                 UShort_t sNN, 
794                                                 Bool_t   sat)
795 {
796   // 
797   // Read in the event selection efficiency 
798   // 
799   // Parameters:
800   //    sys   Collision system
801   //    sNN   Center of mass energy [GeV]
802   //    field Magnetic field in the L3 magnet [kG]
803   // 
804   // Return:
805   //    True on success, false otherwise 
806   //
807   if (fInit) { 
808     AliWarning("Corrections manager initialised, do a forced Init(...)");
809     return kFALSE;
810   }
811
812   TObject* o = GetObject(kAcceptance, runNo, sys, sNN, 0, false, sat);
813   if (!o) return kFALSE;
814
815   CHECK_TYPE(o, fAcceptance, AliFMDCorrAcceptance);
816   if (!fAcceptance) return kFALSE;
817
818   return kTRUE;
819 }
820 //____________________________________________________________________
821 void
822 AliForwardCorrectionManagerOADB::Print(Option_t* option) const
823 {
824   // 
825   // Print stuff 
826   //
827   char ind[gROOT->GetDirLevel()+1];
828   for (Int_t i = 0; i < gROOT->GetDirLevel(); i++) ind[i] = ' ';
829   ind[gROOT->GetDirLevel()] = '\0';
830
831   std::cout << ind << "AliForwardCorrectionManagerOADB:\n"
832             << ind << "  Initialised:      " 
833             << (fInit ? "yes" : "no") << std::endl;
834   if (fInit) 
835     std::cout << ind << "  Run number:       " << fRunNo << "\n"
836               << ind << "  Collision system: " 
837               << AliForwardUtil::CollisionSystemString(fSys) << "\n"
838               << ind << "  Sqrt(s_NN):       "
839               << AliForwardUtil::CenterOfMassEnergyString(fSNN) << "\n"
840               << ind << "  Magnetic field:   " 
841               << AliForwardUtil::MagneticFieldString(fField) << "\n"
842               << ind << "  For simulations:  " << (fMC ? "yes" : "no") << "\n"
843               << ind << "  For satellites:   " << (fSat ? "yes" : "no") << "\n"
844               << std::endl;
845   std::cout << ind << "  Paths:\n" 
846             << ind << "    ELoss Fits:     " << fELossFitsPath << "\n"
847             << ind << "    Merging eff.:   " << fMergingEffPath << "\n"
848             << ind << "    Secondary maps: " << fSecondaryMapPath << "\n"
849             << ind << "    2-hit corr.:    " << fSecondaryMapPath << "\n"
850             << ind << "    Vertex bias:    " << fVertexBiasPath << "\n"
851             << ind << "    Acceptance:     " << fAcceptancePath << std::endl;
852   TString opt(option);
853   opt.ToUpper();
854   if (!opt.Contains("R")) return;
855   
856   gROOT->IncreaseDirLevel();
857   if (fELossFit)          fELossFit->Print(option);
858   else 
859     std::cout << ind << "  Energy loss fits  not initialised" << std::endl;
860   
861   if (fSecondaryMap)      fSecondaryMap->Print(option);
862   else 
863     std::cout << ind << "  Secondary particle correction not initialised" 
864               << std::endl;
865
866   if (fDoubleHit)         fDoubleHit->Print(option);
867   else 
868     std::cout << ind << "  Double hit corr. not initialised" << std::endl;
869
870   if (fVertexBias)        fVertexBias->Print(option);
871   else 
872     std::cout << ind << "  Vertex bias correction not initialised" << std::endl;
873   if (fMergingEfficiency) fMergingEfficiency->Print(option);
874   else 
875     std::cout << ind << "  Merging eff.  not initialised" << std::endl;
876
877   if (fAcceptance)        fAcceptance->Print(option);
878   else 
879     std::cout << ind << "  Acceptance corr.  not initialised" << std::endl;
880   gROOT->DecreaseDirLevel();  
881 }
882
883 //____________________________________________________________________
884 void
885 AliForwardCorrectionManagerOADB::Browse(TBrowser* b)
886 {
887   // 
888   // Browse thos
889   // 
890   if (fELossFit)          b->Add(fELossFit,          "Energy loss fits");
891   if (fSecondaryMap)      b->Add(fSecondaryMap,      "Secondary particle corr");
892   if (fDoubleHit)         b->Add(fDoubleHit,         "Double hit corr");
893   if (fVertexBias)        b->Add(fVertexBias,        "Vertex bias corr");
894   if (fMergingEfficiency) b->Add(fMergingEfficiency, "Merging eff");
895   if (fAcceptance)        b->Add(fAcceptance,        "Acceptance corr");
896 }
897
898 //____________________________________________________________________
899 Bool_t
900 AliForwardCorrectionManagerOADB::StoreObject(ULong_t        runNo,
901                                              UShort_t       sys, 
902                                              UShort_t       sNN, 
903                                              Short_t        fld, 
904                                              Bool_t         mc,
905                                              Bool_t         sat,
906                                              TObject*       obj, 
907                                              Bool_t         full,
908                                              const char*    meth) const
909 {
910   ECorrection what = GetObjectType(obj);
911   if (what == kAll) { 
912     AliErrorF("Cannot deduce the correction type from object of class %s", 
913               obj->ClassName());
914     return false;
915   }
916   return StoreObject(what, runNo, sys, sNN, fld, mc, sat, obj, full, meth);
917 }
918 //____________________________________________________________________
919 Bool_t
920 AliForwardCorrectionManagerOADB::StoreObject(const TString& what, 
921                                              ULong_t        runNo,
922                                              UShort_t       sys, 
923                                              UShort_t       sNN, 
924                                              Short_t        fld, 
925                                              Bool_t         mc,
926                                              Bool_t         sat,
927                                              TObject*       obj, 
928                                              Bool_t         full,
929                                              const char*    meth) const
930 {
931   return StoreObject(GetObjectType(what), runNo, sys, 
932                      sNN, fld, mc, sat, obj, full, meth);
933 }
934
935 //____________________________________________________________________
936 Bool_t
937 AliForwardCorrectionManagerOADB::StoreObject(ECorrection what, 
938                                              ULong_t     runNo,
939                                              UShort_t    sys, 
940                                              UShort_t    sNN, 
941                                              Short_t     fld, 
942                                              Bool_t      mc,
943                                              Bool_t      sat,
944                                              TObject*    obj, 
945                                              Bool_t      full,
946                                              const char* meth) const
947 {
948   // 
949   // Write correction output to (a temporary) file 
950   // 
951   // Parameters: 
952   //   What     What to write 
953   //   sys      Collision system (1: pp, 2: PbPb)
954   //   sNN      Center of mass energy per nucleon (GeV)
955   //   fld      Field (kG)
956   //   mc       MC-only flag 
957   //   obj      Object to write 
958   //   full     if true, write to full path, otherwise locally
959   // 
960   // Return: 
961   //   true on success. 
962   TString         tableName = GetObjectName(what);
963   Bool_t          local     = !(full && fDB);
964   TString         fileName  = (local ? DB_NAME : GetFilePath(what));
965   AliOADBForward* db        = (local ? new AliOADBForward : fDB);
966   ECorrection     otype     = GetObjectType(obj);
967   if (otype != what) { 
968     AliErrorF("Correction type 0x%02x (%s) and type of object 0x%02x (%s) "
969               "does not match", what, tableName.Data(), otype, 
970               obj->ClassName());
971     return false;
972   }
973
974   if (!db->Open(fileName, Form("%s/%s", tableName.Data(), meth), true, true)) {
975     Error("StoreObject", 
976           "Failed to open table %s/%s in %s for read+write (%s)", 
977           tableName.Data(), meth, fileName.Data(), 
978           (local ? "local" : "global"));
979     return false;
980   }
981
982   // Filter out fields we do not need
983   switch (what) {
984   case kELossFits:                      break; // Need all fields
985   case kDoubleHit:         sat = false; break; // Depreacted
986   case kAcceptance:        fld = 0;            // Fall through
987   case kSecondaryMap:                          // Fall through
988   case kVertexBias:                            // Fall through
989   case kMergingEfficiency: mc  = false; break; // No specific for MC
990   default:                              break; // Never here 
991   }
992     
993
994   if (!db->Insert(tableName, obj, runNo, sys, sNN, fld, mc, sat)) { 
995     AliErrorF("Failed to write %s to database for "
996               "table=%s run=%lu, sys=%hu, sNN=%hu, field=%hd, mc=%d, sat=%d", 
997               obj->GetName(), tableName.Data(), runNo, sys, sNN, 
998               fld, mc, sat);
999     return false;
1000   }
1001     
1002   if (local) {
1003     db->Close();
1004     delete db;
1005     
1006     AliInfoF("Correction object %s written to DB in %s - merge this with "
1007              "%s to store for good", obj->GetName(), DB_NAME, 
1008              GetFilePath(what).Data());
1009   }
1010   return true;
1011 }
1012
1013 #ifndef DOXY_INPUT
1014 //______________________________________________________________________________
1015 void AliForwardCorrectionManagerOADB::Streamer(TBuffer &R__b)
1016 {
1017   //
1018   // Stream an object of class AliForwardCorrectionManagerOADB.
1019   //
1020   if (R__b.IsReading()) {
1021      R__b.ReadClassBuffer(AliForwardCorrectionManagerOADB::Class(),this);
1022      if (fgInstance) 
1023        AliWarning(Form("Singleton instance already set (%p) when reading "
1024                        "singleton object (%p).  Read object will be new "
1025                        "singleton object", fgInstance, this));
1026      fgInstance = this;
1027   } else {
1028     R__b.WriteClassBuffer(AliForwardCorrectionManagerOADB::Class(),this);
1029   }
1030 }
1031 #endif
1032 #if 0
1033 //______________________________________________________________________________
1034 void AliForwardCorrectionManagerOADB::Streamer(TBuffer &R__b)
1035 {
1036    // Stream an object of class AliForwardCorrectionManagerOADB.
1037
1038    UInt_t R__s, R__c;
1039    if (R__b.IsReading()) {
1040       Version_t R__v = R__b.ReadVersion(&R__s, &R__c); if (R__v) { }
1041       TObject::Streamer(R__b);
1042       R__b >> fInit;
1043       R__b >> fSys;
1044       R__b >> fSNN;
1045       R__b >> fField;
1046       fELossFitsPath.Streamer(R__b);
1047       fMergingEffPath.Streamer(R__b);
1048       fSecondaryMapPath.Streamer(R__b);
1049       fDoubleHitPath.Streamer(R__b);
1050       fVertexBiasPath.Streamer(R__b);
1051       fAcceptancePath.Streamer(R__b);
1052       R__b >> fELossFit;
1053       R__b >> fSecondaryMap;
1054       R__b >> fDoubleHit;
1055       R__b >> fVertexBias;
1056       R__b >> fMergingEfficiency;
1057       R__b >> fAcceptance;
1058       R__b.CheckByteCount(R__s, R__c, AliForwardCorrectionManagerOADB::IsA());
1059    } else {
1060       R__c = R__b.WriteVersion(AliForwardCorrectionManagerOADB::IsA(), kTRUE);
1061       TObject::Streamer(R__b);
1062       R__b << fInit;
1063       R__b << fSys;
1064       R__b << fSNN;
1065       R__b << fField;
1066       fELossFitsPath.Streamer(R__b);
1067       fMergingEffPath.Streamer(R__b);
1068       fSecondaryMapPath.Streamer(R__b);
1069       fDoubleHitPath.Streamer(R__b);
1070       fVertexBiasPath.Streamer(R__b);
1071       fAcceptancePath.Streamer(R__b);
1072       R__b << fELossFit;
1073       R__b << fSecondaryMap;
1074       R__b << fDoubleHit;
1075       R__b << fVertexBias;
1076       R__b << fMergingEfficiency;
1077       R__b << fAcceptance;
1078       R__b.SetByteCount(R__c, kTRUE);
1079    }
1080 }
1081 #endif
1082
1083 //____________________________________________________________________
1084 //
1085 // EOF
1086 //