]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDParameters.cxx
Coding conventions
[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 /** @file    AliFMDParameters.cxx
17     @author  Christian Holm Christensen <cholm@nbi.dk>
18     @date    Mon Mar 27 12:44:26 2006
19     @brief   Manager of FMD parameters     
20 */
21 //____________________________________________________________________
22 //                                                                          
23 // Forward Multiplicity Detector based on Silicon wafers. 
24 //
25 // This class is a singleton that handles various parameters of
26 // the FMD detectors.  
27 // The manager normally serves the parameters from the Conditions
28 // Database (CDB).  These are retrivied by the member function
29 // `Init'.  Optionally, the class can serve hard-coded constants, if
30 // no CDB is available. 
31 //                                                       
32 #include "AliLog.h"                // ALILOG_H
33 #include "AliFMDParameters.h"      // ALIFMDPARAMETERS_H
34 #include "AliFMDGeometry.h"        // ALIFMDGEOMETRY_H
35 #include "AliFMDRing.h"            // ALIFMDRING_H
36 #include "AliFMDCalibGain.h"       // ALIFMDCALIBGAIN_H
37 #include "AliFMDCalibPedestal.h"   // ALIFMDCALIBPEDESTAL_H
38 #include "AliFMDCalibSampleRate.h" // ALIFMDCALIBSAMPLERATE_H
39 #include "AliFMDCalibStripRange.h" // ALIFMDCALIBSTRIPRANGE_H
40 #include "AliFMDAltroMapping.h"    // ALIFMDALTROMAPPING_H
41 #include <AliCDBManager.h>         // ALICDBMANAGER_H
42 #include <AliCDBEntry.h>           // ALICDBMANAGER_H
43 #include <Riostream.h>
44 #include <sstream>
45 #include <TArrayF.h>
46 #include <TH2D.h>
47
48 //====================================================================
49 ClassImp(AliFMDParameters)
50 #if 0
51   ; // This is here to keep Emacs for indenting the next line
52 #endif
53
54 //____________________________________________________________________
55 AliFMDParameters* AliFMDParameters::fgInstance = 0;
56
57 //____________________________________________________________________
58 const char* AliFMDParameters::fgkPulseGain       = "FMD/Calib/PulseGain";
59 const char* AliFMDParameters::fgkPedestal        = "FMD/Calib/Pedestal";
60 const char* AliFMDParameters::fgkDead            = "FMD/Calib/Dead";
61 const char* AliFMDParameters::fgkSampleRate      = "FMD/Calib/SampleRate";
62 const char* AliFMDParameters::fgkAltroMap        = "FMD/Calib/AltroMap";
63 const char* AliFMDParameters::fgkZeroSuppression = "FMD/Calib/ZeroSuppression";
64 const char* AliFMDParameters::fgkStripRange      = "FMD/Calib/StripRange";
65
66
67 //____________________________________________________________________
68 AliFMDParameters* 
69 AliFMDParameters::Instance() 
70 {
71   // Get static instance 
72   if (!fgInstance) fgInstance = new AliFMDParameters;
73   return fgInstance;
74 }
75
76 //____________________________________________________________________
77 AliFMDParameters::AliFMDParameters() 
78   : fIsInit(kFALSE),
79     fkSiDeDxMip(1.664), 
80     fFixedPulseGain(0), 
81     fEdepMip(0),
82     fZeroSuppression(0), 
83     fSampleRate(0), 
84     fPedestal(0), 
85     fPulseGain(0), 
86     fDeadMap(0), 
87     fAltroMap(0), 
88     fStripRange(0)
89 {
90   // Default constructor 
91   SetVA1MipRange();
92   SetAltroChannelSize();
93   SetChannelsPerAltro();
94   SetZeroSuppression();
95   SetSampleRate();
96   SetPedestal();
97   SetPedestalWidth();
98   SetPedestalFactor();
99   SetThreshold();
100   SetStripRange();
101 }
102
103 //__________________________________________________________________
104 void
105 AliFMDParameters::Init()
106 {
107   // Initialize the parameters manager.  We need to get stuff from the
108   // CDB here. 
109   if (fIsInit) return;
110   InitPulseGain();
111   InitPedestal();
112   InitDeadMap();
113   InitSampleRate();
114   InitZeroSuppression();
115   InitAltroMap();
116   fIsInit = kTRUE;
117   
118 }
119
120 //__________________________________________________________________
121 #define DET2IDX(det,ring,sec,str) \
122   (det * 10000 + (ring == 'I' ? 0 : 1000) + str)  
123   
124 //__________________________________________________________________
125 void
126 AliFMDParameters::Draw(Option_t* option)
127 {
128   TString opt(option);
129   enum {
130     kPulseGain,       // Path to PulseGain calib object
131     kThreshold,       // Path to PulseGain calib object
132     kPedestal,        // Path to Pedestal calib object
133     kPedestalWidth,   // Path to Pedestal calib object
134     kDead,            // Path to Dead calib object
135     kSampleRate,      // Path to SampleRate calib object
136     kAltroMap,        // Path to AltroMap calib object
137     kZeroSuppression, // Path to ZeroSuppression cal object
138     kMinStripRange,   // Path to strip range cal object
139     kMaxStripRange    // Path to strip range cal object
140   } what;
141   
142     
143   if      (opt.Contains("dead", TString::kIgnoreCase)) 
144     what = kDead;
145   else if (opt.Contains("threshold",TString::kIgnoreCase)) 
146     what = kThreshold;
147   else if (opt.Contains("gain",TString::kIgnoreCase)) 
148     what = kPulseGain;
149   else if (opt.Contains("pedestal",TString::kIgnoreCase)) 
150     what = kPedestal;
151   else if (opt.Contains("noise",TString::kIgnoreCase)) 
152     what = kPedestalWidth;
153   else if (opt.Contains("zero",TString::kIgnoreCase)) 
154     what = kZeroSuppression;
155   else if (opt.Contains("rate",TString::kIgnoreCase)) 
156     what = kSampleRate;
157   else if (opt.Contains("min",TString::kIgnoreCase)) 
158     what = kMinStripRange;
159   else if (opt.Contains("max",TString::kIgnoreCase)) 
160     what = kMaxStripRange;
161   else if (opt.Contains("map",TString::kIgnoreCase)) 
162     what = kAltroMap;
163   else {
164     Warning("Draw", "unknown parameter: %s\n\tShould be one of\n\t"
165             "dead, threshold, gain, pedestal, noise, zero, rate, "
166             "min, max, map",  
167             option); 
168     return;
169   }
170
171   TArrayD xbins(3 * 512 + 2 * 256 + 5);
172   Int_t i = 1;
173   Bool_t skip = kTRUE;
174   for (UShort_t det = 1; det <= 3; det++) {
175     UShort_t nRings = (det == 1 ? 1 : 2);
176     for (UShort_t iring = 0; iring < nRings; iring++) {
177       UShort_t nStrip  = (iring == 0 ? 512 : 256);
178       Char_t   ring    = (iring == 0 ? 'I' : 'O');
179       for (UShort_t str = 0; str < nStrip; str++) {
180         Int_t idx = DET2IDX(det, ring, 0, str);
181         if (skip) {
182           xbins[i-1] = idx - .5;
183           skip  = kFALSE;
184         }
185         xbins[i] = idx + .5;
186         i++;
187       }
188       skip = kTRUE;
189       i++;
190     }
191   }
192   TArrayD ybins(41);
193   for (Int_t i = 0; i < 41; i++) ybins[i] = Float_t(i - .5);
194   TH2D* hist = new TH2D("calib", Form("Calibration %s", option), 
195                         xbins.fN-1, xbins.fArray,  
196                         ybins.fN-1, ybins.fArray);
197
198   // hist->Draw("Lego");
199   // return;
200   
201   for (UShort_t det = 1; det <= 3; det++) {
202     UShort_t nRings = (det == 1 ? 1 : 2);
203     for (UShort_t iring = 0; iring < nRings; iring++) {
204       UShort_t nSector = (iring == 0 ?  20 : 40);
205       UShort_t nStrip  = (iring == 0 ? 512 : 256);
206       Char_t   ring    = (iring == 0 ? 'I' : 'O');
207       for (UShort_t sec = 0; sec < nSector; sec++) {
208         for (UShort_t str = 0; str < nStrip; str++) {
209           Int_t idx = DET2IDX(det, ring, sec, str);
210           UInt_t ddl, addr;
211           Double_t val = 0;
212           switch (what) {
213           case kPulseGain:       // Path to PulseGain calib object
214             val = GetPulseGain(det,ring,sec,str); break;
215           case kThreshold:       // Path to PulseGain calib object
216             val = GetThreshold(); break;
217           case kPedestal:        // Path to Pedestal calib object
218             val = GetPedestal(det,ring,sec,str); break;
219           case kPedestalWidth:   // Path to Pedestal calib object
220             val = GetPedestalWidth(det,ring,sec,str); break;
221           case kDead:            // Path to Dead calib object
222             val = IsDead(det,ring,sec,str); break;
223           case kSampleRate:      // Path to SampleRate calib object
224             val = GetSampleRate(det,ring,sec,str); break;
225           case kAltroMap:        // Path to AltroMap calib object
226             Detector2Hardware(det,ring,sec,str, ddl, addr); 
227             val = addr; break;
228           case kZeroSuppression: // Path to ZeroSuppression cal object
229             val = GetZeroSuppression(det,ring,sec,str); break;
230           case kMinStripRange:   // Path to strip range cal object
231             val = GetMinStrip(det,ring,sec,str); break;
232           case kMaxStripRange:    // Path to strip range cal object
233             val = GetMaxStrip(det,ring,sec,str); break;
234           }
235           hist->Fill(idx,sec,val);
236         }
237       }
238     }
239   }
240   hist->Draw("lego");
241 }
242
243 //__________________________________________________________________
244 void
245 AliFMDParameters::Print(Option_t* option) const
246 {
247   // Print information. 
248   // If option contains an 'A' then everything is printed. 
249   TString opt(option);
250   Bool_t showStrips = opt.Contains("a", TString::kIgnoreCase);
251   if (opt.Contains("fmd",TString::kIgnoreCase)) {
252     size_t   i   = opt.Index("fmd",TString::kIgnoreCase);
253     size_t   j   = opt.Index("]",TString::kIgnoreCase);
254     UShort_t det, sec, str;
255     Char_t ring, lbrack, rbrack, comma;
256     UInt_t ddl, addr;
257     std::stringstream s(opt(i+4, j-i-3).Data());
258     s >> det >> ring >> lbrack >> sec >> comma >> str >> rbrack;
259     Detector2Hardware(det, ring, sec, str, ddl, addr);
260     std::cout 
261       << "     Strip    |     Pedestal      |    Gain    | ZS thr. | Address\n"
262       << "--------------+-------------------+------------+---------+---------" 
263       << "\nFMD" << det << ring << "[" << std::setw(2) << sec << "," 
264       << std::setw(3) << str << "] | " 
265       << std::setw(7) << GetPedestal(det, ring, sec, str) 
266       << "+/-" << std::setw(7) 
267       << GetPedestalWidth(det, ring, sec, str) 
268       << " | " << std::setw(10) 
269       << GetPulseGain(det, ring, sec, str) 
270       << " | " << std::setw(7) 
271       << GetZeroSuppression(det, ring, sec, str) 
272       << " | 0x" << std::hex << std::setw(4) 
273       << std::setfill('0') << ddl << ",0x" << std::setw(3) 
274       << addr << std::dec << std::setfill(' ') << std::endl;
275     return;
276   }
277   for (UShort_t det=1 ; det <= 3; det++) {
278     std::cout << "FMD" << det << std::endl;
279     Char_t rings[] = { 'I', (det == 1 ? '\0' : 'O'), '\0' };
280     for (Char_t* ring = rings; *ring != '\0'; ring++) {
281       std::cout << " Ring " << *ring << std::endl;
282       UShort_t nSec = ( *ring == 'I' ? 20  :  40 );
283       UShort_t nStr = ( *ring == 'I' ? 512 : 256 );
284       for (UShort_t sec = 0; sec < nSec; sec++) {
285         UShort_t min  = GetMinStrip(det, *ring, sec, 0);
286         UShort_t max  = GetMaxStrip(det, *ring, sec, 0);
287         UShort_t rate = GetSampleRate(det, *ring, sec, 0);
288         std::cout << "  Sector " << std::setw(2) << sec 
289                   << "  Strip range: " << std::setw(3) << min << "," 
290                   << std::setw(3) << max << "  Rate: " << std::setw(2) 
291                   << rate << std::endl;
292         if (!showStrips) continue;
293         std::cout 
294           << "  Strip |     Pedestal      |    Gain    | ZS thr. | Address\n" 
295           << "--------+-------------------+------------+---------+---------" 
296           << std::endl;
297         for (UShort_t str = 0; str < nStr; str++) {
298           std::cout << "    " << std::setw(3) << str << " | ";
299           if (IsDead(det, *ring, sec, str)) {
300             std::cout << "dead" << std::endl;
301             continue;
302           }
303           UInt_t ddl, addr;
304           Detector2Hardware(det, *ring, sec, str, ddl, addr);
305           std::cout << std::setw(7) << GetPedestal(det, *ring, sec, str) 
306                     << "+/-" << std::setw(7) 
307                     << GetPedestalWidth(det, *ring, sec, str) 
308                     << " | " << std::setw(10) 
309                     << GetPulseGain(det, *ring, sec, str) 
310                     << " | " << std::setw(5) 
311                     << GetZeroSuppression(det, *ring, sec, str) 
312                     << " | 0x" << std::hex << std::setw(4) 
313                     << std::setfill('0') << ddl << ",0x" << std::setw(3) 
314                     << addr << std::dec << std::setfill(' ') << std::endl;
315         }
316       }
317     }
318   }
319 }
320
321 //__________________________________________________________________
322 void
323 AliFMDParameters::SetStripRange(UShort_t min, UShort_t max) 
324 {
325   // Set fixed strip range 
326   fFixedMinStrip = min;
327   fFixedMaxStrip = max;
328 }
329
330 //__________________________________________________________________
331 void
332 AliFMDParameters::InitPulseGain()
333 {
334   // Get pulse gain from CDB or used fixed 
335   AliCDBManager* cdb      = AliCDBManager::Instance();
336   AliCDBEntry*   gain     = cdb->Get(fgkPulseGain);
337   if (!gain) {
338     AliWarning(Form("No %s found in CDB, perhaps you need to "
339                     "use AliFMDCalibFaker?", fgkPulseGain));
340     return;
341   }
342   
343   AliDebug(1, Form("Got gain from CDB"));
344   fPulseGain = dynamic_cast<AliFMDCalibGain*>(gain->GetObject());
345   if (!fPulseGain) AliWarning("Invalid pulser gain object from CDB");
346 }
347 //__________________________________________________________________
348 void
349 AliFMDParameters::InitPedestal()
350 {
351   // Initialize the pedestals from CDB 
352   AliCDBManager* cdb      = AliCDBManager::Instance();
353   AliCDBEntry*   pedestal = cdb->Get(fgkPedestal);
354   if (!pedestal) {
355     AliWarning(Form("No %s found in CDB, perhaps you need to "
356                     "use AliFMDCalibFaker?", fgkPedestal));
357     return;
358   }
359   AliDebug(1, Form("Got pedestal from CDB"));
360   fPedestal = dynamic_cast<AliFMDCalibPedestal*>(pedestal->GetObject());
361   if (!fPedestal) AliWarning("Invalid pedestal object from CDB");
362 }
363
364 //__________________________________________________________________
365 void
366 AliFMDParameters::InitDeadMap()
367 {
368   // Get Dead-channel-map from CDB 
369   AliCDBManager* cdb      = AliCDBManager::Instance();
370   AliCDBEntry*   deadMap  = cdb->Get(fgkDead);
371   if (!deadMap) {
372     AliWarning(Form("No %s found in CDB, perhaps you need to "
373                     "use AliFMDCalibFaker?", fgkDead));
374     return;
375   }
376   AliDebug(1, Form("Got dead map from CDB"));
377   fDeadMap = dynamic_cast<AliFMDCalibDeadMap*>(deadMap->GetObject());
378   if (!fDeadMap) AliWarning("Invalid dead map object from CDB");
379 }
380
381 //__________________________________________________________________
382 void
383 AliFMDParameters::InitZeroSuppression()
384 {
385   // Get 0-suppression from CDB 
386   AliCDBManager* cdb      = AliCDBManager::Instance();
387   AliCDBEntry*   zeroSup  = cdb->Get(fgkZeroSuppression);
388   if (!zeroSup) {
389     AliWarning(Form("No %s found in CDB, perhaps you need to "
390                     "use AliFMDCalibFaker?", fgkZeroSuppression));
391     return;
392   }
393   AliDebug(1, Form("Got zero suppression from CDB"));
394   fZeroSuppression = 
395     dynamic_cast<AliFMDCalibZeroSuppression*>(zeroSup->GetObject());
396   if (!fZeroSuppression)AliWarning("Invalid zero suppression object from CDB");
397 }
398
399 //__________________________________________________________________
400 void
401 AliFMDParameters::InitSampleRate()
402 {
403   // get Sample rate from CDB
404   AliCDBManager* cdb      = AliCDBManager::Instance();
405   AliCDBEntry*   sampRat  = cdb->Get(fgkSampleRate);
406   if (!sampRat) {
407     AliWarning(Form("No %s found in CDB, perhaps you need to "
408                     "use AliFMDCalibFaker?", fgkSampleRate));
409     return;
410   }
411   AliDebug(1, Form("Got zero suppression from CDB"));
412   fSampleRate = dynamic_cast<AliFMDCalibSampleRate*>(sampRat->GetObject());
413   if (!fSampleRate) AliWarning("Invalid zero suppression object from CDB");
414 }
415
416 //__________________________________________________________________
417 void
418 AliFMDParameters::InitAltroMap()
419 {
420   // Get hardware mapping from CDB
421   AliCDBManager* cdb      = AliCDBManager::Instance();
422   AliCDBEntry*   hwMap    = cdb->Get(fgkAltroMap);       
423   if (!hwMap) {
424     AliWarning(Form("No %s found in CDB, perhaps you need to "
425                     "use AliFMDCalibFaker?", fgkAltroMap));
426     fAltroMap = new AliFMDAltroMapping;
427     return;
428   }
429   AliDebug(1, Form("Got ALTRO map from CDB"));
430   fAltroMap = dynamic_cast<AliFMDAltroMapping*>(hwMap->GetObject());
431   if (!fAltroMap) {
432     AliWarning("Invalid ALTRO map object from CDB");
433     fAltroMap = new AliFMDAltroMapping;
434   }
435 }
436
437 //__________________________________________________________________
438 void
439 AliFMDParameters::InitStripRange()
440 {
441   // Get strips read-out from CDB
442   AliCDBManager* cdb      = AliCDBManager::Instance();
443   AliCDBEntry*   range    = cdb->Get(fgkStripRange);
444   if (!range) {
445     AliWarning(Form("No %s found in CDB, perhaps you need to "
446                     "use AliFMDCalibFaker?", fgkStripRange));
447     return;
448   }
449   AliDebug(1, Form("Got strip range from CDB"));
450   fStripRange = dynamic_cast<AliFMDCalibStripRange*>(range->GetObject());
451   if (!fStripRange) AliWarning("Invalid strip range object from CDB");
452 }
453
454
455 //__________________________________________________________________
456 Float_t
457 AliFMDParameters::GetThreshold() const
458 {
459   // Get threshold from CDB
460   if (!fPulseGain) return fFixedThreshold;
461   return fPulseGain->Threshold();
462 }
463
464 //__________________________________________________________________
465 Float_t
466 AliFMDParameters::GetPulseGain(UShort_t detector, Char_t ring, 
467                                UShort_t sector, UShort_t strip) const
468 {
469   // Returns the pulser calibrated gain for strip # strip in sector #
470   // sector or ring id ring of detector # detector. 
471   // 
472   // For simulation, this is normally set to 
473   // 
474   //       VA1_MIP_Range 
475   //    ------------------ * MIP_Energy_Loss
476   //    ALTRO_channel_size
477   // 
478   if (!fPulseGain) { 
479     if (fFixedPulseGain <= 0)
480       fFixedPulseGain = fVA1MipRange * GetEdepMip() / fAltroChannelSize;
481     return fFixedPulseGain;
482   }  
483   AliDebug(50, Form("pulse gain for FMD%d%c[%2d,%3d]=%f",
484                     detector, ring, sector, strip,
485                     fPulseGain->Value(detector, ring, sector, strip)));
486   return fPulseGain->Value(detector, ring, sector, strip);
487 }
488
489 //__________________________________________________________________
490 Bool_t
491 AliFMDParameters::IsDead(UShort_t detector, Char_t ring, 
492                          UShort_t sector, UShort_t strip) const
493 {
494   // Check if the channel is dead 
495   if (!fDeadMap) return kFALSE;
496   AliDebug(50, Form("Dead for FMD%d%c[%2d,%3d]=%s",
497                     detector, ring, sector, strip,
498                     fDeadMap->operator()(detector, ring, sector, strip) ? 
499                     "no" : "yes"));
500   return fDeadMap->operator()(detector, ring, sector, strip);
501 }
502
503 //__________________________________________________________________
504 UShort_t
505 AliFMDParameters::GetZeroSuppression(UShort_t detector, Char_t ring, 
506                                      UShort_t sector, UShort_t strip) const
507 {
508   // Get zero suppression threshold 
509   if (!fZeroSuppression) return fFixedZeroSuppression;
510   // Need to map strip to ALTRO chip. 
511   AliDebug(50, Form("zero sup. for FMD%d%c[%2d,%3d]=%f",
512                     detector, ring, sector, strip,
513                     fZeroSuppression->operator()(detector, ring, 
514                                                  sector, strip)));
515   return fZeroSuppression->operator()(detector, ring, sector, strip/128);
516 }
517
518 //__________________________________________________________________
519 UShort_t
520 AliFMDParameters::GetSampleRate(UShort_t det, Char_t ring, UShort_t sector, 
521                                 UShort_t str) const
522 {
523   // Get sampl rate 
524   if (!fSampleRate) return fFixedSampleRate;
525   // Need to map sector to digitizier card. 
526   UInt_t ret = fSampleRate->Rate(det, ring, sector, str);
527   AliDebug(50, Form("Sample rate for FMD%d%c[%2d,%3d]=%d", 
528                     det, ring, sector, str, ret));
529   return ret;
530 }
531
532 //__________________________________________________________________
533 UShort_t
534 AliFMDParameters::GetMinStrip(UShort_t det, Char_t ring, UShort_t sector, 
535                               UShort_t str) const
536 {
537   // Get strip range read out 
538   if (!fStripRange) return fFixedMinStrip;
539   // Need to map sector to digitizier card. 
540   UInt_t ret = fStripRange->Min(det, ring, sector, str);
541   AliDebug(50, Form("Min strip # for FMD%d%c[%2d,%3d]=%d", 
542                     det, ring, sector, str, ret));
543   return ret;
544 }
545
546 //__________________________________________________________________
547 UShort_t
548 AliFMDParameters::GetMaxStrip(UShort_t det, Char_t ring, UShort_t sector, 
549                               UShort_t str) const
550 {
551   // Get strip range read out 
552   if (!fStripRange) return fFixedMaxStrip;
553   // Need to map sector to digitizier card. 
554   UInt_t ret = fStripRange->Max(det, ring, sector, str);
555   AliDebug(50, Form("Max strip # for FMD%d%c[%2d,%3d]=%d", 
556                     det, ring, sector, str, ret));
557   return ret;
558 }
559
560 //__________________________________________________________________
561 Float_t
562 AliFMDParameters::GetPedestal(UShort_t detector, Char_t ring, 
563                               UShort_t sector, UShort_t strip) const
564 {
565   // Get the pedesal 
566   if (!fPedestal) return fFixedPedestal;
567   AliDebug(50, Form("pedestal for FMD%d%c[%2d,%3d]=%f",
568                     detector, ring, sector, strip,
569                     fPedestal->Value(detector, ring, sector, strip)));
570   return fPedestal->Value(detector, ring, sector, strip);
571 }
572
573 //__________________________________________________________________
574 Float_t
575 AliFMDParameters::GetPedestalWidth(UShort_t detector, Char_t ring, 
576                                    UShort_t sector, UShort_t strip) const
577 {
578   // Get the pedesal 
579   if (!fPedestal) return fFixedPedestalWidth;
580   AliDebug(50, Form("pedetal width for FMD%d%c[%2d,%3d]=%f",
581                     detector, ring, sector, strip,
582                     fPedestal->Width(detector, ring, sector, strip)));
583   return fPedestal->Width(detector, ring, sector, strip);
584 }
585   
586 //__________________________________________________________________
587 AliFMDAltroMapping*
588 AliFMDParameters::GetAltroMap() const
589 {
590   // Get the hardware address to detector index map 
591   return fAltroMap;
592 }
593
594
595 //__________________________________________________________________
596 Bool_t
597 AliFMDParameters::Hardware2Detector(UInt_t ddl, UInt_t addr, UShort_t& det,
598                                     Char_t& ring, UShort_t& sec, 
599                                     UShort_t& str) const
600 {
601   // Map hardware address to detector index
602   if (!fAltroMap) return kFALSE;
603   return fAltroMap->Hardware2Detector(ddl, addr, det, ring, sec, str);
604 }
605
606 //__________________________________________________________________
607 Bool_t
608 AliFMDParameters::Detector2Hardware(UShort_t det, Char_t ring, UShort_t sec, 
609                                     UShort_t str, UInt_t& ddl, 
610                                     UInt_t& addr) const                       
611 {
612   // Map detector index to hardware address
613   if (!fAltroMap) return kFALSE;
614   return fAltroMap->Detector2Hardware(det, ring, sec, str, ddl, addr);
615 }
616
617
618 //__________________________________________________________________
619 Float_t
620 AliFMDParameters::GetEdepMip() const 
621
622   // Get energy deposited by a MIP in the silicon sensors
623   if (fEdepMip <= 0){
624     AliFMDGeometry* fmd = AliFMDGeometry::Instance();
625     fEdepMip = (fkSiDeDxMip 
626                 * fmd->GetRing('I')->GetSiThickness() 
627                 * fmd->GetSiDensity());
628   }
629   return fEdepMip;
630 }
631
632
633   
634   
635   
636 //____________________________________________________________________
637 //
638 // EOF
639 //