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);
175 AliMpCDB::LoadMpSegmentation();
176 AliMpCDB::LoadDDLStore();
177 fManuList = AliMpManuList::ManuList();
178 AliInfo("Manu List generated.");
183 //_____________________________________________________________________________
185 AliMUONCDB::Diff(AliMUONVStore& store1, AliMUONVStore& store2,
188 /// creates a store which contains store1-store2
189 /// if opt="abs" the difference is absolute one,
190 /// if opt="rel" then what is stored is (store1-store2)/store1
191 /// if opt="percent" then what is stored is rel*100
193 /// WARNING Works only for stores which holds AliMUONVCalibParam objects
198 if ( !sopt.Contains("ABS") && !sopt.Contains("REL") && !sopt.Contains("PERCENT") )
200 AliErrorClass(Form("opt %s not supported. Only ABS, REL, PERCENT are",opt));
204 AliMUONVStore* d = static_cast<AliMUONVStore*>(store1.Clone());
206 TIter next(d->CreateIterator());
208 AliMUONVCalibParam* param;
210 while ( ( param = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
212 Int_t detElemId = param->ID0();
213 Int_t manuId = param->ID1();
215 AliMUONVCalibParam* param2 = dynamic_cast<AliMUONVCalibParam*>(store2.FindObject(detElemId,manuId));
216 //FIXME: this might happen. Handle it.
219 cerr << "param2 is null : FIXME : this might happen !" << endl;
224 for ( Int_t i = 0; i < param->Size(); ++i )
226 for ( Int_t j = 0; j < param->Dimension(); ++j )
229 if ( sopt.Contains("ABS") )
231 value = param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j);
233 else if ( sopt.Contains("REL") || sopt.Contains("PERCENT") )
235 if ( param->ValueAsFloat(i,j) )
237 value = (param->ValueAsFloat(i,j) - param2->ValueAsFloat(i,j))/param->ValueAsFloat(i,j);
243 if ( sopt.Contains("PERCENT") ) value *= 100.0;
245 param->SetValueAsFloat(i,j,value);
252 //_____________________________________________________________________________
254 AliMUONCDB::Plot(const AliMUONVStore& store, const char* name, Int_t nbins)
256 /// Make histograms of each dimension of the AliMUONVCalibParam
257 /// contained inside store.
258 /// It produces histograms named name_0, name_1, etc...
260 TIter next(store.CreateIterator());
261 AliMUONVCalibParam* param;
263 const Int_t kNStations = AliMpConstants::NofTrackingChambers()/2;
264 Int_t* nPerStation = new Int_t[kNStations];
267 for ( Int_t i = 0; i < kNStations; ++i ) nPerStation[i]=0;
269 while ( ( param = static_cast<AliMUONVCalibParam*>(next()) ) )
273 Int_t dim = param->Dimension();
275 Float_t* xmin = new Float_t[dim];
276 Float_t* xmax = new Float_t[dim];
277 getBoundaries(store,dim,xmin,xmax);
279 for ( Int_t i = 0; i < dim; ++i )
281 h[i] = new TH1F(Form("%s_%d",name,i),Form("%s_%d",name,i),
282 nbins,xmin[i],xmax[i]);
283 AliInfo(Form("Created histogram %s",h[i]->GetName()));
287 Int_t detElemId = param->ID0();
288 Int_t manuId = param->ID1();
289 Int_t station = AliMpDEManager::GetChamberId(detElemId)/2;
291 const AliMpVSegmentation* seg =
292 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
294 for ( Int_t manuChannel = 0; manuChannel < param->Size(); ++manuChannel )
296 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
297 if (!pad.IsValid()) continue;
300 ++nPerStation[station];
302 for ( Int_t dim = 0; dim < param->Dimension(); ++dim )
304 h[dim]->Fill(param->ValueAsFloat(manuChannel,dim));
309 for ( Int_t i = 0; i < kNStations; ++i )
311 AliInfo(Form("Station %d %d ",(i+1),nPerStation[i]));
314 AliInfo(Form("Number of channels = %d",n));
316 delete[] nPerStation;
319 //_____________________________________________________________________________
321 AliMUONCDB::MakeHVStore(TMap& aliasMap, Bool_t defaultValues)
323 /// Create a HV store
325 AliMUONHVNamer hvNamer;
327 TObjArray* aliases = hvNamer.GenerateAliases();
332 for ( Int_t i = 0; i < aliases->GetEntries(); ++i )
334 TObjString* alias = static_cast<TObjString*>(aliases->At(i));
335 TString& aliasName = alias->String();
336 if ( aliasName.Contains("sw") )
338 // HV Switch (St345 only)
339 TObjArray* valueSet = new TObjArray;
340 valueSet->SetOwner(kTRUE);
342 Bool_t value = kTRUE;
346 Float_t r = gRandom->Uniform();
347 if ( r < 0.007 ) value = kFALSE;
350 for ( UInt_t timeStamp = 0; timeStamp < 60*3; timeStamp += 60 )
352 AliDCSValue* dcsValue = new AliDCSValue(value,timeStamp);
353 valueSet->Add(dcsValue);
355 aliasMap.Add(new TObjString(*alias),valueSet);
360 TObjArray* valueSet = new TObjArray;
361 valueSet->SetOwner(kTRUE);
362 for ( UInt_t timeStamp = 0; timeStamp < 60*15; timeStamp += 120 )
364 Float_t value = 1500;
365 if (!defaultValues) value = GetRandom(1750,62.5,true);
366 AliDCSValue* dcsValue = new AliDCSValue(value,timeStamp);
367 valueSet->Add(dcsValue);
369 aliasMap.Add(new TObjString(*alias),valueSet);
376 AliInfo(Form("%d HV channels and %d switches",nChannels,nSwitch));
378 return nChannels+nSwitch;
381 //_____________________________________________________________________________
383 AliMUONCDB::MakePedestalStore(AliMUONVStore& pedestalStore, Bool_t defaultValues)
385 /// Create a pedestal store. if defaultValues=true, ped.mean=ped.sigma=1,
386 /// otherwise mean and sigma are from a gaussian (with parameters
387 /// defined below by the kPedestal* constants)
389 TIter next(ManuList());
396 const Int_t kChannels(AliMpConstants::ManuNofChannels());
397 const Float_t kPedestalMeanMean(200);
398 const Float_t kPedestalMeanSigma(10);
399 const Float_t kPedestalSigmaMean(1.0);
400 const Float_t kPedestalSigmaSigma(0.2);
402 while ( ( p = (AliMpIntPair*)next() ) )
406 Int_t detElemId = p->GetFirst();
407 Int_t manuId = p->GetSecond();
409 AliMUONVCalibParam* ped =
410 new AliMUONCalibParamNF(2,kChannels,detElemId,manuId,AliMUONVCalibParam::InvalidFloatValue());
412 const AliMpVSegmentation* seg =
413 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
415 for ( Int_t manuChannel = 0; manuChannel < kChannels; ++manuChannel )
417 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
418 if (!pad.IsValid()) continue;
422 Float_t meanPedestal;
423 Float_t sigmaPedestal;
432 Bool_t positive(kTRUE);
434 while ( meanPedestal == 0.0 ) // avoid strict zero
436 meanPedestal = GetRandom(kPedestalMeanMean,kPedestalMeanSigma,positive);
438 sigmaPedestal = GetRandom(kPedestalSigmaMean,kPedestalSigmaSigma,positive);
440 ped->SetValueAsFloat(manuChannel,0,meanPedestal);
441 ped->SetValueAsFloat(manuChannel,1,sigmaPedestal);
444 Bool_t ok = pedestalStore.Add(ped);
447 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
449 if ( fMaxNofChannelsToGenerate > 0 && nchannels >= fMaxNofChannelsToGenerate ) break;
452 AliInfo(Form("%d Manus and %d channels.",nmanus,nchannels));
457 //_____________________________________________________________________________
459 AliMUONCDB::MakeCapacitanceStore(AliMUONVStore& capaStore, const char* file)
461 /// Read the capacitance values from file and append them to the capaStore
463 ifstream in(gSystem->ExpandPathName(file));
464 if (in.bad()) return 0;
469 Int_t serialNumber(-1);
470 AliMUONVCalibParam* param(0x0);
472 while ( in.getline(line,1024,'\n') )
474 if ( isdigit(line[0]) )
476 serialNumber = atoi(line);
477 param = static_cast<AliMUONVCalibParam*>(capaStore.FindObject(serialNumber));
480 AliError(Form("serialNumber %d appears several times !",serialNumber));
484 param = new AliMUONCalibParamNF(2,AliMpConstants::ManuNofChannels(),serialNumber,0,1.0);
485 Bool_t ok = capaStore.Add(param);
488 AliError(Form("Could not set serialNumber=%d",serialNumber));
495 Float_t injectionGain;
496 sscanf(line,"%d %f %f",&channel,&capaValue,&injectionGain);
497 AliDebug(1,Form("SerialNumber %10d Channel %3d Capa %f injectionGain %f",
498 serialNumber,channel,capaValue,injectionGain));
499 param->SetValueAsFloat(channel,0,capaValue);
500 param->SetValueAsFloat(channel,1,injectionGain);
509 //_____________________________________________________________________________
511 AliMUONCDB::MakeCapacitanceStore(AliMUONVStore& capaStore, Bool_t defaultValues)
513 /// Create a capacitance store. if defaultValues=true, all capa are 1.0,
514 /// otherwise they are from a gaussian with parameters defined in the
515 /// kCapa* constants below.
517 TIter next(ManuList());
523 Int_t nmanusOK(0); // manus for which we got the serial number
525 const Float_t kCapaMean(0.3);
526 const Float_t kCapaSigma(0.1);
527 const Float_t kInjectionGainMean(3);
528 const Float_t kInjectionGainSigma(1);
530 while ( ( p = (AliMpIntPair*)next() ) )
534 Int_t detElemId = p->GetFirst();
535 Int_t manuId = p->GetSecond();
537 AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(detElemId);
538 Int_t serialNumber = de->GetManuSerialFromId(manuId);
540 if ( serialNumber <= 0 ) continue;
544 const AliMpVSegmentation* seg =
545 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
547 AliMUONVCalibParam* capa = static_cast<AliMUONVCalibParam*>(capaStore.FindObject(serialNumber));
551 capa = new AliMUONCalibParamNF(2,AliMpConstants::ManuNofChannels(),serialNumber,0,1.0);
552 Bool_t ok = capaStore.Add(capa);
555 AliError(Form("Could not set serialNumber=%d manuId=%d",serialNumber,manuId));
559 for ( Int_t manuChannel = 0; manuChannel < capa->Size(); ++manuChannel )
561 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
562 if (!pad.IsValid()) continue;
567 Float_t injectionGain;
576 capaValue = GetRandom(kCapaMean,kCapaSigma,kTRUE);
577 injectionGain = GetRandom(kInjectionGainMean,kInjectionGainSigma,kTRUE);
579 capa->SetValueAsFloat(manuChannel,0,capaValue);
580 capa->SetValueAsFloat(manuChannel,1,injectionGain);
585 if ( nmanus ) percent = 100*nmanusOK/nmanus;
586 AliInfo(Form("%5d manus with serial number (out of %5d manus = %3.0f%%)",
587 nmanusOK,nmanus,percent));
588 AliInfo(Form("%5d channels",nchannels));
591 AliWarning("Did not get all serial numbers. capaStore is incomplete !!!!");
597 //_____________________________________________________________________________
599 AliMUONCDB::MakeGainStore(AliMUONVStore& gainStore, Bool_t defaultValues)
601 /// Create a gain store. if defaultValues=true, all gains set so that
602 /// charge = (adc-ped)
604 /// otherwise parameters are taken from gaussians with parameters
605 /// defined in the k* constants below.
607 TIter next(ManuList());
614 const Int_t kSaturation(3000);
615 const Double_t kA0Mean(1.2);
616 const Double_t kA0Sigma(0.1);
617 const Double_t kA1Mean(1E-5);
618 const Double_t kA1Sigma(1E-6);
619 const Double_t kQualMean(0xFF);
620 const Double_t kQualSigma(0x10);
621 const Int_t kThresMean(1600);
622 const Int_t kThresSigma(100);
624 while ( ( p = (AliMpIntPair*)next() ) )
628 Int_t detElemId = p->GetFirst();
629 Int_t manuId = p->GetSecond();
631 AliMUONVCalibParam* gain =
632 new AliMUONCalibParamNF(5,AliMpConstants::ManuNofChannels(),
635 AliMUONVCalibParam::InvalidFloatValue());
638 const AliMpVSegmentation* seg =
639 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
641 for ( Int_t manuChannel = 0; manuChannel < gain->Size(); ++manuChannel )
643 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
644 if (!pad.IsValid()) continue;
650 gain->SetValueAsFloat(manuChannel,0,1.0);
651 gain->SetValueAsFloat(manuChannel,1,0.0);
652 gain->SetValueAsInt(manuChannel,2,4095);
653 gain->SetValueAsInt(manuChannel,3,1);
654 gain->SetValueAsInt(manuChannel,4,kSaturation);
658 Bool_t positive(kTRUE);
659 gain->SetValueAsFloat(manuChannel,0,GetRandom(kA0Mean,kA0Sigma,positive));
660 gain->SetValueAsFloat(manuChannel,1,GetRandom(kA1Mean,kA1Sigma,!positive));
661 gain->SetValueAsInt(manuChannel,2,(Int_t)TMath::Nint(GetRandom(kThresMean,kThresSigma,positive)));
662 gain->SetValueAsInt(manuChannel,3,(Int_t)TMath::Nint(GetRandom(kQualMean,kQualSigma,positive)));
663 gain->SetValueAsInt(manuChannel,4,kSaturation);
667 Bool_t ok = gainStore.Add(gain);
670 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
672 if ( fMaxNofChannelsToGenerate > 0 && nchannels >= fMaxNofChannelsToGenerate ) break;
675 AliInfo(Form("%d Manus and %d channels.",nmanus,nchannels));
679 //_____________________________________________________________________________
681 AliMUONCDB::MakeLocalTriggerMaskStore(AliMUONVStore& localBoardMasks) const
683 /// Generate local trigger masks store. All masks are set to FFFF
686 // Generate fake mask values for 234 localboards and put that into
687 // one single container (localBoardMasks)
688 for ( Int_t i = 1; i <= 234; ++i )
690 AliMUONVCalibParam* localBoard = new AliMUONCalibParamNI(1,8,i,0,0);
691 for ( Int_t x = 0; x < 2; ++x )
693 for ( Int_t y = 0; y < 4; ++y )
696 localBoard->SetValueAsInt(index,0,0xFFFF);
700 localBoardMasks.Add(localBoard);
705 //_____________________________________________________________________________
707 AliMUONCDB::MakeRegionalTriggerMaskStore(AliMUONVStore& rtm) const
709 /// Make a regional trigger masks store. All masks are set to 3F
712 for ( Int_t i = 0; i < 16; ++i )
714 AliMUONVCalibParam* regionalBoard = new AliMUONCalibParamNI(1,16,i,0,0);
715 for ( Int_t j = 0; j < 16; ++j )
717 regionalBoard->SetValueAsInt(j,0,0x3F);
720 rtm.Add(regionalBoard);
726 //_____________________________________________________________________________
728 AliMUONCDB::MakeGlobalTriggerMaskStore(AliMUONVCalibParam& gtm) const
730 /// Make a global trigger masks store. All masks set to FFF
734 for ( Int_t j = 0; j < 16; ++j )
736 gtm.SetValueAsInt(j,0,0xFFF);
742 //_____________________________________________________________________________
744 AliMUONCDB::MakeTriggerLUT(const char* file) const
746 /// Make a triggerlut object, from a file.
748 AliMUONTriggerLut* lut = new AliMUONTriggerLut;
749 lut->ReadFromFile(file);
753 //_____________________________________________________________________________
754 AliMUONTriggerEfficiencyCells*
755 AliMUONCDB::MakeTriggerEfficiency(const char* file) const
757 /// Make a trigger efficiency object from a file.
759 return new AliMUONTriggerEfficiencyCells(file);
762 //_____________________________________________________________________________
764 AliMUONCDB::WriteToCDB(const char* calibpath, TObject* object,
765 Int_t startRun, Int_t endRun,
766 const char* filename)
768 /// Write a given object to OCDB
770 AliCDBId id(calibpath,startRun,endRun);
772 md.SetAliRootVersion(gROOT->GetVersion());
773 md.SetComment(gSystem->ExpandPathName(filename));
774 md.SetResponsible("Uploaded using AliMUONCDB class");
775 AliCDBManager* man = AliCDBManager::Instance();
776 man->SetDefaultStorage(fCDBPath);
777 man->Put(object,id,&md);
780 //_____________________________________________________________________________
782 AliMUONCDB::WriteToCDB(const char* calibpath, TObject* object,
783 Int_t startRun, Int_t endRun, Bool_t defaultValues)
785 /// Write a given object to OCDB
787 AliCDBId id(calibpath,startRun,endRun);
789 md.SetAliRootVersion(gROOT->GetVersion());
792 md.SetComment("Test with default values");
796 md.SetComment("Test with random values");
798 md.SetResponsible("AliMUONCDB tester class");
800 AliCDBManager* man = AliCDBManager::Instance();
801 man->SetDefaultStorage(fCDBPath);
802 man->Put(object,id,&md);
805 //_____________________________________________________________________________
807 AliMUONCDB::MakeNeighbourStore(AliMUONVStore& neighbourStore)
809 /// Fill the neighbours store with, for each channel, a TObjArray of its
810 /// neighbouring pads (including itself)
812 AliInfo("Generating NeighbourStore. This will take a while. Please be patient.");
818 TIter next(ManuList());
826 while ( ( p = (AliMpIntPair*)next() ) )
828 Int_t detElemId = p->GetFirst();
829 Int_t manuId = p->GetSecond();
831 const AliMpVSegmentation* seg =
832 AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId,manuId);
834 AliMUONVCalibParam* calibParam = static_cast<AliMUONVCalibParam*>(neighbourStore.FindObject(detElemId,manuId));
838 Int_t size(AliMpConstants::ManuNofChannels());
839 Int_t defaultValue(-1);
840 Int_t packingFactor(size);
842 calibParam = new AliMUONCalibParamNI(dimension,size,detElemId,manuId,defaultValue,packingFactor);
843 Bool_t ok = neighbourStore.Add(calibParam);
846 AliError(Form("Could not set DetElemId=%d manuId=%d",detElemId,manuId));
851 for ( Int_t manuChannel = 0; manuChannel < AliMpConstants::ManuNofChannels(); ++manuChannel )
853 AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,manuChannel),kFALSE);
859 seg->GetNeighbours(pad,tmp,true,true);
860 Int_t nofPadNeighbours = tmp.GetEntriesFast();
862 for ( Int_t i = 0; i < nofPadNeighbours; ++i )
864 AliMpPad* pad = static_cast<AliMpPad*>(tmp.At(i));
866 Bool_t ok = calibParam->PackValues(pad->GetLocation().GetFirst(),pad->GetLocation().GetSecond(),x);
869 AliError("Could not pack value. Something is seriously wrong. Please check");
870 StdoutToAliError(pad->Print(););
873 calibParam->SetValueAsInt(manuChannel,i,x);
884 //_____________________________________________________________________________
886 AliMUONCDB::SetMaxNofChannelsToGenerate(Int_t n)
888 /// Set the maximum number of channels to generate (used for testing only)
889 /// n < 0 means no limit
890 fMaxNofChannelsToGenerate = n;
893 //_____________________________________________________________________________
895 AliMUONCDB::WriteLocalTriggerMasks(Int_t startRun, Int_t endRun)
897 /// Write local trigger masks to OCDB
899 AliMUONVStore* ltm = new AliMUON1DArray(235);
900 Int_t ngenerated = MakeLocalTriggerMaskStore(*ltm);
901 AliInfo(Form("Ngenerated = %d",ngenerated));
904 WriteToCDB("MUON/Calib/LocalTriggerBoardMasks",ltm,startRun,endRun,true);
909 //_____________________________________________________________________________
911 AliMUONCDB::WriteRegionalTriggerMasks(Int_t startRun, Int_t endRun)
913 /// Write regional trigger masks to OCDB
915 AliMUONVStore* rtm = new AliMUON1DArray(16);
916 Int_t ngenerated = MakeRegionalTriggerMaskStore(*rtm);
917 AliInfo(Form("Ngenerated = %d",ngenerated));
920 WriteToCDB("MUON/Calib/RegionalTriggerBoardMasks",rtm,startRun,endRun,true);
925 //_____________________________________________________________________________
927 AliMUONCDB::WriteGlobalTriggerMasks(Int_t startRun, Int_t endRun)
929 /// Write global trigger masks to OCDB
931 AliMUONVCalibParam* gtm = new AliMUONCalibParamNI(1,16,1,0,0);
933 Int_t ngenerated = MakeGlobalTriggerMaskStore(*gtm);
934 AliInfo(Form("Ngenerated = %d",ngenerated));
937 WriteToCDB("MUON/Calib/GlobalTriggerBoardMasks",gtm,startRun,endRun,true);
942 //_____________________________________________________________________________
944 AliMUONCDB::WriteTriggerLut(Int_t startRun, Int_t endRun)
946 /// Write trigger LUT to OCDB
948 AliMUONTriggerLut* lut = MakeTriggerLUT();
951 WriteToCDB("MUON/Calib/TriggerLut",lut,startRun,endRun,true);
956 //_____________________________________________________________________________
958 AliMUONCDB::WriteTriggerEfficiency(Int_t startRun, Int_t endRun)
960 /// Write trigger efficiency to OCDB
962 AliMUONTriggerEfficiencyCells* eff = MakeTriggerEfficiency();
965 WriteToCDB("MUON/Calib/TriggerEfficiency",eff,startRun,endRun,true);
970 //_____________________________________________________________________________
972 AliMUONCDB::WriteNeighbours(Int_t startRun, Int_t endRun)
974 /// Write neighbours to OCDB
976 AliMUONVStore* neighbours = new AliMUON2DMap(kTRUE);
977 Int_t ngenerated = MakeNeighbourStore(*neighbours);
978 AliInfo(Form("Ngenerated = %d",ngenerated));
981 WriteToCDB("MUON/Calib/Neighbours",neighbours,startRun,endRun,true);
986 //_____________________________________________________________________________
988 AliMUONCDB::WriteHV(Bool_t defaultValues,
989 Int_t startRun, Int_t endRun)
991 /// generate HV values (either cste = 1500 V) if defaultValues=true or random
992 /// if defaultValues=false, see makeHVStore) and
993 /// store them into CDB located at cdbpath, with a validity period
994 /// ranging from startRun to endRun
996 TMap* hvStore = new TMap;
997 Int_t ngenerated = MakeHVStore(*hvStore,defaultValues);
998 AliInfo(Form("Ngenerated = %d",ngenerated));
1001 WriteToCDB("MUON/Calib/HV",hvStore,startRun,endRun,defaultValues);
1006 //_____________________________________________________________________________
1008 AliMUONCDB::WritePedestals(Bool_t defaultValues,
1009 Int_t startRun, Int_t endRun)
1011 /// generate pedestal values (either 0 if defaultValues=true or random
1012 /// if defaultValues=false, see makePedestalStore) and
1013 /// store them into CDB located at cdbpath, with a validity period
1014 /// ranging from startRun to endRun
1016 AliMUONVStore* pedestalStore = new AliMUON2DMap(true);
1017 Int_t ngenerated = MakePedestalStore(*pedestalStore,defaultValues);
1018 AliInfo(Form("Ngenerated = %d",ngenerated));
1019 WriteToCDB("MUON/Calib/Pedestals",pedestalStore,startRun,endRun,defaultValues);
1020 delete pedestalStore;
1024 //_____________________________________________________________________________
1026 AliMUONCDB::WriteGains(Bool_t defaultValues,
1027 Int_t startRun, Int_t endRun)
1029 /// generate gain values (either 1 if defaultValues=true or random
1030 /// if defaultValues=false, see makeGainStore) and
1031 /// store them into CDB located at cdbpath, with a validity period
1032 /// ranging from startRun to endRun
1034 AliMUONVStore* gainStore = new AliMUON2DMap(true);
1035 Int_t ngenerated = MakeGainStore(*gainStore,defaultValues);
1036 AliInfo(Form("Ngenerated = %d",ngenerated));
1037 WriteToCDB("MUON/Calib/Gains",gainStore,startRun,endRun,defaultValues);
1041 //_____________________________________________________________________________
1043 AliMUONCDB::WriteCapacitances(const char* filename,
1044 Int_t startRun, Int_t endRun)
1046 /// read manu capacitance and injection gain values from file
1047 /// and store them into CDB located at cdbpath, with a validity period
1048 /// ranging from startRun to endRun
1050 AliMUONVStore* capaStore = new AliMUON1DMap(16828);
1051 Int_t ngenerated = MakeCapacitanceStore(*capaStore,filename);
1052 AliInfo(Form("Ngenerated = %d",ngenerated));
1053 if ( ngenerated > 0 )
1055 WriteToCDB("MUON/Calib/Capacitances",capaStore,startRun,endRun,filename);
1060 //_____________________________________________________________________________
1062 AliMUONCDB::WriteCapacitances(Bool_t defaultValues,
1063 Int_t startRun, Int_t endRun)
1065 /// generate manu capacitance values (either 1 if defaultValues=true or random
1066 /// if defaultValues=false, see makeCapacitanceStore) and
1067 /// store them into CDB located at cdbpath, with a validity period
1068 /// ranging from startRun to endRun
1070 AliMUONVStore* capaStore = new AliMUON1DMap(16828);
1071 Int_t ngenerated = MakeCapacitanceStore(*capaStore,defaultValues);
1072 AliInfo(Form("Ngenerated = %d",ngenerated));
1073 WriteToCDB("MUON/Calib/Capacitances",capaStore,startRun,endRun,defaultValues);
1077 //_____________________________________________________________________________
1079 AliMUONCDB::WriteTrigger(Int_t startRun, Int_t endRun)
1081 /// Writes all Trigger related calibration to CDB
1082 WriteLocalTriggerMasks(startRun,endRun);
1083 WriteRegionalTriggerMasks(startRun,endRun);
1084 WriteGlobalTriggerMasks(startRun,endRun);
1085 WriteTriggerLut(startRun,endRun);
1086 WriteTriggerEfficiency(startRun,endRun);
1089 //_____________________________________________________________________________
1091 AliMUONCDB::WriteTracker(Bool_t defaultValues, Int_t startRun, Int_t endRun)
1093 /// Writes all Tracker related calibration to CDB
1094 WriteHV(defaultValues,startRun,endRun);
1095 WritePedestals(defaultValues,startRun,endRun);
1096 WriteGains(defaultValues,startRun,endRun);
1097 WriteCapacitances(defaultValues,startRun,endRun);
1098 WriteNeighbours(startRun,endRun);