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