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