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 "AliMUONTriggerEfficiencyCells.h"
41 #include "AliMUONTriggerLut.h"
42 #include "AliMUONVStore.h"
43 #include "AliMUONVCalibParam.h"
44 #include "AliMUONVCalibParam.h"
47 #include "AliMpConstants.h"
48 #include "AliMpDDLStore.h"
49 #include "AliMpDEIterator.h"
50 #include "AliMpDEManager.h"
51 #include "AliMpDetElement.h"
52 #include "AliMpHVNamer.h"
53 #include "AliMpManuIterator.h"
54 #include "AliMpSegmentation.h"
55 #include "AliMpStationType.h"
56 #include "AliMpVSegmentation.h"
58 #include "AliCodeTimer.h"
59 #include "AliCDBEntry.h"
60 #include "AliCDBManager.h"
61 #include "AliDCSValue.h"
64 #include <Riostream.h>
70 #include <TObjString.h>
73 #include <TStopwatch.h>
84 //_____________________________________________________________________________
85 AliMUONVStore* Create2DMap()
87 return new AliMUON2DMap(true);
90 //_____________________________________________________________________________
91 void getBoundaries(const AliMUONVStore& store, Int_t dim,
92 Float_t* xmin, Float_t* xmax)
94 /// Assuming the store contains AliMUONVCalibParam objects, compute the
95 /// limits of the value contained in the VCalibParam, for each of its dimensions
96 /// xmin and xmax must be of dimension dim
98 for ( Int_t i = 0; i < dim; ++i )
104 TIter next(store.CreateIterator());
105 AliMUONVCalibParam* value;
107 while ( ( value = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
109 Int_t detElemId = value->ID0();
110 Int_t manuId = value->ID1();
112 const AliMpVSegmentation* seg =
113 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
115 for ( Int_t manuChannel = 0; manuChannel < value->Size(); ++manuChannel )
117 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
118 if (!pad.IsValid()) continue;
120 for ( Int_t i = 0; i < dim; ++i )
122 Float_t x0 = value->ValueAsFloat(manuChannel,i);
124 xmin[i] = TMath::Min(xmin[i],x0);
125 xmax[i] = TMath::Max(xmax[i],x0);
130 for ( Int_t i = 0; i < dim; ++i )
132 if ( TMath::Abs(xmin[i]-xmax[i]) < 1E-3 )
140 //_____________________________________________________________________________
141 Double_t GetRandom(Double_t mean, Double_t sigma, Bool_t mustBePositive)
144 if ( mustBePositive )
148 x = gRandom->Gaus(mean,sigma);
153 x = gRandom->Gaus(mean,sigma);
160 //_____________________________________________________________________________
161 AliMUONCDB::AliMUONCDB(const char* cdbpath)
164 fMaxNofChannelsToGenerate(-1)
168 if ( ! AliMpCDB::LoadDDLStore() ) {
169 AliFatal("Could not access mapping from OCDB !");
173 //_____________________________________________________________________________
174 AliMUONCDB::~AliMUONCDB()
179 //_____________________________________________________________________________
181 AliMUONCDB::Diff(AliMUONVStore& store1, AliMUONVStore& store2,
184 /// creates a store which contains store1-store2
185 /// if opt="abs" the difference is absolute one,
186 /// if opt="rel" then what is stored is (store1-store2)/store1
187 /// if opt="percent" then what is stored is rel*100
189 /// WARNING Works only for stores which holds AliMUONVCalibParam objects
194 if ( !sopt.Contains("ABS") && !sopt.Contains("REL") && !sopt.Contains("PERCENT") )
196 AliErrorClass(Form("opt %s not supported. Only ABS, REL, PERCENT are",opt));
200 AliMUONVStore* d = static_cast<AliMUONVStore*>(store1.Clone());
202 TIter next(d->CreateIterator());
204 AliMUONVCalibParam* param;
206 while ( ( param = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
208 Int_t detElemId = param->ID0();
209 Int_t manuId = param->ID1();
211 AliMUONVCalibParam* param2 = dynamic_cast<AliMUONVCalibParam*>(store2.FindObject(detElemId,manuId));
212 //FIXME: this might happen. Handle it.
215 cerr << "param2 is null : FIXME : this might happen !" << endl;
220 for ( Int_t i = 0; i < param->Size(); ++i )
222 for ( Int_t j = 0; j < param->Dimension(); ++j )
225 if ( sopt.Contains("ABS") )
227 value = param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j);
229 else if ( sopt.Contains("REL") || sopt.Contains("PERCENT") )
231 if ( param->ValueAsFloat(i,j) )
233 value = (param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j))/param->ValueAsFloat(i,j);
239 if ( sopt.Contains("PERCENT") ) value *= 100.0;
241 param->SetValueAsFloat(i,j,value);
248 //_____________________________________________________________________________
250 AliMUONCDB::Plot(const AliMUONVStore& store, const char* name, Int_t nbins)
252 /// Make histograms of each dimension of the AliMUONVCalibParam
253 /// contained inside store.
254 /// It produces histograms named name_0, name_1, etc...
256 TIter next(store.CreateIterator());
257 AliMUONVCalibParam* param;
259 const Int_t kNStations = AliMpConstants::NofTrackingChambers()/2;
260 Int_t* nPerStation = new Int_t[kNStations];
263 for ( Int_t i = 0; i < kNStations; ++i ) nPerStation[i]=0;
265 while ( ( param = static_cast<AliMUONVCalibParam*>(next()) ) )
269 Int_t dim = param->Dimension();
271 Float_t* xmin = new Float_t[dim];
272 Float_t* xmax = new Float_t[dim];
273 getBoundaries(store,dim,xmin,xmax);
275 for ( Int_t i = 0; i < dim; ++i )
277 h[i] = new TH1F(Form("%s_%d",name,i),Form("%s_%d",name,i),
278 nbins,xmin[i],xmax[i]);
279 AliInfo(Form("Created histogram %s",h[i]->GetName()));
283 Int_t detElemId = param->ID0();
284 Int_t manuId = param->ID1();
285 Int_t station = AliMpDEManager::GetChamberId(detElemId)/2;
287 const AliMpVSegmentation* seg =
288 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
290 for ( Int_t manuChannel = 0; manuChannel < param->Size(); ++manuChannel )
292 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
293 if (!pad.IsValid()) continue;
296 ++nPerStation[station];
298 for ( Int_t dim = 0; dim < param->Dimension(); ++dim )
300 h[dim]->Fill(param->ValueAsFloat(manuChannel,dim));
305 for ( Int_t i = 0; i < kNStations; ++i )
307 AliInfo(Form("Station %d %d ",(i+1),nPerStation[i]));
310 AliInfo(Form("Number of channels = %d",n));
312 delete[] nPerStation;
315 //_____________________________________________________________________________
317 AliMUONCDB::MakeHVStore(TMap& aliasMap, Bool_t defaultValues)
319 /// Create a HV store
321 AliMpHVNamer hvNamer;
323 TObjArray* aliases = hvNamer.GenerateAliases();
328 for ( Int_t i = 0; i < aliases->GetEntries(); ++i )
330 TObjString* alias = static_cast<TObjString*>(aliases->At(i));
331 TString& aliasName = alias->String();
332 if ( aliasName.Contains("sw") )
334 // HV Switch (St345 only)
335 TObjArray* valueSet = new TObjArray;
336 valueSet->SetOwner(kTRUE);
338 Bool_t value = kTRUE;
342 Float_t r = gRandom->Uniform();
343 if ( r < 0.007 ) value = kFALSE;
346 for ( UInt_t timeStamp = 0; timeStamp < 60*3; timeStamp += 60 )
348 AliDCSValue* dcsValue = new AliDCSValue(value,timeStamp);
349 valueSet->Add(dcsValue);
351 aliasMap.Add(new TObjString(*alias),valueSet);
356 TObjArray* valueSet = new TObjArray;
357 valueSet->SetOwner(kTRUE);
358 for ( UInt_t timeStamp = 0; timeStamp < 60*15; timeStamp += 120 )
360 Float_t value = 1500;
361 if (!defaultValues) value = GetRandom(1750,62.5,true);
362 AliDCSValue* dcsValue = new AliDCSValue(value,timeStamp);
363 valueSet->Add(dcsValue);
365 aliasMap.Add(new TObjString(*alias),valueSet);
372 AliInfo(Form("%d HV channels and %d switches",nChannels,nSwitch));
374 return nChannels+nSwitch;
377 //_____________________________________________________________________________
379 AliMUONCDB::MakePedestalStore(AliMUONVStore& pedestalStore, Bool_t defaultValues)
381 /// Create a pedestal store. if defaultValues=true, ped.mean=ped.sigma=1,
382 /// otherwise mean and sigma are from a gaussian (with parameters
383 /// defined below by the kPedestal* constants)
385 AliCodeTimerAuto("");
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);
399 AliMpManuIterator it;
401 while ( it.Next(detElemId,manuId) )
405 AliMUONVCalibParam* ped =
406 new AliMUONCalibParamNF(2,kChannels,detElemId,manuId,AliMUONVCalibParam::InvalidFloatValue());
408 AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(detElemId);
410 for ( Int_t manuChannel = 0; manuChannel < kChannels; ++manuChannel )
412 if ( ! de->IsConnectedChannel(manuId,manuChannel) ) 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));
450 //_____________________________________________________________________________
452 AliMUONCDB::MakeCapacitanceStore(AliMUONVStore& capaStore, const char* file)
454 /// Read the capacitance values from file and append them to the capaStore
456 AliCodeTimerAuto(Form("file=%s",file));
458 ifstream in(gSystem->ExpandPathName(file));
459 if (in.bad()) return 0;
464 Int_t serialNumber(-1);
465 AliMUONVCalibParam* param(0x0);
467 while ( in.getline(line,1024,'\n') )
469 if ( isdigit(line[0]) )
471 serialNumber = atoi(line);
472 param = static_cast<AliMUONVCalibParam*>(capaStore.FindObject(serialNumber));
475 AliError(Form("serialNumber %d appears several times !",serialNumber));
479 param = new AliMUONCalibParamNF(2,AliMpConstants::ManuNofChannels(),serialNumber,0,1.0);
480 Bool_t ok = capaStore.Add(param);
483 AliError(Form("Could not set serialNumber=%d",serialNumber));
490 Float_t injectionGain;
491 sscanf(line,"%d %f %f",&channel,&capaValue,&injectionGain);
492 AliDebug(1,Form("SerialNumber %10d Channel %3d Capa %f injectionGain %f",
493 serialNumber,channel,capaValue,injectionGain));
494 param->SetValueAsFloat(channel,0,capaValue);
495 param->SetValueAsFloat(channel,1,injectionGain);
504 //_____________________________________________________________________________
506 AliMUONCDB::MakeCapacitanceStore(AliMUONVStore& capaStore, Bool_t defaultValues)
508 /// Create a capacitance store. if defaultValues=true, all capa are 1.0,
509 /// otherwise they are from a gaussian with parameters defined in the
510 /// kCapa* constants below.
512 AliCodeTimerAuto("");
516 Int_t nmanusOK(0); // manus for which we got the serial number
518 const Float_t kCapaMean(0.3);
519 const Float_t kCapaSigma(0.1);
520 const Float_t kInjectionGainMean(3);
521 const Float_t kInjectionGainSigma(1);
526 AliMpManuIterator it;
528 while ( it.Next(detElemId,manuId) )
532 AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(detElemId);
533 Int_t serialNumber = de->GetManuSerialFromId(manuId);
535 if ( serialNumber <= 0 ) continue;
539 AliMUONVCalibParam* capa = static_cast<AliMUONVCalibParam*>(capaStore.FindObject(serialNumber));
543 capa = new AliMUONCalibParamNF(2,AliMpConstants::ManuNofChannels(),serialNumber,0,1.0);
544 Bool_t ok = capaStore.Add(capa);
547 AliError(Form("Could not set serialNumber=%d manuId=%d",serialNumber,manuId));
551 for ( Int_t manuChannel = 0; manuChannel < capa->Size(); ++manuChannel )
553 if ( ! de->IsConnectedChannel(manuId,manuChannel) ) continue;
558 Float_t injectionGain;
567 capaValue = GetRandom(kCapaMean,kCapaSigma,kTRUE);
568 injectionGain = GetRandom(kInjectionGainMean,kInjectionGainSigma,kTRUE);
570 capa->SetValueAsFloat(manuChannel,0,capaValue);
571 capa->SetValueAsFloat(manuChannel,1,injectionGain);
576 if ( nmanus ) percent = 100*nmanusOK/nmanus;
577 AliInfo(Form("%5d manus with serial number (out of %5d manus = %3.0f%%)",
578 nmanusOK,nmanus,percent));
579 AliInfo(Form("%5d channels",nchannels));
582 AliWarning("Did not get all serial numbers. capaStore is incomplete !!!!");
588 //_____________________________________________________________________________
590 AliMUONCDB::MakeGainStore(AliMUONVStore& gainStore, Bool_t defaultValues)
592 /// Create a gain store. if defaultValues=true, all gains set so that
593 /// charge = (adc-ped)
595 /// otherwise parameters are taken from gaussians with parameters
596 /// defined in the k* constants below.
598 AliCodeTimerAuto("");
603 const Int_t kSaturation(3000);
604 const Double_t kA0Mean(1.2);
605 const Double_t kA0Sigma(0.1);
606 const Double_t kA1Mean(1E-5);
607 const Double_t kA1Sigma(1E-6);
608 const Double_t kQualMean(0xFF);
609 const Double_t kQualSigma(0x10);
610 const Int_t kThresMean(1600);
611 const Int_t kThresSigma(100);
616 AliMpManuIterator it;
618 while ( it.Next(detElemId,manuId) )
622 AliMUONVCalibParam* gain =
623 new AliMUONCalibParamNF(5,AliMpConstants::ManuNofChannels(),
626 AliMUONVCalibParam::InvalidFloatValue());
628 AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(detElemId);
630 for ( Int_t manuChannel = 0; manuChannel < gain->Size(); ++manuChannel )
632 if ( ! de->IsConnectedChannel(manuId,manuChannel) ) continue;
638 gain->SetValueAsFloat(manuChannel,0,1.0);
639 gain->SetValueAsFloat(manuChannel,1,0.0);
640 gain->SetValueAsInt(manuChannel,2,4095);
641 gain->SetValueAsInt(manuChannel,3,1);
642 gain->SetValueAsInt(manuChannel,4,kSaturation);
646 Bool_t positive(kTRUE);
647 gain->SetValueAsFloat(manuChannel,0,GetRandom(kA0Mean,kA0Sigma,positive));
648 gain->SetValueAsFloat(manuChannel,1,GetRandom(kA1Mean,kA1Sigma,!positive));
649 gain->SetValueAsInt(manuChannel,2,(Int_t)TMath::Nint(GetRandom(kThresMean,kThresSigma,positive)));
650 gain->SetValueAsInt(manuChannel,3,(Int_t)TMath::Nint(GetRandom(kQualMean,kQualSigma,positive)));
651 gain->SetValueAsInt(manuChannel,4,kSaturation);
655 Bool_t ok = gainStore.Add(gain);
658 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
660 if ( fMaxNofChannelsToGenerate > 0 && nchannels >= fMaxNofChannelsToGenerate ) break;
663 AliInfo(Form("%d Manus and %d channels.",nmanus,nchannels));
667 //_____________________________________________________________________________
669 AliMUONCDB::MakeLocalTriggerMaskStore(AliMUONVStore& localBoardMasks) const
671 /// Generate local trigger masks store. All masks are set to FFFF
673 AliCodeTimerAuto("");
676 // Generate fake mask values for 234 localboards and put that into
677 // one single container (localBoardMasks)
678 for ( Int_t i = 1; i <= AliMpConstants::TotalNofLocalBoards(); ++i )
680 AliMUONVCalibParam* localBoard = new AliMUONCalibParamNI(1,8,i,0,0);
681 for ( Int_t x = 0; x < 2; ++x )
683 for ( Int_t y = 0; y < 4; ++y )
686 localBoard->SetValueAsInt(index,0,0xFFFF);
690 localBoardMasks.Add(localBoard);
695 //_____________________________________________________________________________
697 AliMUONCDB::MakeRegionalTriggerMaskStore(AliMUONVStore& rtm) const
699 /// Make a regional trigger masks store. Mask is set to FFFF for each local board (Ch.F.)
701 AliCodeTimerAuto("");
704 for ( Int_t i = 0; i < 16; ++i )
706 AliMUONVCalibParam* regionalBoard = new AliMUONCalibParamNI(1,1,i,0,0);
708 regionalBoard->SetValueAsInt(0,0,0xFFFF);
711 rtm.Add(regionalBoard);
717 //_____________________________________________________________________________
719 AliMUONCDB::MakeGlobalTriggerMaskStore(AliMUONVCalibParam& gtm) const
721 /// Make a global trigger masks store. All masks (disable) set to 0x00 for each Darc board (Ch.F.)
723 AliCodeTimerAuto("");
727 for ( Int_t j = 0; j < 2; ++j )
729 gtm.SetValueAsInt(j,0,0x00);
735 //_____________________________________________________________________________
737 AliMUONCDB::MakeTriggerLUT(const char* file) const
739 /// Make a triggerlut object, from a file.
741 AliCodeTimerAuto("");
743 AliMUONTriggerLut* lut = new AliMUONTriggerLut;
744 lut->ReadFromFile(file);
748 //_____________________________________________________________________________
749 AliMUONTriggerEfficiencyCells*
750 AliMUONCDB::MakeTriggerEfficiency(const char* file) const
752 /// Make a trigger efficiency object from a file.
754 AliCodeTimerAuto("");
756 return new AliMUONTriggerEfficiencyCells(file);
759 //_____________________________________________________________________________
761 AliMUONCDB::WriteToCDB(const char* calibpath, TObject* object,
762 Int_t startRun, Int_t endRun,
763 const char* filename)
765 /// Write a given object to OCDB
767 AliCDBId id(calibpath,startRun,endRun);
769 md.SetAliRootVersion(gROOT->GetVersion());
770 md.SetComment(gSystem->ExpandPathName(filename));
771 md.SetResponsible("Uploaded using AliMUONCDB class");
772 AliCDBManager* man = AliCDBManager::Instance();
773 man->SetDefaultStorage(fCDBPath);
774 man->Put(object,id,&md);
777 //_____________________________________________________________________________
779 AliMUONCDB::WriteToCDB(const char* calibpath, TObject* object,
780 Int_t startRun, Int_t endRun, Bool_t defaultValues)
782 /// Write a given object to OCDB
784 AliCDBId id(calibpath,startRun,endRun);
786 md.SetAliRootVersion(gROOT->GetVersion());
789 md.SetComment("Test with default values");
793 md.SetComment("Test with random values");
795 md.SetResponsible("AliMUONCDB tester class");
797 AliCDBManager* man = AliCDBManager::Instance();
798 man->SetDefaultStorage(fCDBPath);
799 man->Put(object,id,&md);
802 //_____________________________________________________________________________
804 AliMUONCDB::MakeNeighbourStore(AliMUONVStore& neighbourStore)
806 /// Fill the neighbours store with, for each channel, a TObjArray of its
807 /// neighbouring pads (including itself)
809 AliCodeTimerAuto("");
811 AliInfo("Generating NeighbourStore. This will take a while. Please be patient.");
820 AliMpManuIterator it;
822 while ( it.Next(detElemId,manuId) )
824 const AliMpVSegmentation* seg =
825 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
827 AliMUONVCalibParam* calibParam = static_cast<AliMUONVCalibParam*>(neighbourStore.FindObject(detElemId,manuId));
831 Int_t size(AliMpConstants::ManuNofChannels());
832 Int_t defaultValue(-1);
833 Int_t packingFactor(size);
835 calibParam = new AliMUONCalibParamNI(dimension,size,detElemId,manuId,defaultValue,packingFactor);
836 Bool_t ok = neighbourStore.Add(calibParam);
839 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
844 for ( Int_t manuChannel = 0; manuChannel < AliMpConstants::ManuNofChannels(); ++manuChannel )
846 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
852 seg->GetNeighbours(pad,tmp,true,true);
853 Int_t nofPadNeighbours = tmp.GetEntriesFast();
855 for ( Int_t i = 0; i < nofPadNeighbours; ++i )
857 AliMpPad* pad = static_cast<AliMpPad*>(tmp.UncheckedAt(i));
860 calibParam->PackValues(pad->GetLocation().GetFirst(),pad->GetLocation().GetSecond(),x);
863 // AliError("Could not pack value. Something is seriously wrong. Please check");
864 // StdoutToAliError(pad->Print(););
867 calibParam->SetValueAsInt(manuChannel,i,x);
876 //_____________________________________________________________________________
878 AliMUONCDB::SetMaxNofChannelsToGenerate(Int_t n)
880 /// Set the maximum number of channels to generate (used for testing only)
881 /// n < 0 means no limit
882 fMaxNofChannelsToGenerate = n;
885 //_____________________________________________________________________________
887 AliMUONCDB::WriteLocalTriggerMasks(Int_t startRun, Int_t endRun)
889 /// Write local trigger masks to OCDB
891 AliMUONVStore* ltm = new AliMUON1DArray(235);
892 Int_t ngenerated = MakeLocalTriggerMaskStore(*ltm);
893 AliInfo(Form("Ngenerated = %d",ngenerated));
896 WriteToCDB("MUON/Calib/LocalTriggerBoardMasks",ltm,startRun,endRun,true);
901 //_____________________________________________________________________________
903 AliMUONCDB::WriteRegionalTriggerMasks(Int_t startRun, Int_t endRun)
905 /// Write regional trigger masks to OCDB
907 AliMUONVStore* rtm = new AliMUON1DArray(16);
908 Int_t ngenerated = MakeRegionalTriggerMaskStore(*rtm);
909 AliInfo(Form("Ngenerated = %d",ngenerated));
912 WriteToCDB("MUON/Calib/RegionalTriggerBoardMasks",rtm,startRun,endRun,true);
917 //_____________________________________________________________________________
919 AliMUONCDB::WriteGlobalTriggerMasks(Int_t startRun, Int_t endRun)
921 /// Write global trigger masks to OCDB
923 AliMUONVCalibParam* gtm = new AliMUONCalibParamNI(1,2,1,0,0);
925 Int_t ngenerated = MakeGlobalTriggerMaskStore(*gtm);
926 AliInfo(Form("Ngenerated = %d",ngenerated));
929 WriteToCDB("MUON/Calib/GlobalTriggerBoardMasks",gtm,startRun,endRun,true);
934 //_____________________________________________________________________________
936 AliMUONCDB::WriteTriggerLut(Int_t startRun, Int_t endRun)
938 /// Write trigger LUT to OCDB
940 AliMUONTriggerLut* lut = MakeTriggerLUT();
943 WriteToCDB("MUON/Calib/TriggerLut",lut,startRun,endRun,true);
948 //_____________________________________________________________________________
950 AliMUONCDB::WriteTriggerEfficiency(Int_t startRun, Int_t endRun)
952 /// Write trigger efficiency to OCDB
954 AliMUONTriggerEfficiencyCells* eff = MakeTriggerEfficiency();
957 WriteToCDB("MUON/Calib/TriggerEfficiency",eff,startRun,endRun,true);
962 //_____________________________________________________________________________
964 AliMUONCDB::WriteNeighbours(Int_t startRun, Int_t endRun)
966 /// Write neighbours to OCDB
968 AliMUONVStore* neighbours = Create2DMap();
969 Int_t ngenerated = MakeNeighbourStore(*neighbours);
970 AliInfo(Form("Ngenerated = %d",ngenerated));
973 WriteToCDB("MUON/Calib/Neighbours",neighbours,startRun,endRun,true);
978 //_____________________________________________________________________________
980 AliMUONCDB::WriteHV(Bool_t defaultValues,
981 Int_t startRun, Int_t endRun)
983 /// generate HV values (either cste = 1500 V) if defaultValues=true or random
984 /// if defaultValues=false, see makeHVStore) and
985 /// store them into CDB located at cdbpath, with a validity period
986 /// ranging from startRun to endRun
988 TMap* hvStore = new TMap;
989 Int_t ngenerated = MakeHVStore(*hvStore,defaultValues);
990 AliInfo(Form("Ngenerated = %d",ngenerated));
993 WriteToCDB("MUON/Calib/HV",hvStore,startRun,endRun,defaultValues);
998 //_____________________________________________________________________________
1000 AliMUONCDB::WritePedestals(Bool_t defaultValues,
1001 Int_t startRun, Int_t endRun)
1003 /// generate pedestal values (either 0 if defaultValues=true or random
1004 /// if defaultValues=false, see makePedestalStore) and
1005 /// store them into CDB located at cdbpath, with a validity period
1006 /// ranging from startRun to endRun
1008 AliMUONVStore* pedestalStore = Create2DMap();
1009 Int_t ngenerated = MakePedestalStore(*pedestalStore,defaultValues);
1010 AliInfo(Form("Ngenerated = %d",ngenerated));
1011 WriteToCDB("MUON/Calib/Pedestals",pedestalStore,startRun,endRun,defaultValues);
1012 delete pedestalStore;
1016 //_____________________________________________________________________________
1018 AliMUONCDB::WriteGains(Bool_t defaultValues,
1019 Int_t startRun, Int_t endRun)
1021 /// generate gain values (either 1 if defaultValues=true or random
1022 /// if defaultValues=false, see makeGainStore) and
1023 /// store them into CDB located at cdbpath, with a validity period
1024 /// ranging from startRun to endRun
1026 AliMUONVStore* gainStore = Create2DMap();
1027 Int_t ngenerated = MakeGainStore(*gainStore,defaultValues);
1028 AliInfo(Form("Ngenerated = %d",ngenerated));
1029 WriteToCDB("MUON/Calib/Gains",gainStore,startRun,endRun,defaultValues);
1033 //_____________________________________________________________________________
1035 AliMUONCDB::WriteCapacitances(const char* filename,
1036 Int_t startRun, Int_t endRun)
1038 /// read manu capacitance and injection gain values from file
1039 /// and store them into CDB located at cdbpath, with a validity period
1040 /// ranging from startRun to endRun
1042 AliMUONVStore* capaStore = new AliMUON1DMap(16828);
1043 Int_t ngenerated = MakeCapacitanceStore(*capaStore,filename);
1044 AliInfo(Form("Ngenerated = %d",ngenerated));
1045 if ( ngenerated > 0 )
1047 WriteToCDB("MUON/Calib/Capacitances",capaStore,startRun,endRun,filename);
1052 //_____________________________________________________________________________
1054 AliMUONCDB::WriteCapacitances(Bool_t defaultValues,
1055 Int_t startRun, Int_t endRun)
1057 /// generate manu capacitance values (either 1 if defaultValues=true or random
1058 /// if defaultValues=false, see makeCapacitanceStore) and
1059 /// store them into CDB located at cdbpath, with a validity period
1060 /// ranging from startRun to endRun
1062 AliMUONVStore* capaStore = new AliMUON1DMap(16828);
1063 Int_t ngenerated = MakeCapacitanceStore(*capaStore,defaultValues);
1064 AliInfo(Form("Ngenerated = %d",ngenerated));
1065 WriteToCDB("MUON/Calib/Capacitances",capaStore,startRun,endRun,defaultValues);
1069 //_____________________________________________________________________________
1071 AliMUONCDB::WriteTrigger(Int_t startRun, Int_t endRun)
1073 /// Writes all Trigger related calibration to CDB
1074 WriteLocalTriggerMasks(startRun,endRun);
1075 WriteRegionalTriggerMasks(startRun,endRun);
1076 WriteGlobalTriggerMasks(startRun,endRun);
1077 WriteTriggerLut(startRun,endRun);
1078 WriteTriggerEfficiency(startRun,endRun);
1081 //_____________________________________________________________________________
1083 AliMUONCDB::WriteTracker(Bool_t defaultValues, Int_t startRun, Int_t endRun)
1085 /// Writes all Tracker related calibration to CDB
1086 WriteHV(defaultValues,startRun,endRun);
1087 WritePedestals(defaultValues,startRun,endRun);
1088 WriteGains(defaultValues,startRun,endRun);
1089 WriteCapacitances(defaultValues,startRun,endRun);
1090 WriteNeighbours(startRun,endRun);