]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDParameters.cxx
Fix to avoid memory resident tree in the vertexer (Marco)
[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 "AliFMDDebug.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     fVA1MipRange(0),
81     fAltroChannelSize(0),
82     fChannelsPerAltro(0),
83     fPedestalFactor(0),
84     fFixedPedestal(0),
85     fFixedPedestalWidth(0),
86     fFixedZeroSuppression(0),
87     fFixedSampleRate(0),
88     fFixedThreshold(0),
89     fFixedMinStrip(0),
90     fFixedMaxStrip(0),
91     fFixedPulseGain(0), 
92     fEdepMip(0),
93     fZeroSuppression(0), 
94     fSampleRate(0), 
95     fPedestal(0), 
96     fPulseGain(0), 
97     fDeadMap(0), 
98     fAltroMap(0), 
99     fStripRange(0)
100 {
101   // Default constructor 
102   SetVA1MipRange();
103   SetAltroChannelSize();
104   SetChannelsPerAltro();
105   SetZeroSuppression();
106   SetSampleRate();
107   SetPedestal();
108   SetPedestalWidth();
109   SetPedestalFactor();
110   SetThreshold();
111   SetStripRange();
112 }
113
114 //__________________________________________________________________
115 void
116 AliFMDParameters::Init(Bool_t forceReInit)
117 {
118   // Initialize the parameters manager.  We need to get stuff from the
119   // CDB here. 
120   if (forceReInit) fIsInit = kFALSE;
121   if (fIsInit) return;
122   InitPulseGain();
123   InitPedestal();
124   InitDeadMap();
125   InitSampleRate();
126   InitZeroSuppression();
127   InitAltroMap();
128   fIsInit = kTRUE;
129   
130 }
131
132 //__________________________________________________________________
133 #define DET2IDX(det,ring,sec,str) \
134   (det * 1000 + (ring == 'I' ? 0 : 512) + str)  
135   
136 //__________________________________________________________________
137 void
138 AliFMDParameters::Draw(Option_t* option)
139 {
140   TString opt(option);
141   enum {
142     kPulseGain,       // Path to PulseGain calib object
143     kThreshold,       // Path to PulseGain calib object
144     kPedestal,        // Path to Pedestal calib object
145     kPedestalWidth,   // Path to Pedestal calib object
146     kDead,            // Path to Dead calib object
147     kSampleRate,      // Path to SampleRate calib object
148     kAltroMap,        // Path to AltroMap calib object
149     kZeroSuppression, // Path to ZeroSuppression cal object
150     kMinStripRange,   // Path to strip range cal object
151     kMaxStripRange    // Path to strip range cal object
152   } what;
153   
154     
155   if      (opt.Contains("dead", TString::kIgnoreCase)) 
156     what = kDead;
157   else if (opt.Contains("threshold",TString::kIgnoreCase)) 
158     what = kThreshold;
159   else if (opt.Contains("gain",TString::kIgnoreCase)) 
160     what = kPulseGain;
161   else if (opt.Contains("pedestal",TString::kIgnoreCase)) 
162     what = kPedestal;
163   else if (opt.Contains("noise",TString::kIgnoreCase)) 
164     what = kPedestalWidth;
165   else if (opt.Contains("zero",TString::kIgnoreCase)) 
166     what = kZeroSuppression;
167   else if (opt.Contains("rate",TString::kIgnoreCase)) 
168     what = kSampleRate;
169   else if (opt.Contains("min",TString::kIgnoreCase)) 
170     what = kMinStripRange;
171   else if (opt.Contains("max",TString::kIgnoreCase)) 
172     what = kMaxStripRange;
173   else if (opt.Contains("map",TString::kIgnoreCase)) 
174     what = kAltroMap;
175   else {
176     Warning("Draw", "unknown parameter: %s\n\tShould be one of\n\t"
177             "dead, threshold, gain, pedestal, noise, zero, rate, "
178             "min, max, map",  
179             option); 
180     return;
181   }
182
183   TArrayD xbins(3 * 512 + 2 * 256 + 5);
184   Int_t i = 1;
185   Bool_t skip = kTRUE;
186   for (UShort_t det = 1; det <= 3; det++) {
187     UShort_t nRings = (det == 1 ? 1 : 2);
188     for (UShort_t iring = 0; iring < nRings; iring++) {
189       UShort_t nStrip  = (iring == 0 ? 512 : 256);
190       Char_t   ring    = (iring == 0 ? 'I' : 'O');
191       for (UShort_t str = 0; str < nStrip; str++) {
192         // UShort_t nSec    = (iring == 0 ? 20  : 40);
193         // Char_t   ring    = (iring == 0 ? 'I' : 'O');
194         // for (UShort_t sec = 0; sec < nSec; sec++) {
195         Int_t idx = DET2IDX(det, ring, 0, str);
196         // Int_t idx = DET2IDX(det, ring, sec, 0);
197         if (skip) {
198           xbins[i-1] = idx - .5;
199           skip  = kFALSE;
200         }
201         xbins[i] = idx + .5;
202         i++;
203       }
204       skip = kTRUE;
205       i++;
206     }
207   }
208   TArrayD ybins(41);
209   for (Int_t i = 0; i < ybins.fN; i++) ybins[i] = Float_t(i - .5);
210   TH2D* hist = new TH2D("calib", Form("Calibration %s", option), 
211                         xbins.fN-1, xbins.fArray,  
212                         ybins.fN-1, ybins.fArray);
213   hist->GetXaxis()->SetTitle("1000 #times detector + 512 #times ring + strip");
214   hist->GetYaxis()->SetTitle("sector");
215   
216   // hist->Draw("Lego");
217   // return;
218   
219   for (UShort_t det = 1; det <= 3; det++) {
220     UShort_t nRings = (det == 1 ? 1 : 2);
221     for (UShort_t iring = 0; iring < nRings; iring++) {
222       UShort_t nSector = (iring == 0 ?  20 : 40);
223       UShort_t nStrip  = (iring == 0 ? 512 : 256);
224       Char_t   ring    = (iring == 0 ? 'I' : 'O');
225       for (UShort_t sec = 0; sec < nSector; sec++) {
226         for (UShort_t str = 0; str < nStrip; str++) {
227           Int_t idx = DET2IDX(det, ring, sec, str);
228           UInt_t ddl, addr;
229           Double_t val = 0;
230           switch (what) {
231           case kPulseGain:       // Path to PulseGain calib object
232             val = GetPulseGain(det,ring,sec,str); break;
233           case kThreshold:       // Path to PulseGain calib object
234             val = GetThreshold(); break;
235           case kPedestal:        // Path to Pedestal calib object
236             val = GetPedestal(det,ring,sec,str); break;
237           case kPedestalWidth:   // Path to Pedestal calib object
238             val = GetPedestalWidth(det,ring,sec,str); break;
239           case kDead:            // Path to Dead calib object
240             val = IsDead(det,ring,sec,str); break;
241           case kSampleRate:      // Path to SampleRate calib object
242             val = GetSampleRate(det,ring,sec,str); break;
243           case kAltroMap:        // Path to AltroMap calib object
244             Detector2Hardware(det,ring,sec,str, ddl, addr); 
245             val = addr; break;
246           case kZeroSuppression: // Path to ZeroSuppression cal object
247             val = GetZeroSuppression(det,ring,sec,str); break;
248           case kMinStripRange:   // Path to strip range cal object
249             val = GetMinStrip(det,ring,sec,str); break;
250           case kMaxStripRange:    // Path to strip range cal object
251             val = GetMaxStrip(det,ring,sec,str); break;
252           }
253           hist->Fill(idx,sec,val);
254           // hist->Fill(idx,str,val);
255         }
256       }
257     }
258   }
259   hist->Draw("lego");
260 }
261
262 //__________________________________________________________________
263 void
264 AliFMDParameters::Print(Option_t* option) const
265 {
266   // Print information. 
267   // If option contains an 'A' then everything is printed. 
268   // If the option contains the string "FMD" the function will search 
269   // for detector, ring, sector, and strip numbers to print, in the
270   // format 
271   // 
272   //    FMD<detector><ring>[<sector>,<string>] 
273   // 
274   // The wild card '*' means all of <detector>, <ring>, <sector>, or 
275   // <strip>. 
276   TString opt(option);
277   Bool_t showStrips  = opt.Contains("a", TString::kIgnoreCase);
278   UShort_t ds[]      = { 1, 2, 3, 0 };
279   Char_t   rs[]      = { 'I', 'O', '\0' };
280   UShort_t minStrip  = 0;
281   UShort_t maxStrip  = 512;
282   UShort_t minSector = 0;
283   UShort_t maxSector = 40;
284   
285   
286   if (opt.Contains("fmd",TString::kIgnoreCase)) {
287     showStrips    = kTRUE;
288     size_t   i    = opt.Index("fmd",TString::kIgnoreCase);
289     size_t   j    = opt.Index("]",TString::kIgnoreCase);
290     enum {
291       kReadDet, 
292       kReadRing, 
293       kReadLbrack,
294       kReadSector,
295       kReadComma,
296       kReadStrip,
297       kReadRbrack, 
298       kEnd
299     } state = kReadDet;
300     std::stringstream s(opt(i+4, j-i-3).Data());
301     while (state != kEnd) {
302       Char_t tmp = s.peek();
303       if (tmp == ' ' || tmp == '\t') {
304         s.get();
305         continue;
306       }
307       switch (state) {
308       case kReadDet: { // First, try to kRead the detector 
309         if (tmp == '*') s.get();
310         else { 
311           UShort_t det;
312           s >> det;
313           if (!s.bad()) {
314             ds[0] = det;
315             ds[1] = 0;
316           }
317         }
318         state = (s.bad() ? kEnd : kReadRing);
319       } break;
320       case kReadRing: { // Then try to read the ring;
321         Char_t ring;
322         s >> ring;
323         if (ring != '*' && !s.bad()) {
324           rs[0] = ring;
325           rs[1] = '\0';
326         }
327         state = (s.bad() ? kEnd : kReadLbrack);
328       } break;
329       case kReadLbrack: { // Try to read a left bracket 
330         Char_t lbrack;
331         s >> lbrack;
332         state = (s.bad() ? kEnd : kReadSector);
333       } break;
334       case kReadSector: { // Try to read a sector 
335         if (tmp == '*') s.get();
336         else {
337           UShort_t sec;
338           s >> sec;
339           if (!s.bad()) {
340             minSector = sec;
341             maxSector = sec + 1;
342           }
343         }
344         state = (s.bad() ? kEnd : kReadComma);
345       } break;
346       case kReadComma: { // Try to read a left bracket 
347         Char_t comma;
348         s >> comma;
349         state = (s.bad() ? kEnd : kReadStrip);
350       } break;
351       case kReadStrip: { // Try to read a strip 
352         if (tmp == '*') s.get();
353         else {
354           UShort_t str;
355           s >> str;
356           if (!s.bad()) {
357             minStrip = str;
358             maxStrip = str + 1;
359           }
360         }
361         state = (s.bad() ? kEnd : kReadRbrack);
362       } break;
363       case kReadRbrack: { // Try to read a left bracket 
364         Char_t rbrack;
365         s >> rbrack;
366         state = kEnd;
367       } break;
368       case kEnd: 
369         break;
370       }
371     }
372   }
373   UShort_t* dp = ds;
374   UShort_t  det;
375   while ((det = *(dp++))) {
376
377     Char_t* rp = rs;
378     Char_t  ring;
379     while ((ring = *(rp++))) {
380       if (det == 1 && ring == 'O') continue;
381       UShort_t min  = GetMinStrip(det, ring, 0, 0);
382       UShort_t max  = GetMaxStrip(det, ring, 0, 0);
383       UShort_t rate = GetSampleRate(det, ring, 0, 0);
384       std::cout << "FMD" << det << ring 
385                 << "  Strip range: " 
386                 << std::setw(3) << min << "," 
387                 << std::setw(3) << max << "  Rate: " 
388                 << std::setw(2) << rate << std::endl;
389
390       if (!showStrips) continue;
391       UShort_t nSec = ( ring == 'I' ? 20  :  40 );
392       UShort_t nStr = ( ring == 'I' ? 512 : 256 );
393       for (UShort_t sec = minSector; sec < maxSector && sec < nSec; sec++) {
394         std::cout 
395           << "  Strip |     Pedestal      |    Gain    | ZS thr. | Address\n" 
396           << "--------+-------------------+------------+---------+---------" 
397           << std::endl;
398         for (UShort_t str = minStrip; str < nStr && str < maxStrip; str++) {
399           if (str == minStrip) std::cout << std::setw(3) << sec << ",";
400           else std::cout << "    ";
401           std::cout << std::setw(3) << str << " | ";
402           if (IsDead(det, ring, sec, str)) {
403             std::cout << "dead" << std::endl;
404             continue;
405           }
406           UInt_t ddl, addr;
407           Detector2Hardware(det, ring, sec, str, ddl, addr);
408           std::cout << std::setw(7) << GetPedestal(det, ring, sec, str) 
409                     << "+/-" << std::setw(7) 
410                     << GetPedestalWidth(det, ring, sec, str) 
411                     << " | " << std::setw(10) 
412                     << GetPulseGain(det, ring, sec, str) 
413                     << " | " << std::setw(7) 
414                     << GetZeroSuppression(det, ring, sec, str) 
415                     << " | 0x" << std::hex << std::setw(4) 
416                     << std::setfill('0') << ddl << ",0x" << std::setw(3) 
417                     << addr << std::dec << std::setfill(' ') << std::endl;
418         } // for (strip)
419       } // for (sector)
420       std::cout
421         << "=============================================================" 
422         << std::endl;
423     } // while (ring)
424   } // while (det)
425   
426 }
427
428 //__________________________________________________________________
429 void
430 AliFMDParameters::SetStripRange(UShort_t min, UShort_t max) 
431 {
432   // Set fixed strip range 
433   fFixedMinStrip = min;
434   fFixedMaxStrip = max;
435 }
436
437 //__________________________________________________________________
438 void
439 AliFMDParameters::InitPulseGain()
440 {
441   // Get pulse gain from CDB or used fixed 
442   AliCDBManager* cdb      = AliCDBManager::Instance();
443   AliCDBEntry*   gain     = cdb->Get(fgkPulseGain);
444   if (!gain) {
445     AliFatal(Form("No %s found in CDB, perhaps you need to "
446                     "use AliFMDCalibFaker?", fgkPulseGain));
447     return;
448   }
449   
450   AliFMDDebug(1, ("Got gain from CDB"));
451   fPulseGain = dynamic_cast<AliFMDCalibGain*>(gain->GetObject());
452   if (!fPulseGain) AliFatal("Invalid pulser gain object from CDB");
453 }
454 //__________________________________________________________________
455 void
456 AliFMDParameters::InitPedestal()
457 {
458   // Initialize the pedestals from CDB 
459   AliCDBManager* cdb      = AliCDBManager::Instance();
460   AliCDBEntry*   pedestal = cdb->Get(fgkPedestal);
461   if (!pedestal) {
462     AliFatal(Form("No %s found in CDB, perhaps you need to "
463                     "use AliFMDCalibFaker?", fgkPedestal));
464     return;
465   }
466   AliFMDDebug(1, ("Got pedestal from CDB"));
467   fPedestal = dynamic_cast<AliFMDCalibPedestal*>(pedestal->GetObject());
468   if (!fPedestal) AliFatal("Invalid pedestal object from CDB");
469 }
470
471 //__________________________________________________________________
472 void
473 AliFMDParameters::InitDeadMap()
474 {
475   // Get Dead-channel-map from CDB 
476   AliCDBManager* cdb      = AliCDBManager::Instance();
477   AliCDBEntry*   deadMap  = cdb->Get(fgkDead);
478   if (!deadMap) {
479     AliFatal(Form("No %s found in CDB, perhaps you need to "
480                     "use AliFMDCalibFaker?", fgkDead));
481     return;
482   }
483   AliFMDDebug(1, ("Got dead map from CDB"));
484   fDeadMap = dynamic_cast<AliFMDCalibDeadMap*>(deadMap->GetObject());
485   if (!fDeadMap) AliFatal("Invalid dead map object from CDB");
486 }
487
488 //__________________________________________________________________
489 void
490 AliFMDParameters::InitZeroSuppression()
491 {
492   // Get 0-suppression from CDB 
493   AliCDBManager* cdb      = AliCDBManager::Instance();
494   AliCDBEntry*   zeroSup  = cdb->Get(fgkZeroSuppression);
495   if (!zeroSup) {
496     AliFatal(Form("No %s found in CDB, perhaps you need to "
497                     "use AliFMDCalibFaker?", fgkZeroSuppression));
498     return;
499   }
500   AliFMDDebug(1, ("Got zero suppression from CDB"));
501   fZeroSuppression = 
502     dynamic_cast<AliFMDCalibZeroSuppression*>(zeroSup->GetObject());
503   if (!fZeroSuppression)AliFatal("Invalid zero suppression object from CDB");
504 }
505
506 //__________________________________________________________________
507 void
508 AliFMDParameters::InitSampleRate()
509 {
510   // get Sample rate from CDB
511   AliCDBManager* cdb      = AliCDBManager::Instance();
512   AliCDBEntry*   sampRat  = cdb->Get(fgkSampleRate);
513   if (!sampRat) {
514     AliFatal(Form("No %s found in CDB, perhaps you need to "
515                     "use AliFMDCalibFaker?", fgkSampleRate));
516     return;
517   }
518   AliFMDDebug(1, ("Got zero suppression from CDB"));
519   fSampleRate = dynamic_cast<AliFMDCalibSampleRate*>(sampRat->GetObject());
520   if (!fSampleRate) AliFatal("Invalid zero suppression object from CDB");
521 }
522
523 //__________________________________________________________________
524 void
525 AliFMDParameters::InitAltroMap()
526 {
527   // Get hardware mapping from CDB
528   AliCDBManager* cdb      = AliCDBManager::Instance();
529   AliCDBEntry*   hwMap    = cdb->Get(fgkAltroMap);       
530   if (!hwMap) {
531     AliFatal(Form("No %s found in CDB, perhaps you need to "
532                     "use AliFMDCalibFaker?", fgkAltroMap));
533     fAltroMap = new AliFMDAltroMapping;
534     return;
535   }
536   AliFMDDebug(1, ("Got ALTRO map from CDB"));
537   fAltroMap = dynamic_cast<AliFMDAltroMapping*>(hwMap->GetObject());
538   if (!fAltroMap) {
539     AliFatal("Invalid ALTRO map object from CDB");
540     fAltroMap = new AliFMDAltroMapping;
541   }
542 }
543
544 //__________________________________________________________________
545 void
546 AliFMDParameters::InitStripRange()
547 {
548   // Get strips read-out from CDB
549   AliCDBManager* cdb      = AliCDBManager::Instance();
550   AliCDBEntry*   range    = cdb->Get(fgkStripRange);
551   if (!range) {
552     AliFatal(Form("No %s found in CDB, perhaps you need to "
553                     "use AliFMDCalibFaker?", fgkStripRange));
554     return;
555   }
556   AliFMDDebug(1, ("Got strip range from CDB"));
557   fStripRange = dynamic_cast<AliFMDCalibStripRange*>(range->GetObject());
558   if (!fStripRange) AliFatal("Invalid strip range object from CDB");
559 }
560
561
562 //__________________________________________________________________
563 Float_t
564 AliFMDParameters::GetThreshold() const
565 {
566   // Get threshold from CDB
567   if (!fPulseGain) return fFixedThreshold;
568   return fPulseGain->Threshold();
569 }
570
571 //__________________________________________________________________
572 Float_t
573 AliFMDParameters::GetPulseGain(UShort_t detector, Char_t ring, 
574                                UShort_t sector, UShort_t strip) const
575 {
576   // Returns the pulser calibrated gain for strip # strip in sector #
577   // sector or ring id ring of detector # detector. 
578   // 
579   // For simulation, this is normally set to 
580   // 
581   //       VA1_MIP_Range 
582   //    ------------------ * MIP_Energy_Loss
583   //    ALTRO_channel_size
584   // 
585   if (!fPulseGain) { 
586     if (fFixedPulseGain <= 0)
587       fFixedPulseGain = fVA1MipRange * GetEdepMip() / fAltroChannelSize;
588     return fFixedPulseGain;
589   }  
590   AliFMDDebug(50, ("pulse gain for FMD%d%c[%2d,%3d]=%f",
591                     detector, ring, sector, strip,
592                     fPulseGain->Value(detector, ring, sector, strip)));
593   return fPulseGain->Value(detector, ring, sector, strip);
594 }
595
596 //__________________________________________________________________
597 Bool_t
598 AliFMDParameters::IsDead(UShort_t detector, Char_t ring, 
599                          UShort_t sector, UShort_t strip) const
600 {
601   // Check if the channel is dead 
602   if (!fDeadMap) return kFALSE;
603   AliFMDDebug(50, ("Dead for FMD%d%c[%2d,%3d]=%s",
604                     detector, ring, sector, strip,
605                     fDeadMap->operator()(detector, ring, sector, strip) ? 
606                     "no" : "yes"));
607   return fDeadMap->operator()(detector, ring, sector, strip);
608 }
609
610 //__________________________________________________________________
611 UShort_t
612 AliFMDParameters::GetZeroSuppression(UShort_t detector, Char_t ring, 
613                                      UShort_t sector, UShort_t strip) const
614 {
615   // Get zero suppression threshold 
616   if (!fZeroSuppression) return fFixedZeroSuppression;
617   // Need to map strip to ALTRO chip. 
618   AliFMDDebug(50, ("zero sup. for FMD%d%c[%2d,%3d]=%f",
619                     detector, ring, sector, strip,
620                     fZeroSuppression->operator()(detector, ring, 
621                                                  sector, strip)));
622   return fZeroSuppression->operator()(detector, ring, sector, strip/128);
623 }
624
625 //__________________________________________________________________
626 UShort_t
627 AliFMDParameters::GetSampleRate(UShort_t det, Char_t ring, UShort_t sector, 
628                                 UShort_t str) const
629 {
630   // Get sampl rate 
631   if (!fSampleRate) return fFixedSampleRate;
632   // Need to map sector to digitizier card. 
633   UInt_t ret = fSampleRate->Rate(det, ring, sector, str);
634   AliFMDDebug(50, ("Sample rate for FMD%d%c[%2d,%3d]=%d", 
635                     det, ring, sector, str, ret));
636   return ret;
637 }
638
639 //__________________________________________________________________
640 UShort_t
641 AliFMDParameters::GetMinStrip(UShort_t det, Char_t ring, UShort_t sector, 
642                               UShort_t str) const
643 {
644   // Get strip range read out 
645   if (!fStripRange) return fFixedMinStrip;
646   // Need to map sector to digitizier card. 
647   UInt_t ret = fStripRange->Min(det, ring, sector, str);
648   AliFMDDebug(50, ("Min strip # for FMD%d%c[%2d,%3d]=%d", 
649                     det, ring, sector, str, ret));
650   return ret;
651 }
652
653 //__________________________________________________________________
654 UShort_t
655 AliFMDParameters::GetMaxStrip(UShort_t det, Char_t ring, UShort_t sector, 
656                               UShort_t str) const
657 {
658   // Get strip range read out 
659   if (!fStripRange) return fFixedMaxStrip;
660   // Need to map sector to digitizier card. 
661   UInt_t ret = fStripRange->Max(det, ring, sector, str);
662   AliFMDDebug(50, ("Max strip # for FMD%d%c[%2d,%3d]=%d", 
663                     det, ring, sector, str, ret));
664   return ret;
665 }
666
667 //__________________________________________________________________
668 Float_t
669 AliFMDParameters::GetPedestal(UShort_t detector, Char_t ring, 
670                               UShort_t sector, UShort_t strip) const
671 {
672   // Get the pedesal 
673   if (!fPedestal) return fFixedPedestal;
674   AliFMDDebug(50, ("pedestal for FMD%d%c[%2d,%3d]=%f",
675                     detector, ring, sector, strip,
676                     fPedestal->Value(detector, ring, sector, strip)));
677   return fPedestal->Value(detector, ring, sector, strip);
678 }
679
680 //__________________________________________________________________
681 Float_t
682 AliFMDParameters::GetPedestalWidth(UShort_t detector, Char_t ring, 
683                                    UShort_t sector, UShort_t strip) const
684 {
685   // Get the pedesal 
686   if (!fPedestal) return fFixedPedestalWidth;
687   AliFMDDebug(50, ("pedetal width for FMD%d%c[%2d,%3d]=%f",
688                     detector, ring, sector, strip,
689                     fPedestal->Width(detector, ring, sector, strip)));
690   return fPedestal->Width(detector, ring, sector, strip);
691 }
692   
693 //__________________________________________________________________
694 AliFMDAltroMapping*
695 AliFMDParameters::GetAltroMap() const
696 {
697   // Get the hardware address to detector index map 
698   return fAltroMap;
699 }
700
701
702 //__________________________________________________________________
703 Bool_t
704 AliFMDParameters::Hardware2Detector(UInt_t ddl, UInt_t addr, UShort_t& det,
705                                     Char_t& ring, UShort_t& sec, 
706                                     UShort_t& str) const
707 {
708   // Map hardware address to detector index
709   if (!fAltroMap) return kFALSE;
710   return fAltroMap->Hardware2Detector(ddl, addr, det, ring, sec, str);
711 }
712
713 //__________________________________________________________________
714 Bool_t
715 AliFMDParameters::Detector2Hardware(UShort_t det, Char_t ring, UShort_t sec, 
716                                     UShort_t str, UInt_t& ddl, 
717                                     UInt_t& addr) const                       
718 {
719   // Map detector index to hardware address
720   if (!fAltroMap) return kFALSE;
721   return fAltroMap->Detector2Hardware(det, ring, sec, str, ddl, addr);
722 }
723
724
725 //__________________________________________________________________
726 Float_t
727 AliFMDParameters::GetEdepMip() const 
728
729   // Get energy deposited by a MIP in the silicon sensors
730   if (fEdepMip <= 0){
731     AliFMDGeometry* fmd = AliFMDGeometry::Instance();
732     fEdepMip = (fkSiDeDxMip 
733                 * fmd->GetRing('I')->GetSiThickness() 
734                 * fmd->GetSiDensity());
735   }
736   return fEdepMip;
737 }
738
739
740   
741   
742   
743 //____________________________________________________________________
744 //
745 // EOF
746 //