1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //-----------------------------------------------------------------------------
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
26 /// For more information, please see READMECDB
28 // \author Laurent Aphecetche
29 //-----------------------------------------------------------------------------
31 #include "AliMUONCDB.h"
33 #include "AliCDBEntry.h"
34 #include "AliCDBManager.h"
35 #include "AliDCSValue.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>
65 #include <TObjString.h>
68 #include <TStopwatch.h>
78 //_____________________________________________________________________________
79 void getBoundaries(const AliMUONVStore& store, Int_t dim,
80 Float_t* xmin, Float_t* xmax)
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
86 for ( Int_t i = 0; i < dim; ++i )
92 TIter next(store.CreateIterator());
93 AliMUONVCalibParam* value;
95 while ( ( value = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
97 Int_t detElemId = value->ID0();
98 Int_t manuId = value->ID1();
100 const AliMpVSegmentation* seg =
101 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
103 for ( Int_t manuChannel = 0; manuChannel < value->Size(); ++manuChannel )
105 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
106 if (!pad.IsValid()) continue;
108 for ( Int_t i = 0; i < dim; ++i )
110 Float_t x0 = value->ValueAsFloat(manuChannel,i);
112 xmin[i] = TMath::Min(xmin[i],x0);
113 xmax[i] = TMath::Max(xmax[i],x0);
118 for ( Int_t i = 0; i < dim; ++i )
120 if ( TMath::Abs(xmin[i]-xmax[i]) < 1E-3 )
128 //_____________________________________________________________________________
129 Double_t GetRandom(Double_t mean, Double_t sigma, Bool_t mustBePositive)
132 if ( mustBePositive )
136 x = gRandom->Gaus(mean,sigma);
141 x = gRandom->Gaus(mean,sigma);
148 //_____________________________________________________________________________
149 AliMUONCDB::AliMUONCDB(const char* cdbpath)
153 fMaxNofChannelsToGenerate(-1)
158 //_____________________________________________________________________________
159 AliMUONCDB::~AliMUONCDB()
165 //_____________________________________________________________________________
167 AliMUONCDB::ManuList()
169 /// return (and create if necessary) the list of (de,manu) pairs
172 AliInfo("Generating ManuList...");
173 fManuList = AliMpManuList::ManuList();
174 AliInfo("Manu List generated.");
179 //_____________________________________________________________________________
181 AliMUONCDB::Diff(AliMUONVStore& store1, AliMUONVStore& store2,
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
189 /// WARNING Works only for stores which holds AliMUONVCalibParam objects
194 if ( !sopt.Contains("ABS") && !sopt.Contains("REL") && !sopt.Contains("PERCENT") )
196 AliErrorClass(Form("opt %s not supported. Only ABS, REL, PERCENT are",opt));
200 AliMUONVStore* d = static_cast<AliMUONVStore*>(store1.Clone());
202 TIter next(d->CreateIterator());
204 AliMUONVCalibParam* param;
206 while ( ( param = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
208 Int_t detElemId = param->ID0();
209 Int_t manuId = param->ID1();
211 AliMUONVCalibParam* param2 = dynamic_cast<AliMUONVCalibParam*>(store2.FindObject(detElemId,manuId));
212 //FIXME: this might happen. Handle it.
215 cerr << "param2 is null : FIXME : this might happen !" << endl;
220 for ( Int_t i = 0; i < param->Size(); ++i )
222 for ( Int_t j = 0; j < param->Dimension(); ++j )
225 if ( sopt.Contains("ABS") )
227 value = param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j);
229 else if ( sopt.Contains("REL") || sopt.Contains("PERCENT") )
231 if ( param->ValueAsFloat(i,j) )
233 value = (param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j))/param->ValueAsFloat(i,j);
239 if ( sopt.Contains("PERCENT") ) value *= 100.0;
241 param->SetValueAsFloat(i,j,value);
248 //_____________________________________________________________________________
250 AliMUONCDB::Plot(const AliMUONVStore& store, const char* name, Int_t nbins)
252 /// Make histograms of each dimension of the AliMUONVCalibParam
253 /// contained inside store.
254 /// It produces histograms named name_0, name_1, etc...
256 TIter next(store.CreateIterator());
257 AliMUONVCalibParam* param;
259 const Int_t kNStations = AliMpConstants::NofTrackingChambers()/2;
260 Int_t* nPerStation = new Int_t[kNStations];
263 for ( Int_t i = 0; i < kNStations; ++i ) nPerStation[i]=0;
265 while ( ( param = static_cast<AliMUONVCalibParam*>(next()) ) )
269 Int_t dim = param->Dimension();
271 Float_t* xmin = new Float_t[dim];
272 Float_t* xmax = new Float_t[dim];
273 getBoundaries(store,dim,xmin,xmax);
275 for ( Int_t i = 0; i < dim; ++i )
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()));
283 Int_t detElemId = param->ID0();
284 Int_t manuId = param->ID1();
285 Int_t station = AliMpDEManager::GetChamberId(detElemId)/2;
287 const AliMpVSegmentation* seg =
288 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
290 for ( Int_t manuChannel = 0; manuChannel < param->Size(); ++manuChannel )
292 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
293 if (!pad.IsValid()) continue;
296 ++nPerStation[station];
298 for ( Int_t dim = 0; dim < param->Dimension(); ++dim )
300 h[dim]->Fill(param->ValueAsFloat(manuChannel,dim));
305 for ( Int_t i = 0; i < kNStations; ++i )
307 AliInfo(Form("Station %d %d ",(i+1),nPerStation[i]));
310 AliInfo(Form("Number of channels = %d",n));
312 delete[] nPerStation;
315 //_____________________________________________________________________________
317 AliMUONCDB::MakeHVStore(TMap& aliasMap, Bool_t defaultValues)
319 /// Create a HV store
321 AliMUONHVNamer hvNamer;
323 TObjArray* aliases = hvNamer.GenerateAliases();
328 for ( Int_t i = 0; i < aliases->GetEntries(); ++i )
330 TObjString* alias = static_cast<TObjString*>(aliases->At(i));
331 TString& aliasName = alias->String();
332 if ( aliasName.Contains("sw") )
334 // HV Switch (St345 only)
335 TObjArray* valueSet = new TObjArray;
336 valueSet->SetOwner(kTRUE);
338 Bool_t value = kTRUE;
342 Float_t r = gRandom->Uniform();
343 if ( r < 0.007 ) value = kFALSE;
346 for ( UInt_t timeStamp = 0; timeStamp < 60*3; timeStamp += 60 )
348 AliDCSValue* dcsValue = new AliDCSValue(value,timeStamp);
349 valueSet->Add(dcsValue);
351 aliasMap.Add(new TObjString(*alias),valueSet);
356 TObjArray* valueSet = new TObjArray;
357 valueSet->SetOwner(kTRUE);
358 for ( UInt_t timeStamp = 0; timeStamp < 60*15; timeStamp += 120 )
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);
365 aliasMap.Add(new TObjString(*alias),valueSet);
372 AliInfo(Form("%d HV channels and %d switches",nChannels,nSwitch));
374 return nChannels+nSwitch;
377 //_____________________________________________________________________________
379 AliMUONCDB::MakePedestalStore(AliMUONVStore& pedestalStore, Bool_t defaultValues)
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)
385 TIter next(ManuList());
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);
398 while ( ( p = (AliMpIntPair*)next() ) )
402 Int_t detElemId = p->GetFirst();
403 Int_t manuId = p->GetSecond();
405 AliMUONVCalibParam* ped =
406 new AliMUONCalibParamNF(2,kChannels,detElemId,manuId,AliMUONVCalibParam::InvalidFloatValue());
408 const AliMpVSegmentation* seg =
409 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
411 for ( Int_t manuChannel = 0; manuChannel < kChannels; ++manuChannel )
413 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
414 if (!pad.IsValid()) continue;
418 Float_t meanPedestal;
419 Float_t sigmaPedestal;
428 Bool_t positive(kTRUE);
430 while ( meanPedestal == 0.0 ) // avoid strict zero
432 meanPedestal = GetRandom(kPedestalMeanMean,kPedestalMeanSigma,positive);
434 sigmaPedestal = GetRandom(kPedestalSigmaMean,kPedestalSigmaSigma,positive);
436 ped->SetValueAsFloat(manuChannel,0,meanPedestal);
437 ped->SetValueAsFloat(manuChannel,1,sigmaPedestal);
440 Bool_t ok = pedestalStore.Add(ped);
443 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
445 if ( fMaxNofChannelsToGenerate > 0 && nchannels >= fMaxNofChannelsToGenerate ) break;
448 AliInfo(Form("%d Manus and %d channels.",nmanus,nchannels));
453 //_____________________________________________________________________________
455 AliMUONCDB::MakeCapacitanceStore(AliMUONVStore& capaStore, Bool_t defaultValues)
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.
461 TIter next(ManuList());
467 Int_t nmanusOK(0); // manus for which we got the serial number
469 const Float_t kCapaMean(1.0);
470 const Float_t kCapaSigma(0.5);
472 while ( ( p = (AliMpIntPair*)next() ) )
476 Int_t detElemId = p->GetFirst();
477 Int_t manuId = p->GetSecond();
479 AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(detElemId);
480 Int_t serialNumber = de->GetManuSerialFromId(manuId);
482 if ( serialNumber <= 0 ) continue;
486 const AliMpVSegmentation* seg =
487 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
489 AliMUONVCalibParam* capa = static_cast<AliMUONVCalibParam*>(capaStore.FindObject(serialNumber));
493 capa = new AliMUONCalibParamNF(1,AliMpConstants::ManuNofChannels(),serialNumber,0,1.0);
494 Bool_t ok = capaStore.Add(capa);
497 AliError(Form("Could not set serialNumber=%d manuId=%d",serialNumber,manuId));
501 for ( Int_t manuChannel = 0; manuChannel < capa->Size(); ++manuChannel )
503 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
504 if (!pad.IsValid()) continue;
516 capaValue = GetRandom(kCapaMean,kCapaSigma,kTRUE);
518 capa->SetValueAsFloat(manuChannel,0,capaValue);
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));
529 AliWarning("Did not get all serial numbers. capaStore is incomplete !!!!");
535 //_____________________________________________________________________________
537 AliMUONCDB::MakeGainStore(AliMUONVStore& gainStore, Bool_t defaultValues)
539 /// Create a gain store. if defaultValues=true, all gains set so that
540 /// charge = (adc-ped)
542 /// otherwise parameters are taken from gaussians with parameters
543 /// defined in the k* constants below.
545 TIter next(ManuList());
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);
562 while ( ( p = (AliMpIntPair*)next() ) )
566 Int_t detElemId = p->GetFirst();
567 Int_t manuId = p->GetSecond();
569 AliMUONVCalibParam* gain =
570 new AliMUONCalibParamNF(5,AliMpConstants::ManuNofChannels(),
573 AliMUONVCalibParam::InvalidFloatValue());
576 const AliMpVSegmentation* seg =
577 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
579 for ( Int_t manuChannel = 0; manuChannel < gain->Size(); ++manuChannel )
581 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
582 if (!pad.IsValid()) continue;
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);
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);
605 Bool_t ok = gainStore.Add(gain);
608 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
610 if ( fMaxNofChannelsToGenerate > 0 && nchannels >= fMaxNofChannelsToGenerate ) break;
613 AliInfo(Form("%d Manus and %d channels.",nmanus,nchannels));
617 //_____________________________________________________________________________
619 AliMUONCDB::MakeLocalTriggerMaskStore(AliMUONVStore& localBoardMasks) const
621 /// Generate local trigger masks store. All masks are set to FFFF
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 )
628 AliMUONVCalibParam* localBoard = new AliMUONCalibParamNI(1,8,i,0,0);
629 for ( Int_t x = 0; x < 2; ++x )
631 for ( Int_t y = 0; y < 4; ++y )
634 localBoard->SetValueAsInt(index,0,0xFFFF);
638 localBoardMasks.Add(localBoard);
643 //_____________________________________________________________________________
645 AliMUONCDB::MakeRegionalTriggerMaskStore(AliMUONVStore& rtm) const
647 /// Make a regional trigger masks store. All masks are set to 3F
650 for ( Int_t i = 0; i < 16; ++i )
652 AliMUONVCalibParam* regionalBoard = new AliMUONCalibParamNI(1,16,i,0,0);
653 for ( Int_t j = 0; j < 16; ++j )
655 regionalBoard->SetValueAsInt(j,0,0x3F);
658 rtm.Add(regionalBoard);
664 //_____________________________________________________________________________
666 AliMUONCDB::MakeGlobalTriggerMaskStore(AliMUONVCalibParam& gtm) const
668 /// Make a global trigger masks store. All masks set to FFF
672 for ( Int_t j = 0; j < 16; ++j )
674 gtm.SetValueAsInt(j,0,0xFFF);
680 //_____________________________________________________________________________
682 AliMUONCDB::MakeTriggerLUT(const char* file) const
684 /// Make a triggerlut object, from a file.
686 AliMUONTriggerLut* lut = new AliMUONTriggerLut;
687 lut->ReadFromFile(file);
691 //_____________________________________________________________________________
692 AliMUONTriggerEfficiencyCells*
693 AliMUONCDB::MakeTriggerEfficiency(const char* file) const
695 /// Make a trigger efficiency object from a file.
697 return new AliMUONTriggerEfficiencyCells(file);
700 //_____________________________________________________________________________
702 AliMUONCDB::WriteToCDB(const char* calibpath, TObject* object,
703 Int_t startRun, Int_t endRun, Bool_t defaultValues)
705 /// Write a given object to OCDB
707 AliCDBId id(calibpath,startRun,endRun);
709 md.SetAliRootVersion(gROOT->GetVersion());
712 md.SetComment("Test with default values");
716 md.SetComment("Test with random values");
718 md.SetResponsible("AliMUONCDB tester class");
720 AliCDBManager* man = AliCDBManager::Instance();
721 man->SetDefaultStorage(fCDBPath);
722 man->Put(object,id,&md);
725 //_____________________________________________________________________________
727 AliMUONCDB::MakeNeighbourStore(AliMUONVStore& neighbourStore)
729 /// Fill the neighbours store with, for each channel, a TObjArray of its
730 /// neighbouring pads (including itself)
732 AliInfo("Generating NeighbourStore. This will take a while. Please be patient.");
738 TIter next(ManuList());
746 while ( ( p = (AliMpIntPair*)next() ) )
748 Int_t detElemId = p->GetFirst();
749 Int_t manuId = p->GetSecond();
751 const AliMpVSegmentation* seg =
752 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
754 AliMUONVCalibParam* calibParam = static_cast<AliMUONVCalibParam*>(neighbourStore.FindObject(detElemId,manuId));
758 Int_t size(AliMpConstants::ManuNofChannels());
759 Int_t defaultValue(-1);
760 Int_t packingFactor(size);
762 calibParam = new AliMUONCalibParamNI(dimension,size,detElemId,manuId,defaultValue,packingFactor);
763 Bool_t ok = neighbourStore.Add(calibParam);
766 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
771 for ( Int_t manuChannel = 0; manuChannel < AliMpConstants::ManuNofChannels(); ++manuChannel )
773 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
779 seg->GetNeighbours(pad,tmp,true,true);
780 Int_t nofPadNeighbours = tmp.GetEntriesFast();
782 for ( Int_t i = 0; i < nofPadNeighbours; ++i )
784 AliMpPad* pad = static_cast<AliMpPad*>(tmp.At(i));
786 Bool_t ok = calibParam->PackValues(pad->GetLocation().GetFirst(),pad->GetLocation().GetSecond(),x);
789 AliError("Could not pack value. Something is seriously wrong. Please check");
790 StdoutToAliError(pad->Print(););
793 calibParam->SetValueAsInt(manuChannel,i,x);
804 //_____________________________________________________________________________
806 AliMUONCDB::SetMaxNofChannelsToGenerate(Int_t n)
808 /// Set the maximum number of channels to generate (used for testing only)
809 /// n < 0 means no limit
810 fMaxNofChannelsToGenerate = n;
813 //_____________________________________________________________________________
815 AliMUONCDB::WriteLocalTriggerMasks(Int_t startRun, Int_t endRun)
817 /// Write local trigger masks to OCDB
819 AliMUONVStore* ltm = new AliMUON1DArray(235);
820 Int_t ngenerated = MakeLocalTriggerMaskStore(*ltm);
821 AliInfo(Form("Ngenerated = %d",ngenerated));
824 WriteToCDB("MUON/Calib/LocalTriggerBoardMasks",ltm,startRun,endRun,true);
829 //_____________________________________________________________________________
831 AliMUONCDB::WriteRegionalTriggerMasks(Int_t startRun, Int_t endRun)
833 /// Write regional trigger masks to OCDB
835 AliMUONVStore* rtm = new AliMUON1DArray(16);
836 Int_t ngenerated = MakeRegionalTriggerMaskStore(*rtm);
837 AliInfo(Form("Ngenerated = %d",ngenerated));
840 WriteToCDB("MUON/Calib/RegionalTriggerBoardMasks",rtm,startRun,endRun,true);
845 //_____________________________________________________________________________
847 AliMUONCDB::WriteGlobalTriggerMasks(Int_t startRun, Int_t endRun)
849 /// Write global trigger masks to OCDB
851 AliMUONVCalibParam* gtm = new AliMUONCalibParamNI(1,16,1,0,0);
853 Int_t ngenerated = MakeGlobalTriggerMaskStore(*gtm);
854 AliInfo(Form("Ngenerated = %d",ngenerated));
857 WriteToCDB("MUON/Calib/GlobalTriggerBoardMasks",gtm,startRun,endRun,true);
862 //_____________________________________________________________________________
864 AliMUONCDB::WriteTriggerLut(Int_t startRun, Int_t endRun)
866 /// Write trigger LUT to OCDB
868 AliMUONTriggerLut* lut = MakeTriggerLUT();
871 WriteToCDB("MUON/Calib/TriggerLut",lut,startRun,endRun,true);
876 //_____________________________________________________________________________
878 AliMUONCDB::WriteTriggerEfficiency(Int_t startRun, Int_t endRun)
880 /// Write trigger efficiency to OCDB
882 AliMUONTriggerEfficiencyCells* eff = MakeTriggerEfficiency();
885 WriteToCDB("MUON/Calib/TriggerEfficiency",eff,startRun,endRun,true);
890 //_____________________________________________________________________________
892 AliMUONCDB::WriteNeighbours(Int_t startRun, Int_t endRun)
894 /// Write neighbours to OCDB
896 AliMUONVStore* neighbours = new AliMUON2DMap(kTRUE);
897 Int_t ngenerated = MakeNeighbourStore(*neighbours);
898 AliInfo(Form("Ngenerated = %d",ngenerated));
901 WriteToCDB("MUON/Calib/Neighbours",neighbours,startRun,endRun,true);
906 //_____________________________________________________________________________
908 AliMUONCDB::WriteHV(Bool_t defaultValues,
909 Int_t startRun, Int_t endRun)
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
916 TMap* hvStore = new TMap;
917 Int_t ngenerated = MakeHVStore(*hvStore,defaultValues);
918 AliInfo(Form("Ngenerated = %d",ngenerated));
921 WriteToCDB("MUON/Calib/HV",hvStore,startRun,endRun,defaultValues);
926 //_____________________________________________________________________________
928 AliMUONCDB::WritePedestals(Bool_t defaultValues,
929 Int_t startRun, Int_t endRun)
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
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;
944 //_____________________________________________________________________________
946 AliMUONCDB::WriteGains(Bool_t defaultValues,
947 Int_t startRun, Int_t endRun)
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
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);
961 //_____________________________________________________________________________
963 AliMUONCDB::WriteCapacitances(Bool_t defaultValues,
964 Int_t startRun, Int_t endRun)
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
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);
978 //_____________________________________________________________________________
980 AliMUONCDB::WriteTrigger(Int_t startRun, Int_t endRun)
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);
990 //_____________________________________________________________________________
992 AliMUONCDB::WriteTracker(Bool_t defaultValues, Int_t startRun, Int_t endRun)
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);