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