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