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