]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDParameters.cxx
Fixes for ROOT trunk - #error is seen by Cling-based rootcint
[u/mrichter/AliRoot.git] / FMD / AliFMDParameters.cxx
1 /*************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  *************************************************************************
15  * $Id$ */
16 /**
17  * @file    AliFMDParameters.cxx
18  * @author  Christian Holm Christensen <cholm@nbi.dk>
19  * @date    Mon Mar 27 12:44:26 2006
20  * @brief   Manager of FMD parameters     
21  */
22 //____________________________________________________________________
23 //                                                                          
24 // Forward Multiplicity Detector based on Silicon wafers. 
25 //
26 // This class is a singleton that handles various parameters of
27 // the FMD detectors.  
28 // The manager normally serves the parameters from the Conditions
29 // Database (CDB).  These are retrivied by the member function
30 // `Init'.  Optionally, the class can serve hard-coded constants, if
31 // no CDB is available. 
32 //                                                       
33 #include "AliFMDDebug.h"                   // ALILOG_H
34 #include "AliFMDParameters.h"      // ALIFMDPARAMETERS_H
35 #include "AliFMDGeometry.h"        // ALIFMDGEOMETRY_H
36 #include "AliFMDRing.h"            // ALIFMDRING_H
37 #include "AliFMDCalibGain.h"       // ALIFMDCALIBGAIN_H
38 #include "AliFMDCalibPedestal.h"   // ALIFMDCALIBPEDESTAL_H
39 #include "AliFMDCalibSampleRate.h" // ALIFMDCALIBSAMPLERATE_H
40 #include "AliFMDCalibStripRange.h" // ALIFMDCALIBSTRIPRANGE_H
41 #include "AliFMDAltroMapping.h"    // ALIFMDALTROMAPPING_H
42 #include <AliCDBManager.h>         // ALICDBMANAGER_H
43 #include <AliCDBEntry.h>           // ALICDBMANAGER_H
44 #include <AliFMDPreprocessor.h>
45 #include <AliLog.h>
46 #include <Riostream.h>
47 #include <sstream>
48 #include <TSystem.h>
49 #include <TArrayF.h>
50 #include <TH2D.h>
51
52 //====================================================================
53 ClassImp(AliFMDParameters)
54 #if 0
55   ; // This is here to keep Emacs for indenting the next line
56 #endif
57
58 //____________________________________________________________________
59 AliFMDParameters* AliFMDParameters::fgInstance = 0;
60
61 //____________________________________________________________________
62 const char* AliFMDParameters::fgkPulseGain          = "FMD/Calib/PulseGain";
63 const char* AliFMDParameters::fgkPedestal           = "FMD/Calib/Pedestal";
64 const char* AliFMDParameters::fgkDead               = "FMD/Calib/Dead";
65 const char* AliFMDParameters::fgkSampleRate         = "FMD/Calib/SampleRate";
66 const char* AliFMDParameters::fgkAltroMap           = "FMD/Calib/AltroMap";
67 const char* AliFMDParameters::fgkZeroSuppression    = "FMD/Calib/ZeroSuppression";
68 const char* AliFMDParameters::fgkStripRange         = "FMD/Calib/StripRange";
69 const char* AliFMDParameters::fkPedestalShuttleID   = "pedestals";
70 const char* AliFMDParameters::fkGainShuttleID       = "gains";
71 const char* AliFMDParameters::fkConditionsShuttleID = "conditions";
72
73 //____________________________________________________________________
74 AliFMDParameters* 
75 AliFMDParameters::Instance() 
76 {
77   // 
78   // Get static instance 
79   //
80   if (!fgInstance) fgInstance = new AliFMDParameters;
81   return fgInstance;
82 }
83
84 //____________________________________________________________________
85 AliFMDParameters::AliFMDParameters() 
86   : fIsInit(kFALSE),
87     fkSiDeDxMip(1.664), 
88     fVA1MipRange(0),
89     fAltroChannelSize(0),
90     fChannelsPerAltro(0),
91     fPedestalFactor(0),
92     fZSPre(1),
93     fZSPost(1),
94     fZSPedSubtract(kTRUE),
95     fFixedPedestal(100),
96     fFixedPedestalWidth(2),
97     fFixedZeroSuppression(1),
98     fFixedSampleRate(2),
99     fFixedThreshold(0),
100     fFixedMinStrip(0),
101     fFixedMaxStrip(127),
102     fFixedPulseGain(2), 
103     fEdepMip(0),
104     fHasCompleteHeader(kTRUE),
105     fZeroSuppression(0), 
106     fSampleRate(0), 
107     fPedestal(0), 
108     fPulseGain(0), 
109     fDeadMap(0), 
110     fAltroMap(0), 
111   fStripRange(0),
112   fRunNo(-1)
113 {
114   //
115   // Default constructor 
116   //
117   SetVA1MipRange();
118   SetAltroChannelSize();
119   SetChannelsPerAltro();
120   SetZeroSuppression();
121   SetSampleRate();
122   SetPedestal();
123   SetPedestalWidth();
124   SetPedestalFactor();
125   SetThreshold();
126   SetStripRange();
127   SetGain();
128   fAltroMap = new AliFMDAltroMapping;
129 }
130
131 //__________________________________________________________________
132 Bool_t
133 AliFMDParameters::CheckForNewRun()
134 {
135   Int_t run = AliCDBManager::Instance()->GetRun();
136   if (run != fRunNo) { 
137     fIsInit = false;
138     fRunNo = run;
139   }
140   return run != fRunNo;
141 }
142
143 //__________________________________________________________________
144 UShort_t
145 AliFMDParameters::Init(Bool_t forceReInit, UInt_t what)
146 {
147   // 
148   // Initialize the manager.  This tries to read the parameters from
149   // CDB.  If that fails, the class uses the hard-coded parameters.
150   // 
151   // Parameters:
152   //    forceReInit Force (re-)initalize flag
153   //    what        What to initialize 
154   //
155   if (forceReInit) fIsInit = kFALSE;
156   CheckForNewRun();
157
158   if (fIsInit) return 0;
159
160   UShort_t errMask = 0;
161   if (what & kPulseGain)       errMask |= InitPulseGain();
162   if (what & kPedestal)        errMask |= InitPedestal();
163   if (what & kDeadMap)         errMask |= InitDeadMap();
164   if (what & kSampleRate)      errMask |= InitSampleRate();
165   if (what & kZeroSuppression) errMask |= InitZeroSuppression();
166   if (what & kAltroMap)        errMask |= InitAltroMap();
167   if (what & kStripRange)      errMask |= InitStripRange();
168   fIsInit = kTRUE;
169
170   return errMask;
171 }
172 //__________________________________________________________________
173 UShort_t
174 AliFMDParameters::Init(AliFMDPreprocessor* pp, Bool_t forceReInit, UInt_t what)
175 {
176   // 
177   // Initialize the manager.  This tries to read the parameters from
178   // CDB.  If that fails, the class uses the hard-coded parameters.
179   // 
180   // Parameters:
181   //    pp          Preprocessor 
182   //    forceReInit Force (re-)initalize flag
183   //    what        What to initialize 
184   //
185   if (forceReInit) fIsInit = kFALSE;
186   CheckForNewRun();
187
188   if (fIsInit) return 0;
189
190   UShort_t errMask = 0;
191   if (what & kPulseGain)       errMask |= InitPulseGain(pp);
192   if (what & kPedestal)        errMask |= InitPedestal(pp);
193   if (what & kDeadMap)         errMask |= InitDeadMap(pp);
194   if (what & kSampleRate)      errMask |= InitSampleRate(pp);
195   if (what & kZeroSuppression) errMask |= InitZeroSuppression(pp);
196   if (what & kAltroMap)        errMask |= InitAltroMap(pp);
197   if (what & kStripRange)      errMask |= InitStripRange(pp);
198   fIsInit = kTRUE;
199
200   return errMask;
201 }
202
203 //__________________________________________________________________
204 Bool_t
205 AliFMDParameters::CheckFile(const char* prefix, 
206                             const char* path, 
207                             int         number, 
208                             TString&    f) const
209 {
210   // 
211   // Check if the file <i>prefix</i><i>number</i> exists in @a path, 
212   // and write the full path to @a f.  
213   // 
214   // Parameters:
215   //    prefix  File prefix (cond, peds, gains, ...)
216   //    path    Path to files
217   //    number  Detector number (1, 2, or 3)
218   //    f       On return full path to file (if found)
219   // 
220   // Return:
221   //    @c true if file exists and is readable, @c false otherwise
222   //
223   f = (Form("%s%d.csv", prefix, number));
224   AliFMDDebug(5, ("Checking if %s exists in %s ...", f.Data(), path));
225   f = gSystem->Which(path, f.Data());
226   AliFMDDebug(5, ("Got back '%s'", f.Data()));
227   return !f.IsNull();
228 }
229
230 //__________________________________________________________________
231 UShort_t
232 AliFMDParameters::Init(const char* path, Bool_t forceReInit, UInt_t what)
233 {
234   // 
235   // Initialize the manager.  This will try to read some calibrations
236   // (sample rate, strip range, gains, pedestals) from local comma
237   // separated value (CSV) files in the directory pointed at by @a
238   // path.  If they are not found, then they will be retrieved from
239   // OCDB as appropriately.   Other calibrations are always read from
240   // OCDB.  
241   // 
242   // The CSV files should be named as 
243   // 
244   // - Pedestals: <tt>peds</tt><i>det_number</i><tt>.csv</tt>
245   // - Gains: <tt>gains</tt><i>det_number</i><tt>.csv</tt>
246   // - Sample Rate: <tt>conditions</tt><i>det_number</i><tt>.csv</tt>
247   // - Strip Range: <tt>conditions</tt><i>det_number</i><tt>.csv</tt>
248   //
249   // where <i>det_number</i> is the detector number (1, 2, or 3). 
250   // 
251   // Parameters:
252   //    path        Where to look for the CSV files
253   //    forceReInit Always reinitialise 
254   //    what        What calibrations to load. 
255   //  
256   if (forceReInit) fIsInit = kFALSE;
257   CheckForNewRun();
258
259   if (fIsInit) return 0;
260
261   AliFMDCalibStripRange*  range = 0;
262   AliFMDCalibSampleRate*  rate  = 0;
263   AliFMDCalibPedestal*    peds  = 0;
264   AliFMDCalibGain*        gains = 0;
265
266   for (Int_t i = 1; i <= 3; i++) { 
267     TString f;
268     if (((what & kSampleRate) || (what & kStripRange)) && 
269         CheckFile("conditions", path, i, f)) {
270       if (!rate  && (what & kSampleRate)) rate  = new AliFMDCalibSampleRate;
271       if (!range && (what & kStripRange)) range = new AliFMDCalibStripRange;
272       std::ifstream in(f.Data());
273       if (range) range->ReadFromFile(in);
274       if (rate)  rate->ReadFromFile(in);
275       in.close();
276     }
277     if ((what & kPedestal) && CheckFile("peds", path, i, f)) {
278       if (!peds) peds  = new AliFMDCalibPedestal;
279       std::ifstream in(f.Data());
280       peds->ReadFromFile(in);
281       in.close();
282     }
283     if ((what & kPulseGain) && CheckFile("gains", path, i, f)) { 
284       if (!gains) gains = new AliFMDCalibGain;
285       std::ifstream in(f.Data());
286       gains->ReadFromFile(in);
287       in.close();
288     }
289   }
290
291   if (range) what &= ~kStripRange;
292   if (rate)  what &= ~kSampleRate;
293   if (peds)  what &= ~kPedestal;
294   if (gains) what &= ~kPulseGain;
295
296   UShort_t ret = Init(kFALSE, what);
297   
298   if (range) SetStripRange(range);
299   if (rate)  SetSampleRate(rate);
300   if (peds)  SetPedestal(peds);
301   if (gains) SetGain(gains);
302
303   fIsInit = kTRUE;
304
305   return ret;
306 }
307
308 //__________________________________________________________________
309 void
310 AliFMDParameters::MakeDeadMap(Float_t maxNoise, 
311                               Float_t minGain, 
312                               Float_t maxGain)
313 {
314   // 
315   // Automatically generate a dead map from the pedestals and gains.
316   // A channel is marked as dead of the noise is too high (currently
317   // more than 10 ADC counts), or the gain is unreasonable (currently
318   // larger than 10, or smaller than 0.1). 
319   // 
320   // The procedure does not overwrite channels previously marked as
321   // dead - e.g., channels marked as dead in the calibration loaded
322   // from OCDB will continue to be marked as dead.  That is, this
323   // procedure will never make a channel un-dead. 
324   // 
325   // Parameters:
326   //    maxNoise  Maximum noise value before a channel is marked
327   // as dead. 
328   //    minGain   Minimum value of the calibrated gain before a
329   // channel is considered dead. 
330   //    maxGain   Maximum value of the calibrated gain before a
331   // channel is considered dead. 
332   //
333   if (fPedestal)  
334     fDeadMap = fPedestal->MakeDeadMap(maxNoise, fDeadMap);
335   if (fPulseGain) 
336     fDeadMap = fPulseGain->MakeDeadMap(minGain, maxGain, fDeadMap);
337 }
338 //__________________________________________________________________
339 #define DET2IDX(det,ring,sec,str) \
340   (det * 1000 + (ring == 'I' ? 0 : 512) + str)  
341   
342 //__________________________________________________________________
343 void
344 AliFMDParameters::Draw(Option_t* option)
345 {
346   // 
347   // Draw parameters. 
348   // 
349   // Parameters:
350   //    option What to draw. Should be one of 
351   // - dead       Dead channels
352   // - threshold Threshold
353   // - gain       Gain
354   // - pedestal  Pedestal
355   // - noise      Noise (or pedestal width)
356   // - zero       Zero suppression
357   // - rate       Sampling rate (VA1 clock / ALTRO clock)
358   // - min        Minimum strip read out
359   // - max        Maximum strip read out
360   // - map        hardware address
361   //
362   TString opt(option);
363   enum {
364     kLocalPulseGain,       // Path to PulseGain calib object
365     kLocalThreshold,       // Path to PulseGain calib object
366     kLocalPedestal,        // Path to Pedestal calib object
367     kLocalPedestalWidth,   // Path to Pedestal calib object
368     kLocalDead,            // Path to Dead calib object
369     kLocalSampleRate,      // Path to SampleRate calib object
370     kLocalAltroMap,        // Path to AltroMap calib object
371     kLocalZeroSuppression, // Path to ZeroSuppression cal object
372     kLocalMinStripRange,   // Path to strip range cal object
373     kLocalMaxStripRange    // Path to strip range cal object
374   } what;
375     
376   if      (opt.Contains("dead", TString::kIgnoreCase)) 
377     what = kLocalDead;
378   else if (opt.Contains("threshold",TString::kIgnoreCase)) 
379     what = kLocalThreshold;
380   else if (opt.Contains("gain",TString::kIgnoreCase)) 
381     what = kLocalPulseGain;
382   else if (opt.Contains("pedestal",TString::kIgnoreCase)) 
383     what = kLocalPedestal;
384   else if (opt.Contains("noise",TString::kIgnoreCase)) 
385     what = kLocalPedestalWidth;
386   else if (opt.Contains("zero",TString::kIgnoreCase)) 
387     what = kLocalZeroSuppression;
388   else if (opt.Contains("rate",TString::kIgnoreCase)) 
389     what = kLocalSampleRate;
390   else if (opt.Contains("min",TString::kIgnoreCase)) 
391     what = kLocalMinStripRange;
392   else if (opt.Contains("max",TString::kIgnoreCase)) 
393     what = kLocalMaxStripRange;
394   else if (opt.Contains("map",TString::kIgnoreCase)) 
395     what = kLocalAltroMap;
396   else {
397     Warning("Draw", "unknown parameter: %s\n\tShould be one of\n\t"
398             "dead, threshold, gain, pedestal, noise, zero, rate, "
399             "min, max, map",  
400             option); 
401     return;
402   }
403
404   TArrayD xbins(3 * 512 + 2 * 256 + 5);
405   Int_t i = 1;
406   Bool_t skip = kTRUE;
407   for (UShort_t det = 1; det <= 3; det++) {
408     UShort_t nRings = (det == 1 ? 1 : 2);
409     for (UShort_t iring = 0; iring < nRings; iring++) {
410       UShort_t nStrip  = (iring == 0 ? 512 : 256);
411       Char_t   ring    = (iring == 0 ? 'I' : 'O');
412       for (UShort_t str = 0; str < nStrip; str++) {
413         // UShort_t nSec    = (iring == 0 ? 20  : 40);
414         // Char_t   ring    = (iring == 0 ? 'I' : 'O');
415         // for (UShort_t sec = 0; sec < nSec; sec++) {
416         Int_t idx = DET2IDX(det, ring, 0, str);
417         // Int_t idx = DET2IDX(det, ring, sec, 0);
418         if (skip) {
419           xbins[i-1] = idx - .5;
420           skip  = kFALSE;
421         }
422         xbins[i] = idx + .5;
423         i++;
424       }
425       skip = kTRUE;
426       i++;
427     }
428   }
429   TArrayD ybins(41);
430   for (/*Int_t*/ i = 0; i < ybins.fN; i++) ybins[i] = Float_t(i - .5);
431   TH2D* hist = new TH2D("calib", Form("Calibration %s", option), 
432                         xbins.fN-1, xbins.fArray,  
433                         ybins.fN-1, ybins.fArray);
434   hist->GetXaxis()->SetTitle("1000 #times detector + 512 #times ring + strip");
435   hist->GetYaxis()->SetTitle("sector");
436   
437   // hist->Draw("Lego");
438   // return;
439   
440   for (UShort_t det = 1; det <= 3; det++) {
441     UShort_t nRings = (det == 1 ? 1 : 2);
442     for (UShort_t iring = 0; iring < nRings; iring++) {
443       UShort_t nSector = (iring == 0 ?  20 : 40);
444       UShort_t nStrip  = (iring == 0 ? 512 : 256);
445       Char_t   ring    = (iring == 0 ? 'I' : 'O');
446       for (UShort_t sec = 0; sec < nSector; sec++) {
447         for (UShort_t str = 0; str < nStrip; str++) {
448           Int_t idx = DET2IDX(det, ring, sec, str);
449           UShort_t ddl, addr, time, sam=0;
450           Double_t val = 0;
451           switch (what) {
452           case kLocalPulseGain:       // Path to PulseGain calib object
453             val = GetPulseGain(det,ring,sec,str); break;
454           case kLocalThreshold:       // Path to PulseGain calib object
455             val = GetThreshold(); break;
456           case kLocalPedestal:        // Path to Pedestal calib object
457             val = GetPedestal(det,ring,sec,str); break;
458           case kLocalPedestalWidth:   // Path to Pedestal calib object
459             val = GetPedestalWidth(det,ring,sec,str); break;
460           case kLocalDead:            // Path to Dead calib object
461             val = IsDead(det,ring,sec,str); break;
462           case kLocalSampleRate:      // Path to SampleRate calib object
463             val = GetSampleRate(det,ring,sec,str); break;
464           case kLocalAltroMap:        // Path to AltroMap calib object
465             Detector2Hardware(det,ring,sec,str,sam,ddl,addr,time); 
466             val = addr; break;
467           case kLocalZeroSuppression: // Path to ZeroSuppression cal object
468             val = GetZeroSuppression(det,ring,sec,str); break;
469           case kLocalMinStripRange:   // Path to strip range cal object
470             val = GetMinStrip(det,ring,sec,str); break;
471           case kLocalMaxStripRange:    // Path to strip range cal object
472             val = GetMaxStrip(det,ring,sec,str); break;
473           }
474           hist->Fill(idx,sec,val);
475           // hist->Fill(idx,str,val);
476         }
477       }
478     }
479   }
480   hist->Draw("lego");
481 }
482
483 //__________________________________________________________________
484 void
485 AliFMDParameters::Print(Option_t* option) const
486 {
487   // Print information. 
488   // If option contains an 'A' then everything is printed. 
489   // If the option contains the string "FMD" the function will search 
490   // for detector, ring, sector, and strip numbers to print, in the
491   // format 
492   // 
493   //    FMD<detector><ring>[<sector>,<string>] 
494   // 
495   // The wild card '*' means all of <detector>, <ring>, <sector>, or 
496   // <strip>. 
497   TString opt(option);
498   Bool_t showStrips  = opt.Contains("a", TString::kIgnoreCase);
499   UShort_t ds[]      = { 1, 2, 3, 0 };
500   Char_t   rs[]      = { 'I', 'O', '\0' };
501   UShort_t minStrip  = 0;
502   UShort_t maxStrip  = 512;
503   UShort_t minSector = 0;
504   UShort_t maxSector = 40;
505   
506   
507   if (opt.Contains("fmd",TString::kIgnoreCase)) {
508     Int_t   i    = opt.Index("fmd",TString::kIgnoreCase);
509     Int_t   j    = opt.Index("]",TString::kIgnoreCase);
510     if (j != kNPOS)
511       showStrips    = kTRUE;
512     else 
513       j = opt.Length();
514     enum {
515       kReadDet, 
516       kReadRing, 
517       kReadLbrack,
518       kReadSector,
519       kReadComma,
520       kReadStrip,
521       kReadRbrack, 
522       kEnd
523     } state = kReadDet;
524     std::stringstream s(opt(i+4, j-i-3).Data());
525     while (state != kEnd) {
526       Char_t tmp = s.peek();
527       if (tmp == ' ' || tmp == '\t') {
528         s.get();
529         continue;
530       }
531       switch (state) {
532       case kReadDet: { // First, try to kRead the detector 
533         if (tmp == '*') s.get();
534         else { 
535           UShort_t det;
536           s >> det;
537           if (!s.bad()) {
538             ds[0] = det;
539             ds[1] = 0;
540           }
541         }
542         state = (s.bad() ? kEnd : kReadRing);
543       } break;
544       case kReadRing: { // Then try to read the ring;
545         Char_t ring;
546         s >> ring;
547         if (ring != '*' && !s.bad()) {
548           rs[0] = ring;
549           rs[1] = '\0';
550         }
551         state = (s.bad() ? kEnd : kReadLbrack);
552       } break;
553       case kReadLbrack: { // Try to read a left bracket 
554         Char_t lbrack;
555         s >> lbrack;
556         state = (s.bad() ? kEnd : kReadSector);
557       } break;
558       case kReadSector: { // Try to read a sector 
559         if (tmp == '*') s.get();
560         else {
561           UShort_t sec;
562           s >> sec;
563           if (!s.bad()) {
564             minSector = sec;
565             maxSector = sec + 1;
566           }
567         }
568         state = (s.bad() ? kEnd : kReadComma);
569       } break;
570       case kReadComma: { // Try to read a left bracket 
571         Char_t comma;
572         s >> comma;
573         state = (s.bad() ? kEnd : kReadStrip);
574       } break;
575       case kReadStrip: { // Try to read a strip 
576         if (tmp == '*') s.get();
577         else {
578           UShort_t str;
579           s >> str;
580           if (!s.bad()) {
581             minStrip = str;
582             maxStrip = str + 1;
583           }
584         }
585         state = (s.bad() ? kEnd : kReadRbrack);
586       } break;
587       case kReadRbrack: { // Try to read a left bracket 
588         Char_t rbrack;
589         s >> rbrack;
590         state = kEnd;
591       } break;
592       case kEnd: 
593         break;
594       }
595     }
596   }
597   UShort_t* dp = ds;
598   UShort_t  det;
599   while ((det = *(dp++))) {
600
601     Char_t* rp = rs;
602     Char_t  ring;
603     while ((ring = *(rp++))) {
604       if (det == 1 && ring == 'O') continue;
605       UShort_t min  = GetMinStrip(det, ring, 0, 0);
606       UShort_t max  = GetMaxStrip(det, ring, 0, 0);
607       std::cout << "FMD" << det << ring 
608                 << "  Strip range: " 
609                 << std::setw(3) << min << "," 
610                 << std::setw(3) << max << std::endl;
611
612       UShort_t nSec = ( ring == 'I' ? 20  :  40 );
613       UShort_t nStr = ( ring == 'I' ? 512 : 256 );
614       for (UShort_t sec = minSector; sec < maxSector && sec < nSec; sec++) {
615
616         UShort_t rate = GetSampleRate(det, ring, sec, 0);
617         std::cout << "FMD" << det << ring << "[" << std::setw(2) << sec 
618                   << "] sample rate: " << rate << std::endl;
619
620         if (!showStrips) continue;
621         std::cout 
622           << "  Strip |     Pedestal      |    Gain    | ZS thr. | Address\n" 
623           << "--------+-------------------+------------+---------+---------" 
624           << std::endl;
625         for (UShort_t str = minStrip; str < nStr && str < maxStrip; str++) {
626           if (str == minStrip) std::cout << std::setw(3) << sec << ",";
627           else std::cout << "    ";
628           std::cout << std::setw(3) << str << " | ";
629           if (IsDead(det, ring, sec, str)) {
630             std::cout << "dead" << std::endl;
631             continue;
632           }
633           UShort_t ddl, addr, time, sam=0;
634           Detector2Hardware(det, ring, sec, str, sam, ddl, addr, time);
635           std::cout << std::setw(7) << GetPedestal(det, ring, sec, str) 
636                     << "+/-" << std::setw(7) 
637                     << GetPedestalWidth(det, ring, sec, str) 
638                     << " | " << std::setw(10) 
639                     << GetPulseGain(det, ring, sec, str) 
640                     << " | " << std::setw(7) 
641                     << GetZeroSuppression(det, ring, sec, str) 
642                     << " | 0x" << std::hex << std::setw(4) 
643                     << std::setfill('0') << ddl << ",0x" << std::setw(3) 
644                     << addr << std::dec << std::setfill(' ') << std::endl;
645         } // for (strip)
646       } // for (sector)
647       std::cout
648         << "=============================================================" 
649         << std::endl;
650     } // while (ring)
651   } // while (det)
652   
653 }
654
655 //__________________________________________________________________
656 AliCDBEntry*
657 AliFMDParameters::GetEntry(const char* path, AliFMDPreprocessor* pp, 
658                            Bool_t fatal) const
659 {
660   // 
661   // Get an entry from either global AliCDBManager or passed
662   // AliFMDPreprocessor. 
663   // 
664   // Parameters:
665   //    path  Path to CDB object. 
666   //    pp    AliFMDPreprocessor 
667   //    fatal If true, raise a fatal flag if we didn't get the entry.
668   // Return:
669   //    AliCDBEntry if found 
670   // 
671   AliCDBEntry* entry = 0;
672   if (!pp) {
673     AliCDBManager* cdb = AliCDBManager::Instance();
674     entry              = cdb->Get(path);
675   }
676   else {
677     const char* third  = gSystem->BaseName(path);
678     const char* second = gSystem->BaseName(gSystem->DirName(path));
679     entry              = pp->GetFromCDB(second, third);
680   }
681   if (!entry) { 
682     TString msg(Form("No %s found in CDB, perhaps you need to "
683                      "use AliFMDCalibFaker?", path));
684     if (fatal) { AliFatal(msg.Data()); }
685     else       AliLog::Message(AliLog::kWarning, msg.Data(), "FMD", 
686                                "AliFMDParameters", "GetEntry", __FILE__, 
687                                __LINE__);
688     return 0;
689   }
690   if (entry && AliLog::GetDebugLevel("FMD", "") > 0) { 
691     AliInfoF("Got entry %p for %s", entry, path);
692     entry->PrintId();
693     entry->PrintMetaData();                     
694     entry->Print();
695   }
696   return entry;
697 }
698
699     
700 //__________________________________________________________________
701 UShort_t
702 AliFMDParameters::InitPulseGain(AliFMDPreprocessor* pp)
703 {
704   // 
705   // Initialize gains.  Try to get them from CDB 
706   // 
707   // Parameters:
708   //    pp Pre-processor if called from shuttle
709   //
710   AliCDBEntry*   gain     = GetEntry(fgkPulseGain, pp);
711   if (!gain) return kPulseGain;
712   
713   AliFMDDebug(5, ("Got gain from CDB"));
714   fPulseGain = dynamic_cast<AliFMDCalibGain*>(gain->GetObject());
715   if (!fPulseGain) { 
716     AliError("Invalid pulser gain object from CDB");
717     return kPulseGain;
718   }
719   if (!fPulseGain->Values().Ptr()) {
720     AliError("Empty pulser gain object from CDB");
721     return kPulseGain;
722   }
723   return 0;
724 }
725 //__________________________________________________________________
726 UShort_t
727 AliFMDParameters::InitPedestal(AliFMDPreprocessor* pp)
728 {
729   //
730   // Initialize pedestals.  Try to get them from CDB
731   // 
732   // Parameters:
733   //    pp Pre-processor if called from shuttle
734   //
735   AliCDBEntry*   pedestal = GetEntry(fgkPedestal, pp);
736   if (!pedestal) return kPedestal;
737
738   AliFMDDebug(5, ("Got pedestal from CDB"));
739   fPedestal = dynamic_cast<AliFMDCalibPedestal*>(pedestal->GetObject());
740   if (!fPedestal) {
741     AliError("Invalid pedestal object from CDB");
742     return kPedestal;
743   }
744   if (!fPedestal->Values().Ptr()) {
745     AliError("Empty pedestal object from CDB");
746     return kPedestal;
747   }
748   return 0;
749 }
750
751 //__________________________________________________________________
752 UShort_t
753 AliFMDParameters::InitDeadMap(AliFMDPreprocessor* pp)
754 {
755   //
756   // Initialize dead map.  Try to get it from CDB
757   // 
758   // Parameters:
759   //    pp Pre-processor if called from shuttle
760   //
761   AliCDBEntry*   deadMap  = GetEntry(fgkDead, pp);
762   if (!deadMap) return kDeadMap;
763   
764   AliFMDDebug(5, ("Got dead map from CDB"));
765   fDeadMap = dynamic_cast<AliFMDCalibDeadMap*>(deadMap->GetObject());
766   if (!fDeadMap) { 
767     AliError("Invalid dead map object from CDB");
768     return kDeadMap;
769   }
770   if (!fDeadMap->Ptr()) {
771     AliError("Empty dead map object from CDB");
772     return kDeadMap;
773   }
774   return 0;
775 }
776
777 //__________________________________________________________________
778 UShort_t
779 AliFMDParameters::InitZeroSuppression(AliFMDPreprocessor* pp)
780 {
781   //
782   // Initialize zero suppression thresholds.  Try to get them from CDB
783   // 
784   // Parameters:
785   //    pp Pre-processor if called from shuttle
786   //
787   AliCDBEntry*   zeroSup  = GetEntry(fgkZeroSuppression, pp);
788   if (!zeroSup) return kZeroSuppression;
789
790   AliFMDDebug(5, ("Got zero suppression from CDB"));
791   fZeroSuppression = 
792     dynamic_cast<AliFMDCalibZeroSuppression*>(zeroSup->GetObject());
793   if (!fZeroSuppression) {
794     AliError("Invalid zero suppression object from CDB");
795     return kZeroSuppression;
796   }
797   if (!fZeroSuppression->Ptr()) {
798     AliWarningF("Empty zero suppression object from CDB, assuming %d",
799                 fFixedZeroSuppression);
800     AliCDBManager* cdbMan = AliCDBManager::Instance();
801     if(!cdbMan || !cdbMan->GetCacheFlag())
802       delete fZeroSuppression;
803     fZeroSuppression = 0;
804   }
805   return 0;
806 }
807
808 //__________________________________________________________________
809 UShort_t
810 AliFMDParameters::InitSampleRate(AliFMDPreprocessor* pp)
811 {
812   //
813   // Initialize sample rates.  Try to get them from CDB
814   // 
815   // Parameters:
816   //    pp Pre-processor if called from shuttle
817   //
818   AliCDBEntry*   sampRat  = GetEntry(fgkSampleRate, pp);
819   if (!sampRat) return kSampleRate;
820
821   AliFMDDebug(5, ("Got zero suppression from CDB"));
822   fSampleRate = dynamic_cast<AliFMDCalibSampleRate*>(sampRat->GetObject());
823   if (!fSampleRate) {
824     AliError("Invalid sample rate object from CDB");
825     return kSampleRate;
826   }
827   if (!fSampleRate->Rates().Ptr()) { 
828     AliError("empty sample rate object from CDB");
829     return kSampleRate;
830   }
831   return 0;
832 }
833
834 //__________________________________________________________________
835 UShort_t
836 AliFMDParameters::InitAltroMap(AliFMDPreprocessor* pp)
837 {
838   //
839   // Initialize hardware map.  Try to get it from CDB
840   // 
841   // Parameters:
842   //    pp Pre-processor if called from shuttle
843   //
844   if (fAltroMap) { 
845     delete fAltroMap;
846     fAltroMap = 0;
847   }
848   AliCDBEntry*   hwMap    = GetEntry(fgkAltroMap, pp, kFALSE);       
849   if (hwMap) {
850     AliFMDDebug(5, ("Got ALTRO map from CDB"));
851     fAltroMap = dynamic_cast<AliFMDAltroMapping*>(hwMap->GetObject());
852   }
853   if (!fAltroMap) {
854     AliError("Invalid ALTRO map object from CDB");
855     fAltroMap = new AliFMDAltroMapping;
856     // return kAltroMap;
857   }
858   return 0;
859 }
860
861 //__________________________________________________________________
862 UShort_t
863 AliFMDParameters::InitStripRange(AliFMDPreprocessor* pp)
864 {
865   //
866   // Initialize strip range.  Try to get it from CDB
867   // 
868   // Parameters:
869   //    pp Pre-processor if called from shuttle
870   //
871   AliCDBEntry*   range    = GetEntry(fgkStripRange, pp);
872   if (!range) return kStripRange;
873
874   AliFMDDebug(5, ("Got strip range from CDB"));
875   fStripRange = dynamic_cast<AliFMDCalibStripRange*>(range->GetObject());
876
877   if (!fStripRange) {
878     AliError("Invalid strip range object from CDB");
879     return kStripRange;
880   }
881   if (!fStripRange->Ranges().Ptr()) {
882     AliError("Empty strip range object from CDB");
883     return kStripRange;
884   }
885   return 0;
886 }
887
888
889 //__________________________________________________________________
890 Float_t
891 AliFMDParameters::GetThreshold() const
892 {
893   // 
894   // Get the threshold in the pulser gain 
895   // 
896   // 
897   // Return:
898   //    Threshold from pulser 
899   //
900   if (!fPulseGain) return fFixedThreshold;
901   return fPulseGain->Threshold();
902 }
903
904 //__________________________________________________________________
905 Float_t
906 AliFMDParameters::GetPulseGain(UShort_t detector, Char_t ring, 
907                                UShort_t sector, UShort_t strip) const
908 {
909   // 
910   // Gain of pre-amp. for strip, sector, ring, detector 
911   //
912   // For simulations this is normally set to 
913   //
914   // @f[ 
915   //  \frac{\mbox{VA1_MIP_Range}{\mbox{ALTRO_channel_size}}\mbox{MIP_Energy_Loss}
916   // @f]
917   // 
918   // 
919   // Parameters:
920   //    detector Detector # (1-3)
921   //    ring     Ring ID ('I' or 'O')
922   //    sector   Sector number (0-39)
923   //    strip    Strip number (0-511)
924   //
925   // Return:
926   //    Gain of pre-amp.  
927   //
928   if (!fPulseGain) { 
929     if (fFixedPulseGain <= 0)
930       fFixedPulseGain = fVA1MipRange * GetEdepMip() / fAltroChannelSize;
931     return fFixedPulseGain;
932   }  
933   AliFMDDebug(50, ("pulse gain for FMD%d%c[%2d,%3d]=%f",
934                     detector, ring, sector, strip,
935                     fPulseGain->Value(detector, ring, sector, strip)));
936   return fPulseGain->Value(detector, ring, sector, strip);
937 }
938
939 //__________________________________________________________________
940 Bool_t
941 AliFMDParameters::IsDead(UShort_t detector, Char_t ring, 
942                          UShort_t sector, UShort_t strip) const
943 {
944   // 
945   // Whether the strip is considered dead
946   // 
947   // Parameters:
948   //    detector Detector # (1-3)
949   //    ring     Ring ID ('I' or 'O')
950   //    sector   Sector number (0-39)
951   //    strip    Strip number (0-511)
952   //
953   // Return:
954   //    @c true if the strip is considered dead, @c false if it's
955   // OK.
956   //
957   if (!fDeadMap) return kFALSE;
958   AliFMDDebug(50, ("Dead for FMD%d%c[%2d,%3d]=%s",
959                     detector, ring, sector, strip,
960                     fDeadMap->operator()(detector, ring, sector, strip) ? 
961                     "no" : "yes"));
962   return fDeadMap->operator()(detector, ring, sector, strip);
963 }
964
965 //__________________________________________________________________
966 UShort_t
967 AliFMDParameters::GetZeroSuppression(UShort_t detector, Char_t ring, 
968                                      UShort_t sector, UShort_t strip) const
969 {
970   // 
971   // zero suppression threshold (in ADC counts)
972   // 
973   // Parameters:
974   //    detector Detector # (1-3)
975   //    ring     Ring ID ('I' or 'O')
976   //    sector   Sector number (0-39)
977   //    strip    Strip number (0-511)
978   //
979   // Return:
980   //    zero suppression threshold (in ADC counts) 
981   //
982   if (!fZeroSuppression) return fFixedZeroSuppression;
983
984   // In case of empty zero suppression objects. 
985   if (!fZeroSuppression->Ptr() || 
986       fZeroSuppression->MaxIndex() <= 0) return fFixedZeroSuppression;
987
988   // Need to map strip to ALTRO chip. 
989   AliFMDDebug(50, ("zero sup. for FMD%d%c[%2d,%3d]=%d",
990                     detector, ring, sector, strip,
991                     fZeroSuppression->operator()(detector, ring, 
992                                                  sector, strip)));
993   return fZeroSuppression->operator()(detector, ring, sector, strip/128);
994 }
995
996 //__________________________________________________________________
997 UShort_t
998 AliFMDParameters::GetSampleRate(UShort_t det, Char_t ring, UShort_t sector, 
999                                 UShort_t str) const
1000 {
1001   // 
1002   // Get the sampling rate
1003   // 
1004   // Parameters:
1005   //    detector Detector # (1-3)
1006   //    ring     Ring ID ('I' or 'O')
1007   //    sector   Sector number (0-39)
1008   //    strip    Strip number (0-511)
1009   //
1010   // Return:
1011   //    The sampling rate 
1012   //
1013   if (!fSampleRate) return fFixedSampleRate;
1014   // Need to map sector to digitizier card. 
1015   UInt_t ret = fSampleRate->Rate(det, ring, sector, str);
1016   AliFMDDebug(50, ("Sample rate for FMD%d%c[%2d,%3d]=%d", 
1017                     det, ring, sector, str, ret));
1018   return ret;
1019 }
1020
1021 //__________________________________________________________________
1022 UShort_t
1023 AliFMDParameters::GetMinStrip(UShort_t det, Char_t ring, UShort_t sector, 
1024                               UShort_t str) const
1025 {
1026   // 
1027   // Get the minimum strip in the read-out range
1028   // 
1029   // Parameters:
1030   //    detector Detector # (1-3)
1031   //    ring     Ring ID ('I' or 'O')
1032   //    sector   Sector number (0-39)
1033   //    strip    Strip number (0-511)
1034   //
1035   // Return:
1036   //    Minimum strip 
1037   //
1038   if (!fStripRange) return fFixedMinStrip;
1039   // Need to map sector to digitizier card. 
1040   UInt_t ret = fStripRange->Min(det, ring, sector, str);
1041   AliFMDDebug(50, ("Min strip # for FMD%d%c[%2d,%3d]=%d", 
1042                     det, ring, sector, str, ret));
1043   return ret;
1044 }
1045
1046 //__________________________________________________________________
1047 UShort_t
1048 AliFMDParameters::GetMaxStrip(UShort_t det, Char_t ring, UShort_t sector, 
1049                               UShort_t str) const
1050 {
1051   // 
1052   // Get the maximum strip in the read-out range
1053   // 
1054   // Parameters:
1055   //    detector Detector # (1-3)
1056   //    ring     Ring ID ('I' or 'O')
1057   //    sector   Sector number (0-39)
1058   //    strip    Strip number (0-511)
1059   //
1060   // Return:
1061   //    Maximum strip 
1062   //
1063   if (!fStripRange) return fFixedMaxStrip;
1064   // Need to map sector to digitizier card. 
1065   UInt_t ret = fStripRange->Max(det, ring, sector, str);
1066   AliFMDDebug(50, ("Max strip # for FMD%d%c[%2d,%3d]=%d", 
1067                     det, ring, sector, str, ret));
1068   return ret;
1069 }
1070
1071 //__________________________________________________________________
1072 Float_t
1073 AliFMDParameters::GetPedestal(UShort_t detector, Char_t ring, 
1074                               UShort_t sector, UShort_t strip) const
1075 {
1076   // 
1077   // Get mean of pedestal
1078   // 
1079   // Parameters:
1080   //    detector Detector # (1-3)
1081   //    ring     Ring ID ('I' or 'O')
1082   //    sector   Sector number (0-39)
1083   //    strip    Strip number (0-511)
1084   //
1085   // Return:
1086   //    Mean of pedestal 
1087   //
1088   if (!fPedestal) return fFixedPedestal;
1089   AliFMDDebug(50, ("pedestal for FMD%d%c[%2d,%3d]=%f",
1090                     detector, ring, sector, strip,
1091                     fPedestal->Value(detector, ring, sector, strip)));
1092   return fPedestal->Value(detector, ring, sector, strip);
1093 }
1094
1095 //__________________________________________________________________
1096 Float_t
1097 AliFMDParameters::GetPedestalWidth(UShort_t detector, Char_t ring, 
1098                                    UShort_t sector, UShort_t strip) const
1099 {
1100   // 
1101   // Width of pedestal
1102   // 
1103   // Parameters:
1104   //    detector Detector # (1-3)
1105   //    ring     Ring ID ('I' or 'O')
1106   //    sector   Sector number (0-39)
1107   //    strip    Strip number (0-511)
1108   //
1109   // Return:
1110   //    Width of pedestal 
1111   //
1112   if (!fPedestal) return fFixedPedestalWidth;
1113   AliFMDDebug(50, ("pedetal width for FMD%d%c[%2d,%3d]=%f",
1114                     detector, ring, sector, strip,
1115                     fPedestal->Width(detector, ring, sector, strip)));
1116   return fPedestal->Width(detector, ring, sector, strip);
1117 }
1118   
1119 //__________________________________________________________________
1120 AliFMDAltroMapping*
1121 AliFMDParameters::GetAltroMap() const
1122 {
1123   // 
1124   // Get the map that translates hardware to detector coordinates 
1125   //
1126   // Return:
1127   //    Get the map that translates hardware to detector
1128   // coordinates 
1129   // 
1130   return fAltroMap;
1131 }
1132
1133
1134 //____________________________________________________________________
1135 Bool_t 
1136 AliFMDParameters::Hardware2Detector(UShort_t  ddl,       UShort_t addr,
1137                                     UShort_t  timebin,   
1138                                     UShort_t& det,       Char_t&  ring, 
1139                                     UShort_t& sec,       Short_t& str,
1140                                     UShort_t& sam) const
1141 {
1142   // 
1143   // Map a hardware address into a detector index. 
1144   // 
1145   // Parameters:
1146   //    ddl        Hardware DDL number 
1147   //    addr       Hardware address.  
1148   //    timebin    Timebin 
1149   //    det        On return, the detector #
1150   //    ring       On return, the ring ID
1151   //    sec        On return, the sector #
1152   //    str        On return, the base of strip #
1153   //    sam        On return, the sample number for this strip
1154   //
1155   // Return:
1156   //    @c true on success, false otherwise 
1157   //
1158   if (!fAltroMap) return kFALSE;
1159   UShort_t board, chip, chan;
1160   fAltroMap->ChannelAddress(addr, board, chip, chan);
1161   return Hardware2Detector(ddl,board,chip,chan,timebin,det,ring,sec,str,sam);
1162 }
1163 //____________________________________________________________________
1164 Bool_t 
1165 AliFMDParameters::Hardware2Detector(UShort_t  ddl,       UShort_t   board,
1166                                     UShort_t  chip,      UShort_t   chan,
1167                                     UShort_t  timebin,   
1168                                     UShort_t& det,       Char_t&   ring, 
1169                                     UShort_t& sec,       Short_t& str,
1170                                     UShort_t& sam) const
1171 {
1172   // 
1173   // Map a hardware address into a detector index. 
1174   // 
1175   // Parameters:
1176   //    ddl        Hardware DDL number 
1177   //    board      FEC number
1178   //    altro      ALTRO number 
1179   //    channel    Channel number 
1180   //    timebin    Timebin 
1181   //    det        On return, the detector #
1182   //    ring       On return, the ring ID
1183   //    sec        On return, the sector #
1184   //    str        On return, the base of strip #
1185   //    sam        On return, the sample number for this strip
1186   //
1187   // Return:
1188   //    @c true on success, false otherwise 
1189   //
1190   if (!fAltroMap) {
1191     AliFMDDebug(1, ("No ALTRO map available"));
1192     return kFALSE;
1193   }
1194   if (fAltroMap->DDL2Detector(ddl) < 0) { 
1195     AliFMDDebug(1, ("Invalid DDL number %d", ddl));
1196     return kFALSE;
1197   }
1198   det = fAltroMap->DDL2Detector(ddl);
1199   Short_t stripBase = 0;
1200   if (!fAltroMap->Channel2StripBase(board,chip,chan, ring, sec, stripBase)) {
1201     AliFMDDebug(1, ("Failed to translate  "
1202                     "%d/0x%02x/0x%x/0x%x/%04d -> "
1203                     "FMD%d%c[%2d,%3d] to detector", 
1204                     ddl, board, chip, chan, timebin, 
1205                     det, ring, sec, stripBase));
1206     return kFALSE;
1207   }
1208   UShort_t preSamples = GetPreSamples(det, ring, sec, stripBase);
1209   UShort_t sampleRate = GetSampleRate(det, ring, sec, stripBase);
1210   Short_t stripOff = 0;
1211   fAltroMap->Timebin2Strip(sec, timebin, preSamples, sampleRate, stripOff, sam);
1212   str = stripBase + stripOff;
1213   AliFMDDebug(50, ("%d/0x%02x/0x%x/0x%x/%04d -> FMD%d%c[%02d,%03d]-%d"
1214                   " (pre=%2d, rate=%d)", 
1215                    ddl, board, chip, chan, timebin, 
1216                    det, ring, sec, str, sam, preSamples, sampleRate));
1217   return kTRUE;
1218 }
1219
1220
1221 //____________________________________________________________________
1222 Bool_t 
1223 AliFMDParameters::Detector2Hardware(UShort_t  det,        Char_t    ring, 
1224                                     UShort_t  sec,        UShort_t  str,
1225                                     UShort_t  sam, 
1226                                     UShort_t& ddl,        UShort_t& board, 
1227                                     UShort_t& altro,      UShort_t& channel, 
1228                                     UShort_t& timebin) const
1229 {
1230   // 
1231   // Map a detector index into a hardware address. 
1232   // 
1233   // Parameters:
1234   //    det         The detector #
1235   //    ring        The ring ID
1236   //    sec         The sector #
1237   //    str         The strip #
1238   //    sam         The sample number 
1239   //    ddl         On return, hardware DDL number 
1240   //    board       On return, the FEC board address (local to DDL)
1241   //    altro       On return, the ALTRO number (local to FEC)
1242   //    channel     On return, the channel number (local to ALTRO)
1243   //    timebin     On return, the timebin number (local to ALTRO)
1244   //
1245   // Return:
1246   //    @c true on success, false otherwise 
1247   //
1248   if (!fAltroMap) { 
1249     AliFMDDebug(1, ("No ALTRO map available"));
1250     return kFALSE;
1251   }
1252   UShort_t preSamples = GetPreSamples(det, ring, sec, str);
1253   UShort_t sampleRate = GetSampleRate(det, ring, sec, str);
1254   UShort_t strip      = str - GetMinStrip(det,ring,sec,str);
1255   return fAltroMap->Detector2Hardware(det, ring, sec, strip, sam,
1256                                       preSamples, sampleRate,
1257                                       ddl, board, altro, channel, timebin);
1258 }
1259
1260   
1261
1262 //____________________________________________________________________
1263 Bool_t 
1264 AliFMDParameters::Detector2Hardware(UShort_t  det,        Char_t    ring, 
1265                                     UShort_t  sec,        UShort_t  str,
1266                                     UShort_t  sam, 
1267                                     UShort_t&   ddl,        UShort_t&   addr,
1268                                     UShort_t& timebin) const
1269 {
1270   // 
1271   // Map a detector index into a hardware address. 
1272   // 
1273   // Parameters:
1274   //    det         The detector #
1275   //    ring        The ring ID
1276   //    sec         The sector #
1277   //    str         The strip #
1278   //    sam         The sample number 
1279   //    ddl         On return, hardware DDL number 
1280   //    addr      On return, hardware address.  
1281   //    timebin     On return, the timebin number (local to ALTRO)
1282   //
1283   // Return:
1284   //    @c true on success, false otherwise 
1285   //
1286   if (!fAltroMap) return kFALSE;
1287   UShort_t preSamples = GetPreSamples(det, ring, sec, str);
1288   UShort_t sampleRate = GetSampleRate(det, ring, sec, str);
1289   UShort_t strip      = str - GetMinStrip(det,ring,sec,str);
1290   return fAltroMap->Detector2Hardware(det, ring, sec, strip, sam,
1291                                       preSamples, sampleRate,
1292                                       ddl, addr, timebin);
1293 }
1294
1295
1296 //__________________________________________________________________
1297 Float_t
1298 AliFMDParameters::GetEdepMip() const 
1299
1300   // 
1301   // Return:
1302   //    The average energy deposited by one MIP 
1303   //
1304   if (fEdepMip <= 0){
1305     AliFMDGeometry* fmd = AliFMDGeometry::Instance();
1306     fEdepMip = (fkSiDeDxMip 
1307                 * fmd->GetRing('I')->GetSiThickness() 
1308                 * fmd->GetSiDensity());
1309   }
1310   return fEdepMip;
1311 }
1312 //____________________________________________________________________
1313 Float_t  
1314 AliFMDParameters::GetDACPerMIP() const
1315 {
1316   // 
1317   // This is the conversion from Digital-to-Analog-Converter setting
1318   // to the number of MIPs. The number was measured in the NBI lab during
1319   // August 2008.
1320   //
1321   // Return:
1322   //    The conversion factor from DAC to ADC 
1323   //
1324   return 29.67;
1325   
1326 }
1327  
1328 //____________________________________________________________________
1329 //
1330 // EOF
1331 //