]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONCDB.C
Updated comments (Raffaele)
[u/mrichter/AliRoot.git] / MUON / MUONCDB.C
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 /// By Laurent Aphecetche
19
20 #if !defined(__CINT__) || defined(__MAKECINT__)
21
22 #include "MUONCDB.h"
23
24 #include "AliCDBEntry.h"
25 #include "AliCDBManager.h"
26 #include "AliMUON1DArray.h"
27 #include "AliMUON2DMap.h"
28 #include "AliMUONCalibParam1I.h"
29 #include "AliMUONCalibParam2F.h"
30 #include "AliMUONConstants.h"
31 #include "AliMUONTriggerLut.h"
32 #include "AliMUONV2DStore.h"
33 #include "AliMUONVCalibParam.h"
34 #include "AliMpDEIterator.h"
35 #include "AliMpDEManager.h"
36 #include "AliMpSegFactory.h"
37 #include "AliMpStationType.h"
38 #include "AliMpVSegmentation.h"
39 #include "AliMUONTriggerEfficiencyCells.h"
40 #include "Riostream.h"
41 #include "TH1F.h"
42 #include "TList.h"
43 #include "TRandom.h"
44 #include "TStopwatch.h"
45 #include "TSystem.h"
46 #endif
47
48 //_____________________________________________________________________________
49 Int_t countChannels(AliMpVSegmentation& seg)
50 {
51   Int_t n(0);
52   
53   for ( Int_t ix = 0; ix < seg.MaxPadIndexX(); ++ix )
54   {
55     for ( Int_t iy = 0; iy < seg.MaxPadIndexY(); ++iy )
56     {
57       if ( seg.HasPad(AliMpIntPair(ix,iy)) ) ++n;
58     }
59   }
60   return n;
61 }
62
63 //_____________________________________________________________________________
64 AliMpSegFactory* segFactory()
65 {
66   static AliMpSegFactory* sf = new AliMpSegFactory();
67   return sf;
68 }
69
70 //_____________________________________________________________________________
71 void countChannels()
72 {
73   AliMpDEIterator it;
74   Int_t ntotal(0);
75   Int_t ntracker(0);
76   Int_t ntrigger(0);
77   
78   for ( Int_t station = 0; station < AliMUONConstants::NCh(); ++station )
79   {
80     Int_t n(0);
81     it.First(station);
82     while (!it.IsDone())
83     {
84       Int_t de = it.CurrentDE();
85       for ( Int_t cathode = 0; cathode < 2; ++cathode )
86       {
87         AliMpVSegmentation* seg = segFactory()->CreateMpSegmentation(de,cathode);
88         n += countChannels(*seg);
89       }
90       it.Next();
91     }
92     cout << "Station " << station << " has " << n << " channels" << endl;
93     if ( station < AliMUONConstants::NTrackingCh() )
94     {
95       ntracker += n;
96     }
97     else
98     {
99       ntrigger += n;
100     }
101     ntotal += n;
102   }
103   cout << "Tracker channels = " << ntracker << endl;
104   cout << "Trigger channels = " << ntrigger << endl;
105   cout << "Total channels =" << ntotal << endl;
106 }
107
108 //_____________________________________________________________________________
109 AliMUONV2DStore* read2D(const char* calibType)
110 {
111   AliCDBManager* man = AliCDBManager::Instance();
112   man->SetDefaultStorage(CDBPath);
113
114   AliCDBEntry* entry = man->Get(calibType,0);
115
116   if (entry)
117     {
118       return (AliMUONV2DStore*)entry->GetObject();
119     }
120   return 0;
121 }
122
123 //_____________________________________________________________________________
124 AliMUONV1DStore* read1D(const char* calibType)
125 {
126   AliCDBManager* man = AliCDBManager::Instance();
127   man->SetDefaultStorage(CDBPath);
128   
129   AliCDBEntry* entry = man->Get(calibType,0);
130   
131   if (entry)
132   {
133     return (AliMUONV1DStore*)entry->GetObject();
134   }
135   return 0;
136 }
137
138 //_____________________________________________________________________________
139 void checkCDB(const char* calibType)
140 {
141   TString c(calibType);
142   Float_t refValue(0);
143   
144   if ( c == "MUON/Calib/DeadChannels" )
145   {
146     refValue=5;
147   }
148    
149   AliMUONV2DStore* store = read2D(calibType);
150   if (!store) return;
151   
152   TIter next(manuList());
153   AliMpIntPair* p;
154   
155   while ( ( p = (AliMpIntPair*)next() ) )
156   {
157     Int_t detElemId = p->GetFirst();
158     Int_t manuId = p->GetSecond();
159     
160     AliMUONVCalibParam* value = 
161       dynamic_cast<AliMUONVCalibParam*>(store->Get(detElemId,manuId));
162     
163     if (value)
164     {
165       for ( Int_t manuChannel = 0; manuChannel < value->Size(); ++manuChannel )
166       {
167         Float_t testValue = value->ValueAsFloat(manuChannel,0);
168         if ( testValue && testValue != refValue )
169         {
170           cout << "Got a strange value for DE=" << detElemId << " manuId="
171           << manuId << " manuChannel=" << manuChannel << " was expecting "
172           << refValue << " and I got " << testValue << endl;
173         }
174       }
175     }
176     else
177     {
178       cout << "Got a null value for DE=" << detElemId << " manuId="
179       << manuId << endl;
180     }
181   }
182   
183   delete store;
184 }
185
186
187 //_____________________________________________________________________________
188 void plotCDB(const char* calibType)
189 {
190   TString c(calibType);
191   
192   TH1* h = 0;
193   TH1* h2 = 0;
194   
195   if ( c == "MUON/Calib/Gains" )
196   {
197     h = new TH1F("gains_mean","mean gain",100,0,1.5);
198     h2 = new TH1F("saturation","adc saturation",4096,-0.5,4095.5);
199   }
200   else if ( c == "MUON/Calib/Pedestals" )
201   {
202     h = new TH1F("pedestals_mean","pedestals_mean",4096,-0.5,4095.5);
203     h2 = new TH1F("pedestals_sigma","pedestals_sigma",100,0,20);
204   }
205   else if ( c == "MUON/Calib/DeadChannels" )
206   {
207     h = new TH1F("dead_channels","dead channels per DE",1500,-0.5,1499.5);
208   }
209   else
210   {
211     cerr << "Don't know how to deal with " << calibType << endl;
212     return;
213   }
214   
215   AliMUONV2DStore* store = read2D(calibType);
216   if (!store) return;
217
218   TIter next(manuList());
219   AliMpIntPair* p;
220   Int_t n(0);
221   Int_t ndead(0);
222   Int_t nPerStation[7];
223   
224   for ( Int_t i = 0; i < 7; ++i ) nPerStation[i]=0;
225   
226   while ( ( p = (AliMpIntPair*)next() ) )
227   {
228     Int_t detElemId = p->GetFirst();
229     Int_t manuId = p->GetSecond();
230     
231     Int_t station = detElemId/100 - 1;
232     
233     AliMUONVCalibParam* value = 
234     dynamic_cast<AliMUONVCalibParam*>(store->Get(detElemId,manuId));
235     
236     if (value)
237     {
238       for ( Int_t manuChannel = 0; manuChannel < value->Size(); ++manuChannel )
239       {
240         ++n;
241         ++nPerStation[station];
242         if (h2)
243         {
244           h->Fill(value->ValueAsFloat(manuChannel,0));
245           h2->Fill(value->ValueAsFloat(manuChannel,1));
246         }
247         else
248         {
249           if( value->ValueAsInt(manuChannel) )
250           {
251             h->Fill(detElemId);
252             ++ndead;
253           }
254         }
255       }
256     }
257     else
258     {
259       cout << "Got a null value for DE=" << detElemId << " manuId="
260         << manuId << endl;
261     }
262   }
263   
264   cout << "Number of channels = " << n << endl;
265   for ( Int_t i = 0; i < 7; ++i )
266   {
267     cout << "Station " << (i+1) << " " << nPerStation[i] << endl;
268   }
269   
270   if (n && c == "MUON/Calib/DeadChannels")
271   {
272     cout << "Number of dead channels=" << ndead << endl;
273   }
274   delete store;
275 }
276
277 //_____________________________________________________________________________
278 void testReadStore(const AliMUONV2DStore& store, Int_t n)
279 {
280   TIter next(manuList());
281   AliMpIntPair* p;
282   
283   while ( ( p = (AliMpIntPair*)next() ) )
284   {
285     for ( Int_t i = 0; i < n; ++i )
286     {
287       store.Get(p->GetFirst(),p->GetSecond());
288     }
289   }
290
291
292 //_____________________________________________________________________________
293 Int_t makeStores(AliMUONV2DStore& pedestalStore,
294                  AliMUONV2DStore& gainStore,
295                  AliMUONV2DStore& deadStore,
296                  Bool_t defaultValues)
297 {  
298   TIter next(manuList());
299   
300   AliMpIntPair* p;
301   
302   Int_t ngenerated(0);
303   
304   Bool_t replace = kFALSE;
305   
306   const Int_t nChannels(64);
307   
308   while ( ( p = (AliMpIntPair*)next() ) )
309   {
310     AliMUONVCalibParam* ped = new AliMUONCalibParam2F(nChannels);
311     AliMUONVCalibParam* gain = new AliMUONCalibParam2F(nChannels);
312     AliMUONVCalibParam* dead = new AliMUONCalibParam1I(nChannels);
313
314     for ( Int_t manuChannel = 0; manuChannel < nChannels; ++manuChannel )
315     {
316       Float_t meanPedestal;
317       Float_t sigmaPedestal;
318       Float_t meanGain;
319       Float_t saturation(3000);
320     
321       if ( defaultValues ) 
322       {
323         meanPedestal = 0.0;
324         sigmaPedestal = 1.0;
325         meanGain = 1.0;
326       }
327       else
328       {
329         meanPedestal = -1;
330         while ( meanPedestal < 0 )
331         {
332           meanPedestal = gRandom->Gaus(150,10);
333         }
334         sigmaPedestal = -1;
335         while ( sigmaPedestal < 0 )
336         {
337           sigmaPedestal = gRandom->Gaus(1,0.2);
338         }
339         meanGain = -1;
340         while ( meanGain < 0 )
341         {
342           meanGain = gRandom->Gaus(1,0.05);
343         }
344       }
345       ped->SetValueAsFloat(manuChannel,0,meanPedestal);
346       ped->SetValueAsFloat(manuChannel,1,sigmaPedestal);
347       gain->SetValueAsFloat(manuChannel,0,meanGain);
348       gain->SetValueAsFloat(manuChannel,1,saturation);
349       
350       if (!defaultValues)
351       {
352         // probability that this channel is dead ~ 1%
353         if ( gRandom->Uniform(100.0) < 1.0 ) 
354         {
355           Int_t reason = 5; // that value could be used to distinguish
356           // why the channel is dead or how it was flagged bad (online,
357           // offline, by chance...). 5 is of course a fake number.
358           dead->SetValueAsInt(manuChannel,0,reason);
359         }
360       }
361     }
362     Int_t detElemId = p->GetFirst();
363     Int_t manuId = p->GetSecond();
364     Bool_t ok1 = pedestalStore.Set(detElemId,manuId,ped,replace);
365     Bool_t ok2 = gainStore.Set(detElemId,manuId,gain,replace);
366     Bool_t ok3 = deadStore.Set(detElemId,manuId,dead,replace);
367     if (!ok1 || !ok2 || !ok3)
368     {
369       cout << "Could not set DetElemId=" << detElemId << " manuId="
370         << manuId << endl;
371     }
372     else
373     {
374       ++ngenerated;
375     }
376   }
377   
378   return ngenerated;
379 }
380
381 //_____________________________________________________________________________
382 TList* manuList(Bool_t reset)
383 {
384   static TList* fgmanuList = new TList;
385   
386   if (reset) 
387   {
388     fgmanuList->Delete();
389     return fgmanuList;
390   }
391   
392   if (!fgmanuList->IsEmpty()) return fgmanuList;
393   
394   TStopwatch timer;
395   
396   cout << "Generating manu list. Please wait" << endl;
397   
398   fgmanuList->SetOwner(kTRUE);
399   
400   timer.Start();
401   
402   AliMpDEIterator it;
403   
404   it.First();
405   
406   while ( !it.IsDone() )
407   {
408     Int_t detElemId = it.CurrentDE();
409     AliMpStationType stationType = AliMpDEManager::GetStationType(detElemId);
410     if ( stationType != kStationTrigger ) 
411     {
412       for ( Int_t cath = 0; cath <=1 ; ++cath )
413       {
414         AliMpVSegmentation* seg = segFactory()->CreateMpSegmentation(detElemId,cath);
415         
416         TArrayI manus;
417         
418         seg->GetAllElectronicCardIDs(manus);
419         
420         for ( Int_t im = 0; im < manus.GetSize(); ++im )
421         {
422           fgmanuList->Add(new AliMpIntPair(detElemId,manus[im]));
423         }        
424       }
425     }
426     it.Next();
427   }
428   
429   cout << "Time to make the manu list = ";
430   timer.Print();
431   
432   return fgmanuList;
433 }
434
435 //_____________________________________________________________________________
436 void testMakeStores(Int_t readLoop)
437 {
438   manuList();
439   
440   AliMUONV2DStore* pedestalStore = new AliMUON2DMap;
441   AliMUONV2DStore* gainStore = new AliMUON2DMap;
442   AliMUONV2DStore* deadStore = new AliMUON2DMap;
443   
444   TStopwatch timer;
445   
446   cout << "Creating" << endl;
447   
448   timer.Start(kTRUE);
449   makeStores(*pedestalStore,*gainStore,*deadStore,true);
450   timer.Print();
451   
452   cout << "Reading..." << endl;
453   timer.Start(kTRUE);
454   testReadStore(*pedestalStore,readLoop);
455   testReadStore(*gainStore,readLoop);
456   testReadStore(*deadStore,readLoop);
457   cout << timer.CpuTime()/readLoop << " CPUs (mean of " << readLoop 
458     <<" samples." << endl;
459   
460   delete pedestalStore;
461   delete gainStore;
462   delete deadStore;
463 }
464
465 //_____________________________________________________________________________
466 void generateTrigger(const char* cdbpath)
467 {
468   // 
469   // Generate trigger related conditions :
470   //
471   // - trigger masks for board (locals, regionals, global)
472   // - trigger lut
473   // - trigger efficiency
474   // - trigger switches (to be implemented FIXME)
475   //
476   const Int_t nlboards = 234;
477   AliMUONV1DStore* localBoardMasks = new AliMUON1DArray(nlboards+1);
478   
479   // Generate fake mask values for 234 localboards and put that into
480   // one single container (localBoardMasks)
481   for ( Int_t i = 1; i <= nlboards; ++i )
482   {
483     AliMUONVCalibParam* localBoard = new AliMUONCalibParam1I(8);
484     for ( Int_t x = 0; x < 2; ++x )
485     {
486       for ( Int_t y = 0; y < 4; ++y )
487       {
488         Int_t index = x*4+y;
489         localBoard->SetValueAsInt(index,0,0xFFFF);
490       }
491     }
492     localBoardMasks->Set(i,localBoard,kFALSE);
493   }
494   
495   // Generate values for regional boards
496   const Int_t nrboards = 16;
497   AliMUONV1DStore* regionalBoardMasks = new AliMUON1DArray(16);
498   
499   for ( Int_t i = 0; i < nrboards; ++i )
500   {
501     AliMUONVCalibParam* regionalBoard = new AliMUONCalibParam1I(16);
502     for ( Int_t j = 0; j < 16; ++j )
503     {
504       regionalBoard->SetValueAsInt(j,0,0x3F);
505     }
506     regionalBoardMasks->Set(i,regionalBoard,kFALSE);
507   }
508   
509   // Generate values for global board
510   AliMUONVCalibParam* globalBoardMasks = new AliMUONCalibParam1I(16);
511   
512   for ( Int_t j = 0; j < 16; ++j )
513   {
514     globalBoardMasks->SetValueAsInt(j,0,0xFFF);
515   }
516     
517   AliMUONTriggerLut lut;
518   lut.ReadFromFile("$(ALICE_ROOT)/MUON/data/lutAptLpt1Hpt1p7.root");
519   
520   AliMUONTriggerEfficiencyCells cells("$ALICE_ROOT/MUON/data/efficiencyCells.dat");
521   
522   //--------------------------------------------
523   // Store the resulting containers into the CDB
524   Int_t ever = 99999999;
525   
526   AliCDBId id("MUON/Calib/LocalTriggerBoardMasks",0,ever);
527   
528   AliCDBMetaData md;
529   md.SetBeamPeriod(1);
530   md.SetAliRootVersion(gROOT->GetVersion());
531   md.SetComment("Test with default values");
532   md.SetResponsible("Rachid Guernane");
533   
534   AliCDBManager* man = AliCDBManager::Instance();
535   man->SetDefaultStorage(cdbpath);
536   man->Put(localBoardMasks,id,&md);
537   
538   id.SetPath("MUON/Calib/RegionalTriggerBoardMasks");
539   
540   man->Put(regionalBoardMasks,id,(AliCDBMetaData*)md.Clone());
541   
542   id.SetPath("MUON/Calib/GlobalTriggerBoardMasks");
543   
544   man->Put(globalBoardMasks,id,(AliCDBMetaData*)md.Clone());
545
546   id.SetPath("MUON/Calib/TriggerLut");
547   
548   man->Put(&lut,id,(AliCDBMetaData*)md.Clone());
549   
550   id.SetPath("MUON/Calib/TriggerEfficiency");
551   md.SetResponsible("Diego Stocco");
552   
553   man->Put(&cells,id,(AliCDBMetaData*)md.Clone());
554   
555   delete localBoardMasks;
556   delete regionalBoardMasks;
557   delete globalBoardMasks;
558 }
559
560 //_____________________________________________________________________________
561 void generateCalibrations(const char* cdbpath, Bool_t defaultValues)
562 {
563   //
564   //
565   //
566
567   AliMUONV2DStore* pedestalStore = new AliMUON2DMap;
568   AliMUONV2DStore* gainStore = new AliMUON2DMap;
569   AliMUONV2DStore* deadStore = new AliMUON2DMap; 
570   
571   makeStores(*pedestalStore,*gainStore,*deadStore,defaultValues);
572   
573   Int_t ever = 99999999;
574   
575   AliCDBId id("MUON/Calib/Pedestals",0,ever);
576   AliCDBMetaData md;
577   md.SetBeamPeriod(1);
578   md.SetAliRootVersion(gROOT->GetVersion());
579   if ( defaultValues )
580   {
581     md.SetComment("Test with default values");
582   }
583   else
584   {
585     md.SetComment("Test with random values");
586   }
587   md.SetResponsible("Laurent Aphecetche");
588   
589   AliCDBManager* man = AliCDBManager::Instance();
590   man->SetDefaultStorage(cdbpath);
591   man->Put(pedestalStore,id,&md);
592   
593   id.SetPath("MUON/Calib/Gains");
594   
595   man->Put(gainStore,id,(AliCDBMetaData*)md.Clone());
596   
597   id.SetPath("MUON/Calib/DeadChannels");
598   
599   man->Put(deadStore,id,(AliCDBMetaData*)md.Clone());
600   
601   delete deadStore;
602   delete pedestalStore;
603   delete gainStore;
604 }
605
606