]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - TPC/AliTPCcalibDB.cxx
Fix compilation warnings
[u/mrichter/AliRoot.git] / TPC / AliTPCcalibDB.cxx
... / ...
CommitLineData
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// Calibration data:
27// 0.) Altro mapping
28// Simulation - not yet
29// Reconstruction - AliTPCclustererMI::Digits2Clusters(AliRawReader* rawReader)
30//
31// 1.) pad by pad calibration - AliTPCCalPad
32//
33// a.) fPadGainFactor
34// Simulation: AliTPCDigitizer::ExecFast - Multiply by gain
35// Reconstruction : AliTPCclustererMI::Digits2Clusters - Divide by gain
36//
37// b.) fPadNoise -
38// Simulation: AliTPCDigitizer::ExecFast
39// Reconstruction: AliTPCclustererMI::FindClusters(AliTPCCalROC * noiseROC)
40// Noise depending cut on clusters charge (n sigma)
41// c.) fPedestal:
42// Simulation: Not used yet - To be impleneted - Rounding to the nearest integer
43// Reconstruction: Used in AliTPCclustererMI::Digits2Clusters(AliRawReader* rawReader)
44// if data taken without zero suppression
45// Currently switch in fRecoParam->GetCalcPedestal();
46//
47// d.) fPadTime0
48// Simulation: applied in the AliTPC::MakeSector - adding offset
49// Reconstruction: AliTPCTransform::Transform() - remove offset
50// AliTPCTransform::Transform() - to be called
51// in AliTPCtracker::Transform()
52//
53//
54// 2.) Space points transformation:
55//
56// a.) General coordinate tranformation - AliTPCtransform (see $ALICE_ROOT/TPC/AliTPCtransform.cxx)
57// Created on fly - use the other calibration components
58// Unisochronity - (substract time0 - pad by pad)
59// Drift velocity - Currently common drift velocity - functionality of AliTPCParam
60// ExB effect
61// Simulation - Not used directly (the effects are applied one by one (see AliTPC::MakeSector)
62// Reconstruction -
63// AliTPCclustererMI::AddCluster
64// AliTPCtrackerMI::Transform
65// b.) ExB effect calibration -
66// classes (base class AliTPCExB, implementation- AliTPCExBExact.h AliTPCExBFirst.h)
67// a.a) Simulation: applied in the AliTPC::MakeSector -
68// calib->GetExB()->CorrectInverse(dxyz0,dxyz1);
69// a.b) Reconstruction -
70//
71// in AliTPCtransform::Correct() - called calib->GetExB()->Correct(dxyz0,dxyz1)
72//
73// 3.) cluster error, shape and Q parameterization
74//
75//
76//
77///////////////////////////////////////////////////////////////////////////////
78
79#include <iostream>
80#include <fstream>
81
82
83#include <AliCDBManager.h>
84#include <AliCDBEntry.h>
85#include <AliLog.h>
86#include <AliMagF.h>
87#include <AliMagWrapCheb.h>
88
89#include "AliTPCcalibDB.h"
90#include "AliTPCAltroMapping.h"
91#include "AliTPCExB.h"
92
93#include "AliTPCCalROC.h"
94#include "AliTPCCalPad.h"
95#include "AliTPCSensorTempArray.h"
96#include "AliGRPObject.h"
97#include "AliTPCTransform.h"
98
99class AliCDBStorage;
100class AliTPCCalDet;
101//
102//
103
104#include "TFile.h"
105#include "TKey.h"
106
107#include "TObjArray.h"
108#include "TObjString.h"
109#include "TString.h"
110#include "AliTPCCalPad.h"
111#include "AliTPCCalibPulser.h"
112#include "AliTPCCalibPedestal.h"
113#include "AliTPCCalibCE.h"
114#include "AliTPCExBFirst.h"
115#include "AliTPCTempMap.h"
116
117
118
119
120ClassImp(AliTPCcalibDB)
121
122AliTPCcalibDB* AliTPCcalibDB::fgInstance = 0;
123Bool_t AliTPCcalibDB::fgTerminated = kFALSE;
124TObjArray AliTPCcalibDB::fgExBArray; // array of ExB corrections
125
126
127//_ singleton implementation __________________________________________________
128AliTPCcalibDB* AliTPCcalibDB::Instance()
129{
130 //
131 // Singleton implementation
132 // Returns an instance of this class, it is created if neccessary
133 //
134
135 if (fgTerminated != kFALSE)
136 return 0;
137
138 if (fgInstance == 0)
139 fgInstance = new AliTPCcalibDB();
140
141 return fgInstance;
142}
143
144void AliTPCcalibDB::Terminate()
145{
146 //
147 // Singleton implementation
148 // Deletes the instance of this class and sets the terminated flag, instances cannot be requested anymore
149 // This function can be called several times.
150 //
151
152 fgTerminated = kTRUE;
153
154 if (fgInstance != 0)
155 {
156 delete fgInstance;
157 fgInstance = 0;
158 }
159}
160
161//_____________________________________________________________________________
162AliTPCcalibDB::AliTPCcalibDB():
163 TObject(),
164 fRun(-1),
165 fTransform(0),
166 fExB(0),
167 fPadGainFactor(0),
168 fDedxGainFactor(0),
169 fPadTime0(0),
170 fPadNoise(0),
171 fPedestals(0),
172 fTemperature(0),
173 fMapping(0),
174 fParam(0),
175 fClusterParam(0),
176 fGRPArray(100000), //! array of GRPs - per run - JUST for calibration studies
177 fGoofieArray(100000), //! array of GOOFIE values -per run - Just for calibration studies
178 fTemperatureArray(100000), //! array of temperature sensors - per run - Just for calibration studies
179 fRunList(100000) //! run list - indicates try to get the run param
180
181{
182 //
183 // constructor
184 //
185 //
186 Update(); // temporary
187}
188
189AliTPCcalibDB::AliTPCcalibDB(const AliTPCcalibDB& ):
190 TObject(),
191 fRun(-1),
192 fTransform(0),
193 fExB(0),
194 fPadGainFactor(0),
195 fDedxGainFactor(0),
196 fPadTime0(0),
197 fPadNoise(0),
198 fPedestals(0),
199 fTemperature(0),
200 fMapping(0),
201 fParam(0),
202 fClusterParam(0),
203 fGRPArray(0), //! array of GRPs - per run - JUST for calibration studies
204 fGoofieArray(0), //! array of GOOFIE values -per run - Just for calibration studies
205 fTemperatureArray(0), //! array of temperature sensors - per run - Just for calibration studies
206 fRunList(0) //! run list - indicates try to get the run param
207{
208 //
209 // Copy constructor invalid -- singleton implementation
210 //
211 Error("copy constructor","invalid -- singleton implementation");
212}
213
214AliTPCcalibDB& AliTPCcalibDB::operator= (const AliTPCcalibDB& )
215{
216//
217// Singleton implementation - no assignment operator
218//
219 Error("operator =", "assignment operator not implemented");
220 return *this;
221}
222
223
224
225//_____________________________________________________________________________
226AliTPCcalibDB::~AliTPCcalibDB()
227{
228 //
229 // destructor
230 //
231
232 // don't delete anything, CDB cache is active!
233 //if (fPadGainFactor) delete fPadGainFactor;
234 //if (fPadTime0) delete fPadTime0;
235 //if (fPadNoise) delete fPadNoise;
236}
237
238
239//_____________________________________________________________________________
240AliCDBEntry* AliTPCcalibDB::GetCDBEntry(const char* cdbPath)
241{
242 //
243 // Retrieves an entry with path <cdbPath> from the CDB.
244 //
245 char chinfo[1000];
246
247 AliCDBEntry* entry = AliCDBManager::Instance()->Get(cdbPath, fRun);
248 if (!entry)
249 {
250 sprintf(chinfo,"AliTPCcalibDB: Failed to get entry:\t%s ", cdbPath);
251 AliError(chinfo);
252 return 0;
253 }
254 return entry;
255}
256
257
258//_____________________________________________________________________________
259void AliTPCcalibDB::SetRun(Long64_t run)
260{
261 //
262 // Sets current run number. Calibration data is read from the corresponding file.
263 //
264 if (fRun == run)
265 return;
266 fRun = run;
267 Update();
268}
269
270
271
272void AliTPCcalibDB::Update(){
273 //
274 AliCDBEntry * entry=0;
275
276 Bool_t cdbCache = AliCDBManager::Instance()->GetCacheFlag(); // save cache status
277 AliCDBManager::Instance()->SetCacheFlag(kTRUE); // activate CDB cache
278
279 //
280 entry = GetCDBEntry("TPC/Calib/PadGainFactor");
281 if (entry){
282 //if (fPadGainFactor) delete fPadGainFactor;
283 entry->SetOwner(kTRUE);
284 fPadGainFactor = (AliTPCCalPad*)entry->GetObject();
285 }
286 //
287 entry = GetCDBEntry("TPC/Calib/GainFactorDedx");
288 if (entry){
289 entry->SetOwner(kTRUE);
290 fDedxGainFactor = (AliTPCCalPad*)entry->GetObject();
291 }
292 //
293 entry = GetCDBEntry("TPC/Calib/PadTime0");
294 if (entry){
295 //if (fPadTime0) delete fPadTime0;
296 entry->SetOwner(kTRUE);
297 fPadTime0 = (AliTPCCalPad*)entry->GetObject();
298 }
299 //
300 //
301 entry = GetCDBEntry("TPC/Calib/PadNoise");
302 if (entry){
303 //if (fPadNoise) delete fPadNoise;
304 entry->SetOwner(kTRUE);
305 fPadNoise = (AliTPCCalPad*)entry->GetObject();
306 }
307
308 entry = GetCDBEntry("TPC/Calib/Pedestals");
309 if (entry){
310 //if (fPedestals) delete fPedestals;
311 entry->SetOwner(kTRUE);
312 fPedestals = (AliTPCCalPad*)entry->GetObject();
313 }
314
315 entry = GetCDBEntry("TPC/Calib/Temperature");
316 if (entry){
317 //if (fTemperature) delete fTemperature;
318 entry->SetOwner(kTRUE);
319 fTemperature = (AliTPCSensorTempArray*)entry->GetObject();
320 }
321
322 entry = GetCDBEntry("TPC/Calib/Parameters");
323 if (entry){
324 //if (fPadNoise) delete fPadNoise;
325 entry->SetOwner(kTRUE);
326 fParam = (AliTPCParam*)(entry->GetObject()->Clone());
327 }
328
329 entry = GetCDBEntry("TPC/Calib/ClusterParam");
330 if (entry){
331 //if (fPadNoise) delete fPadNoise;
332 entry->SetOwner(kTRUE);
333 fClusterParam = (AliTPCClusterParam*)(entry->GetObject()->Clone());
334 }
335
336 entry = GetCDBEntry("TPC/Calib/Mapping");
337 if (entry){
338 //if (fPadNoise) delete fPadNoise;
339 entry->SetOwner(kTRUE);
340 TObjArray * array = dynamic_cast<TObjArray*>(entry->GetObject());
341 if (array && array->GetEntriesFast()==6){
342 fMapping = new AliTPCAltroMapping*[6];
343 for (Int_t i=0; i<6; i++){
344 fMapping[i] = dynamic_cast<AliTPCAltroMapping*>(array->At(i));
345 }
346 }
347 }
348
349
350
351 //entry = GetCDBEntry("TPC/Calib/ExB");
352 //if (entry) {
353 // entry->SetOwner(kTRUE);
354 // fExB=dynamic_cast<AliTPCExB*>(entry->GetObject()->Clone());
355 //}
356 //
357 // ExB - calculate during initialization
358 // -
359 fExB = GetExB(-5,kTRUE);
360 //
361 if (!fTransform) {
362 fTransform=new AliTPCTransform();
363 }
364
365 //
366 AliCDBManager::Instance()->SetCacheFlag(cdbCache); // reset original CDB cache
367
368}
369
370
371
372void AliTPCcalibDB::CreateObjectList(const Char_t *filename, TObjArray *calibObjects)
373{
374//
375// Create calibration objects and read contents from OCDB
376//
377 if ( calibObjects == 0x0 ) return;
378 ifstream in;
379 in.open(filename);
380 if ( !in.is_open() ){
381 fprintf(stderr,"Error: cannot open list file '%s'", filename);
382 return;
383 }
384
385 AliTPCCalPad *calPad=0x0;
386
387 TString sFile;
388 sFile.ReadFile(in);
389 in.close();
390
391 TObjArray *arrFileLine = sFile.Tokenize("\n");
392
393 TIter nextLine(arrFileLine);
394
395 TObjString *sObjLine=0x0;
396 while ( (sObjLine = (TObjString*)nextLine()) ){
397 TString sLine(sObjLine->GetString());
398
399 TObjArray *arrNextCol = sLine.Tokenize("\t");
400
401 TObjString *sObjType = (TObjString*)(arrNextCol->At(0));
402 TObjString *sObjFileName = (TObjString*)(arrNextCol->At(1));
403
404 if ( !sObjType || ! sObjFileName ) continue;
405 TString sType(sObjType->GetString());
406 TString sFileName(sObjFileName->GetString());
407 printf("%s\t%s\n",sType.Data(),sFileName.Data());
408
409 TFile *fIn = TFile::Open(sFileName);
410 if ( !fIn ){
411 fprintf(stderr,"File not found: '%s'", sFileName.Data());
412 continue;
413 }
414
415 if ( sType == "CE" ){
416 AliTPCCalibCE *ce = (AliTPCCalibCE*)fIn->Get("AliTPCCalibCE");
417
418 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadT0());
419 calPad->SetNameTitle("CETmean","CETmean");
420 calibObjects->Add(calPad);
421
422 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadQ());
423 calPad->SetNameTitle("CEQmean","CEQmean");
424 calibObjects->Add(calPad);
425
426 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadRMS());
427 calPad->SetNameTitle("CETrms","CETrms");
428 calibObjects->Add(calPad);
429
430 } else if ( sType == "Pulser") {
431 AliTPCCalibPulser *sig = (AliTPCCalibPulser*)fIn->Get("AliTPCCalibPulser");
432
433 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadT0());
434 calPad->SetNameTitle("PulserTmean","PulserTmean");
435 calibObjects->Add(calPad);
436
437 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadQ());
438 calPad->SetNameTitle("PulserQmean","PulserQmean");
439 calibObjects->Add(calPad);
440
441 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadRMS());
442 calPad->SetNameTitle("PulserTrms","PulserTrms");
443 calibObjects->Add(calPad);
444
445 } else if ( sType == "Pedestals") {
446 AliTPCCalibPedestal *ped = (AliTPCCalibPedestal*)fIn->Get("AliTPCCalibPedestal");
447
448 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadPedestal());
449 calPad->SetNameTitle("Pedestals","Pedestals");
450 calibObjects->Add(calPad);
451
452 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadRMS());
453 calPad->SetNameTitle("Noise","Noise");
454 calibObjects->Add(calPad);
455
456 } else {
457 fprintf(stderr,"Undefined Type: '%s'",sType.Data());
458
459 }
460 delete fIn;
461 }
462}
463
464
465
466void AliTPCcalibDB::MakeTree(const char * fileName, TObjArray * array, const char * mapFileName, AliTPCCalPad* outlierPad, Float_t ltmFraction) {
467 //
468 // Write a tree with all available information
469 // if mapFileName is specified, the Map information are also written to the tree
470 // pads specified in outlierPad are not used for calculating statistics
471 // - the same function as AliTPCCalPad::MakeTree -
472 //
473 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
474
475 TObjArray* mapIROCs = 0;
476 TObjArray* mapOROCs = 0;
477 TVectorF *mapIROCArray = 0;
478 TVectorF *mapOROCArray = 0;
479 Int_t mapEntries = 0;
480 TString* mapNames = 0;
481
482 if (mapFileName) {
483 TFile mapFile(mapFileName, "read");
484
485 TList* listOfROCs = mapFile.GetListOfKeys();
486 mapEntries = listOfROCs->GetEntries()/2;
487 mapIROCs = new TObjArray(mapEntries*2);
488 mapOROCs = new TObjArray(mapEntries*2);
489 mapIROCArray = new TVectorF[mapEntries];
490 mapOROCArray = new TVectorF[mapEntries];
491
492 mapNames = new TString[mapEntries];
493 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
494 TString nameROC(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
495 nameROC.Remove(nameROC.Length()-4, 4);
496 mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "IROC").Data()), ivalue);
497 mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "OROC").Data()), ivalue);
498 mapNames[ivalue].Append(nameROC);
499 }
500
501 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
502 mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
503 mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
504
505 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
506 (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
507 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
508 (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
509 }
510
511 } // if (mapFileName)
512
513 TTreeSRedirector cstream(fileName);
514 Int_t arrayEntries = array->GetEntries();
515
516 TString* names = new TString[arrayEntries];
517 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
518 names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
519
520 for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
521 //
522 // get statistic for given sector
523 //
524 TVectorF median(arrayEntries);
525 TVectorF mean(arrayEntries);
526 TVectorF rms(arrayEntries);
527 TVectorF ltm(arrayEntries);
528 TVectorF ltmrms(arrayEntries);
529 TVectorF medianWithOut(arrayEntries);
530 TVectorF meanWithOut(arrayEntries);
531 TVectorF rmsWithOut(arrayEntries);
532 TVectorF ltmWithOut(arrayEntries);
533 TVectorF ltmrmsWithOut(arrayEntries);
534
535 TVectorF *vectorArray = new TVectorF[arrayEntries];
536 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
537 vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
538
539 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
540 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
541 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
542 AliTPCCalROC* outlierROC = 0;
543 if (outlierPad) outlierROC = outlierPad->GetCalROC(isector);
544 if (calROC) {
545 median[ivalue] = calROC->GetMedian();
546 mean[ivalue] = calROC->GetMean();
547 rms[ivalue] = calROC->GetRMS();
548 Double_t ltmrmsValue = 0;
549 ltm[ivalue] = calROC->GetLTM(&ltmrmsValue, ltmFraction);
550 ltmrms[ivalue] = ltmrmsValue;
551 if (outlierROC) {
552 medianWithOut[ivalue] = calROC->GetMedian(outlierROC);
553 meanWithOut[ivalue] = calROC->GetMean(outlierROC);
554 rmsWithOut[ivalue] = calROC->GetRMS(outlierROC);
555 ltmrmsValue = 0;
556 ltmWithOut[ivalue] = calROC->GetLTM(&ltmrmsValue, ltmFraction, outlierROC);
557 ltmrmsWithOut[ivalue] = ltmrmsValue;
558 }
559 }
560 else {
561 median[ivalue] = 0.;
562 mean[ivalue] = 0.;
563 rms[ivalue] = 0.;
564 ltm[ivalue] = 0.;
565 ltmrms[ivalue] = 0.;
566 medianWithOut[ivalue] = 0.;
567 meanWithOut[ivalue] = 0.;
568 rmsWithOut[ivalue] = 0.;
569 ltmWithOut[ivalue] = 0.;
570 ltmrmsWithOut[ivalue] = 0.;
571 }
572 }
573
574 //
575 // fill vectors of variable per pad
576 //
577 TVectorF *posArray = new TVectorF[8];
578 for (Int_t ivalue = 0; ivalue < 8; ivalue++)
579 posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
580
581 Float_t posG[3] = {0};
582 Float_t posL[3] = {0};
583 Int_t ichannel = 0;
584 for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
585 for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
586 tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
587 tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
588 posArray[0][ichannel] = irow;
589 posArray[1][ichannel] = ipad;
590 posArray[2][ichannel] = posL[0];
591 posArray[3][ichannel] = posL[1];
592 posArray[4][ichannel] = posG[0];
593 posArray[5][ichannel] = posG[1];
594 posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
595 posArray[7][ichannel] = ichannel;
596
597 // loop over array containing AliTPCCalPads
598 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
599 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
600 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
601 if (calROC)
602 (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
603 else
604 (vectorArray[ivalue])[ichannel] = 0;
605 }
606 ichannel++;
607 }
608 }
609
610 cstream << "calPads" <<
611 "sector=" << isector;
612
613 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
614 cstream << "calPads" <<
615 (Char_t*)((names[ivalue] + "_Median=").Data()) << median[ivalue] <<
616 (Char_t*)((names[ivalue] + "_Mean=").Data()) << mean[ivalue] <<
617 (Char_t*)((names[ivalue] + "_RMS=").Data()) << rms[ivalue] <<
618 (Char_t*)((names[ivalue] + "_LTM=").Data()) << ltm[ivalue] <<
619 (Char_t*)((names[ivalue] + "_RMS_LTM=").Data()) << ltmrms[ivalue];
620 if (outlierPad) {
621 cstream << "calPads" <<
622 (Char_t*)((names[ivalue] + "_Median_OutlierCutted=").Data()) << medianWithOut[ivalue] <<
623 (Char_t*)((names[ivalue] + "_Mean_OutlierCutted=").Data()) << meanWithOut[ivalue] <<
624 (Char_t*)((names[ivalue] + "_RMS_OutlierCutted=").Data()) << rmsWithOut[ivalue] <<
625 (Char_t*)((names[ivalue] + "_LTM_OutlierCutted=").Data()) << ltmWithOut[ivalue] <<
626 (Char_t*)((names[ivalue] + "_RMS_LTM_OutlierCutted=").Data()) << ltmrmsWithOut[ivalue];
627 }
628 }
629
630 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
631 cstream << "calPads" <<
632 (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
633 }
634
635 if (mapFileName) {
636 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
637 if (isector < 36)
638 cstream << "calPads" <<
639 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
640 else
641 cstream << "calPads" <<
642 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
643 }
644 }
645
646 cstream << "calPads" <<
647 "row.=" << &posArray[0] <<
648 "pad.=" << &posArray[1] <<
649 "lx.=" << &posArray[2] <<
650 "ly.=" << &posArray[3] <<
651 "gx.=" << &posArray[4] <<
652 "gy.=" << &posArray[5] <<
653 "rpad.=" << &posArray[6] <<
654 "channel.=" << &posArray[7];
655
656 cstream << "calPads" <<
657 "\n";
658
659 delete[] posArray;
660 delete[] vectorArray;
661 }
662
663
664 delete[] names;
665 if (mapFileName) {
666 delete mapIROCs;
667 delete mapOROCs;
668 delete[] mapIROCArray;
669 delete[] mapOROCArray;
670 delete[] mapNames;
671 }
672}
673
674
675
676void AliTPCcalibDB::RegisterExB(Int_t index, Float_t bz, Bool_t bdelete){
677 //
678 // Register static ExB correction map
679 // index - registration index - used for visualization
680 // bz - bz field in kGaus
681
682 Float_t factor = bz/(-5.); // default b filed in Cheb with minus sign
683
684 AliMagF* bmap = new AliMagWrapCheb("Maps","Maps", 2, factor, 10., AliMagWrapCheb::k5kG,kTRUE,"$(ALICE_ROOT)/data/maps/mfchebKGI_sym.root");
685
686 AliTPCExBFirst *exb = new AliTPCExBFirst(bmap,0.88*2.6400e+04,50,50,50);
687 AliTPCExB::SetInstance(exb);
688
689 if (bdelete){
690 delete bmap;
691 }else{
692 AliTPCExB::RegisterField(index,bmap);
693 }
694 if (index>=fgExBArray.GetEntries()) fgExBArray.Expand((index+1)*2+11);
695 fgExBArray.AddAt(exb,index);
696}
697
698
699AliTPCExB* AliTPCcalibDB::GetExB(Float_t bz, Bool_t deleteB) {
700 //
701 // bz filed in KGaus not in tesla
702 // Get ExB correction map
703 // if doesn't exist - create it
704 //
705 Int_t index = TMath::Nint(5+bz);
706 if (index>fgExBArray.GetEntries()) fgExBArray.Expand((index+1)*2+11);
707 if (!fgExBArray.At(index)) AliTPCcalibDB::RegisterExB(index,bz,deleteB);
708 return (AliTPCExB*)fgExBArray.At(index);
709}
710
711
712void AliTPCcalibDB::SetExBField(Float_t bz){
713 //
714 // Set magnetic filed for ExB correction
715 //
716 printf("Set magnetic field for ExB correction = %f\n",bz);
717 fExB = GetExB(bz,kFALSE);
718}
719
720
721
722void AliTPCcalibDB::GetRunInformations( Int_t run){
723 //
724 // - > Don't use it for reconstruction - Only for Calibration studies
725 //
726 AliCDBEntry * entry = 0;
727 if (run>= fRunList.GetSize()){
728 fRunList.Set(run*2+1);
729 fGRPArray.Expand(run*2+1);fGoofieArray.Expand(run*2+1); fTemperatureArray.Expand(run*2+1);
730 }
731 if (fRunList[run]>0) return;
732 entry = AliCDBManager::Instance()->Get("GRP/GRP/Data",run);
733 if (entry) fGRPArray.AddAt(entry->GetObject(),run);
734 entry = AliCDBManager::Instance()->Get("TPC/Calib/Goofie",run);
735 if (entry) fGoofieArray.AddAt(entry->GetObject(),run);
736 entry = AliCDBManager::Instance()->Get("TPC/Calib/Temperature",run);
737 if (entry) fTemperatureArray.AddAt(entry->GetObject(),run);
738 fRunList[run]=1; // sign as used
739}
740
741
742Float_t AliTPCcalibDB::GetGain(Int_t sector, Int_t row, Int_t pad){
743 //
744 //
745 AliTPCCalPad *calPad = Instance()->fDedxGainFactor;;
746 if (!calPad) return 0;
747 return calPad->GetCalROC(sector)->GetValue(row,pad);
748}
749
750AliDCSSensor * AliTPCcalibDB::GetPressureSensor(Int_t run){
751 //
752 //
753 AliGRPObject * grpRun = dynamic_cast<AliGRPObject *>(fGRPArray.At(run));
754 if (!grpRun) {
755 GetRunInformations(run);
756 grpRun = dynamic_cast<AliGRPObject *>(fGRPArray.At(run));
757 if (!grpRun) return 0;
758 }
759 AliDCSSensor * sensor = grpRun->GetCavernAtmosPressure();
760 return sensor;
761
762}
763
764AliTPCSensorTempArray * AliTPCcalibDB::GetTemperatureSensor(Int_t run){
765 //
766 // Get temperature sensor array
767 //
768 AliTPCSensorTempArray * tempArray = (AliTPCSensorTempArray *)fTemperatureArray.At(run);
769 if (!tempArray) {
770 GetRunInformations(run);
771 tempArray = (AliTPCSensorTempArray *)fTemperatureArray.At(run);
772 }
773 return tempArray;
774}
775
776AliDCSSensorArray * AliTPCcalibDB::GetGoofieSensors(Int_t run){
777 //
778 // Get temperature sensor array
779 //
780 AliDCSSensorArray * goofieArray = (AliDCSSensorArray *)fGoofieArray.At(run);
781 if (!goofieArray) {
782 GetRunInformations(run);
783 goofieArray = (AliDCSSensorArray *)fGoofieArray.At(run);
784 }
785 return goofieArray;
786}
787
788
789
790
791Float_t AliTPCcalibDB::GetPressure(Int_t timeStamp, Int_t run){
792 //
793 // GetPressure for given time stamp and runt
794 //
795 TTimeStamp stamp(timeStamp);
796 AliDCSSensor * sensor = Instance()->GetPressureSensor(run);
797 if (!sensor) return 0;
798 if (!sensor->GetFit()) return 0;
799 return sensor->GetValue(stamp);
800}
801
802Bool_t AliTPCcalibDB::GetTemperatureFit(Int_t timeStamp, Int_t run, Int_t side,TVectorD& fit){
803 //
804 //
805 //
806 TTimeStamp tstamp(timeStamp);
807 AliTPCSensorTempArray* tempArray = Instance()->GetTemperatureSensor(run);
808 if (! tempArray) return kFALSE;
809 AliTPCTempMap * tempMap = new AliTPCTempMap(tempArray);
810 TLinearFitter * fitter = tempMap->GetLinearFitter(3,side,tstamp);
811 if (fitter){
812 fitter->Eval();
813 fitter->GetParameters(fit);
814 }
815 delete fitter;
816 delete tempMap;
817 if (!fitter) return kFALSE;
818 return kTRUE;
819}
820
821Float_t AliTPCcalibDB::GetTemperature(Int_t timeStamp, Int_t run, Int_t side){
822 //
823 //
824 //
825 TVectorD vec;
826 if (side==0) {
827 GetTemperatureFit(timeStamp,run,0,vec);
828 return vec[0];
829 }
830 if (side==1){
831 GetTemperatureFit(timeStamp,run,0,vec);
832 return vec[0];
833 }
834}
835
836
837
838void AliTPCcalibDB::ProcessEnv(const char * runList){
839 //
840 // Example test function - how to use the environment variables
841 // runList - ascii file with run numbers
842 // output - dcsTime.root file with tree
843
844 ifstream in;
845 in.open(runList);
846 Int_t irun=0;
847 TTreeSRedirector *pcstream = new TTreeSRedirector("dcsTime.root");
848 while(in.good()) {
849 in >> irun;
850 if (irun==0) continue;
851 printf("Processing run %d\n",irun);
852 AliDCSSensor * sensorPressure = AliTPCcalibDB::Instance()->GetPressureSensor(irun);
853 if (!sensorPressure) continue;
854 AliTPCSensorTempArray * tempArray = AliTPCcalibDB::Instance()->GetTemperatureSensor(irun);
855 AliTPCTempMap * tempMap = new AliTPCTempMap(tempArray);
856 AliDCSSensorArray* goofieArray = AliTPCcalibDB::Instance()->GetGoofieSensors(irun);
857 //
858 Int_t startTime = sensorPressure->GetStartTime();
859 Int_t endTime = sensorPressure->GetEndTime();
860 Int_t dtime = TMath::Max((endTime-startTime)/20,10*60);
861 for (Int_t itime=startTime; itime<endTime; itime+=dtime){
862 //
863 TTimeStamp tstamp(itime);
864 Float_t valuePressure = sensorPressure->GetValue(tstamp);
865
866 TLinearFitter * fitter = 0;
867 TVectorD vecTemp[10];
868 if (itime<tempArray->GetStartTime().GetSec() || itime>tempArray->GetEndTime().GetSec()){
869 }else{
870 for (Int_t itype=0; itype<5; itype++)
871 for (Int_t iside=0; iside<2; iside++){
872 fitter= tempMap->GetLinearFitter(itype,iside,tstamp);
873 if (!fitter) continue;
874 fitter->Eval(); fitter->GetParameters(vecTemp[itype+iside*5]);
875 delete fitter;
876 }
877 }
878
879 TVectorD vecGoofie, vecEntries, vecMean, vecMedian,vecRMS;
880 if (goofieArray){
881 vecGoofie.ResizeTo(goofieArray->NumSensors());
882 ProcessGoofie(goofieArray, vecEntries ,vecMedian, vecMean, vecRMS);
883 //
884 for (Int_t isensor=0; isensor<goofieArray->NumSensors();isensor++){
885 AliDCSSensor *gsensor = goofieArray->GetSensor(isensor);
886 if (gsensor){
887 vecGoofie[isensor] = gsensor->GetValue(tstamp);
888 }
889 }
890 }
891
892
893 //tempMap->GetLinearFitter(0,0,itime);
894 (*pcstream)<<"dcs"<<
895 "run="<<irun<<
896 "time="<<itime<<
897 "goofie.="<<&vecGoofie<<
898 "goofieE.="<<&vecEntries<<
899 "goofieMean.="<<&vecMean<<
900 "goofieMedian.="<<&vecMedian<<
901 "goofieRMS.="<<&vecRMS<<
902 "press="<<valuePressure<<
903 "temp00.="<<&vecTemp[0]<<
904 "temp10.="<<&vecTemp[1]<<
905 "temp20.="<<&vecTemp[2]<<
906 "temp30.="<<&vecTemp[3]<<
907 "temp40.="<<&vecTemp[4]<<
908 "temp01.="<<&vecTemp[5]<<
909 "temp11.="<<&vecTemp[6]<<
910 "temp21.="<<&vecTemp[7]<<
911 "temp31.="<<&vecTemp[8]<<
912 "temp41.="<<&vecTemp[9]<<
913 "\n";
914 }
915 }
916 delete pcstream;
917}
918
919
920void AliTPCcalibDB::ProcessGoofie( AliDCSSensorArray* goofieArray, TVectorD & vecEntries, TVectorD & vecMedian, TVectorD &vecMean, TVectorD &vecRMS){
921 /*
922
923 1 TPC_ANODE_I_A00_STAT
924 2 TPC_DVM_CO2
925 3 TPC_DVM_DriftVelocity
926 4 TPC_DVM_FCageHV
927 5 TPC_DVM_GainFar
928 6 TPC_DVM_GainNear
929 7 TPC_DVM_N2
930 8 TPC_DVM_NumberOfSparks
931 9 TPC_DVM_PeakAreaFar
932 10 TPC_DVM_PeakAreaNear
933 11 TPC_DVM_PeakPosFar
934 12 TPC_DVM_PeakPosNear
935 13 TPC_DVM_PickupHV
936 14 TPC_DVM_Pressure
937 15 TPC_DVM_T1_Over_P
938 16 TPC_DVM_T2_Over_P
939 17 TPC_DVM_T_Over_P
940 18 TPC_DVM_TemperatureS1
941 */
942 //
943 //
944 // TVectorD vecMedian; TVectorD vecEntries; TVectorD vecMean; TVectorD vecRMS;
945 Double_t kEpsilon=0.0000000001;
946 Double_t kBig=100000000000.;
947 Int_t nsensors = goofieArray->NumSensors();
948 vecEntries.ResizeTo(nsensors);
949 vecMedian.ResizeTo(nsensors);
950 vecMean.ResizeTo(nsensors);
951 vecRMS.ResizeTo(nsensors);
952 TVectorF values;
953 for (Int_t isensor=0; isensor<goofieArray->NumSensors();isensor++){
954 AliDCSSensor *gsensor = goofieArray->GetSensor(isensor);
955 if (gsensor && gsensor->GetGraph()){
956 Int_t npoints = gsensor->GetGraph()->GetN();
957 // filter zeroes
958 values.ResizeTo(npoints);
959 Int_t nused =0;
960 for (Int_t ipoint=0; ipoint<npoints; ipoint++){
961 if (TMath::Abs(gsensor->GetGraph()->GetY()[ipoint])>kEpsilon &&
962 TMath::Abs(gsensor->GetGraph()->GetY()[ipoint])<kBig ){
963 values[nused]=gsensor->GetGraph()->GetY()[ipoint];
964 nused++;
965 }
966 }
967 //
968 vecEntries[isensor]= nused;
969 if (nused>1){
970 vecMedian[isensor] = TMath::Median(nused,values.GetMatrixArray());
971 vecMean[isensor] = TMath::Mean(nused,values.GetMatrixArray());
972 vecRMS[isensor] = TMath::RMS(nused,values.GetMatrixArray());
973 }
974 }
975 }
976}
977