]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FORWARD/analysis2/AliForwardCorrectionManager.cxx
Renamed AliFMDCorrDeadChannels to AliFMDCorrAcceptance
[u/mrichter/AliRoot.git] / PWG2 / FORWARD / analysis2 / AliForwardCorrectionManager.cxx
1 //
2 // Manager (singleton) of corrections 
3 // 
4 #include "AliForwardCorrectionManager.h"
5 #include "AliForwardUtil.h"
6 #include <TString.h>
7 #include <AliLog.h>
8 #include <TFile.h>
9 #include <TSystem.h>
10
11     
12 //____________________________________________________________________
13 AliForwardCorrectionManager* AliForwardCorrectionManager::fgInstance = 0;
14 const char* AliForwardCorrectionManager::fgkSecondaryMapSkel = "secondary";
15 const char* AliForwardCorrectionManager::fgkDoubleHitSkel    = "doublehit";
16 const char* AliForwardCorrectionManager::fgkELossFitsSkel    = "elossfits";
17 const char* AliForwardCorrectionManager::fgkVertexBiasSkel   = "vertexbias";
18 const char* AliForwardCorrectionManager::fgkMergingEffSkel   = "merging";
19
20 #define PREFIX "$(ALICE_ROOT)/PWG2/FORWARD/corrections/"
21
22 //____________________________________________________________________
23 AliForwardCorrectionManager& AliForwardCorrectionManager::Instance()
24 {
25   // 
26   // Access to the singleton object 
27   // 
28   // Return:
29   //    Reference to the singleton object 
30   //
31   if (!fgInstance) fgInstance= new AliForwardCorrectionManager;
32   return *fgInstance;
33 }
34
35 //____________________________________________________________________
36 AliForwardCorrectionManager::AliForwardCorrectionManager()
37   : TObject(), 
38     fInit(kFALSE),
39     fSys(0),
40     fSNN(0),
41     fField(999),
42     fELossFitsPath(PREFIX "ELossFits"),
43     fMergingEffPath(PREFIX "MergingEfficiency"), 
44     fSecondaryMapPath(PREFIX "SecondaryMap"),
45     fDoubleHitPath(PREFIX "DoubleHit"),
46     fVertexBiasPath(PREFIX "VertexBias"),
47     fELossFit(0),
48     fSecondaryMap(0),
49     fDoubleHit(0),
50     fVertexBias(0),
51     fMergingEfficiency(0)
52 {
53   // 
54   // Default constructor 
55   //
56 }
57 //____________________________________________________________________
58 AliForwardCorrectionManager::AliForwardCorrectionManager(const AliForwardCorrectionManager& o)
59   : TObject(o),
60     fInit(o.fInit),
61     fSys(o.fSys),
62     fSNN(o.fSNN),
63     fField(o.fField),
64     fELossFitsPath(o.fELossFitsPath),
65     fMergingEffPath(o.fMergingEffPath), 
66     fSecondaryMapPath(o.fSecondaryMapPath),
67     fDoubleHitPath(o.fDoubleHitPath),
68     fVertexBiasPath(o.fVertexBiasPath),
69     fELossFit(o.fELossFit),
70     fSecondaryMap(o.fSecondaryMap),
71     fDoubleHit(o.fDoubleHit),
72     fVertexBias(o.fVertexBias),
73     fMergingEfficiency(o.fMergingEfficiency)
74
75 {
76   // 
77   // Copy constructor 
78   // 
79   // Parameters:
80   //    o Object to copy from 
81   //
82 }
83 //____________________________________________________________________
84 AliForwardCorrectionManager&
85 AliForwardCorrectionManager::operator=(const AliForwardCorrectionManager& o)
86 {
87   // 
88   // Assignment operator 
89   // 
90   // Parameters:
91   //    o Object to assign from 
92   // 
93   // Return:
94   //    Reference to this object 
95   //
96   fInit             = o.fInit;
97   fSys              = o.fSys;
98   fSNN              = o.fSNN;
99   fField            = o.fField;
100   fELossFitsPath    = o.fELossFitsPath;
101   fMergingEffPath   = o.fMergingEffPath;
102   fSecondaryMapPath = o.fSecondaryMapPath;
103   fDoubleHitPath    = o.fDoubleHitPath;
104   fVertexBiasPath   = o.fVertexBiasPath;
105   fELossFit         = o.fELossFit;
106   fSecondaryMap     = o.fSecondaryMap;
107   fDoubleHit        = o.fDoubleHit;
108   fVertexBias       = o.fVertexBias;
109   fMergingEfficiency= o.fMergingEfficiency;
110   return *this;
111 }
112
113 //____________________________________________________________________
114 Bool_t
115 AliForwardCorrectionManager::Init(const char* cms, 
116                                   Float_t     sNN, 
117                                   Float_t     field,
118                                   Bool_t      mc,
119                                   UInt_t      what,
120                                   Bool_t      force)
121 {
122   // 
123   // Read in correction based on passed parameters
124   // 
125   // Parameters:
126   //    collisionSystem Collision system string 
127   //    cmsNN           Center of mass energy per nucleon pair [GeV]
128   //    field           Magnetic field [kG]
129   //    mc              Monte-carlo switch
130   //    what            What to read in 
131   //    force           Force (re-)reading of specified things
132   // 
133   // Return:
134   //    true on success
135   //
136   UShort_t col = AliForwardUtil::ParseCollisionSystem(cms);
137   return Init(col, 
138               AliForwardUtil::ParseCenterOfMassEnergy(col, sNN),
139               AliForwardUtil::ParseMagneticField(field), 
140               mc, what, force);
141 }
142
143 //____________________________________________________________________
144 Bool_t
145 AliForwardCorrectionManager::Init(UShort_t cms, 
146                                   UShort_t sNN, 
147                                   Short_t  field,
148                                   Bool_t   mc,
149                                   UInt_t   what,
150                                   Bool_t   force)
151 {
152   // 
153   // Read in corrections based on the parameters given 
154   // 
155   // Parameters:
156   //    collisionSystem Collision system
157   //    cmsNN           Center of mass energy per nuclean pair [GeV]
158   //    field           Magnetic field setting [kG]
159   //    mc              Monte-carlo switch
160   //    what            What to read in. 
161   //    force           Force (re-)reading of specified things
162   // 
163   // Return:
164   //    
165   //
166   if (force) fInit = kFALSE;
167   if (fInit) return kTRUE;
168
169   Bool_t ret = kTRUE;
170   if (fSys == cms && TMath::Abs(fSNN - sNN) < 10 && fField == field) { 
171     // We're already initialised for these settings - do nothing and return
172     fInit = kTRUE;
173     return ret;
174   }
175   // Set cached parameters 
176   fSys   = cms;
177   fSNN   = sNN;
178   fField = field;
179                
180   // Read secondary map if requested 
181   if (what & kSecondaryMap) {
182     if (!ReadSecondaryMap(cms, sNN, field)) {
183       AliWarning(Form("Failed to read in secondary map for "
184                       "cms=%d, sNN=%dGeV, field=%dkG", cms, sNN, field));
185       ret = kFALSE;
186     }
187   }
188   // Read double hit if requested 
189   if (what & kDoubleHit) {
190     if (!ReadDoubleHit(cms, sNN, field)) {
191       AliWarning(Form("Failed to read in double hit correction for "
192                       "cms=%d, sNN=%dGeV, field=%dkG", cms, sNN, field));
193       ret = kFALSE;
194     }
195   }
196   // Read energy loss fits if requested 
197   if (what & kELossFits) {
198     if (!ReadELossFits(cms, sNN, field, mc)) {
199       AliWarning(Form("Failed to read in energy loss fits for "
200                       "cms=%d, sNN=%dGeV, field=%dkG, %s", 
201                       cms, sNN, field, mc ? "MC" : "real"));
202       ret = kFALSE;
203     }
204   }
205   // Read event selection efficiencies if requested 
206   if (what & kVertexBias) {
207     if (!ReadVertexBias(cms, sNN, field)) {
208       AliWarning(Form("Failed to read in event selection efficiency for "
209                       "cms=%d, sNN=%dGeV, field=%dkG", cms, sNN, field));
210       ret = kFALSE;
211     }
212   }
213   // Read merging efficiencies if requested 
214   if (what & kMergingEfficiency) {
215     if (!ReadMergingEfficiency(cms, sNN, field)) {
216       AliWarning(Form("Failed to read in hit merging efficiency for "
217                       "cms=%d, sNN=%dGeV, field=%dkG", 
218                       cms, sNN, field));
219       ret = kFALSE;
220     }
221   }
222   fInit = kTRUE;
223   return ret;
224 }
225
226 //____________________________________________________________________
227 TString 
228 AliForwardCorrectionManager::GetFileName(ECorrection what, 
229                                          UShort_t    sys, 
230                                          UShort_t    sNN,
231                                          Short_t     field,
232                                          Bool_t      mc) const
233 {
234   // 
235   // Get the path to the specified object 
236   // 
237   // Parameters:
238   //    what  Which stuff to get the path for 
239   //    sys   Collision system
240   //    sNN   Center of mass energy [GeV]
241   //    field Magnetic field in the L3 magnet [kG]
242   //    mc    Whether the correction objects should be valid for MC
243   // 
244   // Return:
245   //    The full path or null 
246   //
247   TString fname = "";
248   fname = GetObjectName(what);
249   fname.Append(Form("_%s_%04dGeV_%c%1dkG_%s.root", 
250                     AliForwardUtil::CollisionSystemString(sys), 
251                     sNN, (field < 0 ? 'm' : 'p'), TMath::Abs(field), 
252                     (mc ? "MC" : "real")));
253   return fname;
254 }
255 //____________________________________________________________________
256 TString
257 AliForwardCorrectionManager::GetFileName(ECorrection what) const
258 {
259   // 
260   // Get the file name of the specified object
261   // 
262   // Parameters:
263   //    what Which stuff to get the path for 
264   // 
265   // Return:
266   //    The full path or null
267   //
268   if (!fInit) { 
269     AliWarning("Corrections manager initialised, do a forced Init(...)");
270     return "";
271   }
272   return GetFileName(what, fSys, fSNN, fField, false);
273 }
274
275 //____________________________________________________________________
276 const Char_t*
277 AliForwardCorrectionManager::GetFileDir(ECorrection what) const
278 {
279   // 
280   // Get the path to the specified object 
281   // 
282   // Parameters:
283   //    what  Which stuff to get the path for 
284   // 
285   // Return:
286   //    The full path or null 
287   //
288   if      (what & kSecondaryMap)        return fSecondaryMapPath;
289   else if (what & kDoubleHit)           return fDoubleHitPath;
290   else if (what & kELossFits)           return fELossFitsPath;
291   else if (what & kVertexBias)          return fVertexBiasPath;
292   else if (what & kMergingEfficiency)   return fMergingEffPath;
293
294   AliWarning(Form("Unknown correction: %d", what));
295   return 0;
296 }
297
298 //____________________________________________________________________
299 TString 
300 AliForwardCorrectionManager::GetFilePath(ECorrection what, 
301                                          UShort_t    sys, 
302                                          UShort_t    sNN,
303                                          Short_t     field,
304                                          Bool_t      mc) const
305 {
306   // 
307   // Get the path to the specified object 
308   // 
309   // Parameters:
310   //    what  Which stuff to get the path for 
311   //    sys   Collision system
312   //    sNN   Center of mass energy [GeV]
313   //    field Magnetic field in the L3 magnet [kG]
314   //    mc    Whether the correction objects should be valid for MC
315   // 
316   // Return:
317   //    The full path or null 
318   //
319   TString path = "";
320   const Char_t* dir = GetFileDir(what);
321   if (!dir) return path;
322   
323   TString fname(GetFileName(what, sys, sNN, field, mc));
324   if (fname.IsNull()) return path;
325
326   path = gSystem->ConcatFileName(gSystem->ExpandPathName(dir), fname);
327   
328   return path;
329 }
330 //____________________________________________________________________
331 TString
332 AliForwardCorrectionManager::GetFilePath(ECorrection what) const
333 {
334   // 
335   // Get the full path to the object.  Note, the manager must be
336   // initialised for this to work
337   // 
338   // Parameters:
339   //    what Which stuff to get the path for 
340   // 
341   // Return:
342   //    The full path or null
343   //
344   if (!fInit) { 
345     AliWarning("Corrections manager initialised, do a forced Init(...)");
346     return "";
347   }
348   return GetFilePath(what, fSys, fSNN, fField, false);
349 }
350
351 //____________________________________________________________________
352 TFile*
353 AliForwardCorrectionManager::GetFile(ECorrection what, 
354                                      UShort_t    sys, 
355                                      UShort_t    sNN, 
356                                      Short_t     field, 
357                                      Bool_t      mc, 
358                                      Bool_t      rw, 
359                                      Bool_t      newfile) const
360 {
361   // 
362   // Open the file that contains the correction object specified 
363   // 
364   // Parameters:
365   //    what  Which stuff to get the path for 
366   //    sys   Collision system
367   //    sNN   Center of mass energy [GeV]
368   //    field Magnetic field in the L3 magnet [kG]
369   //    mc    Whether the correction objects should be valid for MC
370   //    rw    Whether to open the file in read/write
371   //    newfile Wheter to make the file if it doesn't exist
372   // 
373   // Return:
374   //    The file that contains the correction object or null 
375   //
376   TString path = GetFilePath(what, sys, sNN, field, mc);
377   if (path.IsNull()) return 0;
378   
379   TString opt;
380   if (newfile) opt="RECREATE";
381   else {
382     if (gSystem->AccessPathName(path.Data(), 
383                                 (rw ? kWritePermission : kReadPermission))) {
384       AliWarning(Form("file %s cannot be found or insufficient permissions", 
385                       path.Data()));
386       return 0;
387     }
388     opt=(rw ? "UPDATE" : "READ");
389   }
390   TFile* file = TFile::Open(path.Data(), opt.Data());
391   if (!file) { 
392     AliWarning(Form("file %s cannot be opened in mode %s", 
393                     path.Data(), opt.Data()));
394     return 0;
395   }
396   return file;
397 }
398 //____________________________________________________________________
399 TFile*
400 AliForwardCorrectionManager::GetFile(ECorrection what) const
401 {
402   // 
403   // Get the file that contains the object specifed.  Note, the manager
404   // must be initialised for this to work. 
405   // 
406   // Parameters:
407   //    what Which stuff to get the path for 
408   // 
409   // Return:
410   //    The file that contains the correction object or null
411   //
412   if (!fInit) { 
413     AliWarning("Corrections manager initialised, do a forced Init(...)");
414     return 0;
415   }
416   return GetFile(what, fSys, fSNN, fField, false);
417 }
418
419 //____________________________________________________________________
420 const Char_t*
421 AliForwardCorrectionManager::GetObjectName(ECorrection what) const
422 {
423   // 
424   // Get the object name corresponding to correction type 
425   // 
426   // Parameters:
427   //    what Correction 
428   // 
429   // Return:
430   //    Object name or null
431   //
432   if      (what & kSecondaryMap)       return fgkSecondaryMapSkel;
433   else if (what & kDoubleHit)          return fgkDoubleHitSkel;
434   else if (what & kELossFits)          return fgkELossFitsSkel;
435   else if (what & kVertexBias)         return fgkVertexBiasSkel;
436   else if (what & kMergingEfficiency)  return fgkMergingEffSkel;
437   return 0;
438 }
439
440 //____________________________________________________________________
441 TObject*
442 AliForwardCorrectionManager::CheckObject(TFile* file, ECorrection what) const
443 {
444   // 
445   // Check if the specified objet exists in the file, and return it
446   // 
447   // Parameters:
448   //    file File to query 
449   //    what Correction type 
450   // 
451   // Return:
452   //    Object found, or null
453   //
454   TObject* o = file->Get(GetObjectName(what));
455   if (!o) { 
456     AliWarning(Form("Object %s not found in %s", 
457                     GetObjectName(what), file->GetName()));
458     file->Close();
459     return 0;
460   }
461   return o;
462 }
463   
464 //____________________________________________________________________
465 TObject*
466 AliForwardCorrectionManager::GetObject(ECorrection what, 
467                                        UShort_t    sys, 
468                                        UShort_t    sNN, 
469                                        Short_t     field,
470                                        Bool_t      mc) const
471 {
472   // 
473   // Get the path to the specified object 
474   // 
475   // Parameters:
476   //    what  Which stuff to get the path for 
477   //    sys   Collision system
478   //    sNN   Center of mass energy [GeV]
479   //    field Magnetic field in the L3 magnet [kG]
480   //    mc    Whether the correction objects should be valid for MC
481   // 
482   // Return:
483   //    The full path or null 
484   //
485   TFile* file = GetFile(what, sys, sNN, field, mc, false, false);
486   if (!file) return 0;
487   
488   return CheckObject(file, what);
489 }
490 //____________________________________________________________________
491 TObject*
492 AliForwardCorrectionManager::GetObject(ECorrection what) const
493 {
494   // 
495   // Get the object that contaisn the specified correction
496   // 
497   // Parameters:
498   //    what Which object to get
499   // 
500   // Return:
501   //    The object or null
502   //
503   if (!fInit) { 
504     AliWarning("Corrections manager initialised, do a forced Init(...)");
505     return 0;
506   }
507   return GetObject(what, fSys, fSNN, fField, false);
508 }
509
510
511 //____________________________________________________________________
512 Bool_t 
513 AliForwardCorrectionManager::ReadSecondaryMap(UShort_t sys, UShort_t sNN, 
514                                               Short_t field)
515 {
516   // 
517   // Read in the secondary map 
518   // 
519   // Parameters:
520   //    sys   Collision system
521   //    sNN   Center of mass energy [GeV]
522   //    field Magnetic field in the L3 magnet [kG]
523   // 
524   // Return:
525   //    True on success, false otherwise 
526   //
527   if (fInit) { 
528     AliWarning("Corrections manager initialised, do a forced Init(...)");
529     return kFALSE;
530   }
531
532   TObject* o = GetObject(kSecondaryMap, sys, sNN, field, false);
533   if (!o) return kFALSE;
534
535   fSecondaryMap = dynamic_cast<AliFMDCorrSecondaryMap*>(o);
536   if (!fSecondaryMap) {
537     AliWarning(Form("Object %s (%p) is not an AliFMDCorrSecondaryMap object, "
538                     "but %s", fgkSecondaryMapSkel, o, o->ClassName())); 
539     return kFALSE;
540   }
541
542   // file->Close();
543   return kTRUE;
544 }
545 //____________________________________________________________________
546 Bool_t 
547 AliForwardCorrectionManager::ReadDoubleHit(UShort_t sys, UShort_t sNN, 
548                                            Short_t field)
549 {
550   // 
551   // Read in the double hit correction
552   // 
553   // Parameters:
554   //    sys   Collision system
555   //    sNN   Center of mass energy [GeV]
556   //    field Magnetic field in the L3 magnet [kG]
557   // 
558   // Return:
559   //    True on success, false otherwise 
560   //
561   if (fInit) { 
562     AliWarning("Corrections manager initialised, do a forced Init(...)");
563     return kFALSE;
564   }
565
566   TObject* o = GetObject(kDoubleHit, sys, sNN, field, false);
567   if (!o) return kFALSE;
568
569   fDoubleHit = dynamic_cast<AliFMDCorrDoubleHit*>(o);
570   if (!fDoubleHit) {
571     AliWarning(Form("Object %s (%p) is not an AliFMDCorrDoubleHit object, "
572                     "but %s", fgkDoubleHitSkel, o, o->ClassName())); 
573     return kFALSE;
574   }
575
576   // file->Close();
577   return kTRUE;
578 }
579
580 //____________________________________________________________________
581 Bool_t 
582 AliForwardCorrectionManager::ReadELossFits(UShort_t sys, UShort_t sNN, 
583                                            Short_t field, Bool_t mc)
584 {
585   // 
586   // Read in the energy loss fits 
587   // 
588   // Parameters:
589   //    sys   Collision system
590   //    sNN   Center of mass energy [GeV]
591   //    field Magnetic field in the L3 magnet [kG]
592   //    mc    Whether the correction objects should be valid for MC
593   // 
594   // Return:
595   //    True on success, false otherwise 
596   //
597   if (fInit) { 
598     AliWarning("Corrections manager initialised, do a forced Init(...)");
599     return kFALSE;
600   }
601
602   TObject* o = GetObject(kELossFits, sys, sNN, field, mc);
603   if (!o) return kFALSE;
604
605   fELossFit = dynamic_cast<AliFMDCorrELossFit*>(o);
606   if (!fELossFit) {
607     AliWarning(Form("Object %s (%p) is not an AliFMDCorrELossFit object, "
608                     "but %s", fgkELossFitsSkel, o, o->ClassName()));
609     return kFALSE;
610   }
611
612   // file->Close();
613   return kTRUE;
614 }
615
616 //____________________________________________________________________
617 Bool_t 
618 AliForwardCorrectionManager::ReadVertexBias(UShort_t sys, 
619                                             UShort_t sNN, 
620                                             Short_t field)
621 {
622   // 
623   // Read in the event selection efficiency 
624   // 
625   // Parameters:
626   //    sys   Collision system
627   //    sNN   Center of mass energy [GeV]
628   //    field Magnetic field in the L3 magnet [kG]
629   // 
630   // Return:
631   //    True on success, false otherwise 
632   //
633   if (fInit) { 
634     AliWarning("Corrections manager initialised, do a forced Init(...)");
635     return kFALSE;
636   }
637
638   TObject* o = GetObject(kVertexBias, sys, sNN, field, false);
639   if (!o) return kFALSE;
640
641   fVertexBias = dynamic_cast<AliFMDCorrVertexBias*>(o);
642   if (!fVertexBias) {
643     AliWarning(Form("Object %s (%p) is not an AliFMDCorrVertexBias object, "
644                     "but %s", fgkVertexBiasSkel, o, o->ClassName()));
645     return kFALSE;
646   }
647
648   // file->Close();
649   return kTRUE;
650 }
651
652 //____________________________________________________________________
653 Bool_t 
654 AliForwardCorrectionManager::ReadMergingEfficiency(UShort_t sys, 
655                                                    UShort_t sNN, 
656                                                    Short_t field)
657 {
658   // 
659   // Read in the merging efficiency 
660   // 
661   // Parameters:
662   //    sys   Collision system
663   //    sNN   Center of mass energy [GeV]
664   //    field Magnetic field in the L3 magnet [kG]
665   // 
666   // Return:
667   //    True on success, false otherwise 
668   //
669   if (fInit) { 
670     AliWarning("Corrections manager initialised, do a forced Init(...)");
671     return kFALSE;
672   }
673
674   TObject* o = GetObject(kMergingEfficiency, sys, sNN, field, false);
675   if (!o) return kFALSE;
676
677   fMergingEfficiency = dynamic_cast<AliFMDCorrMergingEfficiency*>(o);
678   if (!fMergingEfficiency) {
679     AliWarning(Form("Object %s (%p) is not an AliFMDCorrMergingEfficiency "
680                     "object, but %s", fgkMergingEffSkel, o, o->ClassName()));
681     return kFALSE;
682   }
683
684   // file->Close();
685   return kTRUE;
686 }
687
688 //____________________________________________________________________
689 //
690 // EOF
691 //