]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONCDB.cxx
Bug fixed in AliMUONCDB::WriteTracker (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     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   ifstream in(gSystem->ExpandPathName(file));
462   if (in.bad()) return 0;
463   
464   Int_t ngenerated(0);
465   
466   char line[1024];
467   Int_t serialNumber(-1);
468   AliMUONVCalibParam* param(0x0);
469   
470   while ( in.getline(line,1024,'\n') )
471   {
472     if ( isdigit(line[0]) ) 
473     {
474       serialNumber = atoi(line);
475       param = static_cast<AliMUONVCalibParam*>(capaStore.FindObject(serialNumber));
476       if (param)
477       {
478         AliError(Form("serialNumber %d appears several times !",serialNumber));
479         capaStore.Clear();
480         break;
481       }
482       param = new AliMUONCalibParamNF(2,AliMpConstants::ManuNofChannels(),serialNumber,0,1.0);
483       Bool_t ok = capaStore.Add(param);
484       if (!ok)
485       {
486         AliError(Form("Could not set serialNumber=%d",serialNumber));
487         continue;
488       }      
489       continue;
490     }
491     Int_t channel;
492     Float_t capaValue;
493     Float_t injectionGain;
494     sscanf(line,"%d %f %f",&channel,&capaValue,&injectionGain);
495     AliDebug(1,Form("SerialNumber %10d Channel %3d Capa %f injectionGain %f",
496                     serialNumber,channel,capaValue,injectionGain));
497     param->SetValueAsFloat(channel,0,capaValue);
498     param->SetValueAsFloat(channel,1,injectionGain);
499     ++ngenerated;
500   }
501   
502   in.close();
503   
504   return ngenerated;
505 }
506
507 //_____________________________________________________________________________
508 Int_t 
509 AliMUONCDB::MakeCapacitanceStore(AliMUONVStore& capaStore, Bool_t defaultValues)
510 {
511   /// Create a capacitance store. if defaultValues=true, all capa are 1.0,
512   /// otherwise they are from a gaussian with parameters defined in the
513   /// kCapa* constants below.
514   
515   TIter next(ManuList());
516   
517   AliMpIntPair* p;
518   
519   Int_t nchannels(0);
520   Int_t nmanus(0);
521   Int_t nmanusOK(0); // manus for which we got the serial number
522     
523   const Float_t kCapaMean(0.3);
524   const Float_t kCapaSigma(0.1);
525   const Float_t kInjectionGainMean(3);
526   const Float_t kInjectionGainSigma(1);
527   
528   while ( ( p = (AliMpIntPair*)next() ) )
529   {
530     ++nmanus;
531     
532     Int_t detElemId = p->GetFirst();
533     Int_t manuId = p->GetSecond();
534     
535     AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(detElemId); 
536     Int_t serialNumber = de->GetManuSerialFromId(manuId);
537       
538     if ( serialNumber <= 0 ) continue;
539     
540     ++nmanusOK;
541     
542     const AliMpVSegmentation* seg = 
543       AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
544
545     AliMUONVCalibParam* capa = static_cast<AliMUONVCalibParam*>(capaStore.FindObject(serialNumber));
546     
547     if (!capa)
548     {
549       capa = new AliMUONCalibParamNF(2,AliMpConstants::ManuNofChannels(),serialNumber,0,1.0);
550       Bool_t ok = capaStore.Add(capa);
551       if (!ok)
552       {
553         AliError(Form("Could not set serialNumber=%d manuId=%d",serialNumber,manuId));
554       }      
555     }
556     
557     for ( Int_t manuChannel = 0; manuChannel < capa->Size(); ++manuChannel )
558     {
559       AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
560       if (!pad.IsValid()) continue;
561       
562       ++nchannels;
563       
564       Float_t capaValue;
565       Float_t injectionGain;
566       
567       if ( defaultValues ) 
568       {
569         capaValue = 1.0;
570         injectionGain = 1.0;
571       }
572       else
573       {
574         capaValue = GetRandom(kCapaMean,kCapaSigma,kTRUE);
575         injectionGain = GetRandom(kInjectionGainMean,kInjectionGainSigma,kTRUE);
576       }
577       capa->SetValueAsFloat(manuChannel,0,capaValue);
578       capa->SetValueAsFloat(manuChannel,1,injectionGain);
579     }
580   }
581   
582   Float_t percent = 0;
583   if ( nmanus ) percent = 100*nmanusOK/nmanus;
584   AliInfo(Form("%5d manus with serial number (out of %5d manus = %3.0f%%)",
585                nmanusOK,nmanus,percent));
586   AliInfo(Form("%5d channels",nchannels));
587   if ( percent < 100 ) 
588   {
589     AliWarning("Did not get all serial numbers. capaStore is incomplete !!!!");
590   }
591   return nchannels;
592   
593 }
594
595 //_____________________________________________________________________________
596 Int_t 
597 AliMUONCDB::MakeGainStore(AliMUONVStore& gainStore, Bool_t defaultValues)
598 {  
599   /// Create a gain store. if defaultValues=true, all gains set so that
600   /// charge = (adc-ped)
601   ///
602   /// otherwise parameters are taken from gaussians with parameters 
603   /// defined in the k* constants below.
604   
605   TIter next(ManuList());
606   
607   AliMpIntPair* p;
608   
609   Int_t nchannels(0);
610   Int_t nmanus(0);
611     
612   const Int_t kSaturation(3000);
613   const Double_t kA0Mean(1.2);
614     const Double_t kA0Sigma(0.1);
615     const Double_t kA1Mean(1E-5);
616       const Double_t kA1Sigma(1E-6);
617   const Double_t kQualMean(0xFF);
618   const Double_t kQualSigma(0x10);
619   const Int_t kThresMean(1600);
620     const Int_t kThresSigma(100);
621   
622   while ( ( p = (AliMpIntPair*)next() ) )
623   {
624     ++nmanus;
625
626     Int_t detElemId = p->GetFirst();
627     Int_t manuId = p->GetSecond();
628
629     AliMUONVCalibParam* gain = 
630       new AliMUONCalibParamNF(5,AliMpConstants::ManuNofChannels(),
631                               detElemId,
632                               manuId,
633                               AliMUONVCalibParam::InvalidFloatValue());
634
635
636     const AliMpVSegmentation* seg = 
637       AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
638
639     for ( Int_t manuChannel = 0; manuChannel < gain->Size(); ++manuChannel )
640     {
641       AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
642       if (!pad.IsValid()) continue;
643       
644       ++nchannels;
645       
646       if ( defaultValues ) 
647       {
648         gain->SetValueAsFloat(manuChannel,0,1.0);
649         gain->SetValueAsFloat(manuChannel,1,0.0);
650         gain->SetValueAsInt(manuChannel,2,4095); 
651         gain->SetValueAsInt(manuChannel,3,1);
652         gain->SetValueAsInt(manuChannel,4,kSaturation);
653       }
654       else
655       {
656         Bool_t positive(kTRUE);
657         gain->SetValueAsFloat(manuChannel,0,GetRandom(kA0Mean,kA0Sigma,positive));
658         gain->SetValueAsFloat(manuChannel,1,GetRandom(kA1Mean,kA1Sigma,!positive));
659         gain->SetValueAsInt(manuChannel,2,(Int_t)TMath::Nint(GetRandom(kThresMean,kThresSigma,positive)));
660         gain->SetValueAsInt(manuChannel,3,(Int_t)TMath::Nint(GetRandom(kQualMean,kQualSigma,positive)));
661         gain->SetValueAsInt(manuChannel,4,kSaturation);        
662       }
663       
664     }
665     Bool_t ok = gainStore.Add(gain);
666     if (!ok)
667     {
668       AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
669     }
670     if ( fMaxNofChannelsToGenerate > 0 && nchannels >= fMaxNofChannelsToGenerate ) break;
671   }
672   
673   AliInfo(Form("%d Manus and %d channels.",nmanus,nchannels));
674   return nchannels;
675 }
676
677 //_____________________________________________________________________________
678 Int_t
679 AliMUONCDB::MakeLocalTriggerMaskStore(AliMUONVStore& localBoardMasks) const
680 {
681   /// Generate local trigger masks store. All masks are set to FFFF
682   
683   Int_t ngenerated(0);
684   // Generate fake mask values for 234 localboards and put that into
685   // one single container (localBoardMasks)
686   for ( Int_t i = 1; i <= 234; ++i )
687   {
688     AliMUONVCalibParam* localBoard = new AliMUONCalibParamNI(1,8,i,0,0);
689     for ( Int_t x = 0; x < 2; ++x )
690     {
691       for ( Int_t y = 0; y < 4; ++y )
692       {
693         Int_t index = x*4+y;
694         localBoard->SetValueAsInt(index,0,0xFFFF);
695         ++ngenerated;
696       }
697     }
698     localBoardMasks.Add(localBoard);
699   }
700   return ngenerated;
701 }
702
703 //_____________________________________________________________________________
704 Int_t
705 AliMUONCDB::MakeRegionalTriggerMaskStore(AliMUONVStore& rtm) const
706 {
707   /// Make a regional trigger masks store. All masks are set to 3F
708   
709   Int_t ngenerated(0);
710   for ( Int_t i = 0; i < 16; ++i )
711   {
712     AliMUONVCalibParam* regionalBoard = new AliMUONCalibParamNI(1,16,i,0,0);
713     for ( Int_t j = 0; j < 16; ++j )
714     {
715       regionalBoard->SetValueAsInt(j,0,0x3F);
716       ++ngenerated;
717     }
718     rtm.Add(regionalBoard);
719   }
720   
721   return ngenerated;
722 }
723
724 //_____________________________________________________________________________
725 Int_t 
726 AliMUONCDB::MakeGlobalTriggerMaskStore(AliMUONVCalibParam& gtm) const
727 {
728   /// Make a global trigger masks store. All masks set to FFF
729   
730   Int_t ngenerated(0);
731   
732   for ( Int_t j = 0; j < 16; ++j )
733   {
734     gtm.SetValueAsInt(j,0,0xFFF);
735     ++ngenerated;
736   }
737   return ngenerated;
738 }
739
740 //_____________________________________________________________________________
741 AliMUONTriggerLut* 
742 AliMUONCDB::MakeTriggerLUT(const char* file) const
743 {
744   /// Make a triggerlut object, from a file.
745   
746   AliMUONTriggerLut* lut = new AliMUONTriggerLut;
747   lut->ReadFromFile(file);
748   return lut;
749 }
750
751 //_____________________________________________________________________________
752 AliMUONTriggerEfficiencyCells*
753 AliMUONCDB::MakeTriggerEfficiency(const char* file) const
754 {
755   /// Make a trigger efficiency object from a file.
756   
757   return new AliMUONTriggerEfficiencyCells(file);
758 }
759
760 //_____________________________________________________________________________
761 void 
762 AliMUONCDB::WriteToCDB(const char* calibpath, TObject* object, 
763                        Int_t startRun, Int_t endRun, 
764                        const char* filename)
765 {
766   /// Write a given object to OCDB
767   
768   AliCDBId id(calibpath,startRun,endRun);
769   AliCDBMetaData md;
770   md.SetAliRootVersion(gROOT->GetVersion());
771   md.SetComment(gSystem->ExpandPathName(filename));
772   md.SetResponsible("Uploaded using AliMUONCDB class");
773   AliCDBManager* man = AliCDBManager::Instance();
774   man->SetDefaultStorage(fCDBPath);
775   man->Put(object,id,&md);
776 }
777
778 //_____________________________________________________________________________
779 void 
780 AliMUONCDB::WriteToCDB(const char* calibpath, TObject* object, 
781                        Int_t startRun, Int_t endRun, Bool_t defaultValues)
782 {
783   /// Write a given object to OCDB
784   
785   AliCDBId id(calibpath,startRun,endRun);
786   AliCDBMetaData md;
787   md.SetAliRootVersion(gROOT->GetVersion());
788   if ( defaultValues )
789   {
790     md.SetComment("Test with default values");
791   }
792   else
793   {
794     md.SetComment("Test with random values");
795   }
796   md.SetResponsible("AliMUONCDB tester class");
797   
798   AliCDBManager* man = AliCDBManager::Instance();
799   man->SetDefaultStorage(fCDBPath);
800   man->Put(object,id,&md);
801 }
802
803 //_____________________________________________________________________________
804 Int_t 
805 AliMUONCDB::MakeNeighbourStore(AliMUONVStore& neighbourStore)
806 {
807   /// Fill the neighbours store with, for each channel, a TObjArray of its
808   /// neighbouring pads (including itself)
809   
810   AliInfo("Generating NeighbourStore. This will take a while. Please be patient.");
811   
812   TStopwatch timer;
813   
814   timer.Start(kTRUE);
815   
816   TIter next(ManuList());
817   
818   AliMpIntPair* p;
819   
820   Int_t nchannels(0);
821   
822   TObjArray tmp;
823   
824   while ( ( p = (AliMpIntPair*)next() ) )
825   {
826     Int_t detElemId = p->GetFirst();
827     Int_t manuId = p->GetSecond();
828     
829     const AliMpVSegmentation* seg = 
830       AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
831     
832     AliMUONVCalibParam* calibParam = static_cast<AliMUONVCalibParam*>(neighbourStore.FindObject(detElemId,manuId));
833     if (!calibParam)
834     {
835       Int_t dimension(11);
836       Int_t size(AliMpConstants::ManuNofChannels());
837       Int_t defaultValue(-1);
838       Int_t packingFactor(size);
839       
840       calibParam = new AliMUONCalibParamNI(dimension,size,detElemId,manuId,defaultValue,packingFactor);
841       Bool_t ok = neighbourStore.Add(calibParam);
842       if (!ok)
843       {
844         AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
845         return -1;
846       }      
847     }
848     
849     for ( Int_t manuChannel = 0; manuChannel < AliMpConstants::ManuNofChannels(); ++manuChannel )
850     {
851       AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
852       
853       if (pad.IsValid()) 
854       {
855         ++nchannels;
856
857         seg->GetNeighbours(pad,tmp,true,true);
858         Int_t nofPadNeighbours = tmp.GetEntriesFast();
859             
860         for ( Int_t i = 0; i < nofPadNeighbours; ++i )
861         {
862           AliMpPad* pad = static_cast<AliMpPad*>(tmp.At(i));
863           Int_t x;
864           Bool_t ok = calibParam->PackValues(pad->GetLocation().GetFirst(),pad->GetLocation().GetSecond(),x);
865           if (!ok)
866           {
867             AliError("Could not pack value. Something is seriously wrong. Please check");
868             StdoutToAliError(pad->Print(););
869             return -1;
870           }
871           calibParam->SetValueAsInt(manuChannel,i,x);
872         }
873       }
874     }
875     }
876   
877   timer.Print();
878   
879   return nchannels;
880 }
881
882 //_____________________________________________________________________________
883 void
884 AliMUONCDB::SetMaxNofChannelsToGenerate(Int_t n)
885 {
886   /// Set the maximum number of channels to generate (used for testing only)
887   /// n < 0 means no limit
888   fMaxNofChannelsToGenerate = n;
889 }
890
891 //_____________________________________________________________________________
892 void
893 AliMUONCDB::WriteLocalTriggerMasks(Int_t startRun, Int_t endRun)
894 {  
895   /// Write local trigger masks to OCDB
896   
897   AliMUONVStore* ltm = new AliMUON1DArray(235);
898   Int_t ngenerated = MakeLocalTriggerMaskStore(*ltm);
899   AliInfo(Form("Ngenerated = %d",ngenerated));
900   if (ngenerated>0)
901   {
902     WriteToCDB("MUON/Calib/LocalTriggerBoardMasks",ltm,startRun,endRun,true);
903   }
904   delete ltm;
905 }
906
907 //_____________________________________________________________________________
908 void
909 AliMUONCDB::WriteRegionalTriggerMasks(Int_t startRun, Int_t endRun)
910 {  
911   /// Write regional trigger masks to OCDB
912   
913   AliMUONVStore* rtm = new AliMUON1DArray(16);
914   Int_t ngenerated = MakeRegionalTriggerMaskStore(*rtm);
915   AliInfo(Form("Ngenerated = %d",ngenerated));
916   if (ngenerated>0)
917   {
918     WriteToCDB("MUON/Calib/RegionalTriggerBoardMasks",rtm,startRun,endRun,true);
919   }
920   delete rtm;
921 }
922
923 //_____________________________________________________________________________
924 void
925 AliMUONCDB::WriteGlobalTriggerMasks(Int_t startRun, Int_t endRun)
926 {  
927   /// Write global trigger masks to OCDB
928   
929   AliMUONVCalibParam* gtm = new AliMUONCalibParamNI(1,16,1,0,0);
930
931   Int_t ngenerated = MakeGlobalTriggerMaskStore(*gtm);
932   AliInfo(Form("Ngenerated = %d",ngenerated));
933   if (ngenerated>0)
934   {
935     WriteToCDB("MUON/Calib/GlobalTriggerBoardMasks",gtm,startRun,endRun,true);
936   }
937   delete gtm;
938 }
939
940 //_____________________________________________________________________________
941 void
942 AliMUONCDB::WriteTriggerLut(Int_t startRun, Int_t endRun)
943 {  
944   /// Write trigger LUT to OCDB
945   
946   AliMUONTriggerLut* lut = MakeTriggerLUT();
947   if (lut)
948   {
949     WriteToCDB("MUON/Calib/TriggerLut",lut,startRun,endRun,true);
950   }
951   delete lut;
952 }
953
954 //_____________________________________________________________________________
955 void
956 AliMUONCDB::WriteTriggerEfficiency(Int_t startRun, Int_t endRun)
957 {  
958   /// Write trigger efficiency to OCDB
959   
960   AliMUONTriggerEfficiencyCells* eff = MakeTriggerEfficiency();
961   if (eff)
962   {
963     WriteToCDB("MUON/Calib/TriggerEfficiency",eff,startRun,endRun,true);
964   }
965   delete eff;
966 }
967
968 //_____________________________________________________________________________
969 void 
970 AliMUONCDB::WriteNeighbours(Int_t startRun, Int_t endRun)
971 {
972   /// Write neighbours to OCDB
973   
974   AliMUONVStore* neighbours = new AliMUON2DMap(kTRUE);
975   Int_t ngenerated = MakeNeighbourStore(*neighbours);
976   AliInfo(Form("Ngenerated = %d",ngenerated));
977   if (ngenerated>0)
978   {
979     WriteToCDB("MUON/Calib/Neighbours",neighbours,startRun,endRun,true);
980   }
981   delete neighbours;
982 }
983
984 //_____________________________________________________________________________
985 void 
986 AliMUONCDB::WriteHV(Bool_t defaultValues,
987                     Int_t startRun, Int_t endRun)
988 {
989   /// generate HV values (either cste = 1500 V) if defaultValues=true or random
990   /// if defaultValues=false, see makeHVStore) and
991   /// store them into CDB located at cdbpath, with a validity period
992   /// ranging from startRun to endRun
993   
994   TMap* hvStore = new TMap;
995   Int_t ngenerated = MakeHVStore(*hvStore,defaultValues);
996   AliInfo(Form("Ngenerated = %d",ngenerated));
997   if (ngenerated>0)
998   {
999     WriteToCDB("MUON/Calib/HV",hvStore,startRun,endRun,defaultValues);
1000   }
1001   delete hvStore;
1002 }
1003
1004 //_____________________________________________________________________________
1005 void 
1006 AliMUONCDB::WritePedestals(Bool_t defaultValues,
1007                            Int_t startRun, Int_t endRun)
1008 {
1009   /// generate pedestal values (either 0 if defaultValues=true or random
1010   /// if defaultValues=false, see makePedestalStore) and
1011   /// store them into CDB located at cdbpath, with a validity period
1012   /// ranging from startRun to endRun
1013   
1014   AliMUONVStore* pedestalStore = new AliMUON2DMap(true);
1015   Int_t ngenerated = MakePedestalStore(*pedestalStore,defaultValues);
1016   AliInfo(Form("Ngenerated = %d",ngenerated));
1017   WriteToCDB("MUON/Calib/Pedestals",pedestalStore,startRun,endRun,defaultValues);
1018   delete pedestalStore;
1019 }
1020
1021
1022 //_____________________________________________________________________________
1023 void 
1024 AliMUONCDB::WriteGains(Bool_t defaultValues,
1025                        Int_t startRun, Int_t endRun)
1026 {
1027   /// generate gain values (either 1 if defaultValues=true or random
1028   /// if defaultValues=false, see makeGainStore) and
1029   /// store them into CDB located at cdbpath, with a validity period
1030   /// ranging from startRun to endRun
1031   
1032   AliMUONVStore* gainStore = new AliMUON2DMap(true);
1033   Int_t ngenerated = MakeGainStore(*gainStore,defaultValues);
1034   AliInfo(Form("Ngenerated = %d",ngenerated));  
1035   WriteToCDB("MUON/Calib/Gains",gainStore,startRun,endRun,defaultValues);
1036   delete gainStore;
1037 }
1038
1039 //_____________________________________________________________________________
1040 void 
1041 AliMUONCDB::WriteCapacitances(const char* filename,
1042                               Int_t startRun, Int_t endRun)
1043 {
1044   /// read manu capacitance and injection gain values from file 
1045   /// and store them into CDB located at cdbpath, with a validity period
1046   /// ranging from startRun to endRun
1047   
1048   AliMUONVStore* capaStore = new AliMUON1DMap(16828);
1049   Int_t ngenerated = MakeCapacitanceStore(*capaStore,filename);
1050   AliInfo(Form("Ngenerated = %d",ngenerated));
1051   if ( ngenerated > 0 ) 
1052   {
1053     WriteToCDB("MUON/Calib/Capacitances",capaStore,startRun,endRun,filename);
1054   }
1055   delete capaStore;
1056 }
1057
1058 //_____________________________________________________________________________
1059 void 
1060 AliMUONCDB::WriteCapacitances(Bool_t defaultValues,
1061                               Int_t startRun, Int_t endRun)
1062 {
1063   /// generate manu capacitance values (either 1 if defaultValues=true or random
1064   /// if defaultValues=false, see makeCapacitanceStore) and
1065   /// store them into CDB located at cdbpath, with a validity period
1066   /// ranging from startRun to endRun
1067   
1068   AliMUONVStore* capaStore = new AliMUON1DMap(16828);
1069   Int_t ngenerated = MakeCapacitanceStore(*capaStore,defaultValues);
1070   AliInfo(Form("Ngenerated = %d",ngenerated));
1071   WriteToCDB("MUON/Calib/Capacitances",capaStore,startRun,endRun,defaultValues);
1072   delete capaStore;
1073 }
1074
1075 //_____________________________________________________________________________
1076 void
1077 AliMUONCDB::WriteTrigger(Int_t startRun, Int_t endRun)
1078 {
1079   /// Writes all Trigger related calibration to CDB
1080   WriteLocalTriggerMasks(startRun,endRun);
1081   WriteRegionalTriggerMasks(startRun,endRun);
1082   WriteGlobalTriggerMasks(startRun,endRun);
1083   WriteTriggerLut(startRun,endRun);
1084   WriteTriggerEfficiency(startRun,endRun);
1085 }
1086
1087 //_____________________________________________________________________________
1088 void
1089 AliMUONCDB::WriteTracker(Bool_t defaultValues, Int_t startRun, Int_t endRun)
1090 {
1091   /// Writes all Tracker related calibration to CDB
1092   WriteHV(defaultValues,startRun,endRun);
1093   WritePedestals(defaultValues,startRun,endRun);
1094   WriteGains(defaultValues,startRun,endRun);
1095   WriteCapacitances(defaultValues,startRun,endRun);
1096   WriteNeighbours(startRun,endRun);
1097 }