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