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