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