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 "AliMUONObjectPair.h"
44 #include "AliMUONTriggerEfficiencyCells.h"
45 #include "AliMUONTriggerLut.h"
46 #include "AliMUONV2DStore.h"
47 #include "AliMUONVCalibParam.h"
48 #include "AliMUONVCalibParam.h"
49 #include "AliMUONVDataIterator.h"
50 #include "AliMpDDLStore.h"
51 #include "AliMpDEIterator.h"
52 #include "AliMpDEManager.h"
53 #include "AliMpDetElement.h"
54 #include "AliMpManuList.h"
55 #include "AliMpSegmentation.h"
56 #include "AliMpStationType.h"
57 #include "AliMpVSegmentation.h"
58 #include <Riostream.h>
64 #include <TObjString.h>
67 #include <TStopwatch.h>
74 const Int_t AliMUONCDB::fgkMaxNofChannelsPerManu=64;
78 //_____________________________________________________________________________
79 void getBoundaries(const AliMUONV2DStore& store,
80 Float_t& x0min, Float_t& x0max,
81 Float_t& x1min, Float_t& x1max)
88 AliMUONVDataIterator* it = store.Iterator();
92 while ( ( p = dynamic_cast<AliMUONObjectPair*>(it->Next() ) ) )
94 AliMpIntPair* dm = dynamic_cast<AliMpIntPair*>(p->Key());
95 AliMUONVCalibParam* value = dynamic_cast<AliMUONVCalibParam*>(p->Value());
97 Int_t detElemId = dm->GetFirst();
98 Int_t manuId = dm->GetSecond();
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 Float_t x0 = value->ValueAsFloat(manuChannel,0);
110 x0min = TMath::Min(x0min,x0);
111 x0max = TMath::Max(x0max,x0);
112 if ( value->Dimension()>1 )
114 Float_t x1 = value->ValueAsFloat(manuChannel,1);
115 x1min = TMath::Min(x1min,x1);
116 x1max = TMath::Max(x1max,x1);
119 if (it->IsOwner()) delete p;
125 //_____________________________________________________________________________
126 AliMUONCDB::AliMUONCDB(const char* cdbpath)
129 fManuList(AliMpManuList::ManuList())
134 //_____________________________________________________________________________
135 AliMUONCDB::~AliMUONCDB()
141 //_____________________________________________________________________________
143 AliMUONCDB::Diff(AliMUONV2DStore& store1, AliMUONV2DStore& store2,
146 /// creates a store which contains store1-store2
147 /// if opt="abs" the difference is absolute one,
148 /// if opt="rel" then what is stored is (store1-store2)/store1
149 /// WARNING Works only for stores which holds AliMUONVCalibParam objects
154 if ( !sopt.Contains("ABS") && !sopt.Contains("REL") )
156 AliErrorClass(Form("opt %s not supported. Only ABS or REL are",opt));
160 AliMUONV2DStore* d = static_cast<AliMUONV2DStore*>(store1.Clone());
162 AliMUONVDataIterator* it = d->Iterator();
164 AliMUONObjectPair* p;
166 while ( ( p = dynamic_cast<AliMUONObjectPair*>(it->Next() ) ) )
168 AliMpIntPair* dm = dynamic_cast<AliMpIntPair*>(p->Key());
169 //FIXME: this might happen (if a full manu is missing, for instance)
173 cerr << "dm is null : FIXME: this might happen !" << endl;
177 AliMUONVCalibParam* param = dynamic_cast<AliMUONVCalibParam*>(p->Value());
178 //FIXMENOT: this should *not* happen
181 cerr << "param is null" << endl;
186 AliMUONVCalibParam* param2 = dynamic_cast<AliMUONVCalibParam*>(store2.Get(dm->GetFirst(),dm->GetSecond()));
187 //FIXME: this might happen. Handle it.
190 cerr << "param2 is null : FIXME : this might happen !" << endl;
195 for ( Int_t i = 0; i < param->Size(); ++i )
197 for ( Int_t j = 0; j < param->Dimension(); ++j )
200 if ( sopt.Contains("ABS") )
202 value = param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j);
204 else if ( sopt.Contains("REL") )
206 if ( param->ValueAsFloat(i,j) )
208 value = (param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j))/param->ValueAsFloat(i,j);
215 param->SetValueAsFloat(i,j,value);
218 if (it->IsOwner()) delete p;
223 //_____________________________________________________________________________
225 AliMUONCDB::Plot(const AliMUONV2DStore& store, const char* name, Int_t nbins)
227 /// Make a plot of the first 1 or 2 dimensions of the AliMUONVCalibParam
228 /// contained inside store.
229 /// It produces histograms named name_0, name_1, etc...
231 Float_t x0min, x0max, x1min, x1max;
233 getBoundaries(store,x0min,x0max,x1min,x1max);
237 cerr << Form("Something is wrong with boundaries : x0(min,max)=%e,%e",
238 x0min,x0max) << endl;
242 if ( TMath::Abs(x0min-x0max) < 1E-3 )
248 TH1* h0 = new TH1F(Form("%s_0",name),Form("%s_0",name),
255 h1 = new TH1F(Form("%s_1",name),Form("%s_1",name),
259 TIter next(fManuList);
262 Int_t nPerStation[7];
264 for ( Int_t i = 0; i < 7; ++i ) nPerStation[i]=0;
266 while ( ( p = (AliMpIntPair*)next() ) )
268 Int_t detElemId = p->GetFirst();
269 Int_t manuId = p->GetSecond();
270 Int_t station = AliMpDEManager::GetChamberId(detElemId);
272 const AliMpVSegmentation* seg =
273 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
275 AliMUONVCalibParam* value =
276 dynamic_cast<AliMUONVCalibParam*>(store.Get(detElemId,manuId));
280 for ( Int_t manuChannel = 0; manuChannel < value->Size(); ++manuChannel )
282 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
283 if (!pad.IsValid()) continue;
286 ++nPerStation[station];
287 Float_t x = value->ValueAsFloat(manuChannel,0);
290 AliInfo(Form("DE %d Manu %d Ch %d x=%e",detElemId,manuId,manuChannel,x));
295 h1->Fill(value->ValueAsFloat(manuChannel,1));
301 AliWarning(Form("Got a null value for DE=%d manuId=%d",detElemId,manuId));
305 AliInfo(Form("Number of channels = %d",n));
306 for ( Int_t i = 0; i < 7; ++i )
308 AliInfo(Form("Station %d %d ",(i+1),nPerStation[i]));
312 //_____________________________________________________________________________
314 AliMUONCDB::MakeHVStore(TMap& aliasMap, Bool_t defaultValues)
316 /// Create a HV store
318 AliMUONHVNamer hvNamer;
321 TObjArray* aliases = hvNamer.GenerateAliases();
326 for ( Int_t i = 0; i < aliases->GetEntries(); ++i )
328 TObjString* alias = static_cast<TObjString*>(aliases->At(i));
329 TString& aliasName = alias->String();
330 if ( aliasName.Contains("sw") )
332 // HV Switch (St345 only)
333 TObjArray* valueSet = new TObjArray;
334 valueSet->SetOwner(kTRUE);
336 Bool_t value = kTRUE;
340 Float_t r = random.Uniform();
341 if ( r < 0.007 ) value = kFALSE;
344 for ( UInt_t timeStamp = 0; timeStamp < 60*3; timeStamp += 60 )
346 AliDCSValue* dcsValue = new AliDCSValue(value,timeStamp);
347 valueSet->Add(dcsValue);
349 aliasMap.Add(new TObjString(*alias),valueSet);
354 TObjArray* valueSet = new TObjArray;
355 valueSet->SetOwner(kTRUE);
356 for ( UInt_t timeStamp = 0; timeStamp < 60*15; timeStamp += 120 )
358 Float_t value = 1500;
359 if (!defaultValues) value = random.Gaus(1750,62.5);
360 AliDCSValue* dcsValue = new AliDCSValue(value,timeStamp);
361 valueSet->Add(dcsValue);
363 aliasMap.Add(new TObjString(*alias),valueSet);
370 AliInfo(Form("%d HV channels and %d switches",nChannels,nSwitch));
372 return nChannels+nSwitch;
375 //_____________________________________________________________________________
377 AliMUONCDB::MakePedestalStore(AliMUONV2DStore& pedestalStore, Bool_t defaultValues)
379 /// Create a pedestal store. if defaultValues=true, ped.mean=ped.sigma=1,
380 /// otherwise mean and sigma are from a gaussian (with parameters
381 /// defined below by the kPedestal* constants)
383 TIter next(fManuList);
390 Bool_t replace = kFALSE;
392 const Int_t kChannels(64);
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() ) )
401 AliMUONVCalibParam* ped = new AliMUONCalibParamNF(2,kChannels,AliMUONVCalibParam::InvalidFloatValue());
403 Int_t detElemId = p->GetFirst();
405 Int_t manuId = p->GetSecond();
407 const AliMpVSegmentation* seg =
408 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
410 for ( Int_t manuChannel = 0; manuChannel < kChannels; ++manuChannel )
412 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
413 if (!pad.IsValid()) continue;
417 Float_t meanPedestal;
418 Float_t sigmaPedestal;
428 while ( meanPedestal < 0 )
430 meanPedestal = gRandom->Gaus(kPedestalMeanMean,kPedestalMeanSigma);
433 while ( sigmaPedestal < 0 )
435 sigmaPedestal = gRandom->Gaus(kPedestalSigmaMean,kPedestalSigmaSigma);
438 ped->SetValueAsFloat(manuChannel,0,meanPedestal);
439 ped->SetValueAsFloat(manuChannel,1,sigmaPedestal);
442 Bool_t ok = pedestalStore.Set(detElemId,manuId,ped,replace);
445 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
449 AliInfo(Form("%d Manus and %d channels.",nmanus,nchannels));
454 //_____________________________________________________________________________
456 AliMUONCDB::MakeCapacitanceStore(AliMUONV1DStore& capaStore, Bool_t defaultValues)
458 /// Create a capacitance store. if defaultValues=true, all capa are 1.0,
459 /// otherwise they are from a gaussian with parameters defined in the
460 /// kCapa* constants below.
462 TIter next(fManuList);
468 Int_t nmanusOK(0); // manus for which we got the serial number
470 Bool_t replace = kFALSE;
472 const Float_t kCapaMean(1.0);
473 const Float_t kCapaSigma(0.5);
475 while ( ( p = (AliMpIntPair*)next() ) )
479 Int_t detElemId = p->GetFirst();
480 Int_t manuId = p->GetSecond();
482 AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(detElemId);
483 Int_t serialNumber = de->GetManuSerialFromId(manuId);
485 if ( serialNumber > 0 ) ++nmanusOK;
487 const AliMpVSegmentation* seg =
488 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
490 AliMUONVCalibParam* capa = static_cast<AliMUONVCalibParam*>(capaStore.Get(serialNumber));
494 capa = new AliMUONCalibParamNF(1,fgkMaxNofChannelsPerManu,1.0);
495 Bool_t ok = capaStore.Set(serialNumber,capa,replace);
498 AliError(Form("Could not set serialNumber=%d manuId=%d",serialNumber,manuId));
502 for ( Int_t manuChannel = 0; manuChannel < capa->Size(); ++manuChannel )
504 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
505 if (!pad.IsValid()) continue;
518 while ( capaValue < 0 )
520 capaValue = gRandom->Gaus(kCapaMean,kCapaSigma);
523 capa->SetValueAsFloat(manuChannel,0,capaValue);
528 if ( nmanus ) percent = 100*nmanusOK/nmanus;
529 AliInfo(Form("%5d manus with serial number (out of %5d manus = %3.0f%%)",
530 nmanusOK,nmanus,percent));
531 AliInfo(Form("%5d channels",nchannels));
534 AliWarning("Did not get all serial numbers. capaStore is incomplete !!!!");
540 //_____________________________________________________________________________
542 AliMUONCDB::MakeGainStore(AliMUONV2DStore& gainStore, Bool_t defaultValues)
544 /// Create a gain store. if defaultValues=true, all gain are 1.0,
545 /// otherwise they are from a gaussian with parameters defined in the
546 /// kGain* constants below.
548 TIter next(fManuList);
555 Bool_t replace = kFALSE;
557 const Double_t kSaturation(3000);
558 const Double_t kGainMean(1.0);
559 const Double_t kGainSigma(0.05);
561 while ( ( p = (AliMpIntPair*)next() ) )
564 AliMUONVCalibParam* gain =
565 new AliMUONCalibParamNF(2,fgkMaxNofChannelsPerManu,
566 AliMUONVCalibParam::InvalidFloatValue());
568 Int_t detElemId = p->GetFirst();
569 Int_t manuId = p->GetSecond();
571 const AliMpVSegmentation* seg =
572 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
574 for ( Int_t manuChannel = 0; manuChannel < gain->Size(); ++manuChannel )
576 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
577 if (!pad.IsValid()) continue;
582 Float_t saturation(kSaturation);
591 while ( meanGain < 0 )
593 meanGain = gRandom->Gaus(kGainMean,kGainSigma);
596 gain->SetValueAsFloat(manuChannel,0,meanGain);
597 gain->SetValueAsFloat(manuChannel,1,saturation);
600 Bool_t ok = gainStore.Set(detElemId,manuId,gain,replace);
603 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
607 AliInfo(Form("%d Manus and %d channels.",nmanus,nchannels));
611 //_____________________________________________________________________________
613 AliMUONCDB::MakeLocalTriggerMaskStore(AliMUONV1DStore& localBoardMasks) const
615 /// Generate local trigger masks store. All masks are set to FFFF
618 // Generate fake mask values for 234 localboards and put that into
619 // one single container (localBoardMasks)
620 for ( Int_t i = 1; i <= 234; ++i )
622 AliMUONVCalibParam* localBoard = new AliMUONCalibParamNI(1,8);
623 for ( Int_t x = 0; x < 2; ++x )
625 for ( Int_t y = 0; y < 4; ++y )
628 localBoard->SetValueAsInt(index,0,0xFFFF);
632 localBoardMasks.Set(i,localBoard,kFALSE);
637 //_____________________________________________________________________________
639 AliMUONCDB::MakeRegionalTriggerMaskStore(AliMUONV1DStore& rtm) const
641 /// Make a regional trigger masks store. All masks are set to 3F
644 for ( Int_t i = 0; i < 16; ++i )
646 AliMUONVCalibParam* regionalBoard = new AliMUONCalibParamNI(1,16);
647 for ( Int_t j = 0; j < 16; ++j )
649 regionalBoard->SetValueAsInt(j,0,0x3F);
652 rtm.Set(i,regionalBoard,kFALSE);
658 //_____________________________________________________________________________
660 AliMUONCDB::MakeGlobalTriggerMaskStore(AliMUONVCalibParam& gtm) const
662 /// Make a global trigger masks store. All masks set to FFF
666 for ( Int_t j = 0; j < 16; ++j )
668 gtm.SetValueAsInt(j,0,0xFFF);
674 //_____________________________________________________________________________
676 AliMUONCDB::MakeTriggerLUT(const char* file) const
678 /// Make a triggerlut object, from a file.
680 AliMUONTriggerLut* lut = new AliMUONTriggerLut;
681 lut->ReadFromFile(file);
685 //_____________________________________________________________________________
686 AliMUONTriggerEfficiencyCells*
687 AliMUONCDB::MakeTriggerEfficiency(const char* file) const
689 /// Make a trigger efficiency object from a file.
691 return new AliMUONTriggerEfficiencyCells(file);
694 //_____________________________________________________________________________
696 AliMUONCDB::WriteToCDB(const char* calibpath, TObject* object,
697 Int_t startRun, Int_t endRun, Bool_t defaultValues)
699 /// Write a given object to OCDB
701 AliCDBId id(calibpath,startRun,endRun);
703 md.SetAliRootVersion(gROOT->GetVersion());
706 md.SetComment("Test with default values");
710 md.SetComment("Test with random values");
712 md.SetResponsible("AliMUONCDB tester class");
714 AliCDBManager* man = AliCDBManager::Instance();
715 man->SetDefaultStorage(fCDBPath);
716 man->Put(object,id,&md);
719 //_____________________________________________________________________________
721 AliMUONCDB::MakeNeighbourStore(AliMUONV2DStore& neighbourStore)
723 /// Fill the neighbours store with, for each channel, a TObjArray of its
724 /// neighbouring pads (including itself)
730 TIter next(fManuList);
738 while ( ( p = (AliMpIntPair*)next() ) )
740 Int_t detElemId = p->GetFirst();
741 Int_t manuId = p->GetSecond();
743 const AliMpVSegmentation* seg =
744 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
746 AliMUONVCalibParam* calibParam = static_cast<AliMUONVCalibParam*>(neighbourStore.Get(detElemId,manuId));
751 Int_t defaultValue(-1);
752 Int_t packingFactor(64);
754 calibParam = new AliMUONCalibParamNI(dimension,size,defaultValue,packingFactor);
755 Bool_t ok = neighbourStore.Set(detElemId,manuId,calibParam,kFALSE);
758 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
763 for ( Int_t manuChannel = 0; manuChannel < 64; ++manuChannel )
765 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
771 seg->GetNeighbours(pad,tmp,true,true);
772 Int_t nofPadNeighbours = tmp.GetEntriesFast();
774 for ( Int_t i = 0; i < nofPadNeighbours; ++i )
776 AliMpPad* pad = static_cast<AliMpPad*>(tmp.At(i));
778 Bool_t ok = calibParam->PackValues(pad->GetLocation().GetFirst(),pad->GetLocation().GetSecond(),x);
781 AliError("Could not pack value. Something is seriously wrong. Please check");
784 calibParam->SetValueAsInt(manuChannel,i,x);
795 //_____________________________________________________________________________
797 AliMUONCDB::WriteLocalTriggerMasks(Int_t startRun, Int_t endRun)
799 /// Write local trigger masks to OCDB
801 AliMUONV1DStore* ltm = new AliMUON1DArray(235);
802 Int_t ngenerated = MakeLocalTriggerMaskStore(*ltm);
803 AliInfo(Form("Ngenerated = %d",ngenerated));
806 WriteToCDB("MUON/Calib/LocalTriggerBoardMasks",ltm,startRun,endRun,true);
811 //_____________________________________________________________________________
813 AliMUONCDB::WriteRegionalTriggerMasks(Int_t startRun, Int_t endRun)
815 /// Write regional trigger masks to OCDB
817 AliMUONV1DStore* rtm = new AliMUON1DArray(16);
818 Int_t ngenerated = MakeRegionalTriggerMaskStore(*rtm);
819 AliInfo(Form("Ngenerated = %d",ngenerated));
822 WriteToCDB("MUON/Calib/RegionalTriggerBoardMasks",rtm,startRun,endRun,true);
827 //_____________________________________________________________________________
829 AliMUONCDB::WriteGlobalTriggerMasks(Int_t startRun, Int_t endRun)
831 /// Write global trigger masks to OCDB
833 AliMUONVCalibParam* gtm = new AliMUONCalibParamNI(1,16);
835 Int_t ngenerated = MakeGlobalTriggerMaskStore(*gtm);
836 AliInfo(Form("Ngenerated = %d",ngenerated));
839 WriteToCDB("MUON/Calib/GlobalTriggerBoardMasks",gtm,startRun,endRun,true);
844 //_____________________________________________________________________________
846 AliMUONCDB::WriteTriggerLut(Int_t startRun, Int_t endRun)
848 /// Write trigger LUT to OCDB
850 AliMUONTriggerLut* lut = MakeTriggerLUT();
853 WriteToCDB("MUON/Calib/TriggerLut",lut,startRun,endRun,true);
858 //_____________________________________________________________________________
860 AliMUONCDB::WriteTriggerEfficiency(Int_t startRun, Int_t endRun)
862 /// Write trigger efficiency to OCDB
864 AliMUONTriggerEfficiencyCells* eff = MakeTriggerEfficiency();
867 WriteToCDB("MUON/Calib/TriggerEfficiency",eff,startRun,endRun,true);
872 //_____________________________________________________________________________
874 AliMUONCDB::WriteNeighbours(Int_t startRun, Int_t endRun)
876 /// Write neighbours to OCDB
878 AliMUONV2DStore* neighbours = new AliMUON2DMap(kTRUE);
879 Int_t ngenerated = MakeNeighbourStore(*neighbours);
880 AliInfo(Form("Ngenerated = %d",ngenerated));
883 WriteToCDB("MUON/Calib/Neighbours",neighbours,startRun,endRun,true);
888 //_____________________________________________________________________________
890 AliMUONCDB::WriteHV(Bool_t defaultValues,
891 Int_t startRun, Int_t endRun)
893 /// generate HV values (either cste = 1500 V) if defaultValues=true or random
894 /// if defaultValues=false, see makeHVStore) and
895 /// store them into CDB located at cdbpath, with a validity period
896 /// ranging from startRun to endRun
898 TMap* hvStore = new TMap;
899 Int_t ngenerated = MakeHVStore(*hvStore,defaultValues);
900 AliInfo(Form("Ngenerated = %d",ngenerated));
903 WriteToCDB("MUON/Calib/HV",hvStore,startRun,endRun,defaultValues);
908 //_____________________________________________________________________________
910 AliMUONCDB::WritePedestals(Bool_t defaultValues,
911 Int_t startRun, Int_t endRun)
913 /// generate pedestal values (either 0 if defaultValues=true or random
914 /// if defaultValues=false, see makePedestalStore) and
915 /// store them into CDB located at cdbpath, with a validity period
916 /// ranging from startRun to endRun
918 AliMUONV2DStore* pedestalStore = new AliMUON2DMap(true);
919 Int_t ngenerated = MakePedestalStore(*pedestalStore,defaultValues);
920 AliInfo(Form("Ngenerated = %d",ngenerated));
921 WriteToCDB("MUON/Calib/Pedestals",pedestalStore,startRun,endRun,defaultValues);
922 delete pedestalStore;
926 //_____________________________________________________________________________
928 AliMUONCDB::WriteGains(Bool_t defaultValues,
929 Int_t startRun, Int_t endRun)
931 /// generate gain values (either 1 if defaultValues=true or random
932 /// if defaultValues=false, see makeGainStore) and
933 /// store them into CDB located at cdbpath, with a validity period
934 /// ranging from startRun to endRun
936 AliMUONV2DStore* gainStore = new AliMUON2DMap(true);
937 Int_t ngenerated = MakeGainStore(*gainStore,defaultValues);
938 AliInfo(Form("Ngenerated = %d",ngenerated));
939 WriteToCDB("MUON/Calib/Gains",gainStore,startRun,endRun,defaultValues);
943 //_____________________________________________________________________________
945 AliMUONCDB::WriteCapacitances(Bool_t defaultValues,
946 Int_t startRun, Int_t endRun)
948 /// generate manu capacitance values (either 1 if defaultValues=true or random
949 /// if defaultValues=false, see makeCapacitanceStore) and
950 /// store them into CDB located at cdbpath, with a validity period
951 /// ranging from startRun to endRun
953 AliMUONV1DStore* capaStore = new AliMUON1DMap(16828);
954 Int_t ngenerated = MakeCapacitanceStore(*capaStore,defaultValues);
955 AliInfo(Form("Ngenerated = %d",ngenerated));
956 WriteToCDB("MUON/Calib/Capacitances",capaStore,startRun,endRun,defaultValues);
960 //_____________________________________________________________________________
962 AliMUONCDB::WriteTrigger(Int_t startRun, Int_t endRun)
964 /// Writes all Trigger related calibration to CDB
965 WriteLocalTriggerMasks(startRun,endRun);
966 WriteRegionalTriggerMasks(startRun,endRun);
967 WriteGlobalTriggerMasks(startRun,endRun);
968 WriteTriggerLut(startRun,endRun);
969 WriteTriggerEfficiency(startRun,endRun);
972 //_____________________________________________________________________________
974 AliMUONCDB::WriteTracker(Bool_t defaultValues, Int_t startRun, Int_t endRun)
976 /// Writes all Tracker related calibration to CDB
977 WriteHV(defaultValues,startRun,endRun);
978 WritePedestals(defaultValues,startRun,endRun);
979 WriteGains(defaultValues,startRun,endRun);
980 WriteCapacitances(defaultValues,startRun,endRun);
981 WriteNeighbours(startRun,endRun);