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