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