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