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