]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCcalibDB.cxx
Apply pad by pad noise fluctution in the digits simulation (Marian)
[u/mrichter/AliRoot.git] / TPC / AliTPCcalibDB.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16
17 ///////////////////////////////////////////////////////////////////////////////
18 //                                                                           //
19 // Class providing the calibration parameters by accessing the CDB           //
20 //                                                                           //
21 // Request an instance with AliTPCcalibDB::Instance()                        //
22 // If a new event is processed set the event number with SetRun              //
23 // Then request the calibration data                                         ////
24 // Calibration data:
25 // 1.)  pad by pad calibration -  AliTPCCalPad
26 //      a.) fPadGainFactor
27 //          Simulation: AliTPCDigitizer::ExecFast - Multiply by gain
28 //          Reconstruction : AliTPCclustererMI::Digits2Clusters - Divide by gain  
29 //      b.) fPadNoise -
30 //          Simulation:        AliTPCDigitizer::ExecFast
31 //          Reconstruction:    AliTPCclustererMI::FindClusters(AliTPCCalROC * noiseROC)
32 //                             Noise depending cut on clusters (n sigma)
33 //          
34 ///////////////////////////////////////////////////////////////////////////////
35
36 #include <iostream>
37 #include <fstream>
38
39
40 #include <AliCDBManager.h>
41 #include <AliCDBEntry.h>
42 #include <AliLog.h>
43
44 #include "AliTPCcalibDB.h"
45 #include "AliTPCExB.h"
46
47 #include "AliTPCCalROC.h"
48 #include "AliTPCCalPad.h"
49 #include "AliTPCSensorTempArray.h"
50 #include "AliTPCTransform.h"
51 class AliCDBStorage;
52 class AliTPCCalDet;
53 //
54 //
55
56 #include "TFile.h"
57 #include "TKey.h"
58
59 #include "TObjArray.h"
60 #include "TObjString.h"
61 #include "TString.h"
62 #include "AliTPCCalPad.h"
63 #include "AliTPCCalibPulser.h"
64 #include "AliTPCCalibPedestal.h"
65 #include "AliTPCCalibCE.h"
66
67
68
69
70
71 ClassImp(AliTPCcalibDB)
72
73 AliTPCcalibDB* AliTPCcalibDB::fgInstance = 0;
74 Bool_t AliTPCcalibDB::fgTerminated = kFALSE;
75
76
77 //_ singleton implementation __________________________________________________
78 AliTPCcalibDB* AliTPCcalibDB::Instance()
79 {
80   //
81   // Singleton implementation
82   // Returns an instance of this class, it is created if neccessary
83   //
84   
85   if (fgTerminated != kFALSE)
86     return 0;
87
88   if (fgInstance == 0)
89     fgInstance = new AliTPCcalibDB();
90   
91   return fgInstance;
92 }
93
94 void AliTPCcalibDB::Terminate()
95 {
96   //
97   // Singleton implementation
98   // Deletes the instance of this class and sets the terminated flag, instances cannot be requested anymore
99   // This function can be called several times.
100   //
101   
102   fgTerminated = kTRUE;
103   
104   if (fgInstance != 0)
105   {
106     delete fgInstance;
107     fgInstance = 0;
108   }
109 }
110
111 //_____________________________________________________________________________
112 AliTPCcalibDB::AliTPCcalibDB():
113   fRun(-1),
114   fTrafo(0),
115   fExB(0),
116   fPadGainFactor(0),
117   fPadTime0(0),
118   fPadNoise(0),
119   fPedestals(0),
120   fTemperature(0),
121   fParam(0)
122 {
123   //
124   // constructor
125   //  
126   //
127   Update();    // temporary
128 }
129
130 //_____________________________________________________________________________
131 AliTPCcalibDB::~AliTPCcalibDB() 
132 {
133   //
134   // destructor
135   //
136   
137   // don't delete anything, CDB cache is active!
138   //if (fPadGainFactor) delete fPadGainFactor;
139   //if (fPadTime0) delete fPadTime0;
140   //if (fPadNoise) delete fPadNoise;
141 }
142
143
144 //_____________________________________________________________________________
145 AliCDBEntry* AliTPCcalibDB::GetCDBEntry(const char* cdbPath)
146 {
147   // 
148   // Retrieves an entry with path <cdbPath> from the CDB.
149   //
150   char chinfo[1000];
151     
152   AliCDBEntry* entry = AliCDBManager::Instance()->Get(cdbPath, fRun); 
153   if (!entry) 
154   { 
155     sprintf(chinfo,"AliTPCcalibDB: Failed to get entry:\t%s ", cdbPath);
156     AliError(chinfo); 
157     return 0; 
158   }
159   return entry;
160 }
161
162
163 //_____________________________________________________________________________
164 void AliTPCcalibDB::SetRun(Long64_t run)
165 {
166   //
167   // Sets current run number. Calibration data is read from the corresponding file. 
168   //  
169   if (fRun == run)
170     return;  
171   fRun = run;
172   Update();
173 }
174   
175
176
177 void AliTPCcalibDB::Update(){
178   //
179   AliCDBEntry * entry=0;
180   
181   Bool_t cdbCache = AliCDBManager::Instance()->GetCacheFlag(); // save cache status
182   AliCDBManager::Instance()->SetCacheFlag(kTRUE); // activate CDB cache
183   
184   //
185   entry          = GetCDBEntry("TPC/Calib/PadGainFactor");
186   if (entry){
187     //if (fPadGainFactor) delete fPadGainFactor;
188     entry->SetOwner(kTRUE);
189     fPadGainFactor = (AliTPCCalPad*)entry->GetObject();
190   }
191   //
192   entry          = GetCDBEntry("TPC/Calib/PadTime0");
193   if (entry){
194     //if (fPadTime0) delete fPadTime0;
195     entry->SetOwner(kTRUE);
196     fPadTime0 = (AliTPCCalPad*)entry->GetObject();
197   }
198   //
199   //
200   entry          = GetCDBEntry("TPC/Calib/PadNoise");
201   if (entry){
202     //if (fPadNoise) delete fPadNoise;
203     entry->SetOwner(kTRUE);
204     fPadNoise = (AliTPCCalPad*)entry->GetObject();
205   }
206
207   entry          = GetCDBEntry("TPC/Calib/Pedestals");
208   if (entry){
209     //if (fPedestals) delete fPedestals;
210     entry->SetOwner(kTRUE);
211     fPedestals = (AliTPCCalPad*)entry->GetObject();
212   }
213
214   entry          = GetCDBEntry("TPC/Calib/Temperature");
215   if (entry){
216     //if (fTemperature) delete fTemperature;
217     entry->SetOwner(kTRUE);
218     fTemperature = (AliTPCSensorTempArray*)entry->GetObject();
219   }
220
221   entry          = GetCDBEntry("TPC/Calib/Parameters");
222   if (entry){
223     //if (fPadNoise) delete fPadNoise;
224     entry->SetOwner(kTRUE);
225     fParam = (AliTPCParam*)(entry->GetObject()->Clone());
226   }
227
228   entry          = GetCDBEntry("TPC/Calib/ExB");
229   if (entry) {
230     entry->SetOwner(kTRUE);
231     fExB=dynamic_cast<AliTPCExB*>(entry->GetObject()->Clone());
232   }
233
234   delete fTrafo;
235   fTrafo=new AliTPCTransform(); 
236
237   //
238   AliCDBManager::Instance()->SetCacheFlag(cdbCache); // reset original CDB cache
239   
240 }
241 AliTPCcalibDB::AliTPCcalibDB(const AliTPCcalibDB& org)
242 {
243   //
244   // Copy constructor invalid -- singleton implementation
245   //
246    Error("copy constructor","invalid -- singleton implementation");
247 }
248
249 AliTPCcalibDB& AliTPCcalibDB::operator= (const AliTPCcalibDB& rhs)
250 {
251 //
252 // Singleton implementation - no assignment operator
253 //
254   Error("operator =", "assignment operator not implemented");
255   return *this;
256 }
257
258
259
260 void AliTPCcalibDB::CreateObjectList(const Char_t *filename, TObjArray *calibObjects)
261 {
262 //
263 // Create calibration objects and read contents from OCDB
264 //
265    if ( calibObjects == 0x0 ) return;
266    ifstream in;
267    in.open(filename);
268    if ( !in.is_open() ){
269       fprintf(stderr,"Error: cannot open list file '%s'", filename);
270       return;
271    }
272    
273    AliTPCCalPad *calPad=0x0;
274    
275    TString sFile;
276    sFile.ReadFile(in);
277    in.close();
278    
279    TObjArray *arrFileLine = sFile.Tokenize("\n");
280    
281    TIter nextLine(arrFileLine);
282    
283    TObjString *sObjLine=0x0;
284    while ( sObjLine = (TObjString*)nextLine() ){
285       TString sLine(sObjLine->GetString());
286       
287       TObjArray *arrNextCol = sLine.Tokenize("\t");
288       
289       TObjString *sObjType     = (TObjString*)(arrNextCol->At(0));
290       TObjString *sObjFileName = (TObjString*)(arrNextCol->At(1));
291       
292       if ( !sObjType || ! sObjFileName ) continue;
293       TString sType(sObjType->GetString());
294       TString sFileName(sObjFileName->GetString());
295       printf("%s\t%s\n",sType.Data(),sFileName.Data());
296       
297       TFile *fIn = TFile::Open(sFileName);
298       if ( !fIn ){
299          fprintf(stderr,"File not found: '%s'", sFileName.Data());
300          continue;
301       }
302       
303       if ( sType == "CE" ){
304          AliTPCCalibCE *ce = (AliTPCCalibCE*)fIn->Get("AliTPCCalibCE");
305          
306          calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadT0());         
307          calPad->SetNameTitle("CETmean","CETmean");
308          calibObjects->Add(calPad);
309          
310          calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadQ());         
311          calPad->SetNameTitle("CEQmean","CEQmean");
312          calibObjects->Add(calPad);        
313          
314          calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadRMS());
315          calPad->SetNameTitle("CETrms","CETrms");
316          calibObjects->Add(calPad);         
317                   
318       } else if ( sType == "Pulser") {
319          AliTPCCalibPulser *sig = (AliTPCCalibPulser*)fIn->Get("AliTPCCalibPulser");
320          
321          calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadT0());         
322          calPad->SetNameTitle("PulserTmean","PulserTmean");
323          calibObjects->Add(calPad);
324          
325          calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadQ());         
326          calPad->SetNameTitle("PulserQmean","PulserQmean");
327          calibObjects->Add(calPad);        
328          
329          calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadRMS());
330          calPad->SetNameTitle("PulserTrms","PulserTrms");
331          calibObjects->Add(calPad);         
332       
333       } else if ( sType == "Pedestals") {
334          AliTPCCalibPedestal *ped = (AliTPCCalibPedestal*)fIn->Get("AliTPCCalibPedestal");
335          
336          calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadPedestal());         
337          calPad->SetNameTitle("Pedestals","Pedestals");
338          calibObjects->Add(calPad);
339          
340          calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadRMS());         
341          calPad->SetNameTitle("Noise","Noise");
342          calibObjects->Add(calPad);        
343      
344       } else {
345          fprintf(stderr,"Undefined Type: '%s'",sType.Data());
346          
347       }
348       delete fIn;
349    }
350 }
351
352
353
354 void AliTPCcalibDB::MakeTree(const char * fileName, TObjArray * array, const char * mapFileName, AliTPCCalPad* outlierPad, Float_t ltmFraction) {
355   //
356   // Write a tree with all available information
357   // if mapFileName is specified, the Map information are also written to the tree
358   // pads specified in outlierPad are not used for calculating statistics
359   //  - the same function as AliTPCCalPad::MakeTree - 
360   //
361    AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
362
363    TObjArray* mapIROCs = 0;
364    TObjArray* mapOROCs = 0;
365    TVectorF *mapIROCArray = 0;
366    TVectorF *mapOROCArray = 0;
367    Int_t mapEntries = 0;
368    TString* mapNames = 0;
369    
370    if (mapFileName) {
371       TFile mapFile(mapFileName, "read");
372       
373       TList* listOfROCs = mapFile.GetListOfKeys();
374       mapEntries = listOfROCs->GetEntries()/2;
375       mapIROCs = new TObjArray(mapEntries*2);
376       mapOROCs = new TObjArray(mapEntries*2);
377       mapIROCArray = new TVectorF[mapEntries];
378       mapOROCArray = new TVectorF[mapEntries];
379       
380       mapNames = new TString[mapEntries];
381       for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
382         TString nameROC(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
383          nameROC.Remove(nameROC.Length()-4, 4);
384          mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "IROC").Data()), ivalue);
385          mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "OROC").Data()), ivalue);
386          mapNames[ivalue].Append(nameROC);
387       }
388       
389       for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
390          mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
391          mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
392       
393          for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
394             (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
395          for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
396             (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
397       }
398
399    } //  if (mapFileName)
400   
401    TTreeSRedirector cstream(fileName);
402    Int_t arrayEntries = array->GetEntries();
403    
404    TString* names = new TString[arrayEntries];
405    for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
406       names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
407
408    for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
409       //
410       // get statistic for given sector
411       //
412       TVectorF median(arrayEntries);
413       TVectorF mean(arrayEntries);
414       TVectorF rms(arrayEntries);
415       TVectorF ltm(arrayEntries);
416       TVectorF ltmrms(arrayEntries);
417       TVectorF medianWithOut(arrayEntries);
418       TVectorF meanWithOut(arrayEntries);
419       TVectorF rmsWithOut(arrayEntries);
420       TVectorF ltmWithOut(arrayEntries);
421       TVectorF ltmrmsWithOut(arrayEntries);
422       
423       TVectorF *vectorArray = new TVectorF[arrayEntries];
424       for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
425          vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
426       
427       for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
428          AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
429          AliTPCCalROC* calROC = calPad->GetCalROC(isector);
430          AliTPCCalROC* outlierROC = 0;
431          if (outlierPad) outlierROC = outlierPad->GetCalROC(isector);
432          if (calROC) {
433             median[ivalue] = calROC->GetMedian();
434             mean[ivalue] = calROC->GetMean();
435             rms[ivalue] = calROC->GetRMS();
436             Double_t ltmrmsValue = 0;
437             ltm[ivalue] = calROC->GetLTM(&ltmrmsValue, ltmFraction);
438             ltmrms[ivalue] = ltmrmsValue;
439             if (outlierROC) {
440                medianWithOut[ivalue] = calROC->GetMedian(outlierROC);
441                meanWithOut[ivalue] = calROC->GetMean(outlierROC);
442                rmsWithOut[ivalue] = calROC->GetRMS(outlierROC);
443                ltmrmsValue = 0;
444                ltmWithOut[ivalue] = calROC->GetLTM(&ltmrmsValue, ltmFraction, outlierROC);
445                ltmrmsWithOut[ivalue] = ltmrmsValue;
446             }
447          }
448          else {
449             median[ivalue] = 0.;
450             mean[ivalue] = 0.;
451             rms[ivalue] = 0.;
452             ltm[ivalue] = 0.;
453             ltmrms[ivalue] = 0.;
454             medianWithOut[ivalue] = 0.;
455             meanWithOut[ivalue] = 0.;
456             rmsWithOut[ivalue] = 0.;
457             ltmWithOut[ivalue] = 0.;
458             ltmrmsWithOut[ivalue] = 0.;
459          }
460       }
461       
462       //
463       // fill vectors of variable per pad
464       //
465       TVectorF *posArray = new TVectorF[8];
466       for (Int_t ivalue = 0; ivalue < 8; ivalue++)
467          posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
468
469       Float_t posG[3] = {0};
470       Float_t posL[3] = {0};
471       Int_t ichannel = 0;
472       for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
473          for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
474             tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
475             tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
476             posArray[0][ichannel] = irow;
477             posArray[1][ichannel] = ipad;
478             posArray[2][ichannel] = posL[0];
479             posArray[3][ichannel] = posL[1];
480             posArray[4][ichannel] = posG[0];
481             posArray[5][ichannel] = posG[1];
482             posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
483             posArray[7][ichannel] = ichannel;
484             
485             // loop over array containing AliTPCCalPads
486             for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
487                AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
488                AliTPCCalROC* calROC = calPad->GetCalROC(isector);
489                if (calROC)
490                   (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
491                else
492                   (vectorArray[ivalue])[ichannel] = 0;
493             }
494             ichannel++;
495          }
496       }
497       
498       cstream << "calPads" <<
499          "sector=" << isector;
500       
501       for  (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
502          cstream << "calPads" <<
503             (Char_t*)((names[ivalue] + "_Median=").Data()) << median[ivalue] <<
504             (Char_t*)((names[ivalue] + "_Mean=").Data()) << mean[ivalue] <<
505             (Char_t*)((names[ivalue] + "_RMS=").Data()) << rms[ivalue] <<
506             (Char_t*)((names[ivalue] + "_LTM=").Data()) << ltm[ivalue] <<
507             (Char_t*)((names[ivalue] + "_RMS_LTM=").Data()) << ltmrms[ivalue];
508          if (outlierPad) {
509             cstream << "calPads" <<
510                (Char_t*)((names[ivalue] + "_Median_OutlierCutted=").Data()) << medianWithOut[ivalue] <<
511                (Char_t*)((names[ivalue] + "_Mean_OutlierCutted=").Data()) << meanWithOut[ivalue] <<
512                (Char_t*)((names[ivalue] + "_RMS_OutlierCutted=").Data()) << rmsWithOut[ivalue] <<
513                (Char_t*)((names[ivalue] + "_LTM_OutlierCutted=").Data()) << ltmWithOut[ivalue] <<
514                (Char_t*)((names[ivalue] + "_RMS_LTM_OutlierCutted=").Data()) << ltmrmsWithOut[ivalue];
515          }
516       }
517
518       for  (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
519          cstream << "calPads" <<
520             (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
521       }
522
523       if (mapFileName) {
524          for  (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
525             if (isector < 36)
526                cstream << "calPads" <<
527                   (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
528             else
529                cstream << "calPads" <<
530                   (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
531          }
532       }
533
534       cstream << "calPads" <<
535          "row.=" << &posArray[0] <<
536          "pad.=" << &posArray[1] <<
537          "lx.=" << &posArray[2] <<
538          "ly.=" << &posArray[3] <<
539          "gx.=" << &posArray[4] <<
540          "gy.=" << &posArray[5] <<
541          "rpad.=" << &posArray[6] <<
542          "channel.=" << &posArray[7];
543          
544       cstream << "calPads" <<
545          "\n";
546
547       delete[] posArray;
548       delete[] vectorArray;
549    }
550    
551
552    delete[] names;
553    if (mapFileName) {
554       delete mapIROCs;
555       delete mapOROCs;
556       delete[] mapIROCArray;
557       delete[] mapOROCArray;
558       delete[] mapNames;
559    }
560 }
561