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 "AliCDBEntry.h"
34 #include "AliCDBManager.h"
35 #include "AliDCSValue.h"
37 #include "AliMUON1DArray.h"
38 #include "AliMUON1DMap.h"
39 #include "AliMUON2DMap.h"
40 #include "AliMUON2DStoreValidator.h"
41 #include "AliMUONCalibParamNF.h"
42 #include "AliMUONCalibParamNI.h"
43 #include "AliMUONConstants.h"
44 #include "AliMUONHVNamer.h"
45 #include "AliMUONTriggerEfficiencyCells.h"
46 #include "AliMUONTriggerLut.h"
47 #include "AliMUONVStore.h"
48 #include "AliMUONVCalibParam.h"
49 #include "AliMUONVCalibParam.h"
51 #include "AliMpConstants.h"
52 #include "AliMpDDLStore.h"
53 #include "AliMpDEIterator.h"
54 #include "AliMpDEManager.h"
55 #include "AliMpDetElement.h"
56 #include "AliMpManuList.h"
57 #include "AliMpSegmentation.h"
58 #include "AliMpStationType.h"
59 #include "AliMpVSegmentation.h"
60 #include <Riostream.h>
66 #include <TObjString.h>
69 #include <TStopwatch.h>
79 //_____________________________________________________________________________
80 void getBoundaries(const AliMUONVStore& store, Int_t dim,
81 Float_t* xmin, Float_t* xmax)
83 /// Assuming the store contains AliMUONVCalibParam objects, compute the
84 /// limits of the value contained in the VCalibParam, for each of its dimensions
85 /// xmin and xmax must be of dimension dim
87 for ( Int_t i = 0; i < dim; ++i )
93 TIter next(store.CreateIterator());
94 AliMUONVCalibParam* value;
96 while ( ( value = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
98 Int_t detElemId = value->ID0();
99 Int_t manuId = value->ID1();
101 const AliMpVSegmentation* seg =
102 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
104 for ( Int_t manuChannel = 0; manuChannel < value->Size(); ++manuChannel )
106 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
107 if (!pad.IsValid()) continue;
109 for ( Int_t i = 0; i < dim; ++i )
111 Float_t x0 = value->ValueAsFloat(manuChannel,i);
113 xmin[i] = TMath::Min(xmin[i],x0);
114 xmax[i] = TMath::Max(xmax[i],x0);
119 for ( Int_t i = 0; i < dim; ++i )
121 if ( TMath::Abs(xmin[i]-xmax[i]) < 1E-3 )
129 //_____________________________________________________________________________
130 Double_t GetRandom(Double_t mean, Double_t sigma, Bool_t mustBePositive)
133 if ( mustBePositive )
137 x = gRandom->Gaus(mean,sigma);
142 x = gRandom->Gaus(mean,sigma);
149 //_____________________________________________________________________________
150 AliMUONCDB::AliMUONCDB(const char* cdbpath)
154 fMaxNofChannelsToGenerate(-1)
159 //_____________________________________________________________________________
160 AliMUONCDB::~AliMUONCDB()
166 //_____________________________________________________________________________
168 AliMUONCDB::ManuList()
170 /// return (and create if necessary) the list of (de,manu) pairs
173 AliInfo("Generating ManuList...");
174 AliCDBManager::Instance()->SetDefaultStorage(fCDBPath);
177 if ( ! AliMpCDB::LoadDDLStore() ) {
178 AliFatal("Could not access mapping from OCDB !");
181 fManuList = AliMpManuList::ManuList();
182 AliInfo("Manu List generated.");
187 //_____________________________________________________________________________
189 AliMUONCDB::Diff(AliMUONVStore& store1, AliMUONVStore& store2,
192 /// creates a store which contains store1-store2
193 /// if opt="abs" the difference is absolute one,
194 /// if opt="rel" then what is stored is (store1-store2)/store1
195 /// if opt="percent" then what is stored is rel*100
197 /// WARNING Works only for stores which holds AliMUONVCalibParam objects
202 if ( !sopt.Contains("ABS") && !sopt.Contains("REL") && !sopt.Contains("PERCENT") )
204 AliErrorClass(Form("opt %s not supported. Only ABS, REL, PERCENT are",opt));
208 AliMUONVStore* d = static_cast<AliMUONVStore*>(store1.Clone());
210 TIter next(d->CreateIterator());
212 AliMUONVCalibParam* param;
214 while ( ( param = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
216 Int_t detElemId = param->ID0();
217 Int_t manuId = param->ID1();
219 AliMUONVCalibParam* param2 = dynamic_cast<AliMUONVCalibParam*>(store2.FindObject(detElemId,manuId));
220 //FIXME: this might happen. Handle it.
223 cerr << "param2 is null : FIXME : this might happen !" << endl;
228 for ( Int_t i = 0; i < param->Size(); ++i )
230 for ( Int_t j = 0; j < param->Dimension(); ++j )
233 if ( sopt.Contains("ABS") )
235 value = param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j);
237 else if ( sopt.Contains("REL") || sopt.Contains("PERCENT") )
239 if ( param->ValueAsFloat(i,j) )
241 value = (param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j))/param->ValueAsFloat(i,j);
247 if ( sopt.Contains("PERCENT") ) value *= 100.0;
249 param->SetValueAsFloat(i,j,value);
256 //_____________________________________________________________________________
258 AliMUONCDB::Plot(const AliMUONVStore& store, const char* name, Int_t nbins)
260 /// Make histograms of each dimension of the AliMUONVCalibParam
261 /// contained inside store.
262 /// It produces histograms named name_0, name_1, etc...
264 TIter next(store.CreateIterator());
265 AliMUONVCalibParam* param;
267 const Int_t kNStations = AliMpConstants::NofTrackingChambers()/2;
268 Int_t* nPerStation = new Int_t[kNStations];
271 for ( Int_t i = 0; i < kNStations; ++i ) nPerStation[i]=0;
273 while ( ( param = static_cast<AliMUONVCalibParam*>(next()) ) )
277 Int_t dim = param->Dimension();
279 Float_t* xmin = new Float_t[dim];
280 Float_t* xmax = new Float_t[dim];
281 getBoundaries(store,dim,xmin,xmax);
283 for ( Int_t i = 0; i < dim; ++i )
285 h[i] = new TH1F(Form("%s_%d",name,i),Form("%s_%d",name,i),
286 nbins,xmin[i],xmax[i]);
287 AliInfo(Form("Created histogram %s",h[i]->GetName()));
291 Int_t detElemId = param->ID0();
292 Int_t manuId = param->ID1();
293 Int_t station = AliMpDEManager::GetChamberId(detElemId)/2;
295 const AliMpVSegmentation* seg =
296 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
298 for ( Int_t manuChannel = 0; manuChannel < param->Size(); ++manuChannel )
300 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
301 if (!pad.IsValid()) continue;
304 ++nPerStation[station];
306 for ( Int_t dim = 0; dim < param->Dimension(); ++dim )
308 h[dim]->Fill(param->ValueAsFloat(manuChannel,dim));
313 for ( Int_t i = 0; i < kNStations; ++i )
315 AliInfo(Form("Station %d %d ",(i+1),nPerStation[i]));
318 AliInfo(Form("Number of channels = %d",n));
320 delete[] nPerStation;
323 //_____________________________________________________________________________
325 AliMUONCDB::MakeHVStore(TMap& aliasMap, Bool_t defaultValues)
327 /// Create a HV store
329 AliMUONHVNamer hvNamer;
331 TObjArray* aliases = hvNamer.GenerateAliases();
336 for ( Int_t i = 0; i < aliases->GetEntries(); ++i )
338 TObjString* alias = static_cast<TObjString*>(aliases->At(i));
339 TString& aliasName = alias->String();
340 if ( aliasName.Contains("sw") )
342 // HV Switch (St345 only)
343 TObjArray* valueSet = new TObjArray;
344 valueSet->SetOwner(kTRUE);
346 Bool_t value = kTRUE;
350 Float_t r = gRandom->Uniform();
351 if ( r < 0.007 ) value = kFALSE;
354 for ( UInt_t timeStamp = 0; timeStamp < 60*3; timeStamp += 60 )
356 AliDCSValue* dcsValue = new AliDCSValue(value,timeStamp);
357 valueSet->Add(dcsValue);
359 aliasMap.Add(new TObjString(*alias),valueSet);
364 TObjArray* valueSet = new TObjArray;
365 valueSet->SetOwner(kTRUE);
366 for ( UInt_t timeStamp = 0; timeStamp < 60*15; timeStamp += 120 )
368 Float_t value = 1500;
369 if (!defaultValues) value = GetRandom(1750,62.5,true);
370 AliDCSValue* dcsValue = new AliDCSValue(value,timeStamp);
371 valueSet->Add(dcsValue);
373 aliasMap.Add(new TObjString(*alias),valueSet);
380 AliInfo(Form("%d HV channels and %d switches",nChannels,nSwitch));
382 return nChannels+nSwitch;
385 //_____________________________________________________________________________
387 AliMUONCDB::MakePedestalStore(AliMUONVStore& pedestalStore, Bool_t defaultValues)
389 /// Create a pedestal store. if defaultValues=true, ped.mean=ped.sigma=1,
390 /// otherwise mean and sigma are from a gaussian (with parameters
391 /// defined below by the kPedestal* constants)
393 TIter next(ManuList());
400 const Int_t kChannels(AliMpConstants::ManuNofChannels());
401 const Float_t kPedestalMeanMean(200);
402 const Float_t kPedestalMeanSigma(10);
403 const Float_t kPedestalSigmaMean(1.0);
404 const Float_t kPedestalSigmaSigma(0.2);
406 while ( ( p = (AliMpIntPair*)next() ) )
410 Int_t detElemId = p->GetFirst();
411 Int_t manuId = p->GetSecond();
413 AliMUONVCalibParam* ped =
414 new AliMUONCalibParamNF(2,kChannels,detElemId,manuId,AliMUONVCalibParam::InvalidFloatValue());
416 const AliMpVSegmentation* seg =
417 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
419 for ( Int_t manuChannel = 0; manuChannel < kChannels; ++manuChannel )
421 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
422 if (!pad.IsValid()) continue;
426 Float_t meanPedestal;
427 Float_t sigmaPedestal;
436 Bool_t positive(kTRUE);
438 while ( meanPedestal == 0.0 ) // avoid strict zero
440 meanPedestal = GetRandom(kPedestalMeanMean,kPedestalMeanSigma,positive);
442 sigmaPedestal = GetRandom(kPedestalSigmaMean,kPedestalSigmaSigma,positive);
444 ped->SetValueAsFloat(manuChannel,0,meanPedestal);
445 ped->SetValueAsFloat(manuChannel,1,sigmaPedestal);
448 Bool_t ok = pedestalStore.Add(ped);
451 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
453 if ( fMaxNofChannelsToGenerate > 0 && nchannels >= fMaxNofChannelsToGenerate ) break;
456 AliInfo(Form("%d Manus and %d channels.",nmanus,nchannels));
461 //_____________________________________________________________________________
463 AliMUONCDB::MakeCapacitanceStore(AliMUONVStore& capaStore, const char* file)
465 /// Read the capacitance values from file and append them to the capaStore
467 ifstream in(gSystem->ExpandPathName(file));
468 if (in.bad()) return 0;
473 Int_t serialNumber(-1);
474 AliMUONVCalibParam* param(0x0);
476 while ( in.getline(line,1024,'\n') )
478 if ( isdigit(line[0]) )
480 serialNumber = atoi(line);
481 param = static_cast<AliMUONVCalibParam*>(capaStore.FindObject(serialNumber));
484 AliError(Form("serialNumber %d appears several times !",serialNumber));
488 param = new AliMUONCalibParamNF(2,AliMpConstants::ManuNofChannels(),serialNumber,0,1.0);
489 Bool_t ok = capaStore.Add(param);
492 AliError(Form("Could not set serialNumber=%d",serialNumber));
499 Float_t injectionGain;
500 sscanf(line,"%d %f %f",&channel,&capaValue,&injectionGain);
501 AliDebug(1,Form("SerialNumber %10d Channel %3d Capa %f injectionGain %f",
502 serialNumber,channel,capaValue,injectionGain));
503 param->SetValueAsFloat(channel,0,capaValue);
504 param->SetValueAsFloat(channel,1,injectionGain);
513 //_____________________________________________________________________________
515 AliMUONCDB::MakeCapacitanceStore(AliMUONVStore& capaStore, Bool_t defaultValues)
517 /// Create a capacitance store. if defaultValues=true, all capa are 1.0,
518 /// otherwise they are from a gaussian with parameters defined in the
519 /// kCapa* constants below.
521 TIter next(ManuList());
527 Int_t nmanusOK(0); // manus for which we got the serial number
529 const Float_t kCapaMean(0.3);
530 const Float_t kCapaSigma(0.1);
531 const Float_t kInjectionGainMean(3);
532 const Float_t kInjectionGainSigma(1);
534 while ( ( p = (AliMpIntPair*)next() ) )
538 Int_t detElemId = p->GetFirst();
539 Int_t manuId = p->GetSecond();
541 AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(detElemId);
542 Int_t serialNumber = de->GetManuSerialFromId(manuId);
544 if ( serialNumber <= 0 ) continue;
548 const AliMpVSegmentation* seg =
549 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
551 AliMUONVCalibParam* capa = static_cast<AliMUONVCalibParam*>(capaStore.FindObject(serialNumber));
555 capa = new AliMUONCalibParamNF(2,AliMpConstants::ManuNofChannels(),serialNumber,0,1.0);
556 Bool_t ok = capaStore.Add(capa);
559 AliError(Form("Could not set serialNumber=%d manuId=%d",serialNumber,manuId));
563 for ( Int_t manuChannel = 0; manuChannel < capa->Size(); ++manuChannel )
565 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
566 if (!pad.IsValid()) continue;
571 Float_t injectionGain;
580 capaValue = GetRandom(kCapaMean,kCapaSigma,kTRUE);
581 injectionGain = GetRandom(kInjectionGainMean,kInjectionGainSigma,kTRUE);
583 capa->SetValueAsFloat(manuChannel,0,capaValue);
584 capa->SetValueAsFloat(manuChannel,1,injectionGain);
589 if ( nmanus ) percent = 100*nmanusOK/nmanus;
590 AliInfo(Form("%5d manus with serial number (out of %5d manus = %3.0f%%)",
591 nmanusOK,nmanus,percent));
592 AliInfo(Form("%5d channels",nchannels));
595 AliWarning("Did not get all serial numbers. capaStore is incomplete !!!!");
601 //_____________________________________________________________________________
603 AliMUONCDB::MakeGainStore(AliMUONVStore& gainStore, Bool_t defaultValues)
605 /// Create a gain store. if defaultValues=true, all gains set so that
606 /// charge = (adc-ped)
608 /// otherwise parameters are taken from gaussians with parameters
609 /// defined in the k* constants below.
611 TIter next(ManuList());
618 const Int_t kSaturation(3000);
619 const Double_t kA0Mean(1.2);
620 const Double_t kA0Sigma(0.1);
621 const Double_t kA1Mean(1E-5);
622 const Double_t kA1Sigma(1E-6);
623 const Double_t kQualMean(0xFF);
624 const Double_t kQualSigma(0x10);
625 const Int_t kThresMean(1600);
626 const Int_t kThresSigma(100);
628 while ( ( p = (AliMpIntPair*)next() ) )
632 Int_t detElemId = p->GetFirst();
633 Int_t manuId = p->GetSecond();
635 AliMUONVCalibParam* gain =
636 new AliMUONCalibParamNF(5,AliMpConstants::ManuNofChannels(),
639 AliMUONVCalibParam::InvalidFloatValue());
642 const AliMpVSegmentation* seg =
643 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
645 for ( Int_t manuChannel = 0; manuChannel < gain->Size(); ++manuChannel )
647 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
648 if (!pad.IsValid()) continue;
654 gain->SetValueAsFloat(manuChannel,0,1.0);
655 gain->SetValueAsFloat(manuChannel,1,0.0);
656 gain->SetValueAsInt(manuChannel,2,4095);
657 gain->SetValueAsInt(manuChannel,3,1);
658 gain->SetValueAsInt(manuChannel,4,kSaturation);
662 Bool_t positive(kTRUE);
663 gain->SetValueAsFloat(manuChannel,0,GetRandom(kA0Mean,kA0Sigma,positive));
664 gain->SetValueAsFloat(manuChannel,1,GetRandom(kA1Mean,kA1Sigma,!positive));
665 gain->SetValueAsInt(manuChannel,2,(Int_t)TMath::Nint(GetRandom(kThresMean,kThresSigma,positive)));
666 gain->SetValueAsInt(manuChannel,3,(Int_t)TMath::Nint(GetRandom(kQualMean,kQualSigma,positive)));
667 gain->SetValueAsInt(manuChannel,4,kSaturation);
671 Bool_t ok = gainStore.Add(gain);
674 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
676 if ( fMaxNofChannelsToGenerate > 0 && nchannels >= fMaxNofChannelsToGenerate ) break;
679 AliInfo(Form("%d Manus and %d channels.",nmanus,nchannels));
683 //_____________________________________________________________________________
685 AliMUONCDB::MakeLocalTriggerMaskStore(AliMUONVStore& localBoardMasks) const
687 /// Generate local trigger masks store. All masks are set to FFFF
690 // Generate fake mask values for 234 localboards and put that into
691 // one single container (localBoardMasks)
692 for ( Int_t i = 1; i <= 234; ++i )
694 AliMUONVCalibParam* localBoard = new AliMUONCalibParamNI(1,8,i,0,0);
695 for ( Int_t x = 0; x < 2; ++x )
697 for ( Int_t y = 0; y < 4; ++y )
700 localBoard->SetValueAsInt(index,0,0xFFFF);
704 localBoardMasks.Add(localBoard);
709 //_____________________________________________________________________________
711 AliMUONCDB::MakeRegionalTriggerMaskStore(AliMUONVStore& rtm) const
713 /// Make a regional trigger masks store. All masks are set to 3F
716 for ( Int_t i = 0; i < 16; ++i )
718 AliMUONVCalibParam* regionalBoard = new AliMUONCalibParamNI(1,16,i,0,0);
719 for ( Int_t j = 0; j < 16; ++j )
721 regionalBoard->SetValueAsInt(j,0,0x3F);
724 rtm.Add(regionalBoard);
730 //_____________________________________________________________________________
732 AliMUONCDB::MakeGlobalTriggerMaskStore(AliMUONVCalibParam& gtm) const
734 /// Make a global trigger masks store. All masks set to FFF
738 for ( Int_t j = 0; j < 16; ++j )
740 gtm.SetValueAsInt(j,0,0xFFF);
746 //_____________________________________________________________________________
748 AliMUONCDB::MakeTriggerLUT(const char* file) const
750 /// Make a triggerlut object, from a file.
752 AliMUONTriggerLut* lut = new AliMUONTriggerLut;
753 lut->ReadFromFile(file);
757 //_____________________________________________________________________________
758 AliMUONTriggerEfficiencyCells*
759 AliMUONCDB::MakeTriggerEfficiency(const char* file) const
761 /// Make a trigger efficiency object from a file.
763 return new AliMUONTriggerEfficiencyCells(file);
766 //_____________________________________________________________________________
768 AliMUONCDB::WriteToCDB(const char* calibpath, TObject* object,
769 Int_t startRun, Int_t endRun,
770 const char* filename)
772 /// Write a given object to OCDB
774 AliCDBId id(calibpath,startRun,endRun);
776 md.SetAliRootVersion(gROOT->GetVersion());
777 md.SetComment(gSystem->ExpandPathName(filename));
778 md.SetResponsible("Uploaded using AliMUONCDB class");
779 AliCDBManager* man = AliCDBManager::Instance();
780 man->SetDefaultStorage(fCDBPath);
781 man->Put(object,id,&md);
784 //_____________________________________________________________________________
786 AliMUONCDB::WriteToCDB(const char* calibpath, TObject* object,
787 Int_t startRun, Int_t endRun, Bool_t defaultValues)
789 /// Write a given object to OCDB
791 AliCDBId id(calibpath,startRun,endRun);
793 md.SetAliRootVersion(gROOT->GetVersion());
796 md.SetComment("Test with default values");
800 md.SetComment("Test with random values");
802 md.SetResponsible("AliMUONCDB tester class");
804 AliCDBManager* man = AliCDBManager::Instance();
805 man->SetDefaultStorage(fCDBPath);
806 man->Put(object,id,&md);
809 //_____________________________________________________________________________
811 AliMUONCDB::MakeNeighbourStore(AliMUONVStore& neighbourStore)
813 /// Fill the neighbours store with, for each channel, a TObjArray of its
814 /// neighbouring pads (including itself)
816 AliInfo("Generating NeighbourStore. This will take a while. Please be patient.");
822 TIter next(ManuList());
830 while ( ( p = (AliMpIntPair*)next() ) )
832 Int_t detElemId = p->GetFirst();
833 Int_t manuId = p->GetSecond();
835 const AliMpVSegmentation* seg =
836 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
838 AliMUONVCalibParam* calibParam = static_cast<AliMUONVCalibParam*>(neighbourStore.FindObject(detElemId,manuId));
842 Int_t size(AliMpConstants::ManuNofChannels());
843 Int_t defaultValue(-1);
844 Int_t packingFactor(size);
846 calibParam = new AliMUONCalibParamNI(dimension,size,detElemId,manuId,defaultValue,packingFactor);
847 Bool_t ok = neighbourStore.Add(calibParam);
850 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
855 for ( Int_t manuChannel = 0; manuChannel < AliMpConstants::ManuNofChannels(); ++manuChannel )
857 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
863 seg->GetNeighbours(pad,tmp,true,true);
864 Int_t nofPadNeighbours = tmp.GetEntriesFast();
866 for ( Int_t i = 0; i < nofPadNeighbours; ++i )
868 AliMpPad* pad = static_cast<AliMpPad*>(tmp.At(i));
870 Bool_t ok = calibParam->PackValues(pad->GetLocation().GetFirst(),pad->GetLocation().GetSecond(),x);
873 AliError("Could not pack value. Something is seriously wrong. Please check");
874 StdoutToAliError(pad->Print(););
877 calibParam->SetValueAsInt(manuChannel,i,x);
888 //_____________________________________________________________________________
890 AliMUONCDB::SetMaxNofChannelsToGenerate(Int_t n)
892 /// Set the maximum number of channels to generate (used for testing only)
893 /// n < 0 means no limit
894 fMaxNofChannelsToGenerate = n;
897 //_____________________________________________________________________________
899 AliMUONCDB::WriteLocalTriggerMasks(Int_t startRun, Int_t endRun)
901 /// Write local trigger masks to OCDB
903 AliMUONVStore* ltm = new AliMUON1DArray(235);
904 Int_t ngenerated = MakeLocalTriggerMaskStore(*ltm);
905 AliInfo(Form("Ngenerated = %d",ngenerated));
908 WriteToCDB("MUON/Calib/LocalTriggerBoardMasks",ltm,startRun,endRun,true);
913 //_____________________________________________________________________________
915 AliMUONCDB::WriteRegionalTriggerMasks(Int_t startRun, Int_t endRun)
917 /// Write regional trigger masks to OCDB
919 AliMUONVStore* rtm = new AliMUON1DArray(16);
920 Int_t ngenerated = MakeRegionalTriggerMaskStore(*rtm);
921 AliInfo(Form("Ngenerated = %d",ngenerated));
924 WriteToCDB("MUON/Calib/RegionalTriggerBoardMasks",rtm,startRun,endRun,true);
929 //_____________________________________________________________________________
931 AliMUONCDB::WriteGlobalTriggerMasks(Int_t startRun, Int_t endRun)
933 /// Write global trigger masks to OCDB
935 AliMUONVCalibParam* gtm = new AliMUONCalibParamNI(1,16,1,0,0);
937 Int_t ngenerated = MakeGlobalTriggerMaskStore(*gtm);
938 AliInfo(Form("Ngenerated = %d",ngenerated));
941 WriteToCDB("MUON/Calib/GlobalTriggerBoardMasks",gtm,startRun,endRun,true);
946 //_____________________________________________________________________________
948 AliMUONCDB::WriteTriggerLut(Int_t startRun, Int_t endRun)
950 /// Write trigger LUT to OCDB
952 AliMUONTriggerLut* lut = MakeTriggerLUT();
955 WriteToCDB("MUON/Calib/TriggerLut",lut,startRun,endRun,true);
960 //_____________________________________________________________________________
962 AliMUONCDB::WriteTriggerEfficiency(Int_t startRun, Int_t endRun)
964 /// Write trigger efficiency to OCDB
966 AliMUONTriggerEfficiencyCells* eff = MakeTriggerEfficiency();
969 WriteToCDB("MUON/Calib/TriggerEfficiency",eff,startRun,endRun,true);
974 //_____________________________________________________________________________
976 AliMUONCDB::WriteNeighbours(Int_t startRun, Int_t endRun)
978 /// Write neighbours to OCDB
980 AliMUONVStore* neighbours = new AliMUON2DMap(kTRUE);
981 Int_t ngenerated = MakeNeighbourStore(*neighbours);
982 AliInfo(Form("Ngenerated = %d",ngenerated));
985 WriteToCDB("MUON/Calib/Neighbours",neighbours,startRun,endRun,true);
990 //_____________________________________________________________________________
992 AliMUONCDB::WriteHV(Bool_t defaultValues,
993 Int_t startRun, Int_t endRun)
995 /// generate HV values (either cste = 1500 V) if defaultValues=true or random
996 /// if defaultValues=false, see makeHVStore) and
997 /// store them into CDB located at cdbpath, with a validity period
998 /// ranging from startRun to endRun
1000 TMap* hvStore = new TMap;
1001 Int_t ngenerated = MakeHVStore(*hvStore,defaultValues);
1002 AliInfo(Form("Ngenerated = %d",ngenerated));
1005 WriteToCDB("MUON/Calib/HV",hvStore,startRun,endRun,defaultValues);
1010 //_____________________________________________________________________________
1012 AliMUONCDB::WritePedestals(Bool_t defaultValues,
1013 Int_t startRun, Int_t endRun)
1015 /// generate pedestal values (either 0 if defaultValues=true or random
1016 /// if defaultValues=false, see makePedestalStore) and
1017 /// store them into CDB located at cdbpath, with a validity period
1018 /// ranging from startRun to endRun
1020 AliMUONVStore* pedestalStore = new AliMUON2DMap(true);
1021 Int_t ngenerated = MakePedestalStore(*pedestalStore,defaultValues);
1022 AliInfo(Form("Ngenerated = %d",ngenerated));
1023 WriteToCDB("MUON/Calib/Pedestals",pedestalStore,startRun,endRun,defaultValues);
1024 delete pedestalStore;
1028 //_____________________________________________________________________________
1030 AliMUONCDB::WriteGains(Bool_t defaultValues,
1031 Int_t startRun, Int_t endRun)
1033 /// generate gain values (either 1 if defaultValues=true or random
1034 /// if defaultValues=false, see makeGainStore) and
1035 /// store them into CDB located at cdbpath, with a validity period
1036 /// ranging from startRun to endRun
1038 AliMUONVStore* gainStore = new AliMUON2DMap(true);
1039 Int_t ngenerated = MakeGainStore(*gainStore,defaultValues);
1040 AliInfo(Form("Ngenerated = %d",ngenerated));
1041 WriteToCDB("MUON/Calib/Gains",gainStore,startRun,endRun,defaultValues);
1045 //_____________________________________________________________________________
1047 AliMUONCDB::WriteCapacitances(const char* filename,
1048 Int_t startRun, Int_t endRun)
1050 /// read manu capacitance and injection gain values from file
1051 /// and store them into CDB located at cdbpath, with a validity period
1052 /// ranging from startRun to endRun
1054 AliMUONVStore* capaStore = new AliMUON1DMap(16828);
1055 Int_t ngenerated = MakeCapacitanceStore(*capaStore,filename);
1056 AliInfo(Form("Ngenerated = %d",ngenerated));
1057 if ( ngenerated > 0 )
1059 WriteToCDB("MUON/Calib/Capacitances",capaStore,startRun,endRun,filename);
1064 //_____________________________________________________________________________
1066 AliMUONCDB::WriteCapacitances(Bool_t defaultValues,
1067 Int_t startRun, Int_t endRun)
1069 /// generate manu capacitance values (either 1 if defaultValues=true or random
1070 /// if defaultValues=false, see makeCapacitanceStore) and
1071 /// store them into CDB located at cdbpath, with a validity period
1072 /// ranging from startRun to endRun
1074 AliMUONVStore* capaStore = new AliMUON1DMap(16828);
1075 Int_t ngenerated = MakeCapacitanceStore(*capaStore,defaultValues);
1076 AliInfo(Form("Ngenerated = %d",ngenerated));
1077 WriteToCDB("MUON/Calib/Capacitances",capaStore,startRun,endRun,defaultValues);
1081 //_____________________________________________________________________________
1083 AliMUONCDB::WriteTrigger(Int_t startRun, Int_t endRun)
1085 /// Writes all Trigger related calibration to CDB
1086 WriteLocalTriggerMasks(startRun,endRun);
1087 WriteRegionalTriggerMasks(startRun,endRun);
1088 WriteGlobalTriggerMasks(startRun,endRun);
1089 WriteTriggerLut(startRun,endRun);
1090 WriteTriggerEfficiency(startRun,endRun);
1093 //_____________________________________________________________________________
1095 AliMUONCDB::WriteTracker(Bool_t defaultValues, Int_t startRun, Int_t endRun)
1097 /// Writes all Tracker related calibration to CDB
1098 WriteHV(defaultValues,startRun,endRun);
1099 WritePedestals(defaultValues,startRun,endRun);
1100 WriteGains(defaultValues,startRun,endRun);
1101 WriteCapacitances(defaultValues,startRun,endRun);
1102 WriteNeighbours(startRun,endRun);