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"
48 #include "AliMpConstants.h"
49 #include "AliMpDDLStore.h"
50 #include "AliMpDEIterator.h"
51 #include "AliMpDEManager.h"
52 #include "AliMpDetElement.h"
53 #include "AliMpHVNamer.h"
54 #include "AliMpManuIterator.h"
55 #include "AliMpSegmentation.h"
56 #include "AliMpStationType.h"
57 #include "AliMpVSegmentation.h"
59 #include "AliCodeTimer.h"
60 #include "AliCDBEntry.h"
61 #include "AliCDBManager.h"
62 #include "AliDCSValue.h"
65 #include <Riostream.h>
71 #include <TObjString.h>
74 #include <TStopwatch.h>
85 //_____________________________________________________________________________
86 AliMUONVStore* Create2DMap()
88 return new AliMUON2DMap(true);
91 //_____________________________________________________________________________
92 void getBoundaries(const AliMUONVStore& store, Int_t dim,
93 Float_t* xmin, Float_t* xmax)
95 /// Assuming the store contains AliMUONVCalibParam objects, compute the
96 /// limits of the value contained in the VCalibParam, for each of its dimensions
97 /// xmin and xmax must be of dimension dim
99 for ( Int_t i = 0; i < dim; ++i )
105 TIter next(store.CreateIterator());
106 AliMUONVCalibParam* value;
108 while ( ( value = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
110 Int_t detElemId = value->ID0();
111 Int_t manuId = value->ID1();
113 const AliMpVSegmentation* seg =
114 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
116 for ( Int_t manuChannel = 0; manuChannel < value->Size(); ++manuChannel )
118 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
119 if (!pad.IsValid()) continue;
121 for ( Int_t i = 0; i < dim; ++i )
123 Float_t x0 = value->ValueAsFloat(manuChannel,i);
125 xmin[i] = TMath::Min(xmin[i],x0);
126 xmax[i] = TMath::Max(xmax[i],x0);
131 for ( Int_t i = 0; i < dim; ++i )
133 if ( TMath::Abs(xmin[i]-xmax[i]) < 1E-3 )
141 //_____________________________________________________________________________
142 Double_t GetRandom(Double_t mean, Double_t sigma, Bool_t mustBePositive)
145 if ( mustBePositive )
149 x = gRandom->Gaus(mean,sigma);
154 x = gRandom->Gaus(mean,sigma);
161 //_____________________________________________________________________________
162 AliMUONCDB::AliMUONCDB(const char* cdbpath)
165 fMaxNofChannelsToGenerate(-1)
169 if ( ! AliMpCDB::LoadDDLStore() ) {
170 AliFatal("Could not access mapping from OCDB !");
174 //_____________________________________________________________________________
175 AliMUONCDB::~AliMUONCDB()
180 //_____________________________________________________________________________
182 AliMUONCDB::Diff(AliMUONVStore& store1, AliMUONVStore& store2,
185 /// creates a store which contains store1-store2
186 /// if opt="abs" the difference is absolute one,
187 /// if opt="rel" then what is stored is (store1-store2)/store1
188 /// if opt="percent" then what is stored is rel*100
190 /// WARNING Works only for stores which holds AliMUONVCalibParam objects
195 if ( !sopt.Contains("ABS") && !sopt.Contains("REL") && !sopt.Contains("PERCENT") )
197 AliErrorClass(Form("opt %s not supported. Only ABS, REL, PERCENT are",opt));
201 AliMUONVStore* d = static_cast<AliMUONVStore*>(store1.Clone());
203 TIter next(d->CreateIterator());
205 AliMUONVCalibParam* param;
207 while ( ( param = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
209 Int_t detElemId = param->ID0();
210 Int_t manuId = param->ID1();
212 AliMUONVCalibParam* param2 = dynamic_cast<AliMUONVCalibParam*>(store2.FindObject(detElemId,manuId));
213 //FIXME: this might happen. Handle it.
216 cerr << "param2 is null : FIXME : this might happen !" << endl;
221 for ( Int_t i = 0; i < param->Size(); ++i )
223 for ( Int_t j = 0; j < param->Dimension(); ++j )
226 if ( sopt.Contains("ABS") )
228 value = param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j);
230 else if ( sopt.Contains("REL") || sopt.Contains("PERCENT") )
232 if ( param->ValueAsFloat(i,j) )
234 value = (param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j))/param->ValueAsFloat(i,j);
240 if ( sopt.Contains("PERCENT") ) value *= 100.0;
242 param->SetValueAsFloat(i,j,value);
249 //_____________________________________________________________________________
251 AliMUONCDB::Plot(const AliMUONVStore& store, const char* name, Int_t nbins)
253 /// Make histograms of each dimension of the AliMUONVCalibParam
254 /// contained inside store.
255 /// It produces histograms named name_0, name_1, etc...
257 TIter next(store.CreateIterator());
258 AliMUONVCalibParam* param;
260 const Int_t kNStations = AliMpConstants::NofTrackingChambers()/2;
261 Int_t* nPerStation = new Int_t[kNStations];
264 for ( Int_t i = 0; i < kNStations; ++i ) nPerStation[i]=0;
266 while ( ( param = static_cast<AliMUONVCalibParam*>(next()) ) )
270 Int_t dim = param->Dimension();
272 Float_t* xmin = new Float_t[dim];
273 Float_t* xmax = new Float_t[dim];
274 getBoundaries(store,dim,xmin,xmax);
276 for ( Int_t i = 0; i < dim; ++i )
278 h[i] = new TH1F(Form("%s_%d",name,i),Form("%s_%d",name,i),
279 nbins,xmin[i],xmax[i]);
280 AliInfo(Form("Created histogram %s",h[i]->GetName()));
284 Int_t detElemId = param->ID0();
285 Int_t manuId = param->ID1();
286 Int_t station = AliMpDEManager::GetChamberId(detElemId)/2;
288 const AliMpVSegmentation* seg =
289 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
291 for ( Int_t manuChannel = 0; manuChannel < param->Size(); ++manuChannel )
293 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
294 if (!pad.IsValid()) continue;
297 ++nPerStation[station];
299 for ( Int_t dim = 0; dim < param->Dimension(); ++dim )
301 h[dim]->Fill(param->ValueAsFloat(manuChannel,dim));
306 for ( Int_t i = 0; i < kNStations; ++i )
308 AliInfo(Form("Station %d %d ",(i+1),nPerStation[i]));
311 AliInfo(Form("Number of channels = %d",n));
313 delete[] nPerStation;
316 //_____________________________________________________________________________
318 AliMUONCDB::MakeHVStore(TMap& aliasMap, Bool_t defaultValues)
320 /// Create a HV store
322 AliMpHVNamer hvNamer;
324 TObjArray* aliases = hvNamer.GenerateAliases();
329 for ( Int_t i = 0; i < aliases->GetEntries(); ++i )
331 TObjString* alias = static_cast<TObjString*>(aliases->At(i));
332 TString& aliasName = alias->String();
333 if ( aliasName.Contains("sw") )
335 // HV Switch (St345 only)
336 TObjArray* valueSet = new TObjArray;
337 valueSet->SetOwner(kTRUE);
339 Bool_t value = kTRUE;
343 Float_t r = gRandom->Uniform();
344 if ( r < 0.007 ) value = kFALSE;
347 for ( UInt_t timeStamp = 0; timeStamp < 60*3; timeStamp += 60 )
349 AliDCSValue* dcsValue = new AliDCSValue(value,timeStamp);
350 valueSet->Add(dcsValue);
352 aliasMap.Add(new TObjString(*alias),valueSet);
357 TObjArray* valueSet = new TObjArray;
358 valueSet->SetOwner(kTRUE);
359 for ( UInt_t timeStamp = 0; timeStamp < 60*15; timeStamp += 120 )
361 Float_t value = 1500;
362 if (!defaultValues) value = GetRandom(1750,62.5,true);
363 AliDCSValue* dcsValue = new AliDCSValue(value,timeStamp);
364 valueSet->Add(dcsValue);
366 aliasMap.Add(new TObjString(*alias),valueSet);
373 AliInfo(Form("%d HV channels and %d switches",nChannels,nSwitch));
375 return nChannels+nSwitch;
378 //_____________________________________________________________________________
380 AliMUONCDB::MakePedestalStore(AliMUONVStore& pedestalStore, Bool_t defaultValues)
382 /// Create a pedestal store. if defaultValues=true, ped.mean=ped.sigma=1,
383 /// otherwise mean and sigma are from a gaussian (with parameters
384 /// defined below by the kPedestal* constants)
386 AliCodeTimerAuto("");
391 const Int_t kChannels(AliMpConstants::ManuNofChannels());
392 const Float_t kPedestalMeanMean(200);
393 const Float_t kPedestalMeanSigma(10);
394 const Float_t kPedestalSigmaMean(1.0);
395 const Float_t kPedestalSigmaSigma(0.2);
400 AliMpManuIterator it;
402 while ( it.Next(detElemId,manuId) )
406 AliMUONVCalibParam* ped =
407 new AliMUONCalibParamNF(2,kChannels,detElemId,manuId,AliMUONVCalibParam::InvalidFloatValue());
409 AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(detElemId);
411 for ( Int_t manuChannel = 0; manuChannel < kChannels; ++manuChannel )
413 if ( ! de->IsConnectedChannel(manuId,manuChannel) ) continue;
417 Float_t meanPedestal;
418 Float_t sigmaPedestal;
427 Bool_t positive(kTRUE);
429 while ( meanPedestal == 0.0 ) // avoid strict zero
431 meanPedestal = GetRandom(kPedestalMeanMean,kPedestalMeanSigma,positive);
433 sigmaPedestal = GetRandom(kPedestalSigmaMean,kPedestalSigmaSigma,positive);
435 ped->SetValueAsFloat(manuChannel,0,meanPedestal);
436 ped->SetValueAsFloat(manuChannel,1,sigmaPedestal);
439 Bool_t ok = pedestalStore.Add(ped);
442 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
444 if ( fMaxNofChannelsToGenerate > 0 && nchannels >= fMaxNofChannelsToGenerate ) break;
447 AliInfo(Form("%d Manus and %d channels.",nmanus,nchannels));
451 //_____________________________________________________________________________
453 AliMUONCDB::MakeCapacitanceStore(AliMUONVStore& capaStore, const char* file)
455 /// Read the capacitance values from file and append them to the capaStore
457 return AliMUONTrackerIO::ReadCapacitances(file,capaStore);
460 //_____________________________________________________________________________
462 AliMUONCDB::MakeCapacitanceStore(AliMUONVStore& capaStore, Bool_t defaultValues)
464 /// Create a capacitance store. if defaultValues=true, all capa are 1.0,
465 /// otherwise they are from a gaussian with parameters defined in the
466 /// kCapa* constants below.
468 AliCodeTimerAuto("");
472 Int_t nmanusOK(0); // manus for which we got the serial number
474 const Float_t kCapaMean(0.3);
475 const Float_t kCapaSigma(0.1);
476 const Float_t kInjectionGainMean(3);
477 const Float_t kInjectionGainSigma(1);
482 AliMpManuIterator it;
484 while ( it.Next(detElemId,manuId) )
488 AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(detElemId);
489 Int_t serialNumber = de->GetManuSerialFromId(manuId);
491 if ( serialNumber <= 0 ) continue;
495 AliMUONVCalibParam* capa = static_cast<AliMUONVCalibParam*>(capaStore.FindObject(serialNumber));
499 capa = new AliMUONCalibParamNF(2,AliMpConstants::ManuNofChannels(),serialNumber,0,1.0);
500 Bool_t ok = capaStore.Add(capa);
503 AliError(Form("Could not set serialNumber=%d manuId=%d",serialNumber,manuId));
507 for ( Int_t manuChannel = 0; manuChannel < capa->Size(); ++manuChannel )
509 if ( ! de->IsConnectedChannel(manuId,manuChannel) ) continue;
514 Float_t injectionGain;
523 capaValue = GetRandom(kCapaMean,kCapaSigma,kTRUE);
524 injectionGain = GetRandom(kInjectionGainMean,kInjectionGainSigma,kTRUE);
526 capa->SetValueAsFloat(manuChannel,0,capaValue);
527 capa->SetValueAsFloat(manuChannel,1,injectionGain);
532 if ( nmanus ) percent = 100*nmanusOK/nmanus;
533 AliInfo(Form("%5d manus with serial number (out of %5d manus = %3.0f%%)",
534 nmanusOK,nmanus,percent));
535 AliInfo(Form("%5d channels",nchannels));
538 AliWarning("Did not get all serial numbers. capaStore is incomplete !!!!");
544 //_____________________________________________________________________________
546 AliMUONCDB::MakeGainStore(AliMUONVStore& gainStore, Bool_t defaultValues)
548 /// Create a gain store. if defaultValues=true, all gains set so that
549 /// charge = (adc-ped)
551 /// otherwise parameters are taken from gaussians with parameters
552 /// defined in the k* constants below.
554 AliCodeTimerAuto("");
559 const Int_t kSaturation(3000);
560 const Double_t kA0Mean(1.2);
561 const Double_t kA0Sigma(0.1);
562 const Double_t kA1Mean(1E-5);
563 const Double_t kA1Sigma(1E-6);
564 const Double_t kQualMean(0xFF);
565 const Double_t kQualSigma(0x10);
566 const Int_t kThresMean(1600);
567 const Int_t kThresSigma(100);
572 AliMpManuIterator it;
574 while ( it.Next(detElemId,manuId) )
578 AliMUONVCalibParam* gain =
579 new AliMUONCalibParamNF(5,AliMpConstants::ManuNofChannels(),
582 AliMUONVCalibParam::InvalidFloatValue());
584 AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(detElemId);
586 for ( Int_t manuChannel = 0; manuChannel < gain->Size(); ++manuChannel )
588 if ( ! de->IsConnectedChannel(manuId,manuChannel) ) continue;
594 gain->SetValueAsFloat(manuChannel,0,1.0);
595 gain->SetValueAsFloat(manuChannel,1,0.0);
596 gain->SetValueAsInt(manuChannel,2,4095);
597 gain->SetValueAsInt(manuChannel,3,1);
598 gain->SetValueAsInt(manuChannel,4,kSaturation);
602 Bool_t positive(kTRUE);
603 gain->SetValueAsFloat(manuChannel,0,GetRandom(kA0Mean,kA0Sigma,positive));
604 gain->SetValueAsFloat(manuChannel,1,GetRandom(kA1Mean,kA1Sigma,!positive));
605 gain->SetValueAsInt(manuChannel,2,(Int_t)TMath::Nint(GetRandom(kThresMean,kThresSigma,positive)));
606 gain->SetValueAsInt(manuChannel,3,(Int_t)TMath::Nint(GetRandom(kQualMean,kQualSigma,positive)));
607 gain->SetValueAsInt(manuChannel,4,kSaturation);
611 Bool_t ok = gainStore.Add(gain);
614 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
616 if ( fMaxNofChannelsToGenerate > 0 && nchannels >= fMaxNofChannelsToGenerate ) break;
619 AliInfo(Form("%d Manus and %d channels.",nmanus,nchannels));
623 //_____________________________________________________________________________
625 AliMUONCDB::MakeLocalTriggerMaskStore(AliMUONVStore& localBoardMasks) const
627 /// Generate local trigger masks store. All masks are set to FFFF
629 AliCodeTimerAuto("");
632 // Generate fake mask values for 234 localboards and put that into
633 // one single container (localBoardMasks)
634 for ( Int_t i = 1; i <= AliMpConstants::TotalNofLocalBoards(); ++i )
636 AliMUONVCalibParam* localBoard = new AliMUONCalibParamNI(1,8,i,0,0);
637 for ( Int_t x = 0; x < 2; ++x )
639 for ( Int_t y = 0; y < 4; ++y )
642 localBoard->SetValueAsInt(index,0,0xFFFF);
646 localBoardMasks.Add(localBoard);
651 //_____________________________________________________________________________
653 AliMUONCDB::MakeRegionalTriggerMaskStore(AliMUONVStore& rtm) const
655 /// Make a regional trigger masks store. Mask is set to FFFF for each local board (Ch.F.)
657 AliCodeTimerAuto("");
660 for ( Int_t i = 0; i < 16; ++i )
662 AliMUONVCalibParam* regionalBoard = new AliMUONCalibParamNI(1,1,i,0,0);
664 regionalBoard->SetValueAsInt(0,0,0xFFFF);
667 rtm.Add(regionalBoard);
673 //_____________________________________________________________________________
675 AliMUONCDB::MakeGlobalTriggerMaskStore(AliMUONVCalibParam& gtm) const
677 /// Make a global trigger masks store. All masks (disable) set to 0x00 for each Darc board (Ch.F.)
679 AliCodeTimerAuto("");
683 for ( Int_t j = 0; j < 2; ++j )
685 gtm.SetValueAsInt(j,0,0x00);
691 //_____________________________________________________________________________
693 AliMUONCDB::MakeTriggerLUT(const char* file) const
695 /// Make a triggerlut object, from a file.
697 AliCodeTimerAuto("");
699 AliMUONTriggerLut* lut = new AliMUONTriggerLut;
700 lut->ReadFromFile(file);
704 //_____________________________________________________________________________
705 AliMUONTriggerEfficiencyCells*
706 AliMUONCDB::MakeTriggerEfficiency(const char* file) const
708 /// Make a trigger efficiency object from a file.
710 AliCodeTimerAuto("");
712 return new AliMUONTriggerEfficiencyCells(file);
715 //_____________________________________________________________________________
717 AliMUONCDB::WriteToCDB(const char* calibpath, TObject* object,
718 Int_t startRun, Int_t endRun,
719 const char* filename)
721 /// Write a given object to OCDB
723 AliCDBId id(calibpath,startRun,endRun);
725 md.SetAliRootVersion(gROOT->GetVersion());
726 md.SetComment(gSystem->ExpandPathName(filename));
727 md.SetResponsible("Uploaded using AliMUONCDB class");
728 AliCDBManager* man = AliCDBManager::Instance();
729 man->SetDefaultStorage(fCDBPath);
730 man->Put(object,id,&md);
733 //_____________________________________________________________________________
735 AliMUONCDB::WriteToCDB(const char* calibpath, TObject* object,
736 Int_t startRun, Int_t endRun, Bool_t defaultValues)
738 /// Write a given object to OCDB
740 AliCDBId id(calibpath,startRun,endRun);
742 md.SetAliRootVersion(gROOT->GetVersion());
745 md.SetComment("Test with default values");
749 md.SetComment("Test with random values");
751 md.SetResponsible("AliMUONCDB tester class");
753 AliCDBManager* man = AliCDBManager::Instance();
754 man->SetDefaultStorage(fCDBPath);
755 man->Put(object,id,&md);
758 //_____________________________________________________________________________
760 AliMUONCDB::MakeNeighbourStore(AliMUONVStore& neighbourStore)
762 /// Fill the neighbours store with, for each channel, a TObjArray of its
763 /// neighbouring pads (including itself)
765 AliCodeTimerAuto("");
767 AliInfo("Generating NeighbourStore. This will take a while. Please be patient.");
776 AliMpManuIterator it;
778 while ( it.Next(detElemId,manuId) )
780 const AliMpVSegmentation* seg =
781 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
783 AliMUONVCalibParam* calibParam = static_cast<AliMUONVCalibParam*>(neighbourStore.FindObject(detElemId,manuId));
787 Int_t size(AliMpConstants::ManuNofChannels());
788 Int_t defaultValue(-1);
789 Int_t packingFactor(size);
791 calibParam = new AliMUONCalibParamNI(dimension,size,detElemId,manuId,defaultValue,packingFactor);
792 Bool_t ok = neighbourStore.Add(calibParam);
795 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
800 for ( Int_t manuChannel = 0; manuChannel < AliMpConstants::ManuNofChannels(); ++manuChannel )
802 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
808 seg->GetNeighbours(pad,tmp,true,true);
809 Int_t nofPadNeighbours = tmp.GetEntriesFast();
811 for ( Int_t i = 0; i < nofPadNeighbours; ++i )
813 AliMpPad* pad = static_cast<AliMpPad*>(tmp.UncheckedAt(i));
816 calibParam->PackValues(pad->GetLocation().GetFirst(),pad->GetLocation().GetSecond(),x);
819 // AliError("Could not pack value. Something is seriously wrong. Please check");
820 // StdoutToAliError(pad->Print(););
823 calibParam->SetValueAsInt(manuChannel,i,x);
832 //_____________________________________________________________________________
834 AliMUONCDB::SetMaxNofChannelsToGenerate(Int_t n)
836 /// Set the maximum number of channels to generate (used for testing only)
837 /// n < 0 means no limit
838 fMaxNofChannelsToGenerate = n;
841 //_____________________________________________________________________________
843 AliMUONCDB::WriteLocalTriggerMasks(Int_t startRun, Int_t endRun)
845 /// Write local trigger masks to OCDB
847 AliMUONVStore* ltm = new AliMUON1DArray(235);
848 Int_t ngenerated = MakeLocalTriggerMaskStore(*ltm);
849 AliInfo(Form("Ngenerated = %d",ngenerated));
852 WriteToCDB("MUON/Calib/LocalTriggerBoardMasks",ltm,startRun,endRun,true);
857 //_____________________________________________________________________________
859 AliMUONCDB::WriteRegionalTriggerMasks(Int_t startRun, Int_t endRun)
861 /// Write regional trigger masks to OCDB
863 AliMUONVStore* rtm = new AliMUON1DArray(16);
864 Int_t ngenerated = MakeRegionalTriggerMaskStore(*rtm);
865 AliInfo(Form("Ngenerated = %d",ngenerated));
868 WriteToCDB("MUON/Calib/RegionalTriggerBoardMasks",rtm,startRun,endRun,true);
873 //_____________________________________________________________________________
875 AliMUONCDB::WriteGlobalTriggerMasks(Int_t startRun, Int_t endRun)
877 /// Write global trigger masks to OCDB
879 AliMUONVCalibParam* gtm = new AliMUONCalibParamNI(1,2,1,0,0);
881 Int_t ngenerated = MakeGlobalTriggerMaskStore(*gtm);
882 AliInfo(Form("Ngenerated = %d",ngenerated));
885 WriteToCDB("MUON/Calib/GlobalTriggerBoardMasks",gtm,startRun,endRun,true);
890 //_____________________________________________________________________________
892 AliMUONCDB::WriteTriggerLut(Int_t startRun, Int_t endRun)
894 /// Write trigger LUT to OCDB
896 AliMUONTriggerLut* lut = MakeTriggerLUT();
899 WriteToCDB("MUON/Calib/TriggerLut",lut,startRun,endRun,true);
904 //_____________________________________________________________________________
906 AliMUONCDB::WriteTriggerEfficiency(Int_t startRun, Int_t endRun)
908 /// Write trigger efficiency to OCDB
910 AliMUONTriggerEfficiencyCells* eff = MakeTriggerEfficiency();
913 WriteToCDB("MUON/Calib/TriggerEfficiency",eff,startRun,endRun,true);
918 //_____________________________________________________________________________
920 AliMUONCDB::WriteNeighbours(Int_t startRun, Int_t endRun)
922 /// Write neighbours to OCDB
924 AliMUONVStore* neighbours = Create2DMap();
925 Int_t ngenerated = MakeNeighbourStore(*neighbours);
926 AliInfo(Form("Ngenerated = %d",ngenerated));
929 WriteToCDB("MUON/Calib/Neighbours",neighbours,startRun,endRun,true);
934 //_____________________________________________________________________________
936 AliMUONCDB::WriteHV(Bool_t defaultValues,
937 Int_t startRun, Int_t endRun)
939 /// generate HV values (either cste = 1500 V) if defaultValues=true or random
940 /// if defaultValues=false, see makeHVStore) and
941 /// store them into CDB located at cdbpath, with a validity period
942 /// ranging from startRun to endRun
944 TMap* hvStore = new TMap;
945 Int_t ngenerated = MakeHVStore(*hvStore,defaultValues);
946 AliInfo(Form("Ngenerated = %d",ngenerated));
949 WriteToCDB("MUON/Calib/HV",hvStore,startRun,endRun,defaultValues);
954 //_____________________________________________________________________________
956 AliMUONCDB::WritePedestals(Bool_t defaultValues,
957 Int_t startRun, Int_t endRun)
959 /// generate pedestal values (either 0 if defaultValues=true or random
960 /// if defaultValues=false, see makePedestalStore) and
961 /// store them into CDB located at cdbpath, with a validity period
962 /// ranging from startRun to endRun
964 AliMUONVStore* pedestalStore = Create2DMap();
965 Int_t ngenerated = MakePedestalStore(*pedestalStore,defaultValues);
966 AliInfo(Form("Ngenerated = %d",ngenerated));
967 WriteToCDB("MUON/Calib/Pedestals",pedestalStore,startRun,endRun,defaultValues);
968 delete pedestalStore;
972 //_____________________________________________________________________________
974 AliMUONCDB::WriteGains(Bool_t defaultValues,
975 Int_t startRun, Int_t endRun)
977 /// generate gain values (either 1 if defaultValues=true or random
978 /// if defaultValues=false, see makeGainStore) and
979 /// store them into CDB located at cdbpath, with a validity period
980 /// ranging from startRun to endRun
982 AliMUONVStore* gainStore = Create2DMap();
983 Int_t ngenerated = MakeGainStore(*gainStore,defaultValues);
984 AliInfo(Form("Ngenerated = %d",ngenerated));
985 WriteToCDB("MUON/Calib/Gains",gainStore,startRun,endRun,defaultValues);
989 //_____________________________________________________________________________
991 AliMUONCDB::WriteCapacitances(const char* filename,
992 Int_t startRun, Int_t endRun)
994 /// read manu capacitance and injection gain values from file
995 /// and store them into CDB located at cdbpath, with a validity period
996 /// ranging from startRun to endRun
998 AliMUONVStore* capaStore = new AliMUON1DMap(16828);
999 Int_t ngenerated = MakeCapacitanceStore(*capaStore,filename);
1000 AliInfo(Form("Ngenerated = %d",ngenerated));
1001 if ( ngenerated > 0 )
1003 WriteToCDB("MUON/Calib/Capacitances",capaStore,startRun,endRun,filename);
1008 //_____________________________________________________________________________
1010 AliMUONCDB::WriteCapacitances(Bool_t defaultValues,
1011 Int_t startRun, Int_t endRun)
1013 /// generate manu capacitance values (either 1 if defaultValues=true or random
1014 /// if defaultValues=false, see makeCapacitanceStore) and
1015 /// store them into CDB located at cdbpath, with a validity period
1016 /// ranging from startRun to endRun
1018 AliMUONVStore* capaStore = new AliMUON1DMap(16828);
1019 Int_t ngenerated = MakeCapacitanceStore(*capaStore,defaultValues);
1020 AliInfo(Form("Ngenerated = %d",ngenerated));
1021 WriteToCDB("MUON/Calib/Capacitances",capaStore,startRun,endRun,defaultValues);
1025 //_____________________________________________________________________________
1027 AliMUONCDB::WriteTrigger(Int_t startRun, Int_t endRun)
1029 /// Writes all Trigger related calibration to CDB
1030 WriteLocalTriggerMasks(startRun,endRun);
1031 WriteRegionalTriggerMasks(startRun,endRun);
1032 WriteGlobalTriggerMasks(startRun,endRun);
1033 WriteTriggerLut(startRun,endRun);
1034 WriteTriggerEfficiency(startRun,endRun);
1037 //_____________________________________________________________________________
1039 AliMUONCDB::WriteTracker(Bool_t defaultValues, Int_t startRun, Int_t endRun)
1041 /// Writes all Tracker related calibration to CDB
1042 WriteHV(defaultValues,startRun,endRun);
1043 WritePedestals(defaultValues,startRun,endRun);
1044 WriteGains(defaultValues,startRun,endRun);
1045 WriteCapacitances(defaultValues,startRun,endRun);
1046 WriteNeighbours(startRun,endRun);