]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONPadStatusMaker.cxx
Overload the Clone() method of TObject to replace the method CreateCopy().
[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 "AliCDBEntry.h"
30 #include "AliCDBManager.h"
31 #include "AliCodeTimer.h"
32 #include "AliDCSValue.h"
33 #include "AliLog.h"
34 #include "AliMUON2DMap.h"
35 #include "AliMUON2DStoreValidator.h"
36 #include "AliMUONCalibParamNI.h"
37 #include "AliMUONCalibrationData.h"
38 #include "AliMUONStringIntMap.h"
39 #include "AliMUONVCalibParam.h"
40 #include "AliMpArea.h"
41 #include "AliMpArrayI.h"
42 #include "AliMpConstants.h"
43 #include "AliMpDDLStore.h"
44 #include "AliMpDEIterator.h"
45 #include "AliMpDEManager.h"
46 #include "AliMpDetElement.h"
47 #include "AliMpExMap.h"
48 #include "AliMpHVNamer.h"
49 #include "AliMpIntPair.h"
50 #include "AliMpManuUID.h"
51 #include "AliMpMotifMap.h"
52 #include "AliMpMotifPosition.h"
53 #include "AliMpPCB.h"
54 #include "AliMpPad.h"
55 #include "AliMpSegmentation.h"
56 #include "AliMpSlat.h"
57 #include "AliMpSlatSegmentation.h"
58 #include "AliMpStationType.h"
59 #include "AliMpVPadIterator.h"
60 #include "AliMpVSegmentation.h"
61 #include <Riostream.h>
62 #include <TArrayI.h>
63 #include <TExMap.h>
64 #include <TMap.h>
65 #include <TString.h>
66
67 /// \cond CLASSIMP
68 ClassImp(AliMUONPadStatusMaker)
69 /// \endcond
70
71 //_____________________________________________________________________________
72 AliMUONPadStatusMaker::AliMUONPadStatusMaker(const AliMUONCalibrationData& calibData)
73 : fCalibrationData(calibData),
74   fGainA0Limits(0,1E30),
75   fGainA1Limits(-1E-30,1E30),
76   fGainThresLimits(0,4095),
77   fHVSt12Limits(0,5000),
78   fHVSt345Limits(0,5000),
79   fPedMeanLimits(0,4095),
80   fPedSigmaLimits(0,4095),
81   fStatus(new AliMUON2DMap(true)),
82   fHV(new TExMap),
83   fPedestals(calibData.Pedestals()),
84   fGains(calibData.Gains())
85 {
86    /// ctor
87     AliInfo(Form("ped store %s gain store %s",
88                  fPedestals->ClassName(),
89                  fGains->ClassName()));
90 }
91
92 //_____________________________________________________________________________
93 AliMUONPadStatusMaker::~AliMUONPadStatusMaker()
94 {
95   /// dtor.
96   delete fStatus;
97   delete fHV;
98 }
99
100 //_____________________________________________________________________________
101 TString
102 AliMUONPadStatusMaker::AsString(Int_t status)
103 {
104   /// return a human readable version of the integer status
105   Int_t pedStatus;
106   Int_t gainStatus;
107   Int_t hvStatus;
108   
109   DecodeStatus(status,pedStatus,hvStatus,gainStatus);
110   
111 //  /// Gain status
112 //  enum EGainStatus
113 //  {
114 //    kGainOK = 0,
115 //    kGainA0TooLow = (1<<1),
116 //    kGainA0TooHigh = (1<<2),
117 //    kGainA1TooLow = (1<<3),
118 //    kGainA1TooHigh = (1<<4),
119 //    kGainThresTooLow = (1<<5),
120 //    kGainThresTooHigh = (1<<6),
121 //    
122 //    kGainMissing = kMissing // please always use last bit for meaning "missing"
123 //  };
124 //  
125 //  /// Pedestal status
126 //  enum EPedestalStatus
127 //  {
128 //    kPedOK = 0,
129 //    kPedMeanZero = (1<<1),
130 //    kPedMeanTooLow = (1<<2),
131 //    kPedMeanTooHigh = (1<<3),
132 //    kPedSigmaTooLow = (1<<4),
133 //    kPedSigmaTooHigh = (1<<5),
134 //    
135 //    kPedMissing = kMissing // please always use last bit for meaning "missing"
136 //  };
137 //  
138   TString s("PED ");
139   
140   if ( pedStatus == 0 ) s+= " OK";
141   if ( pedStatus & kPedMeanZero ) s += " Mean is Zero. ";
142   if ( pedStatus & kPedMeanTooLow ) s += " Mean Too Low. ";
143   if ( pedStatus & kPedMeanTooHigh ) s += " Mean Too High. ";
144   if ( pedStatus & kPedSigmaTooLow ) s += " Sigma Too Low. ";
145   if ( pedStatus & kPedSigmaTooHigh ) s += " Sigma Too High. ";
146   if ( pedStatus & kPedMissing ) s += " is missing.";
147   
148 //  /// HV Error
149 //  enum EHVError 
150 //  {
151 //    kHVOK = 0,
152 //    kHVError = (1<<0),
153 //    kHVTooLow = (1<<1),
154 //    kHVTooHigh = (1<<2),
155 //    kHVChannelOFF = (1<<3),
156 //    kHVSwitchOFF = (1<<4),
157 //    
158 //    kHVMissing = kMissing // please always use last bit for meaning "missing"
159 //  };
160
161   return s;
162 }
163
164 //_____________________________________________________________________________
165 Int_t
166 AliMUONPadStatusMaker::BuildStatus(Int_t pedStatus, 
167                                    Int_t hvStatus, 
168                                    Int_t gainStatus)
169 {
170   /// Build a complete status from specific parts (ped,hv,gain)
171   
172   return ( hvStatus & 0xFF ) | ( ( pedStatus & 0xFF ) << 8 ) | 
173   ( ( gainStatus & 0xFF ) << 16 );
174 }
175
176 //_____________________________________________________________________________
177 void
178 AliMUONPadStatusMaker::DecodeStatus(Int_t status, 
179                                     Int_t& pedStatus, 
180                                     Int_t& hvStatus, 
181                                     Int_t& gainStatus)
182 {
183   /// Decode complete status into specific parts (ped,hv,gain)
184   
185   gainStatus = ( status & 0xFF0000 ) >> 16;
186   pedStatus = ( status & 0xFF00 ) >> 8;
187   hvStatus = (status & 0xFF);
188 }
189
190 //_____________________________________________________________________________
191 Bool_t 
192 AliMUONPadStatusMaker::HVSt12Status(Int_t detElemId, Int_t sector,
193                                     Bool_t& hvChannelTooLow,
194                                     Bool_t& hvChannelTooHigh,
195                                     Bool_t& hvChannelON) const
196 {
197   /// Get HV status for one HV sector of St12
198   
199   /// For a given PCB in a given DE, get the HV status (both the channel
200   /// and the switch).
201   /// Returns false if hv switch changed during the run.
202   
203   AliCodeTimerAuto("")
204   
205   Bool_t error = kFALSE;
206   hvChannelTooLow = kFALSE;
207   hvChannelTooHigh = kFALSE;
208   hvChannelON = kTRUE;
209
210   AliMpHVNamer hvNamer;
211   
212   TString hvChannel(hvNamer.DCSHVChannelName(detElemId,sector));
213   
214   TMap* hvMap = fCalibrationData.HV();
215   TPair* hvPair = static_cast<TPair*>(hvMap->FindObject(hvChannel.Data()));
216   if (!hvPair)
217   {
218     AliError(Form("Did not find expected alias (%s) for DE %d",
219                   hvChannel.Data(),detElemId));  
220     error = kTRUE;
221   }
222   else
223   {
224     TObjArray* values = static_cast<TObjArray*>(hvPair->Value());
225     if (!values)
226     {
227       AliError(Form("Could not get values for alias %s",hvChannel.Data()));
228       error = kTRUE;
229     }
230     else
231     {
232       // find out min and max value, and makes a cut
233       Float_t hvMin(1E9);
234       Float_t hvMax(0);
235       TIter next(values);
236       AliDCSValue* val;
237       
238       while ( ( val = static_cast<AliDCSValue*>(next()) ) )
239       {
240         Float_t hv = val->GetFloat();
241         hvMin = TMath::Min(hv,hvMin);
242         hvMax = TMath::Max(hv,hvMax);
243       }
244       
245       float lowThreshold = fHVSt12Limits.X();
246       float highThreshold = fHVSt12Limits.Y();
247             
248       if ( hvMin < lowThreshold ) hvChannelTooLow = kTRUE;
249       if ( hvMax > highThreshold ) hvChannelTooHigh = kTRUE;
250       if ( hvMin < 1 ) hvChannelON = kFALSE;
251     }
252   }
253   
254   return error;
255 }
256
257 //_____________________________________________________________________________
258 Bool_t 
259 AliMUONPadStatusMaker::HVSt345Status(Int_t detElemId, Int_t pcbIndex,
260                                      Bool_t& hvChannelTooLow,
261                                      Bool_t& hvChannelTooHigh,
262                                      Bool_t& hvChannelON,
263                                      Bool_t& hvSwitchON) const
264 {
265   /// For a given PCB in a given DE, get the HV status (both the channel
266   /// and the switch).
267   /// Returns false if something goes wrong (in particular if 
268   /// hv switch changed during the run).
269   
270   AliCodeTimerAuto("")
271   
272   Bool_t error = kFALSE;
273   hvChannelTooLow = kFALSE;
274   hvChannelTooHigh = kFALSE;
275   hvSwitchON = kTRUE;
276   hvChannelON = kTRUE;
277   
278   AliMpHVNamer hvNamer;
279   
280   TString hvChannel(hvNamer.DCSHVChannelName(detElemId));
281   
282   TMap* hvMap = fCalibrationData.HV();
283   
284   TPair* hvPair = static_cast<TPair*>(hvMap->FindObject(hvChannel.Data()));
285   if (!hvPair)
286   {
287     AliError(Form("Did not find expected alias (%s) for DE %d",
288                   hvChannel.Data(),detElemId));  
289     error = kTRUE;
290   }
291   else
292   {
293     TObjArray* values = static_cast<TObjArray*>(hvPair->Value());
294     if (!values)
295     {
296       AliError(Form("Could not get values for alias %s",hvChannel.Data()));
297       error = kTRUE;
298     }
299     else
300     {
301       // find out min and max value, and makes a cut
302       Float_t hvMin(1E9);
303       Float_t hvMax(0);
304       TIter next(values);
305       AliDCSValue* val;
306       
307       while ( ( val = static_cast<AliDCSValue*>(next()) ) )
308       {
309         Float_t hv = val->GetFloat();
310         hvMin = TMath::Min(hv,hvMin);
311         hvMax = TMath::Max(hv,hvMax);
312       }
313
314       float lowThreshold = fHVSt345Limits.X();
315       float highThreshold = fHVSt345Limits.Y();
316
317       if ( hvMin < lowThreshold ) hvChannelTooLow = kTRUE;
318       else if ( hvMax > highThreshold ) hvChannelTooHigh = kTRUE;
319       if ( hvMin < 1 ) hvChannelON = kFALSE;
320     }
321   }
322   
323   TString hvSwitch(hvNamer.DCSHVSwitchName(detElemId,pcbIndex));
324   TPair* switchPair = static_cast<TPair*>(hvMap->FindObject(hvSwitch.Data()));
325   if (!switchPair)
326   {
327     AliError(Form("Did not find expected alias (%s) for DE %d PCB %d",
328                   hvSwitch.Data(),detElemId,pcbIndex));
329     error = kTRUE;
330   }
331   else
332   {
333     TObjArray* values = static_cast<TObjArray*>(switchPair->Value());
334     if (!values)
335     {    
336       AliError(Form("Could not get values for alias %s",hvSwitch.Data()));
337       error = kTRUE;
338     }
339     else
340     {
341       // we'll count the number of ON/OFF for this pad, to insure
342       // consistency (i.e. if status changed during the run, we should
343       // at least notify this fact ;-) and hope it's not the norm)
344       Int_t nTrue(0);
345       Int_t nFalse(0);
346       TIter next(values);
347       AliDCSValue* val;
348       
349       while ( ( val = static_cast<AliDCSValue*>(next()) ) )
350       {
351         if ( val->GetBool() )
352         {
353           ++nTrue;
354         }
355         else
356         {
357           ++nFalse;
358         }
359       }
360       
361       if ( (nTrue>0 && nFalse>0) )
362       {
363         AliWarning(Form("Status of HV Switch %s changed during this run nTrue=%d nFalse=%d! Will consider it OFF",
364                         hvSwitch.Data(),nTrue,nFalse));
365         error = kTRUE;
366       }
367       
368       if ( nFalse ) hvSwitchON = kFALSE;
369     }
370   }
371   return error;
372 }
373
374 //_____________________________________________________________________________
375 Int_t
376 AliMUONPadStatusMaker::HVStatus(Int_t detElemId, Int_t manuId) const
377 {
378   /// Get HV status of one manu
379   
380   AliCodeTimerAuto("")
381   
382   if ( !fCalibrationData.HV() ) return kMissing;
383
384   Long_t lint = fHV->GetValue(AliMpManuUID::BuildUniqueID(detElemId,manuId));
385   
386   if ( lint ) 
387   {
388     return (Int_t)(lint - 1);
389   }
390
391   Int_t status(0);
392   
393   AliMpHVNamer hvNamer;
394   
395   switch ( AliMpDEManager::GetStationType(detElemId) )
396   {
397     case AliMp::kStation1:
398     case AliMp::kStation2:
399     {
400       int sector = hvNamer.ManuId2Sector(detElemId,manuId);
401       if ( sector >= 0 ) 
402       {
403         Bool_t hvChannelTooLow, hvChannelTooHigh, hvChannelON;
404         Bool_t error = HVSt12Status(detElemId,sector,
405                                     hvChannelTooLow,
406                                     hvChannelTooHigh,
407                                     hvChannelON);
408         if ( error ) status |= kHVError;
409         if ( hvChannelTooLow ) status |= kHVTooLow;
410         if ( hvChannelTooHigh ) status |= kHVTooHigh; 
411         if ( !hvChannelON ) status |= kHVChannelOFF;
412         // assign this status to all the other manus handled by the same HV channel
413         SetHVStatus(detElemId,sector,status);
414       }
415     }
416       break;
417     case AliMp::kStation345:
418     {
419       int pcbIndex = hvNamer.ManuId2PCBIndex(detElemId,manuId);
420       if ( pcbIndex >= 0 ) 
421       {
422         Bool_t hvChannelTooLow, hvChannelTooHigh, hvChannelON,hvSwitchON;
423         Bool_t error = HVSt345Status(detElemId,pcbIndex,
424                                      hvChannelTooLow,hvChannelTooHigh,
425                                      hvChannelON,hvSwitchON);
426         if ( error ) status |= kHVError;
427         if ( hvChannelTooLow ) status |= kHVTooLow;
428         if ( hvChannelTooHigh ) status |= kHVTooHigh; 
429         if ( !hvSwitchON ) status |= kHVSwitchOFF; 
430         if ( !hvChannelON) status |= kHVChannelOFF;
431         // assign this status to all the other manus handled by the same HV channel
432         SetHVStatus(detElemId,pcbIndex,status);
433       }
434     }
435       break;
436     default:
437       break;
438   }
439   
440   return status;
441 }
442
443 //_____________________________________________________________________________
444 AliMUONVCalibParam* 
445 AliMUONPadStatusMaker::Neighbours(Int_t detElemId, Int_t manuId) const
446 {
447   /// Get the neighbours parameters for a given manu
448   AliMUONVStore* neighbourStore = fCalibrationData.Neighbours();
449   return static_cast<AliMUONVCalibParam*>(neighbourStore->FindObject(detElemId,manuId));
450 }
451
452 //_____________________________________________________________________________
453 AliMUONVStore* 
454 AliMUONPadStatusMaker::NeighboursStore() const
455 {
456   /// Return the store containing all the neighbours
457   return fCalibrationData.Neighbours();
458 }
459
460 //_____________________________________________________________________________
461 AliMUONVCalibParam*
462 AliMUONPadStatusMaker::ComputeStatus(Int_t detElemId, Int_t manuId) const
463 {
464   /// Compute the status of a given manu, using all available information,
465   /// i.e. pedestals, gains, and HV
466   
467 //  AliCodeTimerAuto("")
468   
469 //  AliCodeTimerStart("Param creation");
470   AliMUONVCalibParam* param = new AliMUONCalibParamNI(1,AliMpConstants::ManuNofChannels(),detElemId,manuId,-1);
471   fStatus->Add(param);
472 //  AliCodeTimerStop("Param creation");
473   
474 //  AliCodeTimerStart("FindObject");
475   AliMUONVCalibParam* pedestals = static_cast<AliMUONVCalibParam*>(fPedestals->FindObject(detElemId,manuId));
476
477   AliMUONVCalibParam* gains = static_cast<AliMUONVCalibParam*>(fGains->FindObject(detElemId,manuId));
478 //  AliCodeTimerStop("FindObject");
479   
480   Int_t hvStatus = HVStatus(detElemId,manuId);
481
482 //  AliCodeTimerStart("Loop");
483   
484   for ( Int_t manuChannel = 0; manuChannel < param->Size(); ++manuChannel )
485   {
486     Int_t pedStatus(0);
487     
488     if (pedestals) 
489     {
490       Float_t pedMean = pedestals->ValueAsFloatFast(manuChannel,0);
491       Float_t pedSigma = pedestals->ValueAsFloatFast(manuChannel,1);
492       if ( pedMean < fPedMeanLimits.X() ) pedStatus |= kPedMeanTooLow;
493       else if ( pedMean > fPedMeanLimits.Y() ) pedStatus |= kPedMeanTooHigh;
494       if ( pedSigma < fPedSigmaLimits.X() ) pedStatus |= kPedSigmaTooLow;
495       else if ( pedSigma > fPedSigmaLimits.Y() ) pedStatus |= kPedSigmaTooHigh;
496       if ( pedMean == 0 ) pedStatus |= kPedMeanZero;
497     }
498     else
499     {
500       pedStatus = kPedMissing;
501     }
502     
503     Int_t gainStatus(0);
504   
505     if ( gains ) 
506     {
507       Float_t a0 = gains->ValueAsFloatFast(manuChannel,0);
508       Float_t a1 = gains->ValueAsFloatFast(manuChannel,1);
509       Float_t thres = gains->ValueAsFloatFast(manuChannel,2);
510   
511       if ( a0 < fGainA0Limits.X() ) gainStatus |= kGainA0TooLow;
512       else if ( a0 > fGainA0Limits.Y() ) gainStatus |= kGainA0TooHigh;
513       if ( a1 < fGainA1Limits.X() ) gainStatus |= kGainA1TooLow;
514       else if ( a1 > fGainA1Limits.Y() ) gainStatus |= kGainA1TooHigh;
515       if ( thres < fGainThresLimits.X() ) gainStatus |= kGainThresTooLow;
516       else if ( thres > fGainThresLimits.Y() ) gainStatus |= kGainThresTooHigh;
517     }
518     else
519     {
520       gainStatus = kGainMissing;
521     }
522         
523     Int_t status = BuildStatus(pedStatus,hvStatus,gainStatus);
524       
525     param->SetValueAsIntFast(manuChannel,0,status);
526   }
527   
528 //  AliCodeTimerStop("Loop");
529   
530   return param;
531 }
532
533 //_____________________________________________________________________________
534 AliMUONVCalibParam* 
535 AliMUONPadStatusMaker::PadStatus(Int_t detElemId, Int_t manuId) const
536 {
537   /// Get the status for a given channel
538   
539   AliMUONVCalibParam* param = static_cast<AliMUONVCalibParam*>(fStatus->FindObject(detElemId,manuId));
540   if (!param)
541   {
542     // not already there, so compute it now
543     AliCodeTimerAuto("ComputeStatus");
544     param = ComputeStatus(detElemId,manuId);
545   }
546   return param;
547 }
548
549 //_____________________________________________________________________________
550 Int_t 
551 AliMUONPadStatusMaker::PadStatus(Int_t detElemId, Int_t manuId, Int_t manuChannel) const
552 {
553   /// Get the status for a given channel
554   
555   AliMUONVCalibParam* param = static_cast<AliMUONVCalibParam*>(fStatus->FindObject(detElemId,manuId));
556   if (!param)
557   {
558     // not already there, so compute it now
559     param = ComputeStatus(detElemId,manuId);
560   }
561   return param->ValueAsInt(manuChannel,0);
562 }
563
564 //_____________________________________________________________________________
565 void
566 AliMUONPadStatusMaker::SetHVStatus(Int_t detElemId, Int_t index, Int_t status) const
567 {
568   /// Assign status to all manus in a given HV "zone" (defined by index, meaning
569   /// is different thing from St12 and St345)
570   
571   AliCodeTimerAuto("")
572   
573   AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(detElemId);
574   
575   const AliMpArrayI* manus = de->ManusForHV(index);
576   
577   for ( Int_t i = 0; i < manus->GetSize(); ++ i ) 
578   {
579     Int_t manuId = manus->GetValue(i);
580     fHV->Add(AliMpManuUID::BuildUniqueID(detElemId,manuId),status + 1);
581   }
582 }