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