]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONPadStatusMaker.cxx
Initial version (Laurent)
[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 "AliDCSValue.h"
32 #include "AliLog.h"
33 #include "AliMUON2DMap.h"
34 #include "AliMUON2DStoreValidator.h"
35 #include "AliMUONCalibParamNI.h"
36 #include "AliMUONCalibrationData.h"
37 #include "AliMUONHVNamer.h"
38 #include "AliMUONVCalibParam.h"
39 #include "AliMpArea.h"
40 #include "AliMpConstants.h"
41 #include "AliMpDEIterator.h"
42 #include "AliMpDEManager.h"
43 #include "AliMpIntPair.h"
44 #include "AliMpManuList.h"
45 #include "AliMpMotifMap.h"
46 #include "AliMpMotifPosition.h"
47 #include "AliMpPCB.h"
48 #include "AliMpPad.h"
49 #include "AliMpSector.h"
50 #include "AliMpSectorSegmentation.h"
51 #include "AliMpSegmentation.h"
52 #include "AliMpSlat.h"
53 #include "AliMpSlatSegmentation.h"
54 #include "AliMpStationType.h"
55 #include "AliMpVPadIterator.h"
56 #include <Riostream.h>
57 #include <TMap.h>
58 #include <TStopwatch.h>
59 #include <TString.h>
60
61
62 /// \cond CLASSIMP
63 ClassImp(AliMUONPadStatusMaker)
64 /// \endcond
65
66 //_____________________________________________________________________________
67 AliMUONPadStatusMaker::AliMUONPadStatusMaker(const AliMUONCalibrationData& calibData)
68 : fCalibrationData(calibData),
69   fPedMeanLimits(0,4095),
70   fPedSigmaLimits(0,4095),
71   fHVSt12Limits(0,5000),
72   fHVSt345Limits(0,5000)
73 {
74    /// ctor
75 }
76
77 //_____________________________________________________________________________
78 AliMUONPadStatusMaker::~AliMUONPadStatusMaker()
79 {
80   /// dtor.
81 }
82
83 //_____________________________________________________________________________
84 AliMUONVStore*
85 AliMUONPadStatusMaker::Combine(const AliMUONVStore& store1,
86                                const AliMUONVStore& store2,
87                                Int_t binShift) const
88 {
89   /// Combine two status containers into one, shifting store2 status bits
90   /// to the left by binShift before making an OR with store1.
91   
92   TStopwatch timer;
93   timer.Start(kTRUE);
94   
95   AliMUONVStore* combined = static_cast<AliMUONVStore*>(store1.Clone());
96   
97   TIter next(store1.CreateIterator());
98   AliMUONVCalibParam* param1;
99   
100   while ( ( param1 = static_cast<AliMUONVCalibParam*>(next()) ) )
101   {
102     Int_t detElemId = param1->ID0();
103     Int_t manuId = param1->ID1();
104     AliMUONVCalibParam* param2 = static_cast<AliMUONVCalibParam*>(store2.FindObject(detElemId,manuId));
105     if (!param2)
106     {
107       AliWarning(Form("Could not get statuses for store2 for DE %d ManuId %d. Marking as missing.",
108                     detElemId,manuId));
109       param2 = static_cast<AliMUONVCalibParam*>(param1->Clone());
110       for ( Int_t manuChannel = 0; manuChannel < param2->Size(); ++manuChannel )
111       {
112         param2->SetValueAsInt(manuChannel,0,kMissing);
113       }      
114     }
115     AliMUONVCalibParam* paramCombined = static_cast<AliMUONVCalibParam*>(combined->FindObject(detElemId,manuId));
116     if (!paramCombined)
117     {
118       paramCombined = static_cast<AliMUONVCalibParam*>(param2->Clone());
119       combined->Add(paramCombined);
120     }
121     
122     for ( Int_t manuChannel = 0; manuChannel < param1->Size(); ++manuChannel )
123     {
124       if ( AliMpManuList::DoesChannelExist(detElemId, manuId, manuChannel) )
125       {
126         Int_t status1(param1->ValueAsInt(manuChannel));
127         Int_t status2(param2->ValueAsInt(manuChannel));
128         
129         Int_t status = status1 | (status2 << binShift);
130         
131         paramCombined->SetValueAsInt(manuChannel,0,status);
132       }
133     }
134   }
135   
136   AliInfo("Timer:");
137   StdoutToAliInfo(timer.Print(););
138   
139   return combined;
140 }
141
142 //_____________________________________________________________________________
143 AliMUONVStore* 
144 AliMUONPadStatusMaker::GeneratePadStatus(Int_t value)
145 {
146   /// Generate a "fake" store, with all (detElemId,manuId) present,
147   /// and containing all the same value
148   
149   AliMUONVStore* store = new AliMUON2DMap(true);
150   
151   TList* list = AliMpManuList::ManuList();
152   
153   AliMpIntPair* pair;
154   
155   TIter next(list);
156   
157   while ( ( pair = static_cast<AliMpIntPair*>(next()) ) ) 
158   {
159     Int_t detElemId = pair->GetFirst();
160     Int_t manuId = pair->GetSecond();
161     AliMUONVCalibParam* param = new AliMUONCalibParamNI(1,AliMpConstants::ManuNofChannels(),detElemId,manuId,value);
162     store->Add(param);
163   }
164   
165   delete list;
166   
167   return store;
168 }
169
170 //_____________________________________________________________________________
171 Bool_t 
172 AliMUONPadStatusMaker::GetSt12Status(const TMap& hvMap,
173                                      Int_t detElemId, Int_t sector,
174                                      Bool_t& hvChannelTooLow,
175                                      Bool_t& hvChannelTooHigh,
176                                      Bool_t& hvChannelON) const
177 {
178   /// Get HV status for one HV sector of St12
179   
180   /// For a given PCB in a given DE, get the HV status (both the channel
181   /// and the switch).
182   /// Returns false if hv switch changed during the run.
183   
184   Bool_t error = kFALSE;
185   hvChannelTooLow = kFALSE;
186   hvChannelTooHigh = kFALSE;
187   hvChannelON = kTRUE;
188   
189   AliMUONHVNamer hvNamer;
190   
191   TString hvChannel(hvNamer.DCSHVChannelName(detElemId,sector));
192   
193   TPair* hvPair = static_cast<TPair*>(hvMap.FindObject(hvChannel.Data()));
194   if (!hvPair)
195   {
196     AliError(Form("Did not find expected alias (%s) for DE %d",
197                   hvChannel.Data(),detElemId));  
198     error = kTRUE;
199   }
200   else
201   {
202     TObjArray* values = static_cast<TObjArray*>(hvPair->Value());
203     if (!values)
204     {
205       AliError(Form("Could not get values for alias %s",hvChannel.Data()));
206       error = kTRUE;
207     }
208     else
209     {
210       // find out min and max value, and makes a cut
211       Float_t hvMin(1E9);
212       Float_t hvMax(0);
213       TIter next(values);
214       AliDCSValue* val;
215       
216       while ( ( val = static_cast<AliDCSValue*>(next()) ) )
217       {
218         Float_t hv = val->GetFloat();
219         hvMin = TMath::Min(hv,hvMin);
220         hvMax = TMath::Max(hv,hvMax);
221       }
222       
223       float lowThreshold = fHVSt12Limits.X();
224       float highThreshold = fHVSt12Limits.Y();
225             
226       if ( hvMin < lowThreshold ) hvChannelTooLow = kTRUE;
227       if ( hvMax > highThreshold ) hvChannelTooHigh = kTRUE;
228       if ( hvMin < 1 ) hvChannelON = kFALSE;
229     }
230   }
231   
232   return error;
233 }
234
235 //_____________________________________________________________________________
236 Bool_t 
237 AliMUONPadStatusMaker::GetSt345Status(const TMap& hvMap,
238                                       Int_t detElemId, Int_t pcbIndex,
239                                       Bool_t& hvChannelTooLow,
240                                       Bool_t& hvChannelTooHigh,
241                                       Bool_t& hvChannelON,
242                                       Bool_t& hvSwitchON) const
243 {
244   /// For a given PCB in a given DE, get the HV status (both the channel
245   /// and the switch).
246   /// Returns false if something goes wrong (in particular if 
247   /// hv switch changed during the run).
248   
249   Bool_t error = kFALSE;
250   hvChannelTooLow = kFALSE;
251   hvChannelTooHigh = kFALSE;
252   hvSwitchON = kTRUE;
253   hvChannelON = kTRUE;
254   
255   AliMUONHVNamer hvNamer;
256   
257   TString hvChannel(hvNamer.DCSHVChannelName(detElemId));
258   
259   TPair* hvPair = static_cast<TPair*>(hvMap.FindObject(hvChannel.Data()));
260   if (!hvPair)
261   {
262     AliError(Form("Did not find expected alias (%s) for DE %d",
263                   hvChannel.Data(),detElemId));  
264     error = kTRUE;
265   }
266   else
267   {
268     TObjArray* values = static_cast<TObjArray*>(hvPair->Value());
269     if (!values)
270     {
271       AliError(Form("Could not get values for alias %s",hvChannel.Data()));
272       error = kTRUE;
273     }
274     else
275     {
276       // find out min and max value, and makes a cut
277       Float_t hvMin(1E9);
278       Float_t hvMax(0);
279       TIter next(values);
280       AliDCSValue* val;
281       
282       while ( ( val = static_cast<AliDCSValue*>(next()) ) )
283       {
284         Float_t hv = val->GetFloat();
285         hvMin = TMath::Min(hv,hvMin);
286         hvMax = TMath::Max(hv,hvMax);
287       }
288
289       float lowThreshold = fHVSt345Limits.X();
290       float highThreshold = fHVSt345Limits.Y();
291
292       if ( hvMin < lowThreshold ) hvChannelTooLow = kTRUE;
293       if ( hvMax > highThreshold ) hvChannelTooHigh = kTRUE;
294       if ( hvMin < 1 ) hvChannelON = kFALSE;
295     }
296   }
297   
298   TString hvSwitch(hvNamer.DCSHVSwitchName(detElemId,pcbIndex));
299   TPair* switchPair = static_cast<TPair*>(hvMap.FindObject(hvSwitch.Data()));
300   if (!switchPair)
301   {
302     AliError(Form("Did not find expected alias (%s) for DE %d PCB %d",
303                   hvSwitch.Data(),detElemId,pcbIndex));
304     error = kTRUE;
305   }
306   else
307   {
308     TObjArray* values = static_cast<TObjArray*>(switchPair->Value());
309     if (!values)
310     {    
311       AliError(Form("Could not get values for alias %s",hvSwitch.Data()));
312       error = kTRUE;
313     }
314     else
315     {
316       // we'll count the number of ON/OFF for this pad, to insure
317       // consistency (i.e. if status changed during the run, we should
318       // at least notify this fact ;-) and hope it's not the norm)
319       Int_t nTrue(0);
320       Int_t nFalse(0);
321       TIter next(values);
322       AliDCSValue* val;
323       
324       while ( ( val = static_cast<AliDCSValue*>(next()) ) )
325       {
326         if ( val->GetBool() )
327         {
328           ++nTrue;
329         }
330         else
331         {
332           ++nFalse;
333         }
334       }
335       
336       if ( (nTrue>0 && nFalse>0) )
337       {
338         AliWarning(Form("Status of HV Switch %s changed during this run nTrue=%d nFalse=%d! Will consider it OFF",
339                         hvSwitch.Data(),nTrue,nFalse));
340         error = kTRUE;
341       }
342       
343       if ( nFalse ) hvSwitchON = kFALSE;
344     }
345   }
346   return error;
347 }
348
349 //_____________________________________________________________________________
350 AliMUONVStore* 
351 AliMUONPadStatusMaker::MakeGainStatus(const AliMUONVStore& /*gainValues*/) const
352 {
353   /// FIXME: to be implemented
354   AliWarning("Not implemented yet");
355   return 0x0;
356 }
357
358 //_____________________________________________________________________________
359 AliMUONVStore* 
360 AliMUONPadStatusMaker::MakeHVStatus(const TMap& hvValues) const
361 {
362   /// Scrutinize HV values and deduce an HV status for each pad
363   
364   TStopwatch timerSt12;
365   TStopwatch timerSt345;
366   
367   timerSt12.Start(kTRUE);
368   timerSt12.Stop();
369   timerSt345.Start(kTRUE);
370   timerSt345.Stop();
371   
372   AliMUONHVNamer hvNamer;
373   
374   AliMpDEIterator deIt;
375   
376   deIt.First();
377   
378   AliMUONVStore* hv = new AliMUON2DMap(kTRUE);
379   
380   while ( !deIt.IsDone() )
381   {
382     Int_t detElemId = deIt.CurrentDEId();
383     
384     switch ( AliMpDEManager::GetStationType(detElemId) )
385     {
386       case AliMp::kStation1:
387       case AliMp::kStation2:
388         timerSt12.Start(kFALSE);
389         for ( int sector = 0; sector < 3; ++sector)
390         {
391           AliDebug(1,Form("detElemId %5d sector %d",detElemId,sector));
392
393           Bool_t hvChannelTooLow, hvChannelTooHigh, hvChannelON;
394           Bool_t error = GetSt12Status(hvValues,
395                                        detElemId,sector,
396                                        hvChannelTooLow,hvChannelTooHigh,
397                                        hvChannelON);
398           Int_t status = 0;
399           if ( error ) status |= kHVError;
400           if ( hvChannelTooLow ) status |= kHVTooLow;
401           if ( hvChannelTooHigh ) status |= kHVTooHigh; 
402           if ( !hvChannelON ) status |= kHVChannelOFF;
403           SetStatusSt12(*hv,detElemId,sector,status);
404           
405         }
406           timerSt12.Stop();
407         break;
408       case AliMp::kStation345:
409       {
410         timerSt345.Start(kFALSE);
411         for ( Int_t pcbIndex = 0; pcbIndex < hvNamer.NumberOfPCBs(detElemId); ++pcbIndex)
412         {
413           AliDebug(1,Form("detElemId %5d pcbIndex %d",detElemId,pcbIndex));
414           Bool_t hvChannelTooLow, hvChannelTooHigh, hvChannelON,hvSwitchON;
415           Bool_t error = GetSt345Status(hvValues,
416                                         detElemId,pcbIndex,
417                                         hvChannelTooLow,hvChannelTooHigh,
418                                         hvChannelON,hvSwitchON);
419           Int_t status = 0;
420           if ( error ) status |= kHVError;
421           if ( hvChannelTooLow ) status |= kHVTooLow;
422           if ( hvChannelTooHigh ) status |= kHVTooHigh; 
423           if ( !hvSwitchON ) status |= kHVSwitchOFF; 
424           if ( !hvChannelON) status |= kHVChannelOFF;
425           SetStatusSt345(*hv,detElemId,pcbIndex,status);
426         }
427         timerSt345.Stop();
428       }
429         break;
430       default:
431         break;
432     }
433     deIt.Next();
434   }
435   
436   AliInfo("St12 timer:");
437   StdoutToAliInfo(timerSt12.Print(););
438   AliInfo("St345 timer:");
439   StdoutToAliInfo(timerSt345.Print(););
440   
441   return hv;
442 }
443
444 //_____________________________________________________________________________
445 AliMUONVStore* 
446 AliMUONPadStatusMaker::MakePedestalStatus(const AliMUONVStore& pedValues) const
447 {
448   /// Assign a pedestal status to each pad
449   
450   TStopwatch timer;
451   
452   timer.Start(kTRUE);
453   
454   AliMUONVStore* pedStatuses = new AliMUON2DMap(kTRUE);
455   
456   TIter next(pedValues.CreateIterator());
457   AliMUONVCalibParam* pedestals;
458   Int_t nofManus(0);
459   
460   while ( ( pedestals = static_cast<AliMUONVCalibParam*>(next() ) ) )
461   {
462     Int_t detElemId = pedestals->ID0();
463     Int_t manuId = pedestals->ID1();
464     ++nofManus;
465     for ( Int_t manuChannel = 0; manuChannel < pedestals->Size(); ++manuChannel )
466     {
467       Int_t status(0);
468       if ( AliMpManuList::DoesChannelExist(detElemId, manuId, manuChannel) )
469       {
470         Float_t pedMean = pedestals->ValueAsFloat(manuChannel,0);
471         Float_t pedSigma = pedestals->ValueAsFloat(manuChannel,1);
472         if ( pedMean < fPedMeanLimits.X() ) status |= kPedMeanTooLow;
473         if ( pedMean > fPedMeanLimits.Y() ) status |= kPedMeanTooHigh;
474         if ( pedSigma < fPedSigmaLimits.X() ) status |= kPedSigmaTooLow;
475         if ( pedSigma > fPedSigmaLimits.Y() ) status |= kPedSigmaTooHigh;
476         if ( pedMean == 0 ) status |= kPedMeanZero;
477         
478         AliMUONVCalibParam* vStatus = 
479           static_cast<AliMUONVCalibParam*>(pedStatuses->FindObject(detElemId,manuId));
480         if ( !vStatus ) 
481         {
482           vStatus = new AliMUONCalibParamNI(1,AliMpConstants::ManuNofChannels(),detElemId,manuId,0);
483           pedStatuses->Add(vStatus);
484         }
485         vStatus->SetValueAsInt(manuChannel,0,status);
486       }
487     }
488   }
489   
490   AliInfo(Form("%d manus checked in :",nofManus));
491   StdoutToAliInfo(timer.Print(););
492   return pedStatuses;  
493 }
494
495 //_____________________________________________________________________________
496 AliMUONVStore* 
497 AliMUONPadStatusMaker::MakeStatus() const
498 {
499   /// Read ped, gains and hv values from CDB, apply some Q&A and produces
500   /// a combined status for each pad.
501
502   TMap* hvValues = fCalibrationData.HV();
503   AliMUONVStore* hvStatus(0x0);
504   
505   if (!hvValues)
506   {
507     AliError("Could not get HV values from CDB. Will create dummy ones and mark those as missing");
508     hvStatus = GeneratePadStatus(kHVMissing);
509   }
510   else
511   {
512     hvStatus = MakeHVStatus(*hvValues);
513   }
514   
515   AliMUONVStore* pedValues = fCalibrationData.Pedestals();
516   AliMUONVStore* pedStatus(0x0);
517
518   if (!pedValues)
519   {
520     AliError("Could not get pedestals values from CDB. Will create dummy ones and mark those as missing");
521     pedStatus = GeneratePadStatus(kPedMissing);
522   }
523   else
524   {
525     pedStatus = MakePedestalStatus(*pedValues);
526   }
527   
528   // FIXME: should do the same for gains as for hv and ped.    
529   
530   AliMUONVStore* status = Combine(*hvStatus,*pedStatus,8);
531   
532   delete hvStatus;
533   delete pedStatus;
534   
535   // Insure we get all channels there (some or even all can be bad, but they
536   // must be there somehow).
537   
538   AliMUON2DStoreValidator validator;
539         
540   TObjArray* a = validator.Validate(*status);
541     
542   if (a) 
543   {
544     // this should not happen.
545     AliError("Status store not complete. Crash to follow soon...");
546     StdoutToAliError(a->Print(););
547     AliFatal("this should not happen at all!");
548     delete status;
549     status = 0x0;
550   }
551     
552   return status;
553 }
554
555 //_____________________________________________________________________________
556 void
557 AliMUONPadStatusMaker::SetStatusSt12(AliMUONVStore& hvStatus,
558                                      Int_t detElemId, 
559                                      Int_t isector,
560                                      Int_t status) const
561 {
562   /// Flag all pads of detElemId (for St12) as bad.
563   
564   // FIXME: need a way to iterator on pads over a given HV sector for St12... 
565   // we currently suppose that one sector is about a third of the chamber...
566   // FIXME !! This has to be checked very carefully...
567   
568   const AliMp::CathodType kCathodes[] = { AliMp::kCath0, AliMp::kCath1 };
569   
570   for ( Int_t icathode = 0; icathode < 2; ++icathode )
571   {
572     const AliMpSectorSegmentation* seg = 
573     static_cast<const AliMpSectorSegmentation*>(AliMpSegmentation::Instance()->GetMpSegmentation(detElemId,kCathodes[icathode]));
574     const AliMpSector* sector = seg->GetSector();
575     AliMpMotifMap* mMap = sector->GetMotifMap();
576     TArrayI a;
577     
578     mMap->GetAllMotifPositionsIDs(a);
579     
580     TVector2 dim = seg->Dimensions();
581     Double_t x = dim.X()*2;
582     Double_t xmin = isector*x/3.0;
583     Double_t xmax = xmin + x/3.0;   
584     
585     for ( Int_t i = 0; i < a.GetSize(); ++i ) 
586     {
587       AliMpMotifPosition* pos = mMap->FindMotifPosition(a[i]);
588       Int_t manuId = pos->GetID();
589       TVector2 position = pos->Position();
590       if ( position.X() >= xmin && position.X() <= xmax) 
591       {
592         AliMUONVCalibParam* dead =
593         static_cast<AliMUONVCalibParam*>(hvStatus.FindObject(detElemId,manuId));
594         if (!dead)
595         {
596           dead = new AliMUONCalibParamNI(1,AliMpConstants::ManuNofChannels(),detElemId,manuId,status);
597           hvStatus.Add(dead);
598         }        
599         else
600         {
601           // FIXME: this should really not happen, if we'd know really the
602           // relationship between manuId and HV sector...
603           // For the time being, let's leave it like that, for testing
604           // purposes only. For production, this will have to be fixed.
605           AliWarning("Please fixme.");
606         }
607       }
608     }
609   }  
610 }
611
612 //_____________________________________________________________________________
613 void
614 AliMUONPadStatusMaker::SetStatusSt345(AliMUONVStore& hvStatus,
615                                       Int_t detElemId, Int_t pcbIndex,
616                                       Int_t status) const
617 {
618   /// Flag all pads of pcbIndex-th PCB of detElemId (for St345) as bad.
619   
620   const AliMp::CathodType kCathodes[] = { AliMp::kCath0, AliMp::kCath1 };
621   
622   for ( Int_t icathode = 0; icathode < 2; ++icathode )
623   {
624     const AliMpSlatSegmentation* seg = static_cast<const AliMpSlatSegmentation*>
625     (AliMpSegmentation::Instance()->GetMpSegmentation(detElemId,kCathodes[icathode]));
626     const AliMpSlat* slat = seg->Slat();
627     const AliMpPCB* pcb = slat->GetPCB(pcbIndex);
628
629     for ( Int_t i = 0; i < pcb->GetSize(); ++i ) 
630     {
631       AliMpMotifPosition* pos = pcb->GetMotifPosition(i);
632       Int_t manuId = pos->GetID();
633       AliMUONVCalibParam* dead = 
634         static_cast<AliMUONVCalibParam*>(hvStatus.FindObject(detElemId,manuId));
635       if (dead)
636       {
637         AliError(Form("dead is not null as expected from DE %d manuId %d",
638                       detElemId,manuId));
639       }
640       if (!dead)
641       {
642         dead = new AliMUONCalibParamNI(1,AliMpConstants::ManuNofChannels(),detElemId,manuId,status);
643         hvStatus.Add(dead);
644       }
645     }    
646   }
647 }
648
649
650