]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONCDB.C
hardcoded detector position; bug in alignment pth fixed
[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     Int_t station = AliMpDEmanager::GetChamberId(detElemId);
231     
232     AliMUONVCalibParam* value = 
233     dynamic_cast<AliMUONVCalibParam*>(store->Get(detElemId,manuId));
234     
235     if (value)
236     {
237       for ( Int_t manuChannel = 0; manuChannel < value->Size(); ++manuChannel )
238       {
239         ++n;
240         ++nPerStation[station];
241         if (h2)
242         {
243           h->Fill(value->ValueAsFloat(manuChannel,0));
244           h2->Fill(value->ValueAsFloat(manuChannel,1));
245         }
246         else
247         {
248           if( value->ValueAsInt(manuChannel) )
249           {
250             h->Fill(detElemId);
251             ++ndead;
252           }
253         }
254       }
255     }
256     else
257     {
258       cout << "Got a null value for DE=" << detElemId << " manuId="
259         << manuId << endl;
260     }
261   }
262   
263   cout << "Number of channels = " << n << endl;
264   for ( Int_t i = 0; i < 7; ++i )
265   {
266     cout << "Station " << (i+1) << " " << nPerStation[i] << endl;
267   }
268   
269   if (n && c == "MUON/Calib/DeadChannels")
270   {
271     cout << "Number of dead channels=" << ndead << endl;
272   }
273   delete store;
274 }
275
276 //_____________________________________________________________________________
277 void testReadStore(const AliMUONV2DStore& store, Int_t n)
278 {
279   TIter next(manuList());
280   AliMpIntPair* p;
281   
282   while ( ( p = (AliMpIntPair*)next() ) )
283   {
284     for ( Int_t i = 0; i < n; ++i )
285     {
286       store.Get(p->GetFirst(),p->GetSecond());
287     }
288   }
289
290
291 //_____________________________________________________________________________
292 Int_t makeStores(AliMUONV2DStore& pedestalStore,
293                  AliMUONV2DStore& gainStore,
294                  AliMUONV2DStore& deadStore,
295                  Bool_t defaultValues)
296 {  
297   TIter next(manuList());
298   
299   AliMpIntPair* p;
300   
301   Int_t ngenerated(0);
302   
303   Bool_t replace = kFALSE;
304   
305   const Int_t nChannels(64);
306   
307   while ( ( p = (AliMpIntPair*)next() ) )
308   {
309     AliMUONVCalibParam* ped = new AliMUONCalibParam2F(nChannels);
310     AliMUONVCalibParam* gain = new AliMUONCalibParam2F(nChannels);
311     AliMUONVCalibParam* dead = new AliMUONCalibParam1I(nChannels);
312
313     for ( Int_t manuChannel = 0; manuChannel < nChannels; ++manuChannel )
314     {
315       Float_t meanPedestal;
316       Float_t sigmaPedestal;
317       Float_t meanGain;
318       Float_t saturation(3000);
319     
320       if ( defaultValues ) 
321       {
322         meanPedestal = 0.0;
323         sigmaPedestal = 1.0;
324         meanGain = 1.0;
325       }
326       else
327       {
328         meanPedestal = -1;
329         while ( meanPedestal < 0 )
330         {
331           meanPedestal = gRandom->Gaus(150,10);
332         }
333         sigmaPedestal = -1;
334         while ( sigmaPedestal < 0 )
335         {
336           sigmaPedestal = gRandom->Gaus(1,0.2);
337         }
338         meanGain = -1;
339         while ( meanGain < 0 )
340         {
341           meanGain = gRandom->Gaus(1,0.05);
342         }
343       }
344       ped->SetValueAsFloat(manuChannel,0,meanPedestal);
345       ped->SetValueAsFloat(manuChannel,1,sigmaPedestal);
346       gain->SetValueAsFloat(manuChannel,0,meanGain);
347       gain->SetValueAsFloat(manuChannel,1,saturation);
348       
349       if (!defaultValues)
350       {
351         // probability that this channel is dead ~ 1%
352         if ( gRandom->Uniform(100.0) < 1.0 ) 
353         {
354           Int_t reason = 5; // that value could be used to distinguish
355           // why the channel is dead or how it was flagged bad (online,
356           // offline, by chance...). 5 is of course a fake number.
357           dead->SetValueAsInt(manuChannel,0,reason);
358         }
359       }
360     }
361     Int_t detElemId = p->GetFirst();
362     Int_t manuId = p->GetSecond();
363     Bool_t ok1 = pedestalStore.Set(detElemId,manuId,ped,replace);
364     Bool_t ok2 = gainStore.Set(detElemId,manuId,gain,replace);
365     Bool_t ok3 = deadStore.Set(detElemId,manuId,dead,replace);
366     if (!ok1 || !ok2 || !ok3)
367     {
368       cout << "Could not set DetElemId=" << detElemId << " manuId="
369         << manuId << endl;
370     }
371     else
372     {
373       ++ngenerated;
374     }
375   }
376   
377   return ngenerated;
378 }
379
380 //_____________________________________________________________________________
381 TList* manuList(Bool_t reset)
382 {
383   static TList* fgmanuList = new TList;
384   
385   if (reset) 
386   {
387     fgmanuList->Delete();
388     return fgmanuList;
389   }
390   
391   if (!fgmanuList->IsEmpty()) return fgmanuList;
392   
393   TStopwatch timer;
394   
395   cout << "Generating manu list. Please wait" << endl;
396   
397   fgmanuList->SetOwner(kTRUE);
398   
399   timer.Start();
400   
401   AliMpDEIterator it;
402   
403   it.First();
404   
405   while ( !it.IsDone() )
406   {
407     Int_t detElemId = it.CurrentDE();
408     AliMpStationType stationType = AliMpDEManager::GetStationType(detElemId);
409     if ( stationType != kStationTrigger ) 
410     {
411       for ( Int_t cath = 0; cath <=1 ; ++cath )
412       {
413         AliMpVSegmentation* seg = segFactory()->CreateMpSegmentation(detElemId,cath);
414         
415         TArrayI manus;
416         
417         seg->GetAllElectronicCardIDs(manus);
418         
419         for ( Int_t im = 0; im < manus.GetSize(); ++im )
420         {
421           fgmanuList->Add(new AliMpIntPair(detElemId,manus[im]));
422         }        
423       }
424     }
425     it.Next();
426   }
427   
428   cout << "Time to make the manu list = ";
429   timer.Print();
430   
431   return fgmanuList;
432 }
433
434 //_____________________________________________________________________________
435 void testMakeStores(Int_t readLoop)
436 {
437   manuList();
438   
439   AliMUONV2DStore* pedestalStore = new AliMUON2DMap;
440   AliMUONV2DStore* gainStore = new AliMUON2DMap;
441   AliMUONV2DStore* deadStore = new AliMUON2DMap;
442   
443   TStopwatch timer;
444   
445   cout << "Creating" << endl;
446   
447   timer.Start(kTRUE);
448   makeStores(*pedestalStore,*gainStore,*deadStore,true);
449   timer.Print();
450   
451   cout << "Reading..." << endl;
452   timer.Start(kTRUE);
453   testReadStore(*pedestalStore,readLoop);
454   testReadStore(*gainStore,readLoop);
455   testReadStore(*deadStore,readLoop);
456   cout << timer.CpuTime()/readLoop << " CPUs (mean of " << readLoop 
457     <<" samples." << endl;
458   
459   delete pedestalStore;
460   delete gainStore;
461   delete deadStore;
462 }
463
464 //_____________________________________________________________________________
465 void generateTrigger(const char* cdbpath)
466 {
467   // 
468   // Generate trigger related conditions :
469   //
470   // - trigger masks for board (locals, regionals, global)
471   // - trigger lut
472   // - trigger efficiency
473   // - trigger switches (to be implemented FIXME)
474   //
475   const Int_t nlboards = 234;
476   AliMUONV1DStore* localBoardMasks = new AliMUON1DArray(nlboards+1);
477   
478   // Generate fake mask values for 234 localboards and put that into
479   // one single container (localBoardMasks)
480   for ( Int_t i = 1; i <= nlboards; ++i )
481   {
482     AliMUONVCalibParam* localBoard = new AliMUONCalibParam1I(8);
483     for ( Int_t x = 0; x < 2; ++x )
484     {
485       for ( Int_t y = 0; y < 4; ++y )
486       {
487         Int_t index = x*4+y;
488         localBoard->SetValueAsInt(index,0,0xFFFF);
489       }
490     }
491     localBoardMasks->Set(i,localBoard,kFALSE);
492   }
493   
494   // Generate values for regional boards
495   const Int_t nrboards = 16;
496   AliMUONV1DStore* regionalBoardMasks = new AliMUON1DArray(16);
497   
498   for ( Int_t i = 0; i < nrboards; ++i )
499   {
500     AliMUONVCalibParam* regionalBoard = new AliMUONCalibParam1I(16);
501     for ( Int_t j = 0; j < 16; ++j )
502     {
503       regionalBoard->SetValueAsInt(j,0,0x3F);
504     }
505     regionalBoardMasks->Set(i,regionalBoard,kFALSE);
506   }
507   
508   // Generate values for global board
509   AliMUONVCalibParam* globalBoardMasks = new AliMUONCalibParam1I(16);
510   
511   for ( Int_t j = 0; j < 16; ++j )
512   {
513     globalBoardMasks->SetValueAsInt(j,0,0xFFF);
514   }
515     
516   AliMUONTriggerLut lut;
517   lut.ReadFromFile("$(ALICE_ROOT)/MUON/data/lutAptLpt1Hpt1p7.root");
518   
519   AliMUONTriggerEfficiencyCells cells("$ALICE_ROOT/MUON/data/efficiencyCells.dat");
520   
521   //--------------------------------------------
522   // Store the resulting containers into the CDB
523   Int_t ever = 99999999;
524   
525   AliCDBId id("MUON/Calib/LocalTriggerBoardMasks",0,ever);
526   
527   AliCDBMetaData md;
528   md.SetBeamPeriod(1);
529   md.SetAliRootVersion(gROOT->GetVersion());
530   md.SetComment("Test with default values");
531   md.SetResponsible("Rachid Guernane");
532   
533   AliCDBManager* man = AliCDBManager::Instance();
534   man->SetDefaultStorage(cdbpath);
535   man->Put(localBoardMasks,id,&md);
536   
537   id.SetPath("MUON/Calib/RegionalTriggerBoardMasks");
538   
539   man->Put(regionalBoardMasks,id,(AliCDBMetaData*)md.Clone());
540   
541   id.SetPath("MUON/Calib/GlobalTriggerBoardMasks");
542   
543   man->Put(globalBoardMasks,id,(AliCDBMetaData*)md.Clone());
544
545   id.SetPath("MUON/Calib/TriggerLut");
546   
547   man->Put(&lut,id,(AliCDBMetaData*)md.Clone());
548   
549   id.SetPath("MUON/Calib/TriggerEfficiency");
550   md.SetResponsible("Diego Stocco");
551   
552   man->Put(&cells,id,(AliCDBMetaData*)md.Clone());
553   
554   delete localBoardMasks;
555   delete regionalBoardMasks;
556   delete globalBoardMasks;
557 }
558
559 //_____________________________________________________________________________
560 void generateCalibrations(const char* cdbpath, Bool_t defaultValues)
561 {
562   //
563   //
564   //
565
566   AliMUONV2DStore* pedestalStore = new AliMUON2DMap;
567   AliMUONV2DStore* gainStore = new AliMUON2DMap;
568   AliMUONV2DStore* deadStore = new AliMUON2DMap; 
569   
570   makeStores(*pedestalStore,*gainStore,*deadStore,defaultValues);
571   
572   Int_t ever = 99999999;
573   
574   AliCDBId id("MUON/Calib/Pedestals",0,ever);
575   AliCDBMetaData md;
576   md.SetBeamPeriod(1);
577   md.SetAliRootVersion(gROOT->GetVersion());
578   if ( defaultValues )
579   {
580     md.SetComment("Test with default values");
581   }
582   else
583   {
584     md.SetComment("Test with random values");
585   }
586   md.SetResponsible("Laurent Aphecetche");
587   
588   AliCDBManager* man = AliCDBManager::Instance();
589   man->SetDefaultStorage(cdbpath);
590   man->Put(pedestalStore,id,&md);
591   
592   id.SetPath("MUON/Calib/Gains");
593   
594   man->Put(gainStore,id,(AliCDBMetaData*)md.Clone());
595   
596   id.SetPath("MUON/Calib/DeadChannels");
597   
598   man->Put(deadStore,id,(AliCDBMetaData*)md.Clone());
599   
600   delete deadStore;
601   delete pedestalStore;
602   delete gainStore;
603 }
604
605