]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDParameters.cxx
5cf64f50dde21b69d1d047c3c64045c8692b1171
[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 // Eventually, this class will use the Conditions DB to get the
28 // various parameters, which code can then request from here.
29 //                                                       
30 #include "AliLog.h"                // ALILOG_H
31 #include "AliFMDParameters.h"      // ALIFMDPARAMETERS_H
32 #include "AliFMDGeometry.h"        // ALIFMDGEOMETRY_H
33 #include "AliFMDRing.h"            // ALIFMDRING_H
34 #include "AliFMDCalibGain.h"       // ALIFMDCALIBGAIN_H
35 #include "AliFMDCalibPedestal.h"   // ALIFMDCALIBPEDESTAL_H
36 #include "AliFMDCalibSampleRate.h" // ALIFMDCALIBSAMPLERATE_H
37 #include "AliFMDCalibStripRange.h" // ALIFMDCALIBSTRIPRANGE_H
38 #include "AliFMDAltroMapping.h"    // ALIFMDALTROMAPPING_H
39 #include <AliCDBManager.h>         // ALICDBMANAGER_H
40 #include <AliCDBEntry.h>           // ALICDBMANAGER_H
41 #include <Riostream.h>
42 #include <sstream>
43
44 //====================================================================
45 ClassImp(AliFMDParameters)
46 #if 0
47   ; // This is here to keep Emacs for indenting the next line
48 #endif
49
50 //____________________________________________________________________
51 AliFMDParameters* AliFMDParameters::fgInstance = 0;
52
53 //____________________________________________________________________
54 const char* AliFMDParameters::fgkPulseGain       = "FMD/Calib/PulseGain";
55 const char* AliFMDParameters::fgkPedestal        = "FMD/Calib/Pedestal";
56 const char* AliFMDParameters::fgkDead            = "FMD/Calib/Dead";
57 const char* AliFMDParameters::fgkSampleRate      = "FMD/Calib/SampleRate";
58 const char* AliFMDParameters::fgkAltroMap        = "FMD/Calib/AltroMap";
59 const char* AliFMDParameters::fgkZeroSuppression = "FMD/Calib/ZeroSuppression";
60 const char* AliFMDParameters::fgkStripRange      = "FMD/Calib/StripRange";
61
62
63 //____________________________________________________________________
64 AliFMDParameters* 
65 AliFMDParameters::Instance() 
66 {
67   // Get static instance 
68   if (!fgInstance) fgInstance = new AliFMDParameters;
69   return fgInstance;
70 }
71
72 //____________________________________________________________________
73 AliFMDParameters::AliFMDParameters() 
74   : fIsInit(kFALSE),
75     fkSiDeDxMip(1.664), 
76     fFixedPulseGain(0), 
77     fEdepMip(0),
78     fZeroSuppression(0), 
79     fSampleRate(0), 
80     fPedestal(0), 
81     fPulseGain(0), 
82     fDeadMap(0), 
83     fAltroMap(0), 
84     fStripRange(0)
85 {
86   // Default constructor 
87   SetVA1MipRange();
88   SetAltroChannelSize();
89   SetChannelsPerAltro();
90   SetZeroSuppression();
91   SetSampleRate();
92   SetPedestal();
93   SetPedestalWidth();
94   SetPedestalFactor();
95   SetThreshold();
96   SetStripRange();
97 }
98
99 //__________________________________________________________________
100 void
101 AliFMDParameters::Init()
102 {
103   // Initialize the parameters manager.  We need to get stuff from the
104   // CDB here. 
105   if (fIsInit) return;
106   InitPulseGain();
107   InitPedestal();
108   InitDeadMap();
109   InitSampleRate();
110   InitZeroSuppression();
111   InitAltroMap();
112   fIsInit = kTRUE;
113   
114 }
115
116 //__________________________________________________________________
117 void
118 AliFMDParameters::Print(Option_t* option) const
119 {
120   // Print information. 
121   // If option contains an 'A' then everything is printed. 
122   TString opt(option);
123   Bool_t showStrips = opt.Contains("a", TString::kIgnoreCase);
124   if (opt.Contains("fmd",TString::kIgnoreCase)) {
125     size_t   i   = opt.Index("fmd",TString::kIgnoreCase);
126     size_t   j   = opt.Index("]",TString::kIgnoreCase);
127     UShort_t det, sec, str;
128     Char_t ring, lbrack, rbrack, comma;
129     UInt_t ddl, addr;
130     Detector2Hardware(det, ring, sec, str, ddl, addr);
131     std::stringstream s(opt(i+4, j-i-3).Data());
132     s >> det >> ring >> lbrack >> sec >> comma >> str >> rbrack;
133     std::cout 
134       << "     Strip    |     Pedestal      |    Gain    | ZS thr. | Address\n"
135       << "--------------+-------------------+------------+---------+---------" 
136       << "\nFMD" << det << ring << "[" << std::setw(2) << sec << "," 
137       << std::setw(3) << str << "] | " 
138       << std::setw(7) << GetPedestal(det, ring, sec, str) 
139       << "+/-" << std::setw(7) 
140       << GetPedestalWidth(det, ring, sec, str) 
141       << " | " << std::setw(10) 
142       << GetPulseGain(det, ring, sec, str) 
143       << " | " << std::setw(7) 
144       << GetZeroSuppression(det, ring, sec, str) 
145       << " | 0x" << std::hex << std::setw(4) 
146       << std::setfill('0') << ddl << ",0x" << std::setw(3) 
147       << addr << std::dec << std::setfill(' ') << std::endl;
148     return;
149   }
150   for (UShort_t det=1 ; det <= 3; det++) {
151     std::cout << "FMD" << det << std::endl;
152     Char_t rings[] = { 'I', (det == 1 ? '\0' : 'O'), '\0' };
153     for (Char_t* ring = rings; *ring != '\0'; ring++) {
154       std::cout << " Ring " << *ring << std::endl;
155       UShort_t nSec = ( *ring == 'I' ? 20  :  40 );
156       UShort_t nStr = ( *ring == 'I' ? 512 : 256 );
157       for (UShort_t sec = 0; sec < nSec; sec++) {
158         UShort_t min  = GetMinStrip(det, *ring, sec, 0);
159         UShort_t max  = GetMaxStrip(det, *ring, sec, 0);
160         UShort_t rate = GetSampleRate(det, *ring, sec, 0);
161         std::cout << "  Sector " << std::setw(2) << sec 
162                   << "  Strip range: " << std::setw(3) << min << "," 
163                   << std::setw(3) << max << "  Rate: " << std::setw(2) 
164                   << rate << std::endl;
165         if (!showStrips) continue;
166         std::cout 
167           << "  Strip |     Pedestal      |    Gain    | ZS thr. | Address\n" 
168           << "--------+-------------------+------------+---------+---------" 
169           << std::endl;
170         for (UShort_t str = 0; str < nStr; str++) {
171           std::cout << "    " << std::setw(3) << str << " | ";
172           if (IsDead(det, *ring, sec, str)) {
173             std::cout << "dead" << std::endl;
174             continue;
175           }
176           UInt_t ddl, addr;
177           Detector2Hardware(det, *ring, sec, str, ddl, addr);
178           std::cout << std::setw(7) << GetPedestal(det, *ring, sec, str) 
179                     << "+/-" << std::setw(7) 
180                     << GetPedestalWidth(det, *ring, sec, str) 
181                     << " | " << std::setw(10) 
182                     << GetPulseGain(det, *ring, sec, str) 
183                     << " | " << std::setw(5) 
184                     << GetZeroSuppression(det, *ring, sec, str) 
185                     << " | 0x" << std::hex << std::setw(4) 
186                     << std::setfill('0') << ddl << ",0x" << std::setw(3) 
187                     << addr << std::dec << std::setfill(' ') << std::endl;
188         }
189       }
190     }
191   }
192 }
193
194 //__________________________________________________________________
195 void
196 AliFMDParameters::SetStripRange(UShort_t min, UShort_t max) 
197 {
198   // Set fixed strip range 
199   fFixedMinStrip = min;
200   fFixedMaxStrip = max;
201 }
202
203 //__________________________________________________________________
204 void
205 AliFMDParameters::InitPulseGain()
206 {
207   // Get pulse gain from CDB or used fixed 
208   AliCDBManager* cdb      = AliCDBManager::Instance();
209   AliCDBEntry*   gain     = cdb->Get(fgkPulseGain);
210   if (!gain) {
211     AliWarning(Form("No %s found in CDB, perhaps you need to "
212                     "use AliFMDCalibFaker?", fgkPulseGain));
213     return;
214   }
215   
216   AliDebug(1, Form("Got gain from CDB"));
217   fPulseGain = dynamic_cast<AliFMDCalibGain*>(gain->GetObject());
218   if (!fPulseGain) AliWarning("Invalid pulser gain object from CDB");
219 }
220 //__________________________________________________________________
221 void
222 AliFMDParameters::InitPedestal()
223 {
224   // Initialize the pedestals from CDB 
225   AliCDBManager* cdb      = AliCDBManager::Instance();
226   AliCDBEntry*   pedestal = cdb->Get(fgkPedestal);
227   if (!pedestal) {
228     AliWarning(Form("No %s found in CDB, perhaps you need to "
229                     "use AliFMDCalibFaker?", fgkPedestal));
230     return;
231   }
232   AliDebug(1, Form("Got pedestal from CDB"));
233   fPedestal = dynamic_cast<AliFMDCalibPedestal*>(pedestal->GetObject());
234   if (!fPedestal) AliWarning("Invalid pedestal object from CDB");
235 }
236
237 //__________________________________________________________________
238 void
239 AliFMDParameters::InitDeadMap()
240 {
241   // Get Dead-channel-map from CDB 
242   AliCDBManager* cdb      = AliCDBManager::Instance();
243   AliCDBEntry*   deadMap  = cdb->Get(fgkDead);
244   if (!deadMap) {
245     AliWarning(Form("No %s found in CDB, perhaps you need to "
246                     "use AliFMDCalibFaker?", fgkDead));
247     return;
248   }
249   AliDebug(1, Form("Got dead map from CDB"));
250   fDeadMap = dynamic_cast<AliFMDCalibDeadMap*>(deadMap->GetObject());
251   if (!fDeadMap) AliWarning("Invalid dead map object from CDB");
252 }
253
254 //__________________________________________________________________
255 void
256 AliFMDParameters::InitZeroSuppression()
257 {
258   // Get 0-suppression from CDB 
259   AliCDBManager* cdb      = AliCDBManager::Instance();
260   AliCDBEntry*   zeroSup  = cdb->Get(fgkZeroSuppression);
261   if (!zeroSup) {
262     AliWarning(Form("No %s found in CDB, perhaps you need to "
263                     "use AliFMDCalibFaker?", fgkZeroSuppression));
264     return;
265   }
266   AliDebug(1, Form("Got zero suppression from CDB"));
267   fZeroSuppression = 
268     dynamic_cast<AliFMDCalibZeroSuppression*>(zeroSup->GetObject());
269   if (!fZeroSuppression)AliWarning("Invalid zero suppression object from CDB");
270 }
271
272 //__________________________________________________________________
273 void
274 AliFMDParameters::InitSampleRate()
275 {
276   // get Sample rate from CDB
277   AliCDBManager* cdb      = AliCDBManager::Instance();
278   AliCDBEntry*   sampRat  = cdb->Get(fgkSampleRate);
279   if (!sampRat) {
280     AliWarning(Form("No %s found in CDB, perhaps you need to "
281                     "use AliFMDCalibFaker?", fgkSampleRate));
282     return;
283   }
284   AliDebug(1, Form("Got zero suppression from CDB"));
285   fSampleRate = dynamic_cast<AliFMDCalibSampleRate*>(sampRat->GetObject());
286   if (!fSampleRate) AliWarning("Invalid zero suppression object from CDB");
287 }
288
289 //__________________________________________________________________
290 void
291 AliFMDParameters::InitAltroMap()
292 {
293   // Get hardware mapping from CDB
294   AliCDBManager* cdb      = AliCDBManager::Instance();
295   AliCDBEntry*   hwMap    = cdb->Get(fgkAltroMap);       
296   if (!hwMap) {
297     AliWarning(Form("No %s found in CDB, perhaps you need to "
298                     "use AliFMDCalibFaker?", fgkAltroMap));
299     fAltroMap = new AliFMDAltroMapping;
300     return;
301   }
302   AliDebug(1, Form("Got ALTRO map from CDB"));
303   fAltroMap = dynamic_cast<AliFMDAltroMapping*>(hwMap->GetObject());
304   if (!fAltroMap) {
305     AliWarning("Invalid ALTRO map object from CDB");
306     fAltroMap = new AliFMDAltroMapping;
307   }
308 }
309
310 //__________________________________________________________________
311 void
312 AliFMDParameters::InitStripRange()
313 {
314   // Get strips read-out from CDB
315   AliCDBManager* cdb      = AliCDBManager::Instance();
316   AliCDBEntry*   range    = cdb->Get(fgkStripRange);
317   if (!range) {
318     AliWarning(Form("No %s found in CDB, perhaps you need to "
319                     "use AliFMDCalibFaker?", fgkStripRange));
320     return;
321   }
322   AliDebug(1, Form("Got strip range from CDB"));
323   fStripRange = dynamic_cast<AliFMDCalibStripRange*>(range->GetObject());
324   if (!fStripRange) AliWarning("Invalid strip range object from CDB");
325 }
326
327
328 //__________________________________________________________________
329 Float_t
330 AliFMDParameters::GetThreshold() const
331 {
332   // Get threshold from CDB
333   if (!fPulseGain) return fFixedThreshold;
334   return fPulseGain->Threshold();
335 }
336
337 //__________________________________________________________________
338 Float_t
339 AliFMDParameters::GetPulseGain(UShort_t detector, Char_t ring, 
340                                UShort_t sector, UShort_t strip) const
341 {
342   // Returns the pulser calibrated gain for strip # strip in sector #
343   // sector or ring id ring of detector # detector. 
344   // 
345   // For simulation, this is normally set to 
346   // 
347   //       VA1_MIP_Range 
348   //    ------------------ * MIP_Energy_Loss
349   //    ALTRO_channel_size
350   // 
351   if (!fPulseGain) { 
352     if (fFixedPulseGain <= 0)
353       fFixedPulseGain = fVA1MipRange * GetEdepMip() / fAltroChannelSize;
354     return fFixedPulseGain;
355   }  
356   AliDebug(50, Form("pulse gain for FMD%d%c[%2d,%3d]=%f",
357                     detector, ring, sector, strip,
358                     fPulseGain->Value(detector, ring, sector, strip)));
359   return fPulseGain->Value(detector, ring, sector, strip);
360 }
361
362 //__________________________________________________________________
363 Bool_t
364 AliFMDParameters::IsDead(UShort_t detector, Char_t ring, 
365                          UShort_t sector, UShort_t strip) const
366 {
367   // Check if the channel is dead 
368   if (!fDeadMap) return kFALSE;
369   AliDebug(50, Form("Dead for FMD%d%c[%2d,%3d]=%s",
370                     detector, ring, sector, strip,
371                     fDeadMap->operator()(detector, ring, sector, strip) ? 
372                     "no" : "yes"));
373   return fDeadMap->operator()(detector, ring, sector, strip);
374 }
375
376 //__________________________________________________________________
377 UShort_t
378 AliFMDParameters::GetZeroSuppression(UShort_t detector, Char_t ring, 
379                                      UShort_t sector, UShort_t strip) const
380 {
381   // Get zero suppression threshold 
382   if (!fZeroSuppression) return fFixedZeroSuppression;
383   // Need to map strip to ALTRO chip. 
384   AliDebug(50, Form("zero sup. for FMD%d%c[%2d,%3d]=%f",
385                     detector, ring, sector, strip,
386                     fZeroSuppression->operator()(detector, ring, 
387                                                  sector, strip)));
388   return fZeroSuppression->operator()(detector, ring, sector, strip/128);
389 }
390
391 //__________________________________________________________________
392 UShort_t
393 AliFMDParameters::GetSampleRate(UShort_t det, Char_t ring, UShort_t sector, 
394                                 UShort_t str) const
395 {
396   // Get sampl rate 
397   if (!fSampleRate) return fFixedSampleRate;
398   // Need to map sector to digitizier card. 
399   UInt_t ret = fSampleRate->Rate(det, ring, sector, str);
400   AliDebug(50, Form("Sample rate for FMD%d%c[%2d,%3d]=%d", 
401                     det, ring, sector, str, ret));
402   return ret;
403 }
404
405 //__________________________________________________________________
406 UShort_t
407 AliFMDParameters::GetMinStrip(UShort_t det, Char_t ring, UShort_t sector, 
408                               UShort_t str) const
409 {
410   // Get strip range read out 
411   if (!fStripRange) return fFixedMinStrip;
412   // Need to map sector to digitizier card. 
413   UInt_t ret = fStripRange->Min(det, ring, sector, str);
414   AliDebug(50, Form("Min strip # for FMD%d%c[%2d,%3d]=%d", 
415                     det, ring, sector, str, ret));
416   return ret;
417 }
418
419 //__________________________________________________________________
420 UShort_t
421 AliFMDParameters::GetMaxStrip(UShort_t det, Char_t ring, UShort_t sector, 
422                               UShort_t str) const
423 {
424   // Get strip range read out 
425   if (!fStripRange) return fFixedMaxStrip;
426   // Need to map sector to digitizier card. 
427   UInt_t ret = fStripRange->Max(det, ring, sector, str);
428   AliDebug(50, Form("Max strip # for FMD%d%c[%2d,%3d]=%d", 
429                     det, ring, sector, str, ret));
430   return ret;
431 }
432
433 //__________________________________________________________________
434 Float_t
435 AliFMDParameters::GetPedestal(UShort_t detector, Char_t ring, 
436                               UShort_t sector, UShort_t strip) const
437 {
438   // Get the pedesal 
439   if (!fPedestal) return fFixedPedestal;
440   AliDebug(50, Form("pedestal for FMD%d%c[%2d,%3d]=%f",
441                     detector, ring, sector, strip,
442                     fPedestal->Value(detector, ring, sector, strip)));
443   return fPedestal->Value(detector, ring, sector, strip);
444 }
445
446 //__________________________________________________________________
447 Float_t
448 AliFMDParameters::GetPedestalWidth(UShort_t detector, Char_t ring, 
449                                    UShort_t sector, UShort_t strip) const
450 {
451   // Get the pedesal 
452   if (!fPedestal) return fFixedPedestalWidth;
453   AliDebug(50, Form("pedetal width for FMD%d%c[%2d,%3d]=%f",
454                     detector, ring, sector, strip,
455                     fPedestal->Width(detector, ring, sector, strip)));
456   return fPedestal->Width(detector, ring, sector, strip);
457 }
458   
459 //__________________________________________________________________
460 AliFMDAltroMapping*
461 AliFMDParameters::GetAltroMap() const
462 {
463   // Get the hardware address to detector index map 
464   return fAltroMap;
465 }
466
467
468 //__________________________________________________________________
469 Bool_t
470 AliFMDParameters::Hardware2Detector(UInt_t ddl, UInt_t addr, UShort_t& det,
471                                     Char_t& ring, UShort_t& sec, 
472                                     UShort_t& str) const
473 {
474   // Map hardware address to detector index
475   if (!fAltroMap) return kFALSE;
476   return fAltroMap->Hardware2Detector(ddl, addr, det, ring, sec, str);
477 }
478
479 //__________________________________________________________________
480 Bool_t
481 AliFMDParameters::Detector2Hardware(UShort_t det, Char_t ring, UShort_t sec, 
482                                     UShort_t str, UInt_t& ddl, 
483                                     UInt_t& addr) const                       
484 {
485   // Map detector index to hardware address
486   if (!fAltroMap) return kFALSE;
487   return fAltroMap->Detector2Hardware(det, ring, sec, str, ddl, addr);
488 }
489
490
491 //__________________________________________________________________
492 Float_t
493 AliFMDParameters::GetEdepMip() const 
494
495   // Get energy deposited by a MIP in the silicon sensors
496   if (fEdepMip <= 0){
497     AliFMDGeometry* fmd = AliFMDGeometry::Instance();
498     fEdepMip = (fkSiDeDxMip 
499                 * fmd->GetRing('I')->GetSiThickness() 
500                 * fmd->GetSiDensity());
501   }
502   return fEdepMip;
503 }
504
505
506   
507   
508   
509 //____________________________________________________________________
510 //
511 // EOF
512 //