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