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