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