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