]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONCDB.cxx
Being more realistic with gain calibration (Laurent)
[u/mrichter/AliRoot.git] / MUON / AliMUONCDB.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 AliMUONCDB
19 ///
20 /// Helper class to experience the OCDB
21 /// It allows to generate dummy (but complete) containers for all the
22 /// calibration data types we have for tracker and trigger, and to write
23 /// them into OCDB.
24 ///
25 /// For more information, please see READMECDB
26 ///
27 // \author Laurent Aphecetche
28
29 #include "AliMUONCDB.h"
30
31 #include "AliCDBEntry.h"
32 #include "AliCDBManager.h"
33 #include "AliDCSValue.h"
34 #include "AliLog.h"
35 #include "AliMUON1DArray.h"
36 #include "AliMUON1DMap.h"
37 #include "AliMUON2DMap.h"
38 #include "AliMUON2DStoreValidator.h"
39 #include "AliMUONCalibParamNF.h"
40 #include "AliMUONCalibParamNI.h"
41 #include "AliMUONConstants.h"
42 #include "AliMUONHVNamer.h"
43 #include "AliMUONTriggerEfficiencyCells.h"
44 #include "AliMUONTriggerLut.h"
45 #include "AliMUONVStore.h"
46 #include "AliMUONVCalibParam.h"
47 #include "AliMUONVCalibParam.h"
48 #include "AliMpConstants.h"
49 #include "AliMpDDLStore.h"
50 #include "AliMpDEIterator.h"
51 #include "AliMpDEManager.h"
52 #include "AliMpDetElement.h"
53 #include "AliMpManuList.h"
54 #include "AliMpSegmentation.h"
55 #include "AliMpStationType.h"
56 #include "AliMpVSegmentation.h"
57 #include <Riostream.h>
58 #include <TArrayI.h>
59 #include <TClass.h>
60 #include <TH1F.h>
61 #include <TList.h>
62 #include <TMap.h>
63 #include <TObjString.h>
64 #include <TROOT.h>
65 #include <TRandom.h>
66 #include <TStopwatch.h>
67 #include <TSystem.h>
68 #include <TMath.h>
69
70 /// \cond CLASSIMP
71 ClassImp(AliMUONCDB)
72 /// \endcond
73
74 namespace
75 {
76 //_____________________________________________________________________________
77 void getBoundaries(const AliMUONVStore& store,
78                    Float_t& x0min, Float_t& x0max,
79                    Float_t& x1min, Float_t& x1max)
80 {
81   x0min=1E30;
82   x0max=-1E30;
83   x1min=1E30;
84   x1max=-1E30;
85   
86   TIter next(store.CreateIterator());
87   AliMUONVCalibParam* value;
88   
89   while ( ( value = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
90   {
91     Int_t detElemId = value->ID0();
92     Int_t manuId = value->ID1();
93     
94     const AliMpVSegmentation* seg = 
95       AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
96         
97     for ( Int_t manuChannel = 0; manuChannel < value->Size(); ++manuChannel )
98     {
99       AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
100       if (!pad.IsValid()) continue;
101       
102       Float_t x0 = value->ValueAsFloat(manuChannel,0);
103       
104       x0min = TMath::Min(x0min,x0);
105       x0max = TMath::Max(x0max,x0);
106       if ( value->Dimension()>1 )
107       {
108         Float_t x1 = value->ValueAsFloat(manuChannel,1);
109         x1min = TMath::Min(x1min,x1);
110         x1max = TMath::Max(x1max,x1);
111       }
112     }
113   }  
114 }
115
116 //_____________________________________________________________________________
117 Double_t GetRandom(Double_t mean, Double_t sigma, Bool_t mustBePositive)
118 {
119   Double_t x(-1);
120   if ( mustBePositive ) 
121   {
122     while ( x < 0 ) 
123     {
124       x = gRandom->Gaus(mean,sigma);
125     }
126   }
127   else
128   {
129     x = gRandom->Gaus(mean,sigma);
130   }
131   return x;
132 }
133
134 }
135
136 //_____________________________________________________________________________
137 AliMUONCDB::AliMUONCDB(const char* cdbpath)
138 : TObject(),
139   fCDBPath(cdbpath),
140   fManuList(0x0)
141 {
142     /// ctor
143 }
144
145 //_____________________________________________________________________________
146 AliMUONCDB::~AliMUONCDB()
147 {
148   /// dtor
149   delete fManuList;
150 }
151
152 //_____________________________________________________________________________
153 TList*
154 AliMUONCDB::ManuList()
155 {
156   /// return (and create if necessary) the list of (de,manu) pairs
157   if (!fManuList) 
158   {
159     AliInfo("Generating ManuList...");
160     fManuList = AliMpManuList::ManuList();
161     AliInfo("Manu List generated.");
162   }
163   return fManuList;
164 }
165
166 //_____________________________________________________________________________
167 AliMUONVStore* 
168 AliMUONCDB::Diff(AliMUONVStore& store1, AliMUONVStore& store2, 
169                  const char* opt)
170 {
171   /// creates a store which contains store1-store2
172   /// if opt="abs" the difference is absolute one,
173   /// if opt="rel" then what is stored is (store1-store2)/store1
174   /// WARNING Works only for stores which holds AliMUONVCalibParam objects
175   
176   TString sopt(opt);
177   sopt.ToUpper();
178   
179   if ( !sopt.Contains("ABS") && !sopt.Contains("REL") )
180   {
181     AliErrorClass(Form("opt %s not supported. Only ABS or REL are",opt));
182     return 0x0;
183   }
184   
185   AliMUONVStore* d = static_cast<AliMUONVStore*>(store1.Clone());
186   
187   TIter next(d->CreateIterator());
188   
189   AliMUONVCalibParam* param;
190   
191   while ( ( param = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
192   {
193     Int_t detElemId = param->ID0();
194     Int_t manuId = param->ID1();
195     
196     AliMUONVCalibParam* param2 = dynamic_cast<AliMUONVCalibParam*>(store2.FindObject(detElemId,manuId));
197     //FIXME: this might happen. Handle it.
198     if (!param2) 
199     {
200       cerr << "param2 is null : FIXME : this might happen !" << endl;
201       delete d;
202       return 0;
203     }
204     
205     for ( Int_t i = 0; i < param->Size(); ++i )
206     {
207       for ( Int_t j = 0; j < param->Dimension(); ++j )
208       {
209         Float_t value(0);
210         if ( sopt.Contains("ABS") )
211         {
212           value = param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j);
213         }
214         else if ( sopt.Contains("REL") )
215         {
216           if ( param->ValueAsFloat(i,j) ) 
217           {
218             value = (param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j))/param->ValueAsFloat(i,j);
219           }
220           else 
221           {
222             continue;
223           }
224         }
225         param->SetValueAsFloat(i,j,value);
226       }      
227     }
228   }
229   return d;
230 }
231
232 //_____________________________________________________________________________
233 void 
234 AliMUONCDB::Plot(const AliMUONVStore& store, const char* name, Int_t nbins)
235 {
236   /// Make a plot of the first 1 or 2 dimensions of the AliMUONVCalibParam
237   /// contained inside store.
238   /// It produces histograms named name_0, name_1, etc...
239   
240   Float_t x0min, x0max, x1min, x1max;
241   
242   getBoundaries(store,x0min,x0max,x1min,x1max);
243   
244   if ( x0min > x0max ) 
245   {
246     cerr << Form("Something is wrong with boundaries : x0(min,max)=%e,%e",
247                  x0min,x0max) << endl;
248     return;
249   }
250
251   if ( TMath::Abs(x0min-x0max) < 1E-3 ) 
252   {
253     x0min -= 1;
254     x0max += 1;
255   }
256   
257   TH1* h0 = new TH1F(Form("%s_0",name),Form("%s_0",name),
258                     nbins,x0min,x0max);
259   
260   TH1* h1(0);
261   
262   if ( x1max > x1min )
263   {
264     h1 = new TH1F(Form("%s_1",name),Form("%s_1",name),
265                   nbins,x1min,x1max);
266   }
267   
268   TIter next(ManuList());
269   AliMpIntPair* p;
270   Int_t n(0);
271   Int_t nPerStation[7];
272   
273   for ( Int_t i = 0; i < 7; ++i ) nPerStation[i]=0;
274   
275   while ( ( p = (AliMpIntPair*)next() ) )
276   {
277     Int_t detElemId = p->GetFirst();
278     Int_t manuId = p->GetSecond();
279     Int_t station = AliMpDEManager::GetChamberId(detElemId);
280     
281     const AliMpVSegmentation* seg = 
282       AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
283     
284     AliMUONVCalibParam* value = 
285       dynamic_cast<AliMUONVCalibParam*>(store.FindObject(detElemId,manuId));
286     
287     if (value)
288     {
289       for ( Int_t manuChannel = 0; manuChannel < value->Size(); ++manuChannel )
290       {
291         AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
292         if (!pad.IsValid()) continue;
293
294         ++n;
295         ++nPerStation[station];
296         Float_t x = value->ValueAsFloat(manuChannel,0);
297         if ( x>1E4 ) 
298         {
299           AliInfo(Form("DE %d Manu %d Ch %d x=%e",detElemId,manuId,manuChannel,x));
300         }
301         h0->Fill(x);
302         if (h1)
303         {
304           h1->Fill(value->ValueAsFloat(manuChannel,1));
305         }
306       }
307     }
308     else
309     {
310       AliWarning(Form("Got a null value for DE=%d manuId=%d",detElemId,manuId));
311     }
312   }
313   
314   AliInfo(Form("Number of channels = %d",n));
315   for ( Int_t i = 0; i < 7; ++i )
316   {
317     AliInfo(Form("Station %d %d ",(i+1),nPerStation[i]));
318   }
319 }
320
321 //_____________________________________________________________________________
322 Int_t 
323 AliMUONCDB::MakeHVStore(TMap& aliasMap, Bool_t defaultValues)
324 {
325   /// Create a HV store
326   
327   AliMUONHVNamer hvNamer;
328   
329   TObjArray* aliases = hvNamer.GenerateAliases();
330   
331   Int_t nSwitch(0);
332   Int_t nChannels(0);
333   
334   for ( Int_t i = 0; i < aliases->GetEntries(); ++i ) 
335   {
336     TObjString* alias = static_cast<TObjString*>(aliases->At(i));
337     TString& aliasName = alias->String();
338     if ( aliasName.Contains("sw") ) 
339     {
340       // HV Switch (St345 only)
341       TObjArray* valueSet = new TObjArray;
342       valueSet->SetOwner(kTRUE);
343       
344       Bool_t value = kTRUE;
345       
346       if (!defaultValues)
347       {
348         Float_t r = gRandom->Uniform();
349         if ( r < 0.007 ) value = kFALSE;      
350       } 
351       
352       for ( UInt_t timeStamp = 0; timeStamp < 60*3; timeStamp += 60 )
353       {
354         AliDCSValue* dcsValue = new AliDCSValue(value,timeStamp);
355         valueSet->Add(dcsValue);
356       }
357       aliasMap.Add(new TObjString(*alias),valueSet);
358       ++nSwitch;
359     }
360     else
361     {
362       TObjArray* valueSet = new TObjArray;
363       valueSet->SetOwner(kTRUE);
364       for ( UInt_t timeStamp = 0; timeStamp < 60*15; timeStamp += 120 )
365       {
366         Float_t value = 1500;
367         if (!defaultValues) value = GetRandom(1750,62.5,true);
368         AliDCSValue* dcsValue = new AliDCSValue(value,timeStamp);
369         valueSet->Add(dcsValue);
370       }
371       aliasMap.Add(new TObjString(*alias),valueSet);
372       ++nChannels;
373     }
374   }
375   
376   delete aliases;
377   
378   AliInfo(Form("%d HV channels and %d switches",nChannels,nSwitch));
379   
380   return nChannels+nSwitch;
381 }
382
383 //_____________________________________________________________________________
384 Int_t 
385 AliMUONCDB::MakePedestalStore(AliMUONVStore& pedestalStore, Bool_t defaultValues)
386 {
387   /// Create a pedestal store. if defaultValues=true, ped.mean=ped.sigma=1,
388   /// otherwise mean and sigma are from a gaussian (with parameters
389   /// defined below by the kPedestal* constants)
390   
391   TIter next(ManuList());
392   
393   AliMpIntPair* p;
394   
395   Int_t nchannels(0);
396   Int_t nmanus(0);
397   
398   const Int_t kChannels(AliMpConstants::ManuNofChannels());
399   const Float_t kPedestalMeanMean(200);
400   const Float_t kPedestalMeanSigma(10);
401   const Float_t kPedestalSigmaMean(1.0);
402   const Float_t kPedestalSigmaSigma(0.2);
403   
404   while ( ( p = (AliMpIntPair*)next() ) )
405   {
406     ++nmanus;
407
408     Int_t detElemId = p->GetFirst();    
409     Int_t manuId = p->GetSecond();
410     
411     AliMUONVCalibParam* ped = new AliMUONCalibParamNF(2,kChannels,detElemId,manuId,AliMUONVCalibParam::InvalidFloatValue());
412     
413     const AliMpVSegmentation* seg = 
414       AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
415     
416     for ( Int_t manuChannel = 0; manuChannel < kChannels; ++manuChannel )
417     {
418       AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
419       if (!pad.IsValid()) continue;
420       
421       ++nchannels;
422       
423       Float_t meanPedestal;
424       Float_t sigmaPedestal;
425       
426       if ( defaultValues ) 
427       {
428         meanPedestal = 0.0;
429         sigmaPedestal = 1.0;
430       }
431       else
432       {
433         Bool_t positive(kTRUE);
434         meanPedestal = GetRandom(kPedestalMeanMean,kPedestalMeanSigma,positive);
435         sigmaPedestal = GetRandom(kPedestalSigmaMean,kPedestalSigmaSigma,positive);
436       }
437       ped->SetValueAsFloat(manuChannel,0,meanPedestal);
438       ped->SetValueAsFloat(manuChannel,1,sigmaPedestal);
439       
440     }
441     Bool_t ok = pedestalStore.Add(ped);
442     if (!ok)
443     {
444       AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
445     }
446   }
447   
448   AliInfo(Form("%d Manus and %d channels.",nmanus,nchannels));
449   return nchannels;
450   
451 }
452
453 //_____________________________________________________________________________
454 Int_t 
455 AliMUONCDB::MakeCapacitanceStore(AliMUONVStore& capaStore, Bool_t defaultValues)
456 {
457   /// Create a capacitance store. if defaultValues=true, all capa are 1.0,
458   /// otherwise they are from a gaussian with parameters defined in the
459   /// kCapa* constants below.
460   
461   TIter next(ManuList());
462   
463   AliMpIntPair* p;
464   
465   Int_t nchannels(0);
466   Int_t nmanus(0);
467   Int_t nmanusOK(0); // manus for which we got the serial number
468     
469   const Float_t kCapaMean(1.0);
470   const Float_t kCapaSigma(0.5);
471   
472   while ( ( p = (AliMpIntPair*)next() ) )
473   {
474     ++nmanus;
475     
476     Int_t detElemId = p->GetFirst();
477     Int_t manuId = p->GetSecond();
478     
479     AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(detElemId); 
480     Int_t serialNumber = de->GetManuSerialFromId(manuId);
481       
482     if ( serialNumber <= 0 ) continue;
483     
484     ++nmanusOK;
485     
486     const AliMpVSegmentation* seg = 
487       AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
488
489     AliMUONVCalibParam* capa = static_cast<AliMUONVCalibParam*>(capaStore.FindObject(serialNumber));
490     
491     if (!capa)
492     {
493       capa = new AliMUONCalibParamNF(1,AliMpConstants::ManuNofChannels(),serialNumber,0,1.0);
494       Bool_t ok = capaStore.Add(capa);
495       if (!ok)
496       {
497         AliError(Form("Could not set serialNumber=%d manuId=%d",serialNumber,manuId));
498       }      
499     }
500     
501     for ( Int_t manuChannel = 0; manuChannel < capa->Size(); ++manuChannel )
502     {
503       AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
504       if (!pad.IsValid()) continue;
505       
506       ++nchannels;
507       
508       Float_t capaValue;
509       
510       if ( defaultValues ) 
511       {
512         capaValue = 1.0;
513       }
514       else
515       {
516         capaValue = GetRandom(kCapaMean,kCapaSigma,kTRUE);
517       }
518       capa->SetValueAsFloat(manuChannel,0,capaValue);
519     }
520   }
521   
522   Float_t percent = 0;
523   if ( nmanus ) percent = 100*nmanusOK/nmanus;
524   AliInfo(Form("%5d manus with serial number (out of %5d manus = %3.0f%%)",
525                nmanusOK,nmanus,percent));
526   AliInfo(Form("%5d channels",nchannels));
527   if ( percent < 100 ) 
528   {
529     AliWarning("Did not get all serial numbers. capaStore is incomplete !!!!");
530   }
531   return nchannels;
532   
533 }
534
535 //_____________________________________________________________________________
536 Int_t 
537 AliMUONCDB::MakeGainStore(AliMUONVStore& gainStore, Bool_t defaultValues)
538 {  
539   /// Create a gain store. if defaultValues=true, all gains set so that
540   /// charge = (adc-ped)
541   ///
542   /// otherwise parameters are taken from gaussians with parameters 
543   /// defined in the k* constants below.
544   
545   TIter next(ManuList());
546   
547   AliMpIntPair* p;
548   
549   Int_t nchannels(0);
550   Int_t nmanus(0);
551     
552   const Int_t kSaturation(3000);
553   const Double_t kA0Mean(1.2);
554   const Double_t kA0Sigma(0.1);
555   const Double_t kA1Mean(1E-5);
556   const Double_t kA1Sigma(1E-6);
557   const Double_t kQualMean(0xFF);
558   const Double_t kQualSigma(0x10);
559   const Int_t kThresMean(1600);
560   const Int_t kThresSigma(100);
561   
562   while ( ( p = (AliMpIntPair*)next() ) )
563   {
564     ++nmanus;
565
566     Int_t detElemId = p->GetFirst();
567     Int_t manuId = p->GetSecond();
568
569     AliMUONVCalibParam* gain = 
570       new AliMUONCalibParamNF(5,AliMpConstants::ManuNofChannels(),
571                               detElemId,
572                               manuId,
573                               AliMUONVCalibParam::InvalidFloatValue());
574
575
576     const AliMpVSegmentation* seg = 
577       AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
578
579     for ( Int_t manuChannel = 0; manuChannel < gain->Size(); ++manuChannel )
580     {
581       AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
582       if (!pad.IsValid()) continue;
583       
584       ++nchannels;
585       
586       if ( defaultValues ) 
587       {
588         gain->SetValueAsFloat(manuChannel,0,1.0);
589         gain->SetValueAsFloat(manuChannel,1,0.0);
590         gain->SetValueAsInt(manuChannel,2,4095); 
591         gain->SetValueAsInt(manuChannel,3,1);
592         gain->SetValueAsInt(manuChannel,4,kSaturation);
593       }
594       else
595       {
596         Bool_t positive(kTRUE);
597         gain->SetValueAsFloat(manuChannel,0,GetRandom(kA0Mean,kA0Sigma,positive));
598         gain->SetValueAsFloat(manuChannel,1,GetRandom(kA1Mean,kA1Sigma,!positive));
599         gain->SetValueAsInt(manuChannel,2,(Int_t)TMath::Nint(GetRandom(kThresMean,kThresSigma,positive)));
600         gain->SetValueAsInt(manuChannel,3,(Int_t)TMath::Nint(GetRandom(kQualMean,kQualSigma,positive)));
601         gain->SetValueAsInt(manuChannel,4,kSaturation);        
602       }
603       
604     }
605     Bool_t ok = gainStore.Add(gain);
606     if (!ok)
607     {
608       AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
609     }
610   }
611   
612   AliInfo(Form("%d Manus and %d channels.",nmanus,nchannels));
613   return nchannels;
614 }
615
616 //_____________________________________________________________________________
617 Int_t
618 AliMUONCDB::MakeLocalTriggerMaskStore(AliMUONVStore& localBoardMasks) const
619 {
620   /// Generate local trigger masks store. All masks are set to FFFF
621   
622   Int_t ngenerated(0);
623   // Generate fake mask values for 234 localboards and put that into
624   // one single container (localBoardMasks)
625   for ( Int_t i = 1; i <= 234; ++i )
626   {
627     AliMUONVCalibParam* localBoard = new AliMUONCalibParamNI(1,8,i,0,0);
628     for ( Int_t x = 0; x < 2; ++x )
629     {
630       for ( Int_t y = 0; y < 4; ++y )
631       {
632         Int_t index = x*4+y;
633         localBoard->SetValueAsInt(index,0,0xFFFF);
634         ++ngenerated;
635       }
636     }
637     localBoardMasks.Add(localBoard);
638   }
639   return ngenerated;
640 }
641
642 //_____________________________________________________________________________
643 Int_t
644 AliMUONCDB::MakeRegionalTriggerMaskStore(AliMUONVStore& rtm) const
645 {
646   /// Make a regional trigger masks store. All masks are set to 3F
647   
648   Int_t ngenerated(0);
649   for ( Int_t i = 0; i < 16; ++i )
650   {
651     AliMUONVCalibParam* regionalBoard = new AliMUONCalibParamNI(1,16,i,0,0);
652     for ( Int_t j = 0; j < 16; ++j )
653     {
654       regionalBoard->SetValueAsInt(j,0,0x3F);
655       ++ngenerated;
656     }
657     rtm.Add(regionalBoard);
658   }
659   
660   return ngenerated;
661 }
662
663 //_____________________________________________________________________________
664 Int_t 
665 AliMUONCDB::MakeGlobalTriggerMaskStore(AliMUONVCalibParam& gtm) const
666 {
667   /// Make a global trigger masks store. All masks set to FFF
668   
669   Int_t ngenerated(0);
670   
671   for ( Int_t j = 0; j < 16; ++j )
672   {
673     gtm.SetValueAsInt(j,0,0xFFF);
674     ++ngenerated;
675   }
676   return ngenerated;
677 }
678
679 //_____________________________________________________________________________
680 AliMUONTriggerLut* 
681 AliMUONCDB::MakeTriggerLUT(const char* file) const
682 {
683   /// Make a triggerlut object, from a file.
684   
685   AliMUONTriggerLut* lut = new AliMUONTriggerLut;
686   lut->ReadFromFile(file);
687   return lut;
688 }
689
690 //_____________________________________________________________________________
691 AliMUONTriggerEfficiencyCells*
692 AliMUONCDB::MakeTriggerEfficiency(const char* file) const
693 {
694   /// Make a trigger efficiency object from a file.
695   
696   return new AliMUONTriggerEfficiencyCells(file);
697 }
698
699 //_____________________________________________________________________________
700 void 
701 AliMUONCDB::WriteToCDB(const char* calibpath, TObject* object, 
702                        Int_t startRun, Int_t endRun, Bool_t defaultValues)
703 {
704   /// Write a given object to OCDB
705   
706   AliCDBId id(calibpath,startRun,endRun);
707   AliCDBMetaData md;
708   md.SetAliRootVersion(gROOT->GetVersion());
709   if ( defaultValues )
710   {
711     md.SetComment("Test with default values");
712   }
713   else
714   {
715     md.SetComment("Test with random values");
716   }
717   md.SetResponsible("AliMUONCDB tester class");
718   
719   AliCDBManager* man = AliCDBManager::Instance();
720   man->SetDefaultStorage(fCDBPath);
721   man->Put(object,id,&md);
722 }
723
724 //_____________________________________________________________________________
725 Int_t 
726 AliMUONCDB::MakeNeighbourStore(AliMUONVStore& neighbourStore)
727 {
728   /// Fill the neighbours store with, for each channel, a TObjArray of its
729   /// neighbouring pads (including itself)
730   
731   AliInfo("Generating NeighbourStore. This will take a while. Please be patient.");
732   
733   TStopwatch timer;
734   
735   timer.Start(kTRUE);
736   
737   TIter next(ManuList());
738   
739   AliMpIntPair* p;
740   
741   Int_t nchannels(0);
742   
743   TObjArray tmp;
744   
745   while ( ( p = (AliMpIntPair*)next() ) )
746   {
747     Int_t detElemId = p->GetFirst();
748     Int_t manuId = p->GetSecond();
749     
750     const AliMpVSegmentation* seg = 
751       AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
752     
753     AliMUONVCalibParam* calibParam = static_cast<AliMUONVCalibParam*>(neighbourStore.FindObject(detElemId,manuId));
754     if (!calibParam)
755     {
756       Int_t dimension(11);
757       Int_t size(AliMpConstants::ManuNofChannels());
758       Int_t defaultValue(-1);
759       Int_t packingFactor(size);
760       
761       calibParam = new AliMUONCalibParamNI(dimension,size,detElemId,manuId,defaultValue,packingFactor);
762       Bool_t ok = neighbourStore.Add(calibParam);
763       if (!ok)
764       {
765         AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
766         return -1;
767       }      
768     }
769     
770     for ( Int_t manuChannel = 0; manuChannel < AliMpConstants::ManuNofChannels(); ++manuChannel )
771     {
772       AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
773       
774       if (pad.IsValid()) 
775       {
776         ++nchannels;
777
778         seg->GetNeighbours(pad,tmp,true,true);
779         Int_t nofPadNeighbours = tmp.GetEntriesFast();
780             
781         for ( Int_t i = 0; i < nofPadNeighbours; ++i )
782         {
783           AliMpPad* pad = static_cast<AliMpPad*>(tmp.At(i));
784           Int_t x;
785           Bool_t ok = calibParam->PackValues(pad->GetLocation().GetFirst(),pad->GetLocation().GetSecond(),x);
786           if (!ok)
787           {
788             AliError("Could not pack value. Something is seriously wrong. Please check");
789             StdoutToAliError(pad->Print(););
790             return -1;
791           }
792           calibParam->SetValueAsInt(manuChannel,i,x);
793         }
794       }
795     }
796     }
797   
798   timer.Print();
799   
800   return nchannels;
801 }
802
803 //_____________________________________________________________________________
804 void
805 AliMUONCDB::WriteLocalTriggerMasks(Int_t startRun, Int_t endRun)
806 {  
807   /// Write local trigger masks to OCDB
808   
809   AliMUONVStore* ltm = new AliMUON1DArray(235);
810   Int_t ngenerated = MakeLocalTriggerMaskStore(*ltm);
811   AliInfo(Form("Ngenerated = %d",ngenerated));
812   if (ngenerated>0)
813   {
814     WriteToCDB("MUON/Calib/LocalTriggerBoardMasks",ltm,startRun,endRun,true);
815   }
816   delete ltm;
817 }
818
819 //_____________________________________________________________________________
820 void
821 AliMUONCDB::WriteRegionalTriggerMasks(Int_t startRun, Int_t endRun)
822 {  
823   /// Write regional trigger masks to OCDB
824   
825   AliMUONVStore* rtm = new AliMUON1DArray(16);
826   Int_t ngenerated = MakeRegionalTriggerMaskStore(*rtm);
827   AliInfo(Form("Ngenerated = %d",ngenerated));
828   if (ngenerated>0)
829   {
830     WriteToCDB("MUON/Calib/RegionalTriggerBoardMasks",rtm,startRun,endRun,true);
831   }
832   delete rtm;
833 }
834
835 //_____________________________________________________________________________
836 void
837 AliMUONCDB::WriteGlobalTriggerMasks(Int_t startRun, Int_t endRun)
838 {  
839   /// Write global trigger masks to OCDB
840   
841   AliMUONVCalibParam* gtm = new AliMUONCalibParamNI(1,16,1,0,0);
842
843   Int_t ngenerated = MakeGlobalTriggerMaskStore(*gtm);
844   AliInfo(Form("Ngenerated = %d",ngenerated));
845   if (ngenerated>0)
846   {
847     WriteToCDB("MUON/Calib/GlobalTriggerBoardMasks",gtm,startRun,endRun,true);
848   }
849   delete gtm;
850 }
851
852 //_____________________________________________________________________________
853 void
854 AliMUONCDB::WriteTriggerLut(Int_t startRun, Int_t endRun)
855 {  
856   /// Write trigger LUT to OCDB
857   
858   AliMUONTriggerLut* lut = MakeTriggerLUT();
859   if (lut)
860   {
861     WriteToCDB("MUON/Calib/TriggerLut",lut,startRun,endRun,true);
862   }
863   delete lut;
864 }
865
866 //_____________________________________________________________________________
867 void
868 AliMUONCDB::WriteTriggerEfficiency(Int_t startRun, Int_t endRun)
869 {  
870   /// Write trigger efficiency to OCDB
871   
872   AliMUONTriggerEfficiencyCells* eff = MakeTriggerEfficiency();
873   if (eff)
874   {
875     WriteToCDB("MUON/Calib/TriggerEfficiency",eff,startRun,endRun,true);
876   }
877   delete eff;
878 }
879
880 //_____________________________________________________________________________
881 void 
882 AliMUONCDB::WriteNeighbours(Int_t startRun, Int_t endRun)
883 {
884   /// Write neighbours to OCDB
885   
886   AliMUONVStore* neighbours = new AliMUON2DMap(kTRUE);
887   Int_t ngenerated = MakeNeighbourStore(*neighbours);
888   AliInfo(Form("Ngenerated = %d",ngenerated));
889   if (ngenerated>0)
890   {
891     WriteToCDB("MUON/Calib/Neighbours",neighbours,startRun,endRun,true);
892   }
893   delete neighbours;
894 }
895
896 //_____________________________________________________________________________
897 void 
898 AliMUONCDB::WriteHV(Bool_t defaultValues,
899                     Int_t startRun, Int_t endRun)
900 {
901   /// generate HV values (either cste = 1500 V) if defaultValues=true or random
902   /// if defaultValues=false, see makeHVStore) and
903   /// store them into CDB located at cdbpath, with a validity period
904   /// ranging from startRun to endRun
905   
906   TMap* hvStore = new TMap;
907   Int_t ngenerated = MakeHVStore(*hvStore,defaultValues);
908   AliInfo(Form("Ngenerated = %d",ngenerated));
909   if (ngenerated>0)
910   {
911     WriteToCDB("MUON/Calib/HV",hvStore,startRun,endRun,defaultValues);
912   }
913   delete hvStore;
914 }
915
916 //_____________________________________________________________________________
917 void 
918 AliMUONCDB::WritePedestals(Bool_t defaultValues,
919                            Int_t startRun, Int_t endRun)
920 {
921   /// generate pedestal values (either 0 if defaultValues=true or random
922   /// if defaultValues=false, see makePedestalStore) and
923   /// store them into CDB located at cdbpath, with a validity period
924   /// ranging from startRun to endRun
925   
926   AliMUONVStore* pedestalStore = new AliMUON2DMap(true);
927   Int_t ngenerated = MakePedestalStore(*pedestalStore,defaultValues);
928   AliInfo(Form("Ngenerated = %d",ngenerated));
929   WriteToCDB("MUON/Calib/Pedestals",pedestalStore,startRun,endRun,defaultValues);
930   delete pedestalStore;
931 }
932
933
934 //_____________________________________________________________________________
935 void 
936 AliMUONCDB::WriteGains(Bool_t defaultValues,
937                        Int_t startRun, Int_t endRun)
938 {
939   /// generate gain values (either 1 if defaultValues=true or random
940   /// if defaultValues=false, see makeGainStore) and
941   /// store them into CDB located at cdbpath, with a validity period
942   /// ranging from startRun to endRun
943   
944   AliMUONVStore* gainStore = new AliMUON2DMap(true);
945   Int_t ngenerated = MakeGainStore(*gainStore,defaultValues);
946   AliInfo(Form("Ngenerated = %d",ngenerated));  
947   WriteToCDB("MUON/Calib/Gains",gainStore,startRun,endRun,defaultValues);
948   delete gainStore;
949 }
950
951 //_____________________________________________________________________________
952 void 
953 AliMUONCDB::WriteCapacitances(Bool_t defaultValues,
954                               Int_t startRun, Int_t endRun)
955 {
956   /// generate manu capacitance values (either 1 if defaultValues=true or random
957   /// if defaultValues=false, see makeCapacitanceStore) and
958   /// store them into CDB located at cdbpath, with a validity period
959   /// ranging from startRun to endRun
960   
961   AliMUONVStore* capaStore = new AliMUON1DMap(16828);
962   Int_t ngenerated = MakeCapacitanceStore(*capaStore,defaultValues);
963   AliInfo(Form("Ngenerated = %d",ngenerated));
964   WriteToCDB("MUON/Calib/Capacitances",capaStore,startRun,endRun,defaultValues);
965   delete capaStore;
966 }
967
968 //_____________________________________________________________________________
969 void
970 AliMUONCDB::WriteTrigger(Int_t startRun, Int_t endRun)
971 {
972   /// Writes all Trigger related calibration to CDB
973   WriteLocalTriggerMasks(startRun,endRun);
974   WriteRegionalTriggerMasks(startRun,endRun);
975   WriteGlobalTriggerMasks(startRun,endRun);
976   WriteTriggerLut(startRun,endRun);
977   WriteTriggerEfficiency(startRun,endRun);
978 }
979
980 //_____________________________________________________________________________
981 void
982 AliMUONCDB::WriteTracker(Bool_t defaultValues, Int_t startRun, Int_t endRun)
983 {
984   /// Writes all Tracker related calibration to CDB
985   WriteHV(defaultValues,startRun,endRun);
986   WritePedestals(defaultValues,startRun,endRun);
987   WriteGains(defaultValues,startRun,endRun);
988   WriteCapacitances(defaultValues,startRun,endRun);
989   WriteNeighbours(startRun,endRun);
990 }
991