]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONPadStatusMaker.cxx
Adding 2013 periods
[u/mrichter/AliRoot.git] / MUON / AliMUONPadStatusMaker.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
16 // $Id$
17
18 //-----------------------------------------------------------------------------
19 /// \class AliMUONPadStatusMaker
20 ///
21 /// Make a 2DStore of pad statuses, using different sources of information,
22 /// like pedestal values, gain values, and HV values.
23 ///
24 /// \author Laurent Aphecetche
25 //-----------------------------------------------------------------------------
26
27 #include "AliMUONPadStatusMaker.h"
28
29 #include "AliMUON2DMap.h"
30 #include "AliMUON2DStoreValidator.h"
31 #include "AliMUONCalibParamNI.h"
32 #include "AliMUONCalibrationData.h"
33 #include "AliMUONLogger.h"
34 #include "AliMUONRecoParam.h"
35 #include "AliMUONStringIntMap.h"
36 #include "AliMUONTrackerData.h"
37 #include "AliMUONVCalibParam.h"
38
39 #include "AliMpArea.h"
40 #include "AliMpArrayI.h"
41 #include "AliMpCDB.h"
42 #include "AliMpConstants.h"
43 #include "AliMpDDLStore.h"
44 #include "AliMpDEManager.h"
45 #include "AliMpDetElement.h"
46 #include "AliMpDCSNamer.h"
47 #include "AliMpManuIterator.h"
48 #include "AliMpManuUID.h"
49
50 #include "AliCDBEntry.h"
51 #include "AliCDBManager.h"
52 #include "AliCodeTimer.h"
53 #include "AliDCSValue.h"
54 #include "AliLog.h"
55
56 #include <Riostream.h>
57 #include <TArrayI.h>
58 #include <TExMap.h>
59 #include <TFile.h>
60 #include <TKey.h>
61 #include <TMap.h>
62 #include <TROOT.h>
63 #include <TString.h>
64 #include <TSystem.h>
65
66 using std::cout;
67 using std::endl;
68 /// \cond CLASSIMP
69 ClassImp(AliMUONPadStatusMaker)
70 /// \endcond
71
72 //_____________________________________________________________________________
73 AliMUONPadStatusMaker::AliMUONPadStatusMaker(const AliMUONCalibrationData& calibData)
74 : fkCalibrationData(calibData),
75 fGainA1Limits(0,1E30),
76 fGainA2Limits(-1E-30,1E30),
77 fGainThresLimits(0,4095),
78 fPedMeanLimits(0,4095),
79 fPedSigmaLimits(0,4095),
80 fManuOccupancyLimits(0,1.0),
81 fBuspatchOccupancyLimits(0,1.0),
82 fDEOccupancyLimits(0,1.0),
83 fStatus(new AliMUON2DMap(true)),
84 fHV(0x0),
85 fPedestals(calibData.Pedestals()),
86 fGains(calibData.Gains()),
87 fTrackerData(0x0)
88 {
89   /// ctor
90   if ( calibData.OccupancyMap() )
91   {
92     /// create a tracker data from the occupancy map
93     fTrackerData = new AliMUONTrackerData("OCC","OCC",*(calibData.OccupancyMap()));
94   }     
95   if ( calibData.HV() )
96   {
97     /// Only create the fHV internal store if there are some HV values available
98     fHV = new TExMap;
99   }
100   
101   SetHVLimit(-1,0.0);
102 }
103
104 //_____________________________________________________________________________
105 AliMUONPadStatusMaker::~AliMUONPadStatusMaker()
106 {
107   /// dtor.
108  
109   delete fStatus;
110   delete fHV;
111   delete fTrackerData;
112 }
113
114 //_____________________________________________________________________________
115 TString
116 AliMUONPadStatusMaker::AsString(Int_t status)
117 {
118   /// return a human readable version of the integer status
119   
120   if ( status == 0 ) 
121   {
122     return "Brave New World";
123   }
124   
125   Int_t pedStatus;
126   Int_t gainStatus;
127   Int_t hvStatus;
128   Int_t occStatus;
129   
130   DecodeStatus(status,pedStatus,hvStatus,gainStatus,occStatus);
131   
132   TString s;
133   
134   if ( pedStatus & kPedMeanZero ) s += "& Ped Mean is Zero ";
135   if ( pedStatus & kPedMeanTooLow ) s += "& Ped Mean Too Low ";
136   if ( pedStatus & kPedMeanTooHigh ) s += "& Ped Mean Too High ";
137   if ( pedStatus & kPedSigmaTooLow ) s += "& Ped Sigma Too Low ";
138   if ( pedStatus & kPedSigmaTooHigh ) s += "& Ped Sigma Too High ";
139   if ( pedStatus & kPedMissing ) s += "& Ped is missing ";
140   
141         if ( gainStatus & kGainA1TooLow ) s+="& Gain A1 is Too Low ";
142         if ( gainStatus & kGainA1TooHigh ) s+="& Gain A1 is Too High ";
143         if ( gainStatus & kGainA2TooLow ) s+="& Gain A2 is Too Low ";
144         if ( gainStatus & kGainA2TooHigh ) s+="& Gain A2 is Too High ";
145         if ( gainStatus & kGainThresTooLow ) s+="& Gain Thres is Too Low ";
146         if ( gainStatus & kGainThresTooHigh ) s+="& Gain Thres is Too High ";
147         if ( gainStatus & kGainMissing ) s+="& Gain is missing ";
148         
149         if ( hvStatus & kHVError ) s+="& HV is on error ";
150         if ( hvStatus & kHVTooLow ) s+="& HV is Too Low ";
151         if ( hvStatus & kHVTooHigh ) s+="& HV is Too High ";
152         if ( hvStatus & kHVChannelOFF ) s+="& HV has channel OFF ";
153         if ( hvStatus & kHVSwitchOFF ) s+="& HV has switch OFF ";
154         if ( hvStatus & kHVMissing ) s+="& HV is missing ";
155
156   if ( occStatus & kManuOccupancyTooHigh ) s+="& manu occupancy too high ";
157   if ( occStatus & kManuOccupancyTooLow ) s+="& manu occupancy too low ";
158   if ( occStatus & kBusPatchOccupancyTooHigh ) s+="& bus patch occupancy too high ";
159   if ( occStatus & kBusPatchOccupancyTooLow ) s+="& bus patch occupancy too low ";
160   if ( occStatus & kDEOccupancyTooHigh ) s+="& DE occupancy too high ";
161   if ( occStatus & kDEOccupancyTooLow ) s+="& DE occupancy too low ";
162   
163   if ( s[0] == '&' ) s[0] = ' ';
164   
165   return s;
166 }
167
168 //_____________________________________________________________________________
169 TString
170 AliMUONPadStatusMaker::AsCondition(Int_t mask)
171 {
172   /// return a human readable version of the mask's equivalent condition
173   
174   TString s(AsString(mask));
175   
176   s.ReplaceAll("&","|");
177   
178   return s;
179 }
180
181 //_____________________________________________________________________________
182 Int_t
183 AliMUONPadStatusMaker::BuildStatus(Int_t pedStatus, 
184                                    Int_t hvStatus, 
185                                    Int_t gainStatus,
186                                    Int_t occStatus)
187 {
188   /// Build a complete status from specific parts (ped,hv,gain)
189   
190   return ( hvStatus & 0xFF ) | ( ( pedStatus & 0xFF ) << 8 ) | 
191   ( ( gainStatus & 0xFF ) << 16 ) |
192   ( ( occStatus & 0xFF ) << 24 ) ;
193 }
194
195 //_____________________________________________________________________________
196 void
197 AliMUONPadStatusMaker::DecodeStatus(Int_t status, 
198                                     Int_t& pedStatus, 
199                                     Int_t& hvStatus, 
200                                     Int_t& gainStatus,
201                                     Int_t& occStatus)
202 {
203   /// Decode complete status into specific parts (ped,hv,gain)
204   
205   occStatus = ( status & 0xFF000000 ) >> 24;
206   gainStatus = ( status & 0xFF0000 ) >> 16;
207   pedStatus = ( status & 0xFF00 ) >> 8;
208   hvStatus = (status & 0xFF);
209 }
210
211 //_____________________________________________________________________________
212 Bool_t 
213 AliMUONPadStatusMaker::HVSt12Status(Int_t detElemId, Int_t sector,
214                                     Bool_t& hvChannelTooLow,
215                                     Bool_t& hvChannelTooHigh,
216                                     Bool_t& hvChannelON) const
217 {
218   /// Get HV status for one HV sector of St12
219   
220   /// For a given PCB in a given DE, get the HV status (both the channel
221   /// and the switch).
222   /// Returns false if hv switch changed during the run.
223   
224   AliCodeTimerAuto("",0)
225   
226   if (!fHV) return kFALSE;
227   
228   Bool_t error = kFALSE;
229   hvChannelTooLow = kFALSE;
230   hvChannelTooHigh = kFALSE;
231   hvChannelON = kTRUE;
232
233   Int_t chamberId = AliMpDEManager::GetChamberId(detElemId);
234   
235   AliMpDCSNamer hvNamer("TRACKER");
236   
237   TString hvChannel(hvNamer.DCSAliasName(detElemId,sector));
238   
239   TMap* hvMap = fkCalibrationData.HV();
240   TPair* hvPair = static_cast<TPair*>(hvMap->FindObject(hvChannel.Data()));
241   if (!hvPair)
242   {
243     AliError(Form("Did not find expected alias (%s) for DE %d",
244                   hvChannel.Data(),detElemId));  
245     error = kTRUE;
246   }
247   else
248   {
249     TObjArray* values = static_cast<TObjArray*>(hvPair->Value());
250     if (!values)
251     {
252       AliError(Form("Could not get values for alias %s",hvChannel.Data()));
253       error = kTRUE;
254     }
255     else
256     {
257       // find out min value, and makes a cut
258       Float_t hvMin(1E9);
259       TIter next(values);
260       AliDCSValue* val;
261       
262       while ( ( val = static_cast<AliDCSValue*>(next()) ) )
263       {
264         Float_t hv = val->GetFloat();
265         hvMin = TMath::Min(hv,hvMin);
266       }
267       
268       float lowThreshold = HVLimit(chamberId);
269             
270       if ( hvMin < lowThreshold ) hvChannelTooLow = kTRUE;
271       if ( hvMin < hvNamer.TrackerHVOFF() ) hvChannelON = kFALSE;
272     }
273   }
274   
275   return error;
276 }
277
278 //_____________________________________________________________________________
279 Float_t
280 AliMUONPadStatusMaker::SwitchValue(const TObjArray& dcsArray)
281 {
282   /// Loop over the dcs value for a single switch to decide whether
283   /// we should consider it on or off
284   
285   // we'll count the number of ON/OFF for this pad, to insure
286   // consistency (i.e. if status changed during the run, we should
287   // at least notify this fact ;-) and hope it's not the norm)
288   Int_t nTrue(0);
289   Int_t nFalse(0);
290   TIter next(&dcsArray);
291   AliDCSValue* val;
292   
293   while ( ( val = static_cast<AliDCSValue*>(next()) ) )
294   {
295     if ( val->GetBool() )
296     {
297       ++nTrue;
298     }
299     else
300     {
301       ++nFalse;
302     }
303   }
304   
305   if ( (nTrue>0 && nFalse>0) )
306   {
307     // change of state during the run, consider it off
308     return 0.0;
309   }
310   
311   if ( nFalse ) 
312   {
313     /// switch = FALSE means the HV was flowding up to the PCB.
314     /// i.e. switch = FALSE = ON
315     return 1.0;    
316   }
317   
318   return 0.0;
319 }
320
321 //_____________________________________________________________________________
322 Bool_t 
323 AliMUONPadStatusMaker::HVSt345Status(Int_t detElemId, Int_t pcbIndex,
324                                      Bool_t& hvChannelTooLow,
325                                      Bool_t& hvChannelTooHigh,
326                                      Bool_t& hvChannelON,
327                                      Bool_t& hvSwitchON) const
328 {
329   /// For a given PCB in a given DE, get the HV status (both the channel
330   /// and the switch).
331   /// Returns false if something goes wrong (in particular if 
332   /// hv switch changed during the run).
333   
334   AliCodeTimerAuto("",0)
335   
336   if (!fHV) return kFALSE;
337   
338   Bool_t error = kFALSE;
339   hvChannelTooLow = kFALSE;
340   hvChannelTooHigh = kFALSE;
341   hvSwitchON = kTRUE;
342   hvChannelON = kTRUE;
343   
344   AliMpDCSNamer hvNamer("TRACKER");
345   
346   Int_t chamberId = AliMpDEManager::GetChamberId(detElemId);
347   
348   TString hvChannel(hvNamer.DCSAliasName(detElemId));
349   
350   TMap* hvMap = fkCalibrationData.HV();
351   
352   TPair* hvPair = static_cast<TPair*>(hvMap->FindObject(hvChannel.Data()));
353   if (!hvPair)
354   {
355     AliError(Form("Did not find expected alias (%s) for DE %d",
356                   hvChannel.Data(),detElemId));  
357     error = kTRUE;
358   }
359   else
360   {
361     TObjArray* values = static_cast<TObjArray*>(hvPair->Value());
362     if (!values)
363     {
364       AliError(Form("Could not get values for alias %s",hvChannel.Data()));
365       error = kTRUE;
366     }
367     else
368     {
369       // find out min value, and makes a cut
370       Float_t hvMin(1E9);
371       TIter next(values);
372       AliDCSValue* val;
373       
374       while ( ( val = static_cast<AliDCSValue*>(next()) ) )
375       {
376         Float_t hv = val->GetFloat();
377         hvMin = TMath::Min(hv,hvMin);
378       }
379
380       float lowThreshold = HVLimit(chamberId);
381
382       if ( hvMin < lowThreshold ) hvChannelTooLow = kTRUE;
383       if ( hvMin < hvNamer.TrackerHVOFF() ) hvChannelON = kFALSE;
384     }
385   }
386   
387   TString hvSwitch(hvNamer.DCSSwitchAliasName(detElemId,pcbIndex));
388   TPair* switchPair = static_cast<TPair*>(hvMap->FindObject(hvSwitch.Data()));
389   if (!switchPair)
390   {
391     AliError(Form("Did not find expected alias (%s) for DE %d PCB %d",
392                   hvSwitch.Data(),detElemId,pcbIndex));
393     error = kTRUE;
394   }
395   else
396   {
397     TObjArray* values = static_cast<TObjArray*>(switchPair->Value());
398     if (!values)
399     {    
400       AliError(Form("Could not get values for alias %s",hvSwitch.Data()));
401       error = kTRUE;
402     }
403     else
404     {
405       Float_t sv = SwitchValue(*values);
406       if ( sv < 0.99 ) hvSwitchON = kFALSE;
407     }
408   }
409   return error;
410 }
411
412 //_____________________________________________________________________________
413 Int_t
414 AliMUONPadStatusMaker::HVStatus(Int_t detElemId, Int_t manuId) const
415 {
416   /// Get HV status of one manu
417   
418   AliCodeTimerAuto("",0)
419   
420   if ( !fHV ) return kMissing;
421
422   Long_t lint = fHV->GetValue(AliMpManuUID::BuildUniqueID(detElemId,manuId));
423   
424   if ( lint ) 
425   {
426     return (Int_t)(lint - 1);
427   }
428
429   Int_t status(0);
430   
431   AliMpDCSNamer hvNamer("TRACKER");
432   
433   switch ( AliMpDEManager::GetStationType(detElemId) )
434   {
435     case AliMp::kStation12:
436     {
437       int sector = hvNamer.ManuId2Sector(detElemId,manuId);
438       if ( sector >= 0 ) 
439       {
440         Bool_t hvChannelTooLow, hvChannelTooHigh, hvChannelON;
441         Bool_t error = HVSt12Status(detElemId,sector,
442                                     hvChannelTooLow,
443                                     hvChannelTooHigh,
444                                     hvChannelON);
445         if ( error ) status |= kHVError;
446         if ( hvChannelTooLow ) status |= kHVTooLow;
447         if ( hvChannelTooHigh ) status |= kHVTooHigh; 
448         if ( !hvChannelON ) status |= kHVChannelOFF;
449         // assign this status to all the other manus handled by the same HV channel
450         SetHVStatus(detElemId,sector,status);
451       }
452     }
453       break;
454     case AliMp::kStation345:
455     {
456       int pcbIndex = hvNamer.ManuId2PCBIndex(detElemId,manuId);
457       if ( pcbIndex >= 0 ) 
458       {
459         Bool_t hvChannelTooLow, hvChannelTooHigh, hvChannelON,hvSwitchON;
460         Bool_t error = HVSt345Status(detElemId,pcbIndex,
461                                      hvChannelTooLow,hvChannelTooHigh,
462                                      hvChannelON,hvSwitchON);
463         if ( error ) status |= kHVError;
464         if ( hvChannelTooLow ) status |= kHVTooLow;
465         if ( hvChannelTooHigh ) status |= kHVTooHigh; 
466         if ( !hvSwitchON ) status |= kHVSwitchOFF; 
467         if ( !hvChannelON) status |= kHVChannelOFF;
468         // assign this status to all the other manus handled by the same HV channel
469         SetHVStatus(detElemId,pcbIndex,status);
470       }
471     }
472       break;
473     default:
474       break;
475   }
476   
477   return status;
478 }
479
480 //_____________________________________________________________________________
481 AliMUONVCalibParam* 
482 AliMUONPadStatusMaker::Neighbours(Int_t detElemId, Int_t manuId) const
483 {
484   /// Get the neighbours parameters for a given manu
485   AliMUONVStore* neighbourStore = fkCalibrationData.Neighbours();
486   return static_cast<AliMUONVCalibParam*>(neighbourStore->FindObject(detElemId,manuId));
487 }
488
489 //_____________________________________________________________________________
490 AliMUONVStore* 
491 AliMUONPadStatusMaker::NeighboursStore() const
492 {
493   /// Return the store containing all the neighbours
494   return fkCalibrationData.Neighbours();
495 }
496
497 //_____________________________________________________________________________
498 AliMUONVCalibParam*
499 AliMUONPadStatusMaker::ComputeStatus(Int_t detElemId, Int_t manuId) const
500 {
501   /// Compute the status of a given manu, using all available information,
502   /// i.e. pedestals, gains, and HV
503   
504   AliMUONVCalibParam* param = new AliMUONCalibParamNI(1,AliMpConstants::ManuNofChannels(),detElemId,manuId,-1);
505   fStatus->Add(param);
506
507   AliMUONVCalibParam* pedestals = static_cast<AliMUONVCalibParam*>(fPedestals->FindObject(detElemId,manuId));
508
509   AliMUONVCalibParam* gains = static_cast<AliMUONVCalibParam*>(fGains->FindObject(detElemId,manuId));
510   
511   Int_t hvStatus = HVStatus(detElemId,manuId);
512
513   Int_t occStatus = OccupancyStatus(detElemId,manuId);
514   
515   for ( Int_t manuChannel = 0; manuChannel < param->Size(); ++manuChannel )
516   {
517     Int_t pedStatus(0);
518     
519     if (pedestals) 
520     {
521       Float_t pedMean = pedestals->ValueAsFloatFast(manuChannel,0);
522       Float_t pedSigma = pedestals->ValueAsFloatFast(manuChannel,1);
523       if ( pedMean < fPedMeanLimits.X() ) pedStatus |= kPedMeanTooLow;
524       else if ( pedMean > fPedMeanLimits.Y() ) pedStatus |= kPedMeanTooHigh;
525       if ( pedSigma < fPedSigmaLimits.X() ) pedStatus |= kPedSigmaTooLow;
526       else if ( pedSigma > fPedSigmaLimits.Y() ) pedStatus |= kPedSigmaTooHigh;
527       if ( pedMean == 0 ) pedStatus |= kPedMeanZero;
528     }
529     else
530     {
531       pedStatus = kPedMissing;
532     }
533     
534     Int_t gainStatus(0);
535   
536     if ( gains ) 
537     {
538       Float_t a0 = gains->ValueAsFloatFast(manuChannel,0);
539       Float_t a1 = gains->ValueAsFloatFast(manuChannel,1);
540       Float_t thres = gains->ValueAsFloatFast(manuChannel,2);
541   
542       if ( a0 < fGainA1Limits.X() ) gainStatus |= kGainA1TooLow;
543       else if ( a0 > fGainA1Limits.Y() ) gainStatus |= kGainA1TooHigh;
544       if ( a1 < fGainA2Limits.X() ) gainStatus |= kGainA2TooLow;
545       else if ( a1 > fGainA2Limits.Y() ) gainStatus |= kGainA2TooHigh;
546       if ( thres < fGainThresLimits.X() ) gainStatus |= kGainThresTooLow;
547       else if ( thres > fGainThresLimits.Y() ) gainStatus |= kGainThresTooHigh;
548     }
549     else
550     {
551       gainStatus = kGainMissing;
552     }
553     
554     Int_t status = BuildStatus(pedStatus,hvStatus,gainStatus,occStatus);
555       
556     param->SetValueAsIntFast(manuChannel,0,status);
557   }
558   
559   return param;
560 }
561
562 //_____________________________________________________________________________
563 Int_t 
564 AliMUONPadStatusMaker::OccupancyStatus(Int_t detElemId, Int_t manuId) const
565 {
566   /// Get the "other" status for a given manu
567  
568   Int_t rv(0);
569   
570   if ( fTrackerData ) 
571   {
572     const Int_t occIndex = 2;
573     
574     Double_t occ = fTrackerData->DetectionElement(detElemId,occIndex);
575     
576     if ( occ <= fDEOccupancyLimits.X() )
577     {
578       rv |= kDEOccupancyTooLow;
579     } 
580     else if ( occ > fDEOccupancyLimits.Y() )
581     {
582       rv |= kDEOccupancyTooHigh;
583     }
584     
585     Int_t busPatchId = AliMpDDLStore::Instance()->GetBusPatchId(detElemId,manuId);
586     
587     occ = fTrackerData->BusPatch(busPatchId,occIndex);
588
589     if ( occ <= fBuspatchOccupancyLimits.X() )
590     {
591       rv |= kBusPatchOccupancyTooLow;
592     } 
593     else if ( occ > fBuspatchOccupancyLimits.Y() )
594     {
595       rv |= kBusPatchOccupancyTooHigh;
596     }
597     
598     occ = fTrackerData->Manu(detElemId,manuId,occIndex);
599     
600     if ( occ <= fManuOccupancyLimits.X() )
601     {
602       rv |= kManuOccupancyTooLow;
603     } 
604     else if ( occ > fManuOccupancyLimits.Y() )
605     {
606       rv |= kManuOccupancyTooHigh;
607     }
608   }
609   return rv;
610 }
611
612 //_____________________________________________________________________________
613 AliMUONVCalibParam* 
614 AliMUONPadStatusMaker::PadStatus(Int_t detElemId, Int_t manuId) const
615 {
616   /// Get the status container for a given manu
617   
618   AliMUONVCalibParam* param = static_cast<AliMUONVCalibParam*>(fStatus->FindObject(detElemId,manuId));
619   if (!param)
620   {
621     // not already there, so compute it now
622     AliCodeTimerAuto("ComputeStatus",0);
623     param = ComputeStatus(detElemId,manuId);
624   }
625   return param;
626 }
627
628 //_____________________________________________________________________________
629 Int_t 
630 AliMUONPadStatusMaker::PadStatus(Int_t detElemId, Int_t manuId, Int_t manuChannel) const
631 {
632   /// Get the status for a given channel
633   
634   AliMUONVCalibParam* param = static_cast<AliMUONVCalibParam*>(fStatus->FindObject(detElemId,manuId));
635   if (!param)
636   {
637     // not already there, so compute it now
638     param = ComputeStatus(detElemId,manuId);
639   }
640   return param->ValueAsInt(manuChannel,0);
641 }
642
643 //_____________________________________________________________________________
644 void
645 AliMUONPadStatusMaker::SetHVStatus(Int_t detElemId, Int_t index, Int_t status) const
646 {
647   /// Assign status to all manus in a given HV "zone" (defined by index, meaning
648   /// is different thing from St12 and St345)
649   
650   AliCodeTimerAuto("",0)
651   
652   AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(detElemId);
653   
654   const AliMpArrayI* manus = de->ManusForHV(index);
655   
656   for ( Int_t i = 0; i < manus->GetSize(); ++ i ) 
657   {
658     Int_t manuId = manus->GetValue(i);
659     fHV->Add(AliMpManuUID::BuildUniqueID(detElemId,manuId),status + 1);
660   }
661 }
662
663 //_____________________________________________________________________________
664 Double_t
665 AliMUONPadStatusMaker::HVLimit(Int_t chamberId) const
666 {
667   /// Get HV limit for a given chamber
668   if ( chamberId >=0 && chamberId < 10 ) 
669   {
670     return fHVLimit[chamberId];
671   }
672   return 0.0;
673 }
674
675 //_____________________________________________________________________________
676 void
677 AliMUONPadStatusMaker::SetHVLimit(Int_t chamberId, Double_t hv) 
678 {
679   /// Set hv limit for a given chamber (or all if chamberId==-1)
680   
681   if ( chamberId == -1 ) 
682   {
683     for ( Int_t i = 0; i < 10; ++i ) 
684     {
685       fHVLimit[i] = hv;
686     }
687   }
688   else if ( chamberId >= 0 && chamberId < 10 ) 
689   {
690     fHVLimit[chamberId]=hv;
691   }
692   else
693   {
694     AliError(Form("chamberId=%d is invalid",chamberId));
695   }
696 }
697
698 //_____________________________________________________________________________
699 void
700 AliMUONPadStatusMaker::SetLimits(const AliMUONRecoParam& recoParams) 
701 {
702   /// Set the limits from the recoparam
703   
704   for ( int i = 0; i < 10; ++i )
705   {
706     SetHVLimit(i,recoParams.HVLimit(i));
707   }
708   
709   SetPedMeanLimits(recoParams.PedMeanLowLimit(),recoParams.PedMeanHighLimit());
710   SetPedSigmaLimits(recoParams.PedSigmaLowLimit(),recoParams.PedSigmaHighLimit());
711   
712   SetGainA1Limits(recoParams.GainA1LowLimit(),recoParams.GainA1HighLimit());
713   SetGainA2Limits(recoParams.GainA2LowLimit(),recoParams.GainA2HighLimit());
714   SetGainThresLimits(recoParams.GainThresLowLimit(),recoParams.GainThresHighLimit());
715   
716   SetManuOccupancyLimits(recoParams.ManuOccupancyLowLimit(),recoParams.ManuOccupancyHighLimit());
717   SetBuspatchOccupancyLimits(recoParams.BuspatchOccupancyLowLimit(),recoParams.BuspatchOccupancyHighLimit());
718   SetDEOccupancyLimits(recoParams.DEOccupancyLowLimit(),recoParams.DEOccupancyHighLimit());  
719 }
720
721 //_____________________________________________________________________________
722 void 
723 AliMUONPadStatusMaker::Report(UInt_t mask)
724 {
725   /// Report the number of bad pads, according to the mask,
726   /// and the various reasons why they are bad (with occurence rates)
727   
728   AliInfo("");
729   AliCodeTimerAuto("",0);
730
731   AliMUONLogger log(1064008);
732   
733   Int_t nBadPads(0);
734   Int_t nPads(0);
735   
736   AliMpManuIterator it;
737   
738   Int_t detElemId, manuId;
739   
740   while ( it.Next(detElemId,manuId) )
741   {
742     AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(detElemId);
743     
744     for ( Int_t i = 0; i < AliMpConstants::ManuNofChannels(); ++i )
745     {
746       if ( de->IsConnectedChannel(manuId,i) ) 
747       {
748         ++nPads;
749         
750         Int_t status = PadStatus(detElemId,manuId,i);          
751         
752         if ( mask && ( status & mask) ) // note that if mask == 0, all pads are good...
753         {
754           ++nBadPads;
755           log.Log(AsString(status));
756         }
757       }
758     }
759   }
760   
761   if (!nPads) 
762   {
763     AliError("Got no pad from the iterator ?! That's not normal. Please check !");
764     return;
765   }
766   
767   TString msg;
768   Int_t ntimes;
769   
770   cout << Form("According to mask %x (human readable form below) %6d pads are bad (over a total of %6d, i.e. %7.2f %%)",
771                mask,nBadPads,nPads,nBadPads*100.0/nPads) << endl;
772   cout << AliMUONPadStatusMaker::AsCondition(mask) << endl;
773   cout << "--------" << endl;
774   
775   while ( log.Next(msg,ntimes) )
776   {
777     cout << Form("The message (%120s) occured %15d times (%7.4f %%)",msg.Data(),ntimes,ntimes*100.0/nPads) << endl;
778   }
779   
780 }
781