1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 //-----------------------------------------------------------------------------
21 /// Helper class to experience the OCDB
22 /// It allows to generate dummy (but complete) containers for all the
23 /// calibration data types we have for tracker and trigger, and to write
26 /// For more information, please see READMECDB
28 // \author Laurent Aphecetche
29 //-----------------------------------------------------------------------------
31 #include "AliMUONCDB.h"
33 #include "AliMUON1DArray.h"
34 #include "AliMUON1DMap.h"
35 #include "AliMUON2DMap.h"
36 #include "AliMUON2DStoreValidator.h"
37 #include "AliMUONCalibParamNF.h"
38 #include "AliMUONCalibParamNI.h"
39 #include "AliMUONConstants.h"
40 #include "AliMUONTrackerIO.h"
41 #include "AliMUONTriggerEfficiencyCells.h"
42 #include "AliMUONTriggerLut.h"
43 #include "AliMUONVStore.h"
44 #include "AliMUONVCalibParam.h"
45 #include "AliMUONVCalibParam.h"
46 #include "AliMUONGlobalCrateConfig.h"
47 #include "AliMUONRegionalTriggerConfig.h"
50 #include "AliMpConstants.h"
51 #include "AliMpDDLStore.h"
52 #include "AliMpDEManager.h"
53 #include "AliMpDetElement.h"
54 #include "AliMpFiles.h"
55 #include "AliMpHVNamer.h"
56 #include "AliMpManuIterator.h"
57 #include "AliMpSegmentation.h"
58 #include "AliMpStationType.h"
59 #include "AliMpVSegmentation.h"
61 #include "AliCodeTimer.h"
62 #include "AliCDBEntry.h"
63 #include "AliCDBManager.h"
64 #include "AliDCSValue.h"
67 #include <Riostream.h>
73 #include <TObjString.h>
76 #include <TStopwatch.h>
87 //_____________________________________________________________________________
88 AliMUONVStore* Create2DMap()
90 return new AliMUON2DMap(true);
93 //_____________________________________________________________________________
94 void getBoundaries(const AliMUONVStore& store, Int_t dim,
95 Float_t* xmin, Float_t* xmax)
97 /// Assuming the store contains AliMUONVCalibParam objects, compute the
98 /// limits of the value contained in the VCalibParam, for each of its dimensions
99 /// xmin and xmax must be of dimension dim
101 for ( Int_t i = 0; i < dim; ++i )
107 TIter next(store.CreateIterator());
108 AliMUONVCalibParam* value;
110 while ( ( value = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
112 Int_t detElemId = value->ID0();
113 Int_t manuId = value->ID1();
115 const AliMpVSegmentation* seg =
116 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
118 for ( Int_t manuChannel = 0; manuChannel < value->Size(); ++manuChannel )
120 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
121 if (!pad.IsValid()) continue;
123 for ( Int_t i = 0; i < dim; ++i )
125 Float_t x0 = value->ValueAsFloat(manuChannel,i);
127 xmin[i] = TMath::Min(xmin[i],x0);
128 xmax[i] = TMath::Max(xmax[i],x0);
133 for ( Int_t i = 0; i < dim; ++i )
135 if ( TMath::Abs(xmin[i]-xmax[i]) < 1E-3 )
143 //_____________________________________________________________________________
144 Double_t GetRandom(Double_t mean, Double_t sigma, Bool_t mustBePositive)
147 if ( mustBePositive )
151 x = gRandom->Gaus(mean,sigma);
156 x = gRandom->Gaus(mean,sigma);
163 //_____________________________________________________________________________
164 AliMUONCDB::AliMUONCDB(const char* cdbpath)
167 fMaxNofChannelsToGenerate(-1)
171 if ( ! AliMpCDB::LoadDDLStore() ) {
172 AliFatal("Could not access mapping from OCDB !");
176 //_____________________________________________________________________________
177 AliMUONCDB::~AliMUONCDB()
182 //_____________________________________________________________________________
184 AliMUONCDB::Diff(AliMUONVStore& store1, AliMUONVStore& store2,
187 /// creates a store which contains store1-store2
188 /// if opt="abs" the difference is absolute one,
189 /// if opt="rel" then what is stored is (store1-store2)/store1
190 /// if opt="percent" then what is stored is rel*100
192 /// WARNING Works only for stores which holds AliMUONVCalibParam objects
197 if ( !sopt.Contains("ABS") && !sopt.Contains("REL") && !sopt.Contains("PERCENT") )
199 AliErrorClass(Form("opt %s not supported. Only ABS, REL, PERCENT are",opt));
203 AliMUONVStore* d = static_cast<AliMUONVStore*>(store1.Clone());
205 TIter next(d->CreateIterator());
207 AliMUONVCalibParam* param;
209 while ( ( param = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
211 Int_t detElemId = param->ID0();
212 Int_t manuId = param->ID1();
214 AliMUONVCalibParam* param2 = dynamic_cast<AliMUONVCalibParam*>(store2.FindObject(detElemId,manuId));
215 //FIXME: this might happen. Handle it.
218 cerr << "param2 is null : FIXME : this might happen !" << endl;
223 for ( Int_t i = 0; i < param->Size(); ++i )
225 for ( Int_t j = 0; j < param->Dimension(); ++j )
228 if ( sopt.Contains("ABS") )
230 value = param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j);
232 else if ( sopt.Contains("REL") || sopt.Contains("PERCENT") )
234 if ( param->ValueAsFloat(i,j) )
236 value = (param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j))/param->ValueAsFloat(i,j);
242 if ( sopt.Contains("PERCENT") ) value *= 100.0;
244 param->SetValueAsFloat(i,j,value);
251 //_____________________________________________________________________________
253 AliMUONCDB::Plot(const AliMUONVStore& store, const char* name, Int_t nbins)
255 /// Make histograms of each dimension of the AliMUONVCalibParam
256 /// contained inside store.
257 /// It produces histograms named name_0, name_1, etc...
259 TIter next(store.CreateIterator());
260 AliMUONVCalibParam* param;
262 const Int_t kNStations = AliMpConstants::NofTrackingChambers()/2;
263 Int_t* nPerStation = new Int_t[kNStations];
266 for ( Int_t i = 0; i < kNStations; ++i ) nPerStation[i]=0;
268 while ( ( param = static_cast<AliMUONVCalibParam*>(next()) ) )
272 Int_t dim = param->Dimension();
274 Float_t* xmin = new Float_t[dim];
275 Float_t* xmax = new Float_t[dim];
276 getBoundaries(store,dim,xmin,xmax);
278 for ( Int_t i = 0; i < dim; ++i )
280 h[i] = new TH1F(Form("%s_%d",name,i),Form("%s_%d",name,i),
281 nbins,xmin[i],xmax[i]);
282 AliInfo(Form("Created histogram %s",h[i]->GetName()));
286 Int_t detElemId = param->ID0();
287 Int_t manuId = param->ID1();
288 Int_t station = AliMpDEManager::GetChamberId(detElemId)/2;
290 const AliMpVSegmentation* seg =
291 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
293 for ( Int_t manuChannel = 0; manuChannel < param->Size(); ++manuChannel )
295 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
296 if (!pad.IsValid()) continue;
299 ++nPerStation[station];
301 for ( Int_t dim = 0; dim < param->Dimension(); ++dim )
303 h[dim]->Fill(param->ValueAsFloat(manuChannel,dim));
308 for ( Int_t i = 0; i < kNStations; ++i )
310 AliInfo(Form("Station %d %d ",(i+1),nPerStation[i]));
313 AliInfo(Form("Number of channels = %d",n));
315 delete[] nPerStation;
318 //_____________________________________________________________________________
320 AliMUONCDB::MakeHVStore(TMap& aliasMap, Bool_t defaultValues)
322 /// Create a HV store
324 AliMpHVNamer hvNamer;
326 TObjArray* aliases = hvNamer.GenerateAliases();
331 for ( Int_t i = 0; i < aliases->GetEntries(); ++i )
333 TObjString* alias = static_cast<TObjString*>(aliases->At(i));
334 TString& aliasName = alias->String();
335 if ( aliasName.Contains("sw") )
337 // HV Switch (St345 only)
338 TObjArray* valueSet = new TObjArray;
339 valueSet->SetOwner(kTRUE);
341 Bool_t value = kTRUE;
345 Float_t r = gRandom->Uniform();
346 if ( r < 0.007 ) value = kFALSE;
349 for ( UInt_t timeStamp = 0; timeStamp < 60*3; timeStamp += 60 )
351 AliDCSValue* dcsValue = new AliDCSValue(value,timeStamp);
352 valueSet->Add(dcsValue);
354 aliasMap.Add(new TObjString(*alias),valueSet);
359 TObjArray* valueSet = new TObjArray;
360 valueSet->SetOwner(kTRUE);
361 for ( UInt_t timeStamp = 0; timeStamp < 60*15; timeStamp += 120 )
363 Float_t value = 1500;
364 if (!defaultValues) value = GetRandom(1750,62.5,true);
365 AliDCSValue* dcsValue = new AliDCSValue(value,timeStamp);
366 valueSet->Add(dcsValue);
368 aliasMap.Add(new TObjString(*alias),valueSet);
375 AliInfo(Form("%d HV channels and %d switches",nChannels,nSwitch));
377 return nChannels+nSwitch;
380 //_____________________________________________________________________________
382 AliMUONCDB::MakePedestalStore(AliMUONVStore& pedestalStore, Bool_t defaultValues)
384 /// Create a pedestal store. if defaultValues=true, ped.mean=ped.sigma=1,
385 /// otherwise mean and sigma are from a gaussian (with parameters
386 /// defined below by the kPedestal* constants)
388 AliCodeTimerAuto("");
393 const Int_t kChannels(AliMpConstants::ManuNofChannels());
394 const Float_t kPedestalMeanMean(200);
395 const Float_t kPedestalMeanSigma(10);
396 const Float_t kPedestalSigmaMean(1.0);
397 const Float_t kPedestalSigmaSigma(0.2);
402 AliMpManuIterator it;
404 while ( it.Next(detElemId,manuId) )
408 AliMUONVCalibParam* ped =
409 new AliMUONCalibParamNF(2,kChannels,detElemId,manuId,AliMUONVCalibParam::InvalidFloatValue());
411 AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(detElemId);
413 for ( Int_t manuChannel = 0; manuChannel < kChannels; ++manuChannel )
415 if ( ! de->IsConnectedChannel(manuId,manuChannel) ) continue;
419 Float_t meanPedestal;
420 Float_t sigmaPedestal;
429 Bool_t positive(kTRUE);
431 while ( meanPedestal == 0.0 ) // avoid strict zero
433 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));
446 if ( fMaxNofChannelsToGenerate > 0 && nchannels >= fMaxNofChannelsToGenerate ) break;
449 AliInfo(Form("%d Manus and %d channels.",nmanus,nchannels));
453 //_____________________________________________________________________________
455 AliMUONCDB::MakeCapacitanceStore(AliMUONVStore& capaStore, const char* file)
457 /// Read the capacitance values from file and append them to the capaStore
459 return AliMUONTrackerIO::ReadCapacitances(file,capaStore);
462 //_____________________________________________________________________________
464 AliMUONCDB::MakeCapacitanceStore(AliMUONVStore& capaStore, Bool_t defaultValues)
466 /// Create a capacitance store. if defaultValues=true, all capa are 1.0,
467 /// otherwise they are from a gaussian with parameters defined in the
468 /// kCapa* constants below.
470 AliCodeTimerAuto("");
474 Int_t nmanusOK(0); // manus for which we got the serial number
476 const Float_t kCapaMean(0.3);
477 const Float_t kCapaSigma(0.1);
478 const Float_t kInjectionGainMean(3);
479 const Float_t kInjectionGainSigma(1);
484 AliMpManuIterator it;
486 while ( it.Next(detElemId,manuId) )
490 AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(detElemId);
491 Int_t serialNumber = de->GetManuSerialFromId(manuId);
493 if ( serialNumber <= 0 ) continue;
497 AliMUONVCalibParam* capa = static_cast<AliMUONVCalibParam*>(capaStore.FindObject(serialNumber));
501 capa = new AliMUONCalibParamNF(2,AliMpConstants::ManuNofChannels(),serialNumber,0,1.0);
502 Bool_t ok = capaStore.Add(capa);
505 AliError(Form("Could not set serialNumber=%d manuId=%d",serialNumber,manuId));
509 for ( Int_t manuChannel = 0; manuChannel < capa->Size(); ++manuChannel )
511 if ( ! de->IsConnectedChannel(manuId,manuChannel) ) continue;
516 Float_t injectionGain;
525 capaValue = GetRandom(kCapaMean,kCapaSigma,kTRUE);
526 injectionGain = GetRandom(kInjectionGainMean,kInjectionGainSigma,kTRUE);
528 capa->SetValueAsFloat(manuChannel,0,capaValue);
529 capa->SetValueAsFloat(manuChannel,1,injectionGain);
534 if ( nmanus ) percent = 100*nmanusOK/nmanus;
535 AliInfo(Form("%5d manus with serial number (out of %5d manus = %3.0f%%)",
536 nmanusOK,nmanus,percent));
537 AliInfo(Form("%5d channels",nchannels));
540 AliWarning("Did not get all serial numbers. capaStore is incomplete !!!!");
546 //_____________________________________________________________________________
548 AliMUONCDB::MakeGainStore(AliMUONVStore& gainStore, Bool_t defaultValues)
550 /// Create a gain store. if defaultValues=true, all gains set so that
551 /// charge = (adc-ped)
553 /// otherwise parameters are taken from gaussians with parameters
554 /// defined in the k* constants below.
556 AliCodeTimerAuto("");
561 const Int_t kSaturation(3000);
562 const Double_t kA0Mean(1.2);
563 const Double_t kA0Sigma(0.1);
564 const Double_t kA1Mean(1E-5);
565 const Double_t kA1Sigma(1E-6);
566 const Double_t kQualMean(0xFF);
567 const Double_t kQualSigma(0x10);
568 const Int_t kThresMean(1600);
569 const Int_t kThresSigma(100);
574 AliMpManuIterator it;
576 while ( it.Next(detElemId,manuId) )
580 AliMUONVCalibParam* gain =
581 new AliMUONCalibParamNF(5,AliMpConstants::ManuNofChannels(),
584 AliMUONVCalibParam::InvalidFloatValue());
586 AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(detElemId);
588 for ( Int_t manuChannel = 0; manuChannel < gain->Size(); ++manuChannel )
590 if ( ! de->IsConnectedChannel(manuId,manuChannel) ) continue;
596 gain->SetValueAsFloat(manuChannel,0,1.0);
597 gain->SetValueAsFloat(manuChannel,1,0.0);
598 gain->SetValueAsInt(manuChannel,2,4095);
599 gain->SetValueAsInt(manuChannel,3,1);
600 gain->SetValueAsInt(manuChannel,4,kSaturation);
604 Bool_t positive(kTRUE);
605 gain->SetValueAsFloat(manuChannel,0,GetRandom(kA0Mean,kA0Sigma,positive));
606 gain->SetValueAsFloat(manuChannel,1,GetRandom(kA1Mean,kA1Sigma,!positive));
607 gain->SetValueAsInt(manuChannel,2,(Int_t)TMath::Nint(GetRandom(kThresMean,kThresSigma,positive)));
608 gain->SetValueAsInt(manuChannel,3,(Int_t)TMath::Nint(GetRandom(kQualMean,kQualSigma,positive)));
609 gain->SetValueAsInt(manuChannel,4,kSaturation);
613 Bool_t ok = gainStore.Add(gain);
616 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
618 if ( fMaxNofChannelsToGenerate > 0 && nchannels >= fMaxNofChannelsToGenerate ) break;
621 AliInfo(Form("%d Manus and %d channels.",nmanus,nchannels));
625 //_____________________________________________________________________________
627 AliMUONCDB::MakeLocalTriggerMaskStore(AliMUONVStore& localBoardMasks) const
629 /// Generate local trigger masks store. All masks are set to FFFF
631 AliCodeTimerAuto("");
634 // Generate fake mask values for all localboards and put that into
635 // one single container (localBoardMasks)
636 for ( Int_t i = 1; i <= AliMpConstants::TotalNofLocalBoards(); ++i )
638 AliMUONVCalibParam* localBoard = new AliMUONCalibParamNI(1,8,i,0,0);
639 for ( Int_t x = 0; x < 2; ++x )
641 for ( Int_t y = 0; y < 4; ++y )
644 localBoard->SetValueAsInt(index,0,0xFFFF);
648 localBoardMasks.Add(localBoard);
653 //_____________________________________________________________________________
655 AliMUONCDB::MakeRegionalTriggerConfigStore(AliMUONRegionalTriggerConfig& rtm) const
657 /// Make a regional trigger config store. Mask is set to FFFF for each local board (Ch.F.)
659 AliCodeTimerAuto("");
661 return rtm.ReadData(AliMpFiles::LocalTriggerBoardMapping());
666 //_____________________________________________________________________________
668 AliMUONCDB::MakeGlobalTriggerConfigStore(AliMUONGlobalCrateConfig& gtm) const
670 /// Make a global trigger config store. All masks (disable) set to 0x00 for each Darc board (Ch.F.)
672 AliCodeTimerAuto("");
674 return gtm.ReadData(AliMpFiles::GlobalTriggerBoardMapping());
678 //_____________________________________________________________________________
680 AliMUONCDB::MakeTriggerLUT(const char* file) const
682 /// Make a triggerlut object, from a file.
684 AliCodeTimerAuto("");
686 AliMUONTriggerLut* lut = new AliMUONTriggerLut;
687 lut->ReadFromFile(file);
691 //_____________________________________________________________________________
692 AliMUONTriggerEfficiencyCells*
693 AliMUONCDB::MakeTriggerEfficiency(const char* file) const
695 /// Make a trigger efficiency object from a file.
697 AliCodeTimerAuto("");
699 return new AliMUONTriggerEfficiencyCells(file);
702 //_____________________________________________________________________________
704 AliMUONCDB::WriteToCDB(const char* calibpath, TObject* object,
705 Int_t startRun, Int_t endRun,
706 const char* filename)
708 /// Write a given object to OCDB
710 AliCDBId id(calibpath,startRun,endRun);
712 md.SetAliRootVersion(gROOT->GetVersion());
713 md.SetComment(gSystem->ExpandPathName(filename));
714 md.SetResponsible("Uploaded using AliMUONCDB class");
715 AliCDBManager* man = AliCDBManager::Instance();
716 man->SetDefaultStorage(fCDBPath);
717 man->Put(object,id,&md);
720 //_____________________________________________________________________________
722 AliMUONCDB::WriteToCDB(const char* calibpath, TObject* object,
723 Int_t startRun, Int_t endRun, Bool_t defaultValues)
725 /// Write a given object to OCDB
727 AliCDBId id(calibpath,startRun,endRun);
729 md.SetAliRootVersion(gROOT->GetVersion());
732 md.SetComment("Test with default values");
736 md.SetComment("Test with random values");
738 md.SetResponsible("AliMUONCDB tester class");
740 AliCDBManager* man = AliCDBManager::Instance();
741 man->SetDefaultStorage(fCDBPath);
742 man->Put(object,id,&md);
745 //_____________________________________________________________________________
747 AliMUONCDB::MakeNeighbourStore(AliMUONVStore& neighbourStore)
749 /// Fill the neighbours store with, for each channel, a TObjArray of its
750 /// neighbouring pads (including itself)
752 AliCodeTimerAuto("");
754 AliInfo("Generating NeighbourStore. This will take a while. Please be patient.");
763 AliMpManuIterator it;
765 while ( it.Next(detElemId,manuId) )
767 const AliMpVSegmentation* seg =
768 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
770 AliMUONVCalibParam* calibParam = static_cast<AliMUONVCalibParam*>(neighbourStore.FindObject(detElemId,manuId));
774 Int_t size(AliMpConstants::ManuNofChannels());
775 Int_t defaultValue(-1);
776 Int_t packingFactor(size);
778 calibParam = new AliMUONCalibParamNI(dimension,size,detElemId,manuId,defaultValue,packingFactor);
779 Bool_t ok = neighbourStore.Add(calibParam);
782 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
787 for ( Int_t manuChannel = 0; manuChannel < AliMpConstants::ManuNofChannels(); ++manuChannel )
789 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
795 seg->GetNeighbours(pad,tmp,true,true);
796 Int_t nofPadNeighbours = tmp.GetEntriesFast();
798 for ( Int_t i = 0; i < nofPadNeighbours; ++i )
800 AliMpPad* p = static_cast<AliMpPad*>(tmp.UncheckedAt(i));
803 calibParam->PackValues(p->GetLocation().GetFirst(),p->GetLocation().GetSecond(),x);
806 // AliError("Could not pack value. Something is seriously wrong. Please check");
807 // StdoutToAliError(pad->Print(););
810 calibParam->SetValueAsInt(manuChannel,i,x);
819 //_____________________________________________________________________________
821 AliMUONCDB::SetMaxNofChannelsToGenerate(Int_t n)
823 /// Set the maximum number of channels to generate (used for testing only)
824 /// n < 0 means no limit
825 fMaxNofChannelsToGenerate = n;
828 //_____________________________________________________________________________
830 AliMUONCDB::WriteLocalTriggerMasks(Int_t startRun, Int_t endRun)
832 /// Write local trigger masks to OCDB
834 AliMUONVStore* ltm = new AliMUON1DArray(AliMpConstants::TotalNofLocalBoards()+1);
835 Int_t ngenerated = MakeLocalTriggerMaskStore(*ltm);
836 AliInfo(Form("Ngenerated = %d",ngenerated));
839 WriteToCDB("MUON/Calib/LocalTriggerBoardMasks",ltm,startRun,endRun,true);
844 //_____________________________________________________________________________
846 AliMUONCDB::WriteRegionalTriggerConfig(Int_t startRun, Int_t endRun)
848 /// Write regional trigger masks to OCDB
850 AliMUONRegionalTriggerConfig* rtm = new AliMUONRegionalTriggerConfig();
851 Int_t ngenerated = MakeRegionalTriggerConfigStore(*rtm);
852 AliInfo(Form("Ngenerated = %d",ngenerated));
855 WriteToCDB("MUON/Calib/RegionalTriggerConfig",rtm,startRun,endRun,true);
861 //_____________________________________________________________________________
863 AliMUONCDB::WriteGlobalTriggerConfig(Int_t startRun, Int_t endRun)
865 /// Write global trigger masks to OCDB
867 AliMUONGlobalCrateConfig* gtm = new AliMUONGlobalCrateConfig();
869 Int_t ngenerated = MakeGlobalTriggerConfigStore(*gtm);
870 AliInfo(Form("Ngenerated = %d",ngenerated));
873 WriteToCDB("MUON/Calib/GlobalTriggerCrateConfig",gtm,startRun,endRun,true);
879 //_____________________________________________________________________________
881 AliMUONCDB::WriteTriggerLut(Int_t startRun, Int_t endRun)
883 /// Write trigger LUT to OCDB
885 AliMUONTriggerLut* lut = MakeTriggerLUT();
888 WriteToCDB("MUON/Calib/TriggerLut",lut,startRun,endRun,true);
893 //_____________________________________________________________________________
895 AliMUONCDB::WriteTriggerEfficiency(Int_t startRun, Int_t endRun)
897 /// Write trigger efficiency to OCDB
899 AliMUONTriggerEfficiencyCells* eff = MakeTriggerEfficiency();
902 WriteToCDB("MUON/Calib/TriggerEfficiency",eff,startRun,endRun,true);
907 //_____________________________________________________________________________
909 AliMUONCDB::WriteNeighbours(Int_t startRun, Int_t endRun)
911 /// Write neighbours to OCDB
913 AliMUONVStore* neighbours = Create2DMap();
914 Int_t ngenerated = MakeNeighbourStore(*neighbours);
915 AliInfo(Form("Ngenerated = %d",ngenerated));
918 WriteToCDB("MUON/Calib/Neighbours",neighbours,startRun,endRun,true);
923 //_____________________________________________________________________________
925 AliMUONCDB::WriteHV(Bool_t defaultValues,
926 Int_t startRun, Int_t endRun)
928 /// generate HV values (either cste = 1500 V) if defaultValues=true or random
929 /// if defaultValues=false, see makeHVStore) and
930 /// store them into CDB located at cdbpath, with a validity period
931 /// ranging from startRun to endRun
933 TMap* hvStore = new TMap;
934 Int_t ngenerated = MakeHVStore(*hvStore,defaultValues);
935 AliInfo(Form("Ngenerated = %d",ngenerated));
938 WriteToCDB("MUON/Calib/HV",hvStore,startRun,endRun,defaultValues);
943 //_____________________________________________________________________________
945 AliMUONCDB::WritePedestals(Bool_t defaultValues,
946 Int_t startRun, Int_t endRun)
948 /// generate pedestal values (either 0 if defaultValues=true or random
949 /// if defaultValues=false, see makePedestalStore) and
950 /// store them into CDB located at cdbpath, with a validity period
951 /// ranging from startRun to endRun
953 AliMUONVStore* pedestalStore = Create2DMap();
954 Int_t ngenerated = MakePedestalStore(*pedestalStore,defaultValues);
955 AliInfo(Form("Ngenerated = %d",ngenerated));
956 WriteToCDB("MUON/Calib/Pedestals",pedestalStore,startRun,endRun,defaultValues);
957 delete pedestalStore;
961 //_____________________________________________________________________________
963 AliMUONCDB::WriteGains(Bool_t defaultValues,
964 Int_t startRun, Int_t endRun)
966 /// generate gain values (either 1 if defaultValues=true or random
967 /// if defaultValues=false, see makeGainStore) and
968 /// store them into CDB located at cdbpath, with a validity period
969 /// ranging from startRun to endRun
971 AliMUONVStore* gainStore = Create2DMap();
972 Int_t ngenerated = MakeGainStore(*gainStore,defaultValues);
973 AliInfo(Form("Ngenerated = %d",ngenerated));
974 WriteToCDB("MUON/Calib/Gains",gainStore,startRun,endRun,defaultValues);
978 //_____________________________________________________________________________
980 AliMUONCDB::WriteCapacitances(const char* filename,
981 Int_t startRun, Int_t endRun)
983 /// read manu capacitance and injection gain values from file
984 /// and store them into CDB located at cdbpath, with a validity period
985 /// ranging from startRun to endRun
987 AliMUONVStore* capaStore = new AliMUON1DMap(16828);
988 Int_t ngenerated = MakeCapacitanceStore(*capaStore,filename);
989 AliInfo(Form("Ngenerated = %d",ngenerated));
990 if ( ngenerated > 0 )
992 WriteToCDB("MUON/Calib/Capacitances",capaStore,startRun,endRun,filename);
997 //_____________________________________________________________________________
999 AliMUONCDB::WriteCapacitances(Bool_t defaultValues,
1000 Int_t startRun, Int_t endRun)
1002 /// generate manu capacitance values (either 1 if defaultValues=true or random
1003 /// if defaultValues=false, see makeCapacitanceStore) and
1004 /// store them into CDB located at cdbpath, with a validity period
1005 /// ranging from startRun to endRun
1007 AliMUONVStore* capaStore = new AliMUON1DMap(16828);
1008 Int_t ngenerated = MakeCapacitanceStore(*capaStore,defaultValues);
1009 AliInfo(Form("Ngenerated = %d",ngenerated));
1010 WriteToCDB("MUON/Calib/Capacitances",capaStore,startRun,endRun,defaultValues);
1014 //_____________________________________________________________________________
1016 AliMUONCDB::WriteTrigger(Int_t startRun, Int_t endRun)
1018 /// Writes all Trigger related calibration to CDB
1019 WriteLocalTriggerMasks(startRun,endRun);
1020 WriteRegionalTriggerConfig(startRun,endRun);
1021 WriteGlobalTriggerConfig(startRun,endRun);
1022 WriteTriggerLut(startRun,endRun);
1023 WriteTriggerEfficiency(startRun,endRun);
1026 //_____________________________________________________________________________
1028 AliMUONCDB::WriteTracker(Bool_t defaultValues, Int_t startRun, Int_t endRun)
1030 /// Writes all Tracker related calibration to CDB
1031 WriteHV(defaultValues,startRun,endRun);
1032 WritePedestals(defaultValues,startRun,endRun);
1033 WriteGains(defaultValues,startRun,endRun);
1034 WriteCapacitances(defaultValues,startRun,endRun);
1035 WriteNeighbours(startRun,endRun);