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 **************************************************************************/
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
25 /// For more information, please see READMECDB
27 // \author Laurent Aphecetche
29 #include "AliMUONCDB.h"
31 #include "AliCDBEntry.h"
32 #include "AliCDBManager.h"
33 #include "AliDCSValue.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>
63 #include <TObjString.h>
66 #include <TStopwatch.h>
75 //_____________________________________________________________________________
76 void getBoundaries(const AliMUONVStore& store,
77 Float_t& x0min, Float_t& x0max,
78 Float_t& x1min, Float_t& x1max)
85 TIter next(store.CreateIterator());
86 AliMUONVCalibParam* value;
88 while ( ( value = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
90 Int_t detElemId = value->ID0();
91 Int_t manuId = value->ID1();
93 const AliMpVSegmentation* seg =
94 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
96 for ( Int_t manuChannel = 0; manuChannel < value->Size(); ++manuChannel )
98 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
99 if (!pad.IsValid()) continue;
101 Float_t x0 = value->ValueAsFloat(manuChannel,0);
103 x0min = TMath::Min(x0min,x0);
104 x0max = TMath::Max(x0max,x0);
105 if ( value->Dimension()>1 )
107 Float_t x1 = value->ValueAsFloat(manuChannel,1);
108 x1min = TMath::Min(x1min,x1);
109 x1max = TMath::Max(x1max,x1);
117 //_____________________________________________________________________________
118 AliMUONCDB::AliMUONCDB(const char* cdbpath)
126 //_____________________________________________________________________________
127 AliMUONCDB::~AliMUONCDB()
133 //_____________________________________________________________________________
135 AliMUONCDB::ManuList()
137 /// return (and create if necessary) the list of (de,manu) pairs
140 AliInfo("Generating ManuList...");
141 fManuList = AliMpManuList::ManuList();
142 AliInfo("Manu List generated.");
147 //_____________________________________________________________________________
149 AliMUONCDB::Diff(AliMUONVStore& store1, AliMUONVStore& store2,
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
160 if ( !sopt.Contains("ABS") && !sopt.Contains("REL") )
162 AliErrorClass(Form("opt %s not supported. Only ABS or REL are",opt));
166 AliMUONVStore* d = static_cast<AliMUONVStore*>(store1.Clone());
168 TIter next(d->CreateIterator());
170 AliMUONVCalibParam* param;
172 while ( ( param = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
174 Int_t detElemId = param->ID0();
175 Int_t manuId = param->ID1();
177 AliMUONVCalibParam* param2 = dynamic_cast<AliMUONVCalibParam*>(store2.FindObject(detElemId,manuId));
178 //FIXME: this might happen. Handle it.
181 cerr << "param2 is null : FIXME : this might happen !" << endl;
186 for ( Int_t i = 0; i < param->Size(); ++i )
188 for ( Int_t j = 0; j < param->Dimension(); ++j )
191 if ( sopt.Contains("ABS") )
193 value = param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j);
195 else if ( sopt.Contains("REL") )
197 if ( param->ValueAsFloat(i,j) )
199 value = (param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j))/param->ValueAsFloat(i,j);
206 param->SetValueAsFloat(i,j,value);
213 //_____________________________________________________________________________
215 AliMUONCDB::Plot(const AliMUONVStore& store, const char* name, Int_t nbins)
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...
221 Float_t x0min, x0max, x1min, x1max;
223 getBoundaries(store,x0min,x0max,x1min,x1max);
227 cerr << Form("Something is wrong with boundaries : x0(min,max)=%e,%e",
228 x0min,x0max) << endl;
232 if ( TMath::Abs(x0min-x0max) < 1E-3 )
238 TH1* h0 = new TH1F(Form("%s_0",name),Form("%s_0",name),
245 h1 = new TH1F(Form("%s_1",name),Form("%s_1",name),
249 TIter next(ManuList());
252 Int_t nPerStation[7];
254 for ( Int_t i = 0; i < 7; ++i ) nPerStation[i]=0;
256 while ( ( p = (AliMpIntPair*)next() ) )
258 Int_t detElemId = p->GetFirst();
259 Int_t manuId = p->GetSecond();
260 Int_t station = AliMpDEManager::GetChamberId(detElemId);
262 const AliMpVSegmentation* seg =
263 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
265 AliMUONVCalibParam* value =
266 dynamic_cast<AliMUONVCalibParam*>(store.FindObject(detElemId,manuId));
270 for ( Int_t manuChannel = 0; manuChannel < value->Size(); ++manuChannel )
272 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
273 if (!pad.IsValid()) continue;
276 ++nPerStation[station];
277 Float_t x = value->ValueAsFloat(manuChannel,0);
280 AliInfo(Form("DE %d Manu %d Ch %d x=%e",detElemId,manuId,manuChannel,x));
285 h1->Fill(value->ValueAsFloat(manuChannel,1));
291 AliWarning(Form("Got a null value for DE=%d manuId=%d",detElemId,manuId));
295 AliInfo(Form("Number of channels = %d",n));
296 for ( Int_t i = 0; i < 7; ++i )
298 AliInfo(Form("Station %d %d ",(i+1),nPerStation[i]));
302 //_____________________________________________________________________________
304 AliMUONCDB::MakeHVStore(TMap& aliasMap, Bool_t defaultValues)
306 /// Create a HV store
308 AliMUONHVNamer hvNamer;
310 TObjArray* aliases = hvNamer.GenerateAliases();
315 for ( Int_t i = 0; i < aliases->GetEntries(); ++i )
317 TObjString* alias = static_cast<TObjString*>(aliases->At(i));
318 TString& aliasName = alias->String();
319 if ( aliasName.Contains("sw") )
321 // HV Switch (St345 only)
322 TObjArray* valueSet = new TObjArray;
323 valueSet->SetOwner(kTRUE);
325 Bool_t value = kTRUE;
329 Float_t r = gRandom->Uniform();
330 if ( r < 0.007 ) value = kFALSE;
333 for ( UInt_t timeStamp = 0; timeStamp < 60*3; timeStamp += 60 )
335 AliDCSValue* dcsValue = new AliDCSValue(value,timeStamp);
336 valueSet->Add(dcsValue);
338 aliasMap.Add(new TObjString(*alias),valueSet);
343 TObjArray* valueSet = new TObjArray;
344 valueSet->SetOwner(kTRUE);
345 for ( UInt_t timeStamp = 0; timeStamp < 60*15; timeStamp += 120 )
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);
352 aliasMap.Add(new TObjString(*alias),valueSet);
359 AliInfo(Form("%d HV channels and %d switches",nChannels,nSwitch));
361 return nChannels+nSwitch;
364 //_____________________________________________________________________________
366 AliMUONCDB::MakePedestalStore(AliMUONVStore& pedestalStore, Bool_t defaultValues)
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)
372 TIter next(ManuList());
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);
385 while ( ( p = (AliMpIntPair*)next() ) )
389 Int_t detElemId = p->GetFirst();
390 Int_t manuId = p->GetSecond();
392 AliMUONVCalibParam* ped = new AliMUONCalibParamNF(2,kChannels,detElemId,manuId,AliMUONVCalibParam::InvalidFloatValue());
394 const AliMpVSegmentation* seg =
395 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
397 for ( Int_t manuChannel = 0; manuChannel < kChannels; ++manuChannel )
399 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
400 if (!pad.IsValid()) continue;
404 Float_t meanPedestal;
405 Float_t sigmaPedestal;
415 while ( meanPedestal < 0 )
417 meanPedestal = gRandom->Gaus(kPedestalMeanMean,kPedestalMeanSigma);
420 while ( sigmaPedestal < 0 )
422 sigmaPedestal = gRandom->Gaus(kPedestalSigmaMean,kPedestalSigmaSigma);
425 ped->SetValueAsFloat(manuChannel,0,meanPedestal);
426 ped->SetValueAsFloat(manuChannel,1,sigmaPedestal);
429 Bool_t ok = pedestalStore.Add(ped);
432 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
436 AliInfo(Form("%d Manus and %d channels.",nmanus,nchannels));
441 //_____________________________________________________________________________
443 AliMUONCDB::MakeCapacitanceStore(AliMUONVStore& capaStore, Bool_t defaultValues)
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.
449 TIter next(ManuList());
455 Int_t nmanusOK(0); // manus for which we got the serial number
457 const Float_t kCapaMean(1.0);
458 const Float_t kCapaSigma(0.5);
460 while ( ( p = (AliMpIntPair*)next() ) )
464 Int_t detElemId = p->GetFirst();
465 Int_t manuId = p->GetSecond();
467 AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(detElemId);
468 Int_t serialNumber = de->GetManuSerialFromId(manuId);
470 if ( serialNumber <= 0 ) continue;
474 const AliMpVSegmentation* seg =
475 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
477 AliMUONVCalibParam* capa = static_cast<AliMUONVCalibParam*>(capaStore.FindObject(serialNumber));
481 capa = new AliMUONCalibParamNF(1,AliMpConstants::ManuNofChannels(),serialNumber,0,1.0);
482 Bool_t ok = capaStore.Add(capa);
485 AliError(Form("Could not set serialNumber=%d manuId=%d",serialNumber,manuId));
489 for ( Int_t manuChannel = 0; manuChannel < capa->Size(); ++manuChannel )
491 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
492 if (!pad.IsValid()) continue;
505 while ( capaValue < 0 )
507 capaValue = gRandom->Gaus(kCapaMean,kCapaSigma);
510 capa->SetValueAsFloat(manuChannel,0,capaValue);
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));
521 AliWarning("Did not get all serial numbers. capaStore is incomplete !!!!");
527 //_____________________________________________________________________________
529 AliMUONCDB::MakeGainStore(AliMUONVStore& gainStore, Bool_t defaultValues)
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.
535 TIter next(ManuList());
542 const Double_t kSaturation(3000);
543 const Double_t kGainMean(1.0);
544 const Double_t kGainSigma(0.05);
546 while ( ( p = (AliMpIntPair*)next() ) )
550 Int_t detElemId = p->GetFirst();
551 Int_t manuId = p->GetSecond();
553 AliMUONVCalibParam* gain =
554 new AliMUONCalibParamNF(2,AliMpConstants::ManuNofChannels(),
557 AliMUONVCalibParam::InvalidFloatValue());
560 const AliMpVSegmentation* seg =
561 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
563 for ( Int_t manuChannel = 0; manuChannel < gain->Size(); ++manuChannel )
565 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
566 if (!pad.IsValid()) continue;
571 Float_t saturation(kSaturation);
580 while ( meanGain < 0 )
582 meanGain = gRandom->Gaus(kGainMean,kGainSigma);
585 gain->SetValueAsFloat(manuChannel,0,meanGain);
586 gain->SetValueAsFloat(manuChannel,1,saturation);
589 Bool_t ok = gainStore.Add(gain);
592 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
596 AliInfo(Form("%d Manus and %d channels.",nmanus,nchannels));
600 //_____________________________________________________________________________
602 AliMUONCDB::MakeLocalTriggerMaskStore(AliMUONVStore& localBoardMasks) const
604 /// Generate local trigger masks store. All masks are set to FFFF
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 )
611 AliMUONVCalibParam* localBoard = new AliMUONCalibParamNI(1,8,i,0,0);
612 for ( Int_t x = 0; x < 2; ++x )
614 for ( Int_t y = 0; y < 4; ++y )
617 localBoard->SetValueAsInt(index,0,0xFFFF);
621 localBoardMasks.Add(localBoard);
626 //_____________________________________________________________________________
628 AliMUONCDB::MakeRegionalTriggerMaskStore(AliMUONVStore& rtm) const
630 /// Make a regional trigger masks store. All masks are set to 3F
633 for ( Int_t i = 0; i < 16; ++i )
635 AliMUONVCalibParam* regionalBoard = new AliMUONCalibParamNI(1,16,i,0,0);
636 for ( Int_t j = 0; j < 16; ++j )
638 regionalBoard->SetValueAsInt(j,0,0x3F);
641 rtm.Add(regionalBoard);
647 //_____________________________________________________________________________
649 AliMUONCDB::MakeGlobalTriggerMaskStore(AliMUONVCalibParam& gtm) const
651 /// Make a global trigger masks store. All masks set to FFF
655 for ( Int_t j = 0; j < 16; ++j )
657 gtm.SetValueAsInt(j,0,0xFFF);
663 //_____________________________________________________________________________
665 AliMUONCDB::MakeTriggerLUT(const char* file) const
667 /// Make a triggerlut object, from a file.
669 AliMUONTriggerLut* lut = new AliMUONTriggerLut;
670 lut->ReadFromFile(file);
674 //_____________________________________________________________________________
675 AliMUONTriggerEfficiencyCells*
676 AliMUONCDB::MakeTriggerEfficiency(const char* file) const
678 /// Make a trigger efficiency object from a file.
680 return new AliMUONTriggerEfficiencyCells(file);
683 //_____________________________________________________________________________
685 AliMUONCDB::WriteToCDB(const char* calibpath, TObject* object,
686 Int_t startRun, Int_t endRun, Bool_t defaultValues)
688 /// Write a given object to OCDB
690 AliCDBId id(calibpath,startRun,endRun);
692 md.SetAliRootVersion(gROOT->GetVersion());
695 md.SetComment("Test with default values");
699 md.SetComment("Test with random values");
701 md.SetResponsible("AliMUONCDB tester class");
703 AliCDBManager* man = AliCDBManager::Instance();
704 man->SetDefaultStorage(fCDBPath);
705 man->Put(object,id,&md);
708 //_____________________________________________________________________________
710 AliMUONCDB::MakeNeighbourStore(AliMUONVStore& neighbourStore)
712 /// Fill the neighbours store with, for each channel, a TObjArray of its
713 /// neighbouring pads (including itself)
715 AliInfo("Generating NeighbourStore. This will take a while. Please be patient.");
721 TIter next(ManuList());
729 while ( ( p = (AliMpIntPair*)next() ) )
731 Int_t detElemId = p->GetFirst();
732 Int_t manuId = p->GetSecond();
734 const AliMpVSegmentation* seg =
735 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
737 AliMUONVCalibParam* calibParam = static_cast<AliMUONVCalibParam*>(neighbourStore.FindObject(detElemId,manuId));
741 Int_t size(AliMpConstants::ManuNofChannels());
742 Int_t defaultValue(-1);
743 Int_t packingFactor(size);
745 calibParam = new AliMUONCalibParamNI(dimension,size,detElemId,manuId,defaultValue,packingFactor);
746 Bool_t ok = neighbourStore.Add(calibParam);
749 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
754 for ( Int_t manuChannel = 0; manuChannel < AliMpConstants::ManuNofChannels(); ++manuChannel )
756 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
762 seg->GetNeighbours(pad,tmp,true,true);
763 Int_t nofPadNeighbours = tmp.GetEntriesFast();
765 for ( Int_t i = 0; i < nofPadNeighbours; ++i )
767 AliMpPad* pad = static_cast<AliMpPad*>(tmp.At(i));
769 Bool_t ok = calibParam->PackValues(pad->GetLocation().GetFirst(),pad->GetLocation().GetSecond(),x);
772 AliError("Could not pack value. Something is seriously wrong. Please check");
773 StdoutToAliError(pad->Print(););
776 calibParam->SetValueAsInt(manuChannel,i,x);
787 //_____________________________________________________________________________
789 AliMUONCDB::WriteLocalTriggerMasks(Int_t startRun, Int_t endRun)
791 /// Write local trigger masks to OCDB
793 AliMUONVStore* ltm = new AliMUON1DArray(235);
794 Int_t ngenerated = MakeLocalTriggerMaskStore(*ltm);
795 AliInfo(Form("Ngenerated = %d",ngenerated));
798 WriteToCDB("MUON/Calib/LocalTriggerBoardMasks",ltm,startRun,endRun,true);
803 //_____________________________________________________________________________
805 AliMUONCDB::WriteRegionalTriggerMasks(Int_t startRun, Int_t endRun)
807 /// Write regional trigger masks to OCDB
809 AliMUONVStore* rtm = new AliMUON1DArray(16);
810 Int_t ngenerated = MakeRegionalTriggerMaskStore(*rtm);
811 AliInfo(Form("Ngenerated = %d",ngenerated));
814 WriteToCDB("MUON/Calib/RegionalTriggerBoardMasks",rtm,startRun,endRun,true);
819 //_____________________________________________________________________________
821 AliMUONCDB::WriteGlobalTriggerMasks(Int_t startRun, Int_t endRun)
823 /// Write global trigger masks to OCDB
825 AliMUONVCalibParam* gtm = new AliMUONCalibParamNI(1,16,1,0,0);
827 Int_t ngenerated = MakeGlobalTriggerMaskStore(*gtm);
828 AliInfo(Form("Ngenerated = %d",ngenerated));
831 WriteToCDB("MUON/Calib/GlobalTriggerBoardMasks",gtm,startRun,endRun,true);
836 //_____________________________________________________________________________
838 AliMUONCDB::WriteTriggerLut(Int_t startRun, Int_t endRun)
840 /// Write trigger LUT to OCDB
842 AliMUONTriggerLut* lut = MakeTriggerLUT();
845 WriteToCDB("MUON/Calib/TriggerLut",lut,startRun,endRun,true);
850 //_____________________________________________________________________________
852 AliMUONCDB::WriteTriggerEfficiency(Int_t startRun, Int_t endRun)
854 /// Write trigger efficiency to OCDB
856 AliMUONTriggerEfficiencyCells* eff = MakeTriggerEfficiency();
859 WriteToCDB("MUON/Calib/TriggerEfficiency",eff,startRun,endRun,true);
864 //_____________________________________________________________________________
866 AliMUONCDB::WriteNeighbours(Int_t startRun, Int_t endRun)
868 /// Write neighbours to OCDB
870 AliMUONVStore* neighbours = new AliMUON2DMap(kTRUE);
871 Int_t ngenerated = MakeNeighbourStore(*neighbours);
872 AliInfo(Form("Ngenerated = %d",ngenerated));
875 WriteToCDB("MUON/Calib/Neighbours",neighbours,startRun,endRun,true);
880 //_____________________________________________________________________________
882 AliMUONCDB::WriteHV(Bool_t defaultValues,
883 Int_t startRun, Int_t endRun)
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
890 TMap* hvStore = new TMap;
891 Int_t ngenerated = MakeHVStore(*hvStore,defaultValues);
892 AliInfo(Form("Ngenerated = %d",ngenerated));
895 WriteToCDB("MUON/Calib/HV",hvStore,startRun,endRun,defaultValues);
900 //_____________________________________________________________________________
902 AliMUONCDB::WritePedestals(Bool_t defaultValues,
903 Int_t startRun, Int_t endRun)
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
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;
918 //_____________________________________________________________________________
920 AliMUONCDB::WriteGains(Bool_t defaultValues,
921 Int_t startRun, Int_t endRun)
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
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);
935 //_____________________________________________________________________________
937 AliMUONCDB::WriteCapacitances(Bool_t defaultValues,
938 Int_t startRun, Int_t endRun)
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
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);
952 //_____________________________________________________________________________
954 AliMUONCDB::WriteTrigger(Int_t startRun, Int_t endRun)
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);
964 //_____________________________________________________________________________
966 AliMUONCDB::WriteTracker(Bool_t defaultValues, Int_t startRun, Int_t endRun)
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);