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