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