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