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, Int_t dim,
78 Float_t* xmin, Float_t* xmax)
80 /// Assuming the store contains AliMUONVCalibParam objects, compute the
81 /// limits of the value contained in the VCalibParam, for each of its dimensions
82 /// xmin and xmax must be of dimension dim
84 for ( Int_t i = 0; i < dim; ++i )
90 TIter next(store.CreateIterator());
91 AliMUONVCalibParam* value;
93 while ( ( value = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
95 Int_t detElemId = value->ID0();
96 Int_t manuId = value->ID1();
98 const AliMpVSegmentation* seg =
99 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
101 for ( Int_t manuChannel = 0; manuChannel < value->Size(); ++manuChannel )
103 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
104 if (!pad.IsValid()) continue;
106 for ( Int_t i = 0; i < dim; ++i )
108 Float_t x0 = value->ValueAsFloat(manuChannel,i);
110 xmin[i] = TMath::Min(xmin[i],x0);
111 xmax[i] = TMath::Max(xmax[i],x0);
116 for ( Int_t i = 0; i < dim; ++i )
118 if ( TMath::Abs(xmin[i]-xmax[i]) < 1E-3 )
126 //_____________________________________________________________________________
127 Double_t GetRandom(Double_t mean, Double_t sigma, Bool_t mustBePositive)
130 if ( mustBePositive )
134 x = gRandom->Gaus(mean,sigma);
139 x = gRandom->Gaus(mean,sigma);
146 //_____________________________________________________________________________
147 AliMUONCDB::AliMUONCDB(const char* cdbpath)
151 fMaxNofChannelsToGenerate(-1)
156 //_____________________________________________________________________________
157 AliMUONCDB::~AliMUONCDB()
163 //_____________________________________________________________________________
165 AliMUONCDB::ManuList()
167 /// return (and create if necessary) the list of (de,manu) pairs
170 AliInfo("Generating ManuList...");
171 fManuList = AliMpManuList::ManuList();
172 AliInfo("Manu List generated.");
177 //_____________________________________________________________________________
179 AliMUONCDB::Diff(AliMUONVStore& store1, AliMUONVStore& store2,
182 /// creates a store which contains store1-store2
183 /// if opt="abs" the difference is absolute one,
184 /// if opt="rel" then what is stored is (store1-store2)/store1
185 /// if opt="percent" then what is stored is rel*100
187 /// WARNING Works only for stores which holds AliMUONVCalibParam objects
192 if ( !sopt.Contains("ABS") && !sopt.Contains("REL") && !sopt.Contains("PERCENT") )
194 AliErrorClass(Form("opt %s not supported. Only ABS, REL, PERCENT are",opt));
198 AliMUONVStore* d = static_cast<AliMUONVStore*>(store1.Clone());
200 TIter next(d->CreateIterator());
202 AliMUONVCalibParam* param;
204 while ( ( param = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
206 Int_t detElemId = param->ID0();
207 Int_t manuId = param->ID1();
209 AliMUONVCalibParam* param2 = dynamic_cast<AliMUONVCalibParam*>(store2.FindObject(detElemId,manuId));
210 //FIXME: this might happen. Handle it.
213 cerr << "param2 is null : FIXME : this might happen !" << endl;
218 for ( Int_t i = 0; i < param->Size(); ++i )
220 for ( Int_t j = 0; j < param->Dimension(); ++j )
223 if ( sopt.Contains("ABS") )
225 value = param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j);
227 else if ( sopt.Contains("REL") || sopt.Contains("PERCENT") )
229 if ( param->ValueAsFloat(i,j) )
231 value = (param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j))/param->ValueAsFloat(i,j);
237 if ( sopt.Contains("PERCENT") ) value *= 100.0;
239 param->SetValueAsFloat(i,j,value);
246 //_____________________________________________________________________________
248 AliMUONCDB::Plot(const AliMUONVStore& store, const char* name, Int_t nbins)
250 /// Make histograms of each dimension of the AliMUONVCalibParam
251 /// contained inside store.
252 /// It produces histograms named name_0, name_1, etc...
254 TIter next(store.CreateIterator());
255 AliMUONVCalibParam* param;
257 const Int_t kNStations = AliMpConstants::NofTrackingChambers()/2;
258 Int_t* nPerStation = new Int_t[kNStations];
261 for ( Int_t i = 0; i < kNStations; ++i ) nPerStation[i]=0;
263 while ( ( param = static_cast<AliMUONVCalibParam*>(next()) ) )
267 Int_t dim = param->Dimension();
269 Float_t* xmin = new Float_t[dim];
270 Float_t* xmax = new Float_t[dim];
271 getBoundaries(store,dim,xmin,xmax);
273 for ( Int_t i = 0; i < dim; ++i )
275 h[i] = new TH1F(Form("%s_%d",name,i),Form("%s_%d",name,i),
276 nbins,xmin[i],xmax[i]);
277 AliInfo(Form("Created histogram %s",h[i]->GetName()));
281 Int_t detElemId = param->ID0();
282 Int_t manuId = param->ID1();
283 Int_t station = AliMpDEManager::GetChamberId(detElemId)/2;
285 const AliMpVSegmentation* seg =
286 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
288 for ( Int_t manuChannel = 0; manuChannel < param->Size(); ++manuChannel )
290 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
291 if (!pad.IsValid()) continue;
294 ++nPerStation[station];
296 for ( Int_t dim = 0; dim < param->Dimension(); ++dim )
298 h[dim]->Fill(param->ValueAsFloat(manuChannel,dim));
303 for ( Int_t i = 0; i < kNStations; ++i )
305 AliInfo(Form("Station %d %d ",(i+1),nPerStation[i]));
308 AliInfo(Form("Number of channels = %d",n));
310 delete[] nPerStation;
313 //_____________________________________________________________________________
315 AliMUONCDB::MakeHVStore(TMap& aliasMap, Bool_t defaultValues)
317 /// Create a HV store
319 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 = gRandom->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 = GetRandom(1750,62.5,true);
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(AliMUONVStore& 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(ManuList());
390 const Int_t kChannels(AliMpConstants::ManuNofChannels());
391 const Float_t kPedestalMeanMean(200);
392 const Float_t kPedestalMeanSigma(10);
393 const Float_t kPedestalSigmaMean(1.0);
394 const Float_t kPedestalSigmaSigma(0.2);
396 while ( ( p = (AliMpIntPair*)next() ) )
400 Int_t detElemId = p->GetFirst();
401 Int_t manuId = p->GetSecond();
403 AliMUONVCalibParam* ped =
404 new AliMUONCalibParamNF(2,kChannels,detElemId,manuId,AliMUONVCalibParam::InvalidFloatValue());
406 const AliMpVSegmentation* seg =
407 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
409 for ( Int_t manuChannel = 0; manuChannel < kChannels; ++manuChannel )
411 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
412 if (!pad.IsValid()) continue;
416 Float_t meanPedestal;
417 Float_t sigmaPedestal;
426 Bool_t positive(kTRUE);
428 while ( meanPedestal == 0.0 ) // avoid strict zero
430 meanPedestal = GetRandom(kPedestalMeanMean,kPedestalMeanSigma,positive);
432 sigmaPedestal = GetRandom(kPedestalSigmaMean,kPedestalSigmaSigma,positive);
434 ped->SetValueAsFloat(manuChannel,0,meanPedestal);
435 ped->SetValueAsFloat(manuChannel,1,sigmaPedestal);
438 Bool_t ok = pedestalStore.Add(ped);
441 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
443 if ( fMaxNofChannelsToGenerate > 0 && nchannels >= fMaxNofChannelsToGenerate ) break;
446 AliInfo(Form("%d Manus and %d channels.",nmanus,nchannels));
451 //_____________________________________________________________________________
453 AliMUONCDB::MakeCapacitanceStore(AliMUONVStore& capaStore, Bool_t defaultValues)
455 /// Create a capacitance store. if defaultValues=true, all capa are 1.0,
456 /// otherwise they are from a gaussian with parameters defined in the
457 /// kCapa* constants below.
459 TIter next(ManuList());
465 Int_t nmanusOK(0); // manus for which we got the serial number
467 const Float_t kCapaMean(1.0);
468 const Float_t kCapaSigma(0.5);
470 while ( ( p = (AliMpIntPair*)next() ) )
474 Int_t detElemId = p->GetFirst();
475 Int_t manuId = p->GetSecond();
477 AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(detElemId);
478 Int_t serialNumber = de->GetManuSerialFromId(manuId);
480 if ( serialNumber <= 0 ) continue;
484 const AliMpVSegmentation* seg =
485 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
487 AliMUONVCalibParam* capa = static_cast<AliMUONVCalibParam*>(capaStore.FindObject(serialNumber));
491 capa = new AliMUONCalibParamNF(1,AliMpConstants::ManuNofChannels(),serialNumber,0,1.0);
492 Bool_t ok = capaStore.Add(capa);
495 AliError(Form("Could not set serialNumber=%d manuId=%d",serialNumber,manuId));
499 for ( Int_t manuChannel = 0; manuChannel < capa->Size(); ++manuChannel )
501 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
502 if (!pad.IsValid()) continue;
514 capaValue = GetRandom(kCapaMean,kCapaSigma,kTRUE);
516 capa->SetValueAsFloat(manuChannel,0,capaValue);
521 if ( nmanus ) percent = 100*nmanusOK/nmanus;
522 AliInfo(Form("%5d manus with serial number (out of %5d manus = %3.0f%%)",
523 nmanusOK,nmanus,percent));
524 AliInfo(Form("%5d channels",nchannels));
527 AliWarning("Did not get all serial numbers. capaStore is incomplete !!!!");
533 //_____________________________________________________________________________
535 AliMUONCDB::MakeGainStore(AliMUONVStore& gainStore, Bool_t defaultValues)
537 /// Create a gain store. if defaultValues=true, all gains set so that
538 /// charge = (adc-ped)
540 /// otherwise parameters are taken from gaussians with parameters
541 /// defined in the k* constants below.
543 TIter next(ManuList());
550 const Int_t kSaturation(3000);
551 const Double_t kA0Mean(1.2);
552 const Double_t kA0Sigma(0.1);
553 const Double_t kA1Mean(1E-5);
554 const Double_t kA1Sigma(1E-6);
555 const Double_t kQualMean(0xFF);
556 const Double_t kQualSigma(0x10);
557 const Int_t kThresMean(1600);
558 const Int_t kThresSigma(100);
560 while ( ( p = (AliMpIntPair*)next() ) )
564 Int_t detElemId = p->GetFirst();
565 Int_t manuId = p->GetSecond();
567 AliMUONVCalibParam* gain =
568 new AliMUONCalibParamNF(5,AliMpConstants::ManuNofChannels(),
571 AliMUONVCalibParam::InvalidFloatValue());
574 const AliMpVSegmentation* seg =
575 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
577 for ( Int_t manuChannel = 0; manuChannel < gain->Size(); ++manuChannel )
579 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
580 if (!pad.IsValid()) continue;
586 gain->SetValueAsFloat(manuChannel,0,1.0);
587 gain->SetValueAsFloat(manuChannel,1,0.0);
588 gain->SetValueAsInt(manuChannel,2,4095);
589 gain->SetValueAsInt(manuChannel,3,1);
590 gain->SetValueAsInt(manuChannel,4,kSaturation);
594 Bool_t positive(kTRUE);
595 gain->SetValueAsFloat(manuChannel,0,GetRandom(kA0Mean,kA0Sigma,positive));
596 gain->SetValueAsFloat(manuChannel,1,GetRandom(kA1Mean,kA1Sigma,!positive));
597 gain->SetValueAsInt(manuChannel,2,(Int_t)TMath::Nint(GetRandom(kThresMean,kThresSigma,positive)));
598 gain->SetValueAsInt(manuChannel,3,(Int_t)TMath::Nint(GetRandom(kQualMean,kQualSigma,positive)));
599 gain->SetValueAsInt(manuChannel,4,kSaturation);
603 Bool_t ok = gainStore.Add(gain);
606 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
608 if ( fMaxNofChannelsToGenerate > 0 && nchannels >= fMaxNofChannelsToGenerate ) break;
611 AliInfo(Form("%d Manus and %d channels.",nmanus,nchannels));
615 //_____________________________________________________________________________
617 AliMUONCDB::MakeLocalTriggerMaskStore(AliMUONVStore& localBoardMasks) const
619 /// Generate local trigger masks store. All masks are set to FFFF
622 // Generate fake mask values for 234 localboards and put that into
623 // one single container (localBoardMasks)
624 for ( Int_t i = 1; i <= 234; ++i )
626 AliMUONVCalibParam* localBoard = new AliMUONCalibParamNI(1,8,i,0,0);
627 for ( Int_t x = 0; x < 2; ++x )
629 for ( Int_t y = 0; y < 4; ++y )
632 localBoard->SetValueAsInt(index,0,0xFFFF);
636 localBoardMasks.Add(localBoard);
641 //_____________________________________________________________________________
643 AliMUONCDB::MakeRegionalTriggerMaskStore(AliMUONVStore& rtm) const
645 /// Make a regional trigger masks store. All masks are set to 3F
648 for ( Int_t i = 0; i < 16; ++i )
650 AliMUONVCalibParam* regionalBoard = new AliMUONCalibParamNI(1,16,i,0,0);
651 for ( Int_t j = 0; j < 16; ++j )
653 regionalBoard->SetValueAsInt(j,0,0x3F);
656 rtm.Add(regionalBoard);
662 //_____________________________________________________________________________
664 AliMUONCDB::MakeGlobalTriggerMaskStore(AliMUONVCalibParam& gtm) const
666 /// Make a global trigger masks store. All masks set to FFF
670 for ( Int_t j = 0; j < 16; ++j )
672 gtm.SetValueAsInt(j,0,0xFFF);
678 //_____________________________________________________________________________
680 AliMUONCDB::MakeTriggerLUT(const char* file) const
682 /// Make a triggerlut object, from a file.
684 AliMUONTriggerLut* lut = new AliMUONTriggerLut;
685 lut->ReadFromFile(file);
689 //_____________________________________________________________________________
690 AliMUONTriggerEfficiencyCells*
691 AliMUONCDB::MakeTriggerEfficiency(const char* file) const
693 /// Make a trigger efficiency object from a file.
695 return new AliMUONTriggerEfficiencyCells(file);
698 //_____________________________________________________________________________
700 AliMUONCDB::WriteToCDB(const char* calibpath, TObject* object,
701 Int_t startRun, Int_t endRun, Bool_t defaultValues)
703 /// Write a given object to OCDB
705 AliCDBId id(calibpath,startRun,endRun);
707 md.SetAliRootVersion(gROOT->GetVersion());
710 md.SetComment("Test with default values");
714 md.SetComment("Test with random values");
716 md.SetResponsible("AliMUONCDB tester class");
718 AliCDBManager* man = AliCDBManager::Instance();
719 man->SetDefaultStorage(fCDBPath);
720 man->Put(object,id,&md);
723 //_____________________________________________________________________________
725 AliMUONCDB::MakeNeighbourStore(AliMUONVStore& neighbourStore)
727 /// Fill the neighbours store with, for each channel, a TObjArray of its
728 /// neighbouring pads (including itself)
730 AliInfo("Generating NeighbourStore. This will take a while. Please be patient.");
736 TIter next(ManuList());
744 while ( ( p = (AliMpIntPair*)next() ) )
746 Int_t detElemId = p->GetFirst();
747 Int_t manuId = p->GetSecond();
749 const AliMpVSegmentation* seg =
750 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
752 AliMUONVCalibParam* calibParam = static_cast<AliMUONVCalibParam*>(neighbourStore.FindObject(detElemId,manuId));
756 Int_t size(AliMpConstants::ManuNofChannels());
757 Int_t defaultValue(-1);
758 Int_t packingFactor(size);
760 calibParam = new AliMUONCalibParamNI(dimension,size,detElemId,manuId,defaultValue,packingFactor);
761 Bool_t ok = neighbourStore.Add(calibParam);
764 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
769 for ( Int_t manuChannel = 0; manuChannel < AliMpConstants::ManuNofChannels(); ++manuChannel )
771 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
777 seg->GetNeighbours(pad,tmp,true,true);
778 Int_t nofPadNeighbours = tmp.GetEntriesFast();
780 for ( Int_t i = 0; i < nofPadNeighbours; ++i )
782 AliMpPad* pad = static_cast<AliMpPad*>(tmp.At(i));
784 Bool_t ok = calibParam->PackValues(pad->GetLocation().GetFirst(),pad->GetLocation().GetSecond(),x);
787 AliError("Could not pack value. Something is seriously wrong. Please check");
788 StdoutToAliError(pad->Print(););
791 calibParam->SetValueAsInt(manuChannel,i,x);
802 //_____________________________________________________________________________
804 AliMUONCDB::SetMaxNofChannelsToGenerate(Int_t n)
806 /// Set the maximum number of channels to generate (used for testing only)
807 /// n < 0 means no limit
808 fMaxNofChannelsToGenerate = n;
811 //_____________________________________________________________________________
813 AliMUONCDB::WriteLocalTriggerMasks(Int_t startRun, Int_t endRun)
815 /// Write local trigger masks to OCDB
817 AliMUONVStore* ltm = new AliMUON1DArray(235);
818 Int_t ngenerated = MakeLocalTriggerMaskStore(*ltm);
819 AliInfo(Form("Ngenerated = %d",ngenerated));
822 WriteToCDB("MUON/Calib/LocalTriggerBoardMasks",ltm,startRun,endRun,true);
827 //_____________________________________________________________________________
829 AliMUONCDB::WriteRegionalTriggerMasks(Int_t startRun, Int_t endRun)
831 /// Write regional trigger masks to OCDB
833 AliMUONVStore* rtm = new AliMUON1DArray(16);
834 Int_t ngenerated = MakeRegionalTriggerMaskStore(*rtm);
835 AliInfo(Form("Ngenerated = %d",ngenerated));
838 WriteToCDB("MUON/Calib/RegionalTriggerBoardMasks",rtm,startRun,endRun,true);
843 //_____________________________________________________________________________
845 AliMUONCDB::WriteGlobalTriggerMasks(Int_t startRun, Int_t endRun)
847 /// Write global trigger masks to OCDB
849 AliMUONVCalibParam* gtm = new AliMUONCalibParamNI(1,16,1,0,0);
851 Int_t ngenerated = MakeGlobalTriggerMaskStore(*gtm);
852 AliInfo(Form("Ngenerated = %d",ngenerated));
855 WriteToCDB("MUON/Calib/GlobalTriggerBoardMasks",gtm,startRun,endRun,true);
860 //_____________________________________________________________________________
862 AliMUONCDB::WriteTriggerLut(Int_t startRun, Int_t endRun)
864 /// Write trigger LUT to OCDB
866 AliMUONTriggerLut* lut = MakeTriggerLUT();
869 WriteToCDB("MUON/Calib/TriggerLut",lut,startRun,endRun,true);
874 //_____________________________________________________________________________
876 AliMUONCDB::WriteTriggerEfficiency(Int_t startRun, Int_t endRun)
878 /// Write trigger efficiency to OCDB
880 AliMUONTriggerEfficiencyCells* eff = MakeTriggerEfficiency();
883 WriteToCDB("MUON/Calib/TriggerEfficiency",eff,startRun,endRun,true);
888 //_____________________________________________________________________________
890 AliMUONCDB::WriteNeighbours(Int_t startRun, Int_t endRun)
892 /// Write neighbours to OCDB
894 AliMUONVStore* neighbours = new AliMUON2DMap(kTRUE);
895 Int_t ngenerated = MakeNeighbourStore(*neighbours);
896 AliInfo(Form("Ngenerated = %d",ngenerated));
899 WriteToCDB("MUON/Calib/Neighbours",neighbours,startRun,endRun,true);
904 //_____________________________________________________________________________
906 AliMUONCDB::WriteHV(Bool_t defaultValues,
907 Int_t startRun, Int_t endRun)
909 /// generate HV values (either cste = 1500 V) if defaultValues=true or random
910 /// if defaultValues=false, see makeHVStore) and
911 /// store them into CDB located at cdbpath, with a validity period
912 /// ranging from startRun to endRun
914 TMap* hvStore = new TMap;
915 Int_t ngenerated = MakeHVStore(*hvStore,defaultValues);
916 AliInfo(Form("Ngenerated = %d",ngenerated));
919 WriteToCDB("MUON/Calib/HV",hvStore,startRun,endRun,defaultValues);
924 //_____________________________________________________________________________
926 AliMUONCDB::WritePedestals(Bool_t defaultValues,
927 Int_t startRun, Int_t endRun)
929 /// generate pedestal values (either 0 if defaultValues=true or random
930 /// if defaultValues=false, see makePedestalStore) and
931 /// store them into CDB located at cdbpath, with a validity period
932 /// ranging from startRun to endRun
934 AliMUONVStore* pedestalStore = new AliMUON2DMap(true);
935 Int_t ngenerated = MakePedestalStore(*pedestalStore,defaultValues);
936 AliInfo(Form("Ngenerated = %d",ngenerated));
937 WriteToCDB("MUON/Calib/Pedestals",pedestalStore,startRun,endRun,defaultValues);
938 delete pedestalStore;
942 //_____________________________________________________________________________
944 AliMUONCDB::WriteGains(Bool_t defaultValues,
945 Int_t startRun, Int_t endRun)
947 /// generate gain values (either 1 if defaultValues=true or random
948 /// if defaultValues=false, see makeGainStore) and
949 /// store them into CDB located at cdbpath, with a validity period
950 /// ranging from startRun to endRun
952 AliMUONVStore* gainStore = new AliMUON2DMap(true);
953 Int_t ngenerated = MakeGainStore(*gainStore,defaultValues);
954 AliInfo(Form("Ngenerated = %d",ngenerated));
955 WriteToCDB("MUON/Calib/Gains",gainStore,startRun,endRun,defaultValues);
959 //_____________________________________________________________________________
961 AliMUONCDB::WriteCapacitances(Bool_t defaultValues,
962 Int_t startRun, Int_t endRun)
964 /// generate manu capacitance values (either 1 if defaultValues=true or random
965 /// if defaultValues=false, see makeCapacitanceStore) and
966 /// store them into CDB located at cdbpath, with a validity period
967 /// ranging from startRun to endRun
969 AliMUONVStore* capaStore = new AliMUON1DMap(16828);
970 Int_t ngenerated = MakeCapacitanceStore(*capaStore,defaultValues);
971 AliInfo(Form("Ngenerated = %d",ngenerated));
972 WriteToCDB("MUON/Calib/Capacitances",capaStore,startRun,endRun,defaultValues);
976 //_____________________________________________________________________________
978 AliMUONCDB::WriteTrigger(Int_t startRun, Int_t endRun)
980 /// Writes all Trigger related calibration to CDB
981 WriteLocalTriggerMasks(startRun,endRun);
982 WriteRegionalTriggerMasks(startRun,endRun);
983 WriteGlobalTriggerMasks(startRun,endRun);
984 WriteTriggerLut(startRun,endRun);
985 WriteTriggerEfficiency(startRun,endRun);
988 //_____________________________________________________________________________
990 AliMUONCDB::WriteTracker(Bool_t defaultValues, Int_t startRun, Int_t endRun)
992 /// Writes all Tracker related calibration to CDB
993 WriteHV(defaultValues,startRun,endRun);
994 WritePedestals(defaultValues,startRun,endRun);
995 WriteGains(defaultValues,startRun,endRun);
996 WriteCapacitances(defaultValues,startRun,endRun);
997 WriteNeighbours(startRun,endRun);