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>
76 //_____________________________________________________________________________
77 void getBoundaries(const AliMUONVStore& store,
78 Float_t& x0min, Float_t& x0max,
79 Float_t& x1min, Float_t& x1max)
86 TIter next(store.CreateIterator());
87 AliMUONVCalibParam* value;
89 while ( ( value = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
91 Int_t detElemId = value->ID0();
92 Int_t manuId = value->ID1();
94 const AliMpVSegmentation* seg =
95 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
97 for ( Int_t manuChannel = 0; manuChannel < value->Size(); ++manuChannel )
99 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
100 if (!pad.IsValid()) continue;
102 Float_t x0 = value->ValueAsFloat(manuChannel,0);
104 x0min = TMath::Min(x0min,x0);
105 x0max = TMath::Max(x0max,x0);
106 if ( value->Dimension()>1 )
108 Float_t x1 = value->ValueAsFloat(manuChannel,1);
109 x1min = TMath::Min(x1min,x1);
110 x1max = TMath::Max(x1max,x1);
116 //_____________________________________________________________________________
117 Double_t GetRandom(Double_t mean, Double_t sigma, Bool_t mustBePositive)
120 if ( mustBePositive )
124 x = gRandom->Gaus(mean,sigma);
129 x = gRandom->Gaus(mean,sigma);
136 //_____________________________________________________________________________
137 AliMUONCDB::AliMUONCDB(const char* cdbpath)
145 //_____________________________________________________________________________
146 AliMUONCDB::~AliMUONCDB()
152 //_____________________________________________________________________________
154 AliMUONCDB::ManuList()
156 /// return (and create if necessary) the list of (de,manu) pairs
159 AliInfo("Generating ManuList...");
160 fManuList = AliMpManuList::ManuList();
161 AliInfo("Manu List generated.");
166 //_____________________________________________________________________________
168 AliMUONCDB::Diff(AliMUONVStore& store1, AliMUONVStore& store2,
171 /// creates a store which contains store1-store2
172 /// if opt="abs" the difference is absolute one,
173 /// if opt="rel" then what is stored is (store1-store2)/store1
174 /// WARNING Works only for stores which holds AliMUONVCalibParam objects
179 if ( !sopt.Contains("ABS") && !sopt.Contains("REL") )
181 AliErrorClass(Form("opt %s not supported. Only ABS or REL are",opt));
185 AliMUONVStore* d = static_cast<AliMUONVStore*>(store1.Clone());
187 TIter next(d->CreateIterator());
189 AliMUONVCalibParam* param;
191 while ( ( param = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
193 Int_t detElemId = param->ID0();
194 Int_t manuId = param->ID1();
196 AliMUONVCalibParam* param2 = dynamic_cast<AliMUONVCalibParam*>(store2.FindObject(detElemId,manuId));
197 //FIXME: this might happen. Handle it.
200 cerr << "param2 is null : FIXME : this might happen !" << endl;
205 for ( Int_t i = 0; i < param->Size(); ++i )
207 for ( Int_t j = 0; j < param->Dimension(); ++j )
210 if ( sopt.Contains("ABS") )
212 value = param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j);
214 else if ( sopt.Contains("REL") )
216 if ( param->ValueAsFloat(i,j) )
218 value = (param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j))/param->ValueAsFloat(i,j);
225 param->SetValueAsFloat(i,j,value);
232 //_____________________________________________________________________________
234 AliMUONCDB::Plot(const AliMUONVStore& store, const char* name, Int_t nbins)
236 /// Make a plot of the first 1 or 2 dimensions of the AliMUONVCalibParam
237 /// contained inside store.
238 /// It produces histograms named name_0, name_1, etc...
240 Float_t x0min, x0max, x1min, x1max;
242 getBoundaries(store,x0min,x0max,x1min,x1max);
246 cerr << Form("Something is wrong with boundaries : x0(min,max)=%e,%e",
247 x0min,x0max) << endl;
251 if ( TMath::Abs(x0min-x0max) < 1E-3 )
257 TH1* h0 = new TH1F(Form("%s_0",name),Form("%s_0",name),
264 h1 = new TH1F(Form("%s_1",name),Form("%s_1",name),
268 TIter next(ManuList());
271 Int_t nPerStation[7];
273 for ( Int_t i = 0; i < 7; ++i ) nPerStation[i]=0;
275 while ( ( p = (AliMpIntPair*)next() ) )
277 Int_t detElemId = p->GetFirst();
278 Int_t manuId = p->GetSecond();
279 Int_t station = AliMpDEManager::GetChamberId(detElemId);
281 const AliMpVSegmentation* seg =
282 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
284 AliMUONVCalibParam* value =
285 dynamic_cast<AliMUONVCalibParam*>(store.FindObject(detElemId,manuId));
289 for ( Int_t manuChannel = 0; manuChannel < value->Size(); ++manuChannel )
291 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
292 if (!pad.IsValid()) continue;
295 ++nPerStation[station];
296 Float_t x = value->ValueAsFloat(manuChannel,0);
299 AliInfo(Form("DE %d Manu %d Ch %d x=%e",detElemId,manuId,manuChannel,x));
304 h1->Fill(value->ValueAsFloat(manuChannel,1));
310 AliWarning(Form("Got a null value for DE=%d manuId=%d",detElemId,manuId));
314 AliInfo(Form("Number of channels = %d",n));
315 for ( Int_t i = 0; i < 7; ++i )
317 AliInfo(Form("Station %d %d ",(i+1),nPerStation[i]));
321 //_____________________________________________________________________________
323 AliMUONCDB::MakeHVStore(TMap& aliasMap, Bool_t defaultValues)
325 /// Create a HV store
327 AliMUONHVNamer hvNamer;
329 TObjArray* aliases = hvNamer.GenerateAliases();
334 for ( Int_t i = 0; i < aliases->GetEntries(); ++i )
336 TObjString* alias = static_cast<TObjString*>(aliases->At(i));
337 TString& aliasName = alias->String();
338 if ( aliasName.Contains("sw") )
340 // HV Switch (St345 only)
341 TObjArray* valueSet = new TObjArray;
342 valueSet->SetOwner(kTRUE);
344 Bool_t value = kTRUE;
348 Float_t r = gRandom->Uniform();
349 if ( r < 0.007 ) value = kFALSE;
352 for ( UInt_t timeStamp = 0; timeStamp < 60*3; timeStamp += 60 )
354 AliDCSValue* dcsValue = new AliDCSValue(value,timeStamp);
355 valueSet->Add(dcsValue);
357 aliasMap.Add(new TObjString(*alias),valueSet);
362 TObjArray* valueSet = new TObjArray;
363 valueSet->SetOwner(kTRUE);
364 for ( UInt_t timeStamp = 0; timeStamp < 60*15; timeStamp += 120 )
366 Float_t value = 1500;
367 if (!defaultValues) value = GetRandom(1750,62.5,true);
368 AliDCSValue* dcsValue = new AliDCSValue(value,timeStamp);
369 valueSet->Add(dcsValue);
371 aliasMap.Add(new TObjString(*alias),valueSet);
378 AliInfo(Form("%d HV channels and %d switches",nChannels,nSwitch));
380 return nChannels+nSwitch;
383 //_____________________________________________________________________________
385 AliMUONCDB::MakePedestalStore(AliMUONVStore& pedestalStore, Bool_t defaultValues)
387 /// Create a pedestal store. if defaultValues=true, ped.mean=ped.sigma=1,
388 /// otherwise mean and sigma are from a gaussian (with parameters
389 /// defined below by the kPedestal* constants)
391 TIter next(ManuList());
398 const Int_t kChannels(AliMpConstants::ManuNofChannels());
399 const Float_t kPedestalMeanMean(200);
400 const Float_t kPedestalMeanSigma(10);
401 const Float_t kPedestalSigmaMean(1.0);
402 const Float_t kPedestalSigmaSigma(0.2);
404 while ( ( p = (AliMpIntPair*)next() ) )
408 Int_t detElemId = p->GetFirst();
409 Int_t manuId = p->GetSecond();
411 AliMUONVCalibParam* ped = new AliMUONCalibParamNF(2,kChannels,detElemId,manuId,AliMUONVCalibParam::InvalidFloatValue());
413 const AliMpVSegmentation* seg =
414 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
416 for ( Int_t manuChannel = 0; manuChannel < kChannels; ++manuChannel )
418 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
419 if (!pad.IsValid()) continue;
423 Float_t meanPedestal;
424 Float_t sigmaPedestal;
433 Bool_t positive(kTRUE);
434 meanPedestal = GetRandom(kPedestalMeanMean,kPedestalMeanSigma,positive);
435 sigmaPedestal = GetRandom(kPedestalSigmaMean,kPedestalSigmaSigma,positive);
437 ped->SetValueAsFloat(manuChannel,0,meanPedestal);
438 ped->SetValueAsFloat(manuChannel,1,sigmaPedestal);
441 Bool_t ok = pedestalStore.Add(ped);
444 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
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));
612 AliInfo(Form("%d Manus and %d channels.",nmanus,nchannels));
616 //_____________________________________________________________________________
618 AliMUONCDB::MakeLocalTriggerMaskStore(AliMUONVStore& localBoardMasks) const
620 /// Generate local trigger masks store. All masks are set to FFFF
623 // Generate fake mask values for 234 localboards and put that into
624 // one single container (localBoardMasks)
625 for ( Int_t i = 1; i <= 234; ++i )
627 AliMUONVCalibParam* localBoard = new AliMUONCalibParamNI(1,8,i,0,0);
628 for ( Int_t x = 0; x < 2; ++x )
630 for ( Int_t y = 0; y < 4; ++y )
633 localBoard->SetValueAsInt(index,0,0xFFFF);
637 localBoardMasks.Add(localBoard);
642 //_____________________________________________________________________________
644 AliMUONCDB::MakeRegionalTriggerMaskStore(AliMUONVStore& rtm) const
646 /// Make a regional trigger masks store. All masks are set to 3F
649 for ( Int_t i = 0; i < 16; ++i )
651 AliMUONVCalibParam* regionalBoard = new AliMUONCalibParamNI(1,16,i,0,0);
652 for ( Int_t j = 0; j < 16; ++j )
654 regionalBoard->SetValueAsInt(j,0,0x3F);
657 rtm.Add(regionalBoard);
663 //_____________________________________________________________________________
665 AliMUONCDB::MakeGlobalTriggerMaskStore(AliMUONVCalibParam& gtm) const
667 /// Make a global trigger masks store. All masks set to FFF
671 for ( Int_t j = 0; j < 16; ++j )
673 gtm.SetValueAsInt(j,0,0xFFF);
679 //_____________________________________________________________________________
681 AliMUONCDB::MakeTriggerLUT(const char* file) const
683 /// Make a triggerlut object, from a file.
685 AliMUONTriggerLut* lut = new AliMUONTriggerLut;
686 lut->ReadFromFile(file);
690 //_____________________________________________________________________________
691 AliMUONTriggerEfficiencyCells*
692 AliMUONCDB::MakeTriggerEfficiency(const char* file) const
694 /// Make a trigger efficiency object from a file.
696 return new AliMUONTriggerEfficiencyCells(file);
699 //_____________________________________________________________________________
701 AliMUONCDB::WriteToCDB(const char* calibpath, TObject* object,
702 Int_t startRun, Int_t endRun, Bool_t defaultValues)
704 /// Write a given object to OCDB
706 AliCDBId id(calibpath,startRun,endRun);
708 md.SetAliRootVersion(gROOT->GetVersion());
711 md.SetComment("Test with default values");
715 md.SetComment("Test with random values");
717 md.SetResponsible("AliMUONCDB tester class");
719 AliCDBManager* man = AliCDBManager::Instance();
720 man->SetDefaultStorage(fCDBPath);
721 man->Put(object,id,&md);
724 //_____________________________________________________________________________
726 AliMUONCDB::MakeNeighbourStore(AliMUONVStore& neighbourStore)
728 /// Fill the neighbours store with, for each channel, a TObjArray of its
729 /// neighbouring pads (including itself)
731 AliInfo("Generating NeighbourStore. This will take a while. Please be patient.");
737 TIter next(ManuList());
745 while ( ( p = (AliMpIntPair*)next() ) )
747 Int_t detElemId = p->GetFirst();
748 Int_t manuId = p->GetSecond();
750 const AliMpVSegmentation* seg =
751 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
753 AliMUONVCalibParam* calibParam = static_cast<AliMUONVCalibParam*>(neighbourStore.FindObject(detElemId,manuId));
757 Int_t size(AliMpConstants::ManuNofChannels());
758 Int_t defaultValue(-1);
759 Int_t packingFactor(size);
761 calibParam = new AliMUONCalibParamNI(dimension,size,detElemId,manuId,defaultValue,packingFactor);
762 Bool_t ok = neighbourStore.Add(calibParam);
765 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
770 for ( Int_t manuChannel = 0; manuChannel < AliMpConstants::ManuNofChannels(); ++manuChannel )
772 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
778 seg->GetNeighbours(pad,tmp,true,true);
779 Int_t nofPadNeighbours = tmp.GetEntriesFast();
781 for ( Int_t i = 0; i < nofPadNeighbours; ++i )
783 AliMpPad* pad = static_cast<AliMpPad*>(tmp.At(i));
785 Bool_t ok = calibParam->PackValues(pad->GetLocation().GetFirst(),pad->GetLocation().GetSecond(),x);
788 AliError("Could not pack value. Something is seriously wrong. Please check");
789 StdoutToAliError(pad->Print(););
792 calibParam->SetValueAsInt(manuChannel,i,x);
803 //_____________________________________________________________________________
805 AliMUONCDB::WriteLocalTriggerMasks(Int_t startRun, Int_t endRun)
807 /// Write local trigger masks to OCDB
809 AliMUONVStore* ltm = new AliMUON1DArray(235);
810 Int_t ngenerated = MakeLocalTriggerMaskStore(*ltm);
811 AliInfo(Form("Ngenerated = %d",ngenerated));
814 WriteToCDB("MUON/Calib/LocalTriggerBoardMasks",ltm,startRun,endRun,true);
819 //_____________________________________________________________________________
821 AliMUONCDB::WriteRegionalTriggerMasks(Int_t startRun, Int_t endRun)
823 /// Write regional trigger masks to OCDB
825 AliMUONVStore* rtm = new AliMUON1DArray(16);
826 Int_t ngenerated = MakeRegionalTriggerMaskStore(*rtm);
827 AliInfo(Form("Ngenerated = %d",ngenerated));
830 WriteToCDB("MUON/Calib/RegionalTriggerBoardMasks",rtm,startRun,endRun,true);
835 //_____________________________________________________________________________
837 AliMUONCDB::WriteGlobalTriggerMasks(Int_t startRun, Int_t endRun)
839 /// Write global trigger masks to OCDB
841 AliMUONVCalibParam* gtm = new AliMUONCalibParamNI(1,16,1,0,0);
843 Int_t ngenerated = MakeGlobalTriggerMaskStore(*gtm);
844 AliInfo(Form("Ngenerated = %d",ngenerated));
847 WriteToCDB("MUON/Calib/GlobalTriggerBoardMasks",gtm,startRun,endRun,true);
852 //_____________________________________________________________________________
854 AliMUONCDB::WriteTriggerLut(Int_t startRun, Int_t endRun)
856 /// Write trigger LUT to OCDB
858 AliMUONTriggerLut* lut = MakeTriggerLUT();
861 WriteToCDB("MUON/Calib/TriggerLut",lut,startRun,endRun,true);
866 //_____________________________________________________________________________
868 AliMUONCDB::WriteTriggerEfficiency(Int_t startRun, Int_t endRun)
870 /// Write trigger efficiency to OCDB
872 AliMUONTriggerEfficiencyCells* eff = MakeTriggerEfficiency();
875 WriteToCDB("MUON/Calib/TriggerEfficiency",eff,startRun,endRun,true);
880 //_____________________________________________________________________________
882 AliMUONCDB::WriteNeighbours(Int_t startRun, Int_t endRun)
884 /// Write neighbours to OCDB
886 AliMUONVStore* neighbours = new AliMUON2DMap(kTRUE);
887 Int_t ngenerated = MakeNeighbourStore(*neighbours);
888 AliInfo(Form("Ngenerated = %d",ngenerated));
891 WriteToCDB("MUON/Calib/Neighbours",neighbours,startRun,endRun,true);
896 //_____________________________________________________________________________
898 AliMUONCDB::WriteHV(Bool_t defaultValues,
899 Int_t startRun, Int_t endRun)
901 /// generate HV values (either cste = 1500 V) if defaultValues=true or random
902 /// if defaultValues=false, see makeHVStore) and
903 /// store them into CDB located at cdbpath, with a validity period
904 /// ranging from startRun to endRun
906 TMap* hvStore = new TMap;
907 Int_t ngenerated = MakeHVStore(*hvStore,defaultValues);
908 AliInfo(Form("Ngenerated = %d",ngenerated));
911 WriteToCDB("MUON/Calib/HV",hvStore,startRun,endRun,defaultValues);
916 //_____________________________________________________________________________
918 AliMUONCDB::WritePedestals(Bool_t defaultValues,
919 Int_t startRun, Int_t endRun)
921 /// generate pedestal values (either 0 if defaultValues=true or random
922 /// if defaultValues=false, see makePedestalStore) and
923 /// store them into CDB located at cdbpath, with a validity period
924 /// ranging from startRun to endRun
926 AliMUONVStore* pedestalStore = new AliMUON2DMap(true);
927 Int_t ngenerated = MakePedestalStore(*pedestalStore,defaultValues);
928 AliInfo(Form("Ngenerated = %d",ngenerated));
929 WriteToCDB("MUON/Calib/Pedestals",pedestalStore,startRun,endRun,defaultValues);
930 delete pedestalStore;
934 //_____________________________________________________________________________
936 AliMUONCDB::WriteGains(Bool_t defaultValues,
937 Int_t startRun, Int_t endRun)
939 /// generate gain values (either 1 if defaultValues=true or random
940 /// if defaultValues=false, see makeGainStore) and
941 /// store them into CDB located at cdbpath, with a validity period
942 /// ranging from startRun to endRun
944 AliMUONVStore* gainStore = new AliMUON2DMap(true);
945 Int_t ngenerated = MakeGainStore(*gainStore,defaultValues);
946 AliInfo(Form("Ngenerated = %d",ngenerated));
947 WriteToCDB("MUON/Calib/Gains",gainStore,startRun,endRun,defaultValues);
951 //_____________________________________________________________________________
953 AliMUONCDB::WriteCapacitances(Bool_t defaultValues,
954 Int_t startRun, Int_t endRun)
956 /// generate manu capacitance values (either 1 if defaultValues=true or random
957 /// if defaultValues=false, see makeCapacitanceStore) and
958 /// store them into CDB located at cdbpath, with a validity period
959 /// ranging from startRun to endRun
961 AliMUONVStore* capaStore = new AliMUON1DMap(16828);
962 Int_t ngenerated = MakeCapacitanceStore(*capaStore,defaultValues);
963 AliInfo(Form("Ngenerated = %d",ngenerated));
964 WriteToCDB("MUON/Calib/Capacitances",capaStore,startRun,endRun,defaultValues);
968 //_____________________________________________________________________________
970 AliMUONCDB::WriteTrigger(Int_t startRun, Int_t endRun)
972 /// Writes all Trigger related calibration to CDB
973 WriteLocalTriggerMasks(startRun,endRun);
974 WriteRegionalTriggerMasks(startRun,endRun);
975 WriteGlobalTriggerMasks(startRun,endRun);
976 WriteTriggerLut(startRun,endRun);
977 WriteTriggerEfficiency(startRun,endRun);
980 //_____________________________________________________________________________
982 AliMUONCDB::WriteTracker(Bool_t defaultValues, Int_t startRun, Int_t endRun)
984 /// Writes all Tracker related calibration to CDB
985 WriteHV(defaultValues,startRun,endRun);
986 WritePedestals(defaultValues,startRun,endRun);
987 WriteGains(defaultValues,startRun,endRun);
988 WriteCapacitances(defaultValues,startRun,endRun);
989 WriteNeighbours(startRun,endRun);