1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 ///////////////////////////////////////////////////////////////////////////////
19 // Class providing the calibration parameters by accessing the CDB //
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 ////
28 // Simulation - not yet
29 // Reconstruction - AliTPCclustererMI::Digits2Clusters(AliRawReader* rawReader)
31 // 1.) pad by pad calibration - AliTPCCalPad
34 // Simulation: AliTPCDigitizer::ExecFast - Multiply by gain
35 // Reconstruction : AliTPCclustererMI::Digits2Clusters - Divide by gain
38 // Simulation: AliTPCDigitizer::ExecFast
39 // Reconstruction: AliTPCclustererMI::FindClusters(AliTPCCalROC * noiseROC)
40 // Noise depending cut on clusters charge (n sigma)
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();
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()
54 // 2.) Space points transformation:
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
61 // Simulation - Not used directly (the effects are applied one by one (see AliTPC::MakeSector)
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 -
71 // in AliTPCtransform::Correct() - called calib->GetExB()->Correct(dxyz0,dxyz1)
73 // 3.) cluster error, shape and Q parameterization
77 ///////////////////////////////////////////////////////////////////////////////
83 #include <AliCDBManager.h>
84 #include <AliCDBEntry.h>
87 #include <AliSplineFit.h>
89 #include "AliTPCcalibDB.h"
90 #include "AliTPCAltroMapping.h"
91 #include "AliTPCExB.h"
93 #include "AliTPCCalROC.h"
94 #include "AliTPCCalPad.h"
95 #include "AliTPCSensorTempArray.h"
96 #include "AliGRPObject.h"
97 #include "AliTPCTransform.h"
107 #include "TObjArray.h"
108 #include "TObjString.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 #include "AliTPCCalibVdrift.h"
117 #include "AliTPCCalibRaw.h"
119 #include "AliTPCPreprocessorOnline.h"
122 ClassImp(AliTPCcalibDB)
124 AliTPCcalibDB* AliTPCcalibDB::fgInstance = 0;
125 Bool_t AliTPCcalibDB::fgTerminated = kFALSE;
126 TObjArray AliTPCcalibDB::fgExBArray; // array of ExB corrections
129 //_ singleton implementation __________________________________________________
130 AliTPCcalibDB* AliTPCcalibDB::Instance()
133 // Singleton implementation
134 // Returns an instance of this class, it is created if neccessary
137 if (fgTerminated != kFALSE)
141 fgInstance = new AliTPCcalibDB();
146 void AliTPCcalibDB::Terminate()
149 // Singleton implementation
150 // Deletes the instance of this class and sets the terminated flag, instances cannot be requested anymore
151 // This function can be called several times.
154 fgTerminated = kTRUE;
163 //_____________________________________________________________________________
164 AliTPCcalibDB::AliTPCcalibDB():
183 fTimeGainSplinesArray(100000),
184 fGRPArray(100000), //! array of GRPs - per run - JUST for calibration studies
185 fGRPMaps(100000), //! array of GRPs - per run - JUST for calibration studies
186 fGoofieArray(100000), //! array of GOOFIE values -per run - Just for calibration studies
187 fVoltageArray(100000),
188 fTemperatureArray(100000), //! array of temperature sensors - per run - Just for calibration studies
189 fVdriftArray(100000), //! array of v drift interfaces
190 fDriftCorrectionArray(100000), //! array of drift correction
191 fRunList(100000) //! run list - indicates try to get the run param
198 Update(); // temporary
201 AliTPCcalibDB::AliTPCcalibDB(const AliTPCcalibDB& ):
220 fTimeGainSplinesArray(100000),
221 fGRPArray(0), //! array of GRPs - per run - JUST for calibration studies
222 fGRPMaps(0), //! array of GRPs - per run - JUST for calibration studies
223 fGoofieArray(0), //! array of GOOFIE values -per run - Just for calibration studies
225 fTemperatureArray(0), //! array of temperature sensors - per run - Just for calibration studies
226 fVdriftArray(0), //! array of v drift interfaces
227 fDriftCorrectionArray(0), //! array of v drift interfaces
228 fRunList(0) //! run list - indicates try to get the run param
231 // Copy constructor invalid -- singleton implementation
233 Error("copy constructor","invalid -- singleton implementation");
236 AliTPCcalibDB& AliTPCcalibDB::operator= (const AliTPCcalibDB& )
239 // Singleton implementation - no assignment operator
241 Error("operator =", "assignment operator not implemented");
247 //_____________________________________________________________________________
248 AliTPCcalibDB::~AliTPCcalibDB()
254 // don't delete anything, CDB cache is active!
255 //if (fPadGainFactor) delete fPadGainFactor;
256 //if (fPadTime0) delete fPadTime0;
257 //if (fPadNoise) delete fPadNoise;
261 //_____________________________________________________________________________
262 AliCDBEntry* AliTPCcalibDB::GetCDBEntry(const char* cdbPath)
265 // Retrieves an entry with path <cdbPath> from the CDB.
269 AliCDBEntry* entry = AliCDBManager::Instance()->Get(cdbPath, fRun);
272 sprintf(chinfo,"AliTPCcalibDB: Failed to get entry:\t%s ", cdbPath);
280 //_____________________________________________________________________________
281 void AliTPCcalibDB::SetRun(Long64_t run)
284 // Sets current run number. Calibration data is read from the corresponding file.
294 void AliTPCcalibDB::Update(){
296 AliCDBEntry * entry=0;
298 Bool_t cdbCache = AliCDBManager::Instance()->GetCacheFlag(); // save cache status
299 AliCDBManager::Instance()->SetCacheFlag(kTRUE); // activate CDB cache
302 entry = GetCDBEntry("TPC/Calib/PadGainFactor");
304 //if (fPadGainFactor) delete fPadGainFactor;
305 entry->SetOwner(kTRUE);
306 fPadGainFactor = (AliTPCCalPad*)entry->GetObject();
309 entry = GetCDBEntry("TPC/Calib/TimeGain");
311 //if (fTimeGainSplines) delete fTimeGainSplines;
312 entry->SetOwner(kTRUE);
313 fTimeGainSplines = (TObjArray*)entry->GetObject();
316 entry = GetCDBEntry("TPC/Calib/GainFactorDedx");
318 entry->SetOwner(kTRUE);
319 fDedxGainFactor = (AliTPCCalPad*)entry->GetObject();
322 entry = GetCDBEntry("TPC/Calib/PadTime0");
324 //if (fPadTime0) delete fPadTime0;
325 entry->SetOwner(kTRUE);
326 fPadTime0 = (AliTPCCalPad*)entry->GetObject();
330 entry = GetCDBEntry("TPC/Calib/PadNoise");
332 //if (fPadNoise) delete fPadNoise;
333 entry->SetOwner(kTRUE);
334 fPadNoise = (AliTPCCalPad*)entry->GetObject();
337 entry = GetCDBEntry("TPC/Calib/Pedestals");
339 //if (fPedestals) delete fPedestals;
340 entry->SetOwner(kTRUE);
341 fPedestals = (AliTPCCalPad*)entry->GetObject();
344 entry = GetCDBEntry("TPC/Calib/Temperature");
346 //if (fTemperature) delete fTemperature;
347 entry->SetOwner(kTRUE);
348 fTemperature = (AliTPCSensorTempArray*)entry->GetObject();
351 entry = GetCDBEntry("TPC/Calib/Parameters");
353 //if (fPadNoise) delete fPadNoise;
354 entry->SetOwner(kTRUE);
355 fParam = (AliTPCParam*)(entry->GetObject()->Clone());
358 entry = GetCDBEntry("TPC/Calib/ClusterParam");
360 entry->SetOwner(kTRUE);
361 fClusterParam = (AliTPCClusterParam*)(entry->GetObject()->Clone());
364 //ALTRO configuration data
365 entry = GetCDBEntry("TPC/Calib/AltroConfig");
367 entry->SetOwner(kTRUE);
368 fALTROConfigData=(TObjArray*)(entry->GetObject());
371 //Calibration Pulser data
372 entry = GetCDBEntry("TPC/Calib/Pulser");
374 entry->SetOwner(kTRUE);
375 fPulserData=(TObjArray*)(entry->GetObject());
379 entry = GetCDBEntry("TPC/Calib/CE");
381 entry->SetOwner(kTRUE);
382 fCEData=(TObjArray*)(entry->GetObject());
384 //RAW calibration data
385 entry = GetCDBEntry("TPC/Calib/Raw");
387 entry->SetOwner(kTRUE);
388 TObjArray *arr=(TObjArray*)(entry->GetObject());
389 if (arr) fCalibRaw=(AliTPCCalibRaw*)arr->At(0);
392 entry = GetCDBEntry("TPC/Calib/Mapping");
394 //if (fPadNoise) delete fPadNoise;
395 entry->SetOwner(kTRUE);
396 TObjArray * array = dynamic_cast<TObjArray*>(entry->GetObject());
397 if (array && array->GetEntriesFast()==6){
398 fMapping = new AliTPCAltroMapping*[6];
399 for (Int_t i=0; i<6; i++){
400 fMapping[i] = dynamic_cast<AliTPCAltroMapping*>(array->At(i));
407 //entry = GetCDBEntry("TPC/Calib/ExB");
409 // entry->SetOwner(kTRUE);
410 // fExB=dynamic_cast<AliTPCExB*>(entry->GetObject()->Clone());
413 // ExB - calculate during initialization
415 fExB = GetExB(-5,kTRUE);
418 fTransform=new AliTPCTransform();
419 fTransform->SetCurrentRun(AliCDBManager::Instance()->GetRun());
423 AliCDBManager::Instance()->SetCacheFlag(cdbCache); // reset original CDB cache
429 void AliTPCcalibDB::CreateObjectList(const Char_t *filename, TObjArray *calibObjects)
432 // Create calibration objects and read contents from OCDB
434 if ( calibObjects == 0x0 ) return;
437 if ( !in.is_open() ){
438 fprintf(stderr,"Error: cannot open list file '%s'", filename);
442 AliTPCCalPad *calPad=0x0;
448 TObjArray *arrFileLine = sFile.Tokenize("\n");
450 TIter nextLine(arrFileLine);
452 TObjString *sObjLine=0x0;
453 while ( (sObjLine = (TObjString*)nextLine()) ){
454 TString sLine(sObjLine->GetString());
456 TObjArray *arrNextCol = sLine.Tokenize("\t");
458 TObjString *sObjType = (TObjString*)(arrNextCol->At(0));
459 TObjString *sObjFileName = (TObjString*)(arrNextCol->At(1));
461 if ( !sObjType || ! sObjFileName ) continue;
462 TString sType(sObjType->GetString());
463 TString sFileName(sObjFileName->GetString());
464 printf("%s\t%s\n",sType.Data(),sFileName.Data());
466 TFile *fIn = TFile::Open(sFileName);
468 fprintf(stderr,"File not found: '%s'", sFileName.Data());
472 if ( sType == "CE" ){
473 AliTPCCalibCE *ce = (AliTPCCalibCE*)fIn->Get("AliTPCCalibCE");
475 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadT0());
476 calPad->SetNameTitle("CETmean","CETmean");
477 calibObjects->Add(calPad);
479 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadQ());
480 calPad->SetNameTitle("CEQmean","CEQmean");
481 calibObjects->Add(calPad);
483 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadRMS());
484 calPad->SetNameTitle("CETrms","CETrms");
485 calibObjects->Add(calPad);
487 } else if ( sType == "Pulser") {
488 AliTPCCalibPulser *sig = (AliTPCCalibPulser*)fIn->Get("AliTPCCalibPulser");
490 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadT0());
491 calPad->SetNameTitle("PulserTmean","PulserTmean");
492 calibObjects->Add(calPad);
494 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadQ());
495 calPad->SetNameTitle("PulserQmean","PulserQmean");
496 calibObjects->Add(calPad);
498 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadRMS());
499 calPad->SetNameTitle("PulserTrms","PulserTrms");
500 calibObjects->Add(calPad);
502 } else if ( sType == "Pedestals") {
503 AliTPCCalibPedestal *ped = (AliTPCCalibPedestal*)fIn->Get("AliTPCCalibPedestal");
505 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadPedestal());
506 calPad->SetNameTitle("Pedestals","Pedestals");
507 calibObjects->Add(calPad);
509 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadRMS());
510 calPad->SetNameTitle("Noise","Noise");
511 calibObjects->Add(calPad);
514 fprintf(stderr,"Undefined Type: '%s'",sType.Data());
523 void AliTPCcalibDB::MakeTree(const char * fileName, TObjArray * array, const char * mapFileName, AliTPCCalPad* outlierPad, Float_t ltmFraction) {
525 // Write a tree with all available information
526 // if mapFileName is specified, the Map information are also written to the tree
527 // pads specified in outlierPad are not used for calculating statistics
528 // - the same function as AliTPCCalPad::MakeTree -
530 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
532 TObjArray* mapIROCs = 0;
533 TObjArray* mapOROCs = 0;
534 TVectorF *mapIROCArray = 0;
535 TVectorF *mapOROCArray = 0;
536 Int_t mapEntries = 0;
537 TString* mapNames = 0;
540 TFile mapFile(mapFileName, "read");
542 TList* listOfROCs = mapFile.GetListOfKeys();
543 mapEntries = listOfROCs->GetEntries()/2;
544 mapIROCs = new TObjArray(mapEntries*2);
545 mapOROCs = new TObjArray(mapEntries*2);
546 mapIROCArray = new TVectorF[mapEntries];
547 mapOROCArray = new TVectorF[mapEntries];
549 mapNames = new TString[mapEntries];
550 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
551 TString nameROC(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
552 nameROC.Remove(nameROC.Length()-4, 4);
553 mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "IROC").Data()), ivalue);
554 mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "OROC").Data()), ivalue);
555 mapNames[ivalue].Append(nameROC);
558 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
559 mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
560 mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
562 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
563 (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
564 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
565 (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
568 } // if (mapFileName)
570 TTreeSRedirector cstream(fileName);
571 Int_t arrayEntries = array->GetEntries();
573 TString* names = new TString[arrayEntries];
574 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
575 names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
577 for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
579 // get statistic for given sector
581 TVectorF median(arrayEntries);
582 TVectorF mean(arrayEntries);
583 TVectorF rms(arrayEntries);
584 TVectorF ltm(arrayEntries);
585 TVectorF ltmrms(arrayEntries);
586 TVectorF medianWithOut(arrayEntries);
587 TVectorF meanWithOut(arrayEntries);
588 TVectorF rmsWithOut(arrayEntries);
589 TVectorF ltmWithOut(arrayEntries);
590 TVectorF ltmrmsWithOut(arrayEntries);
592 TVectorF *vectorArray = new TVectorF[arrayEntries];
593 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
594 vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
596 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
597 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
598 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
599 AliTPCCalROC* outlierROC = 0;
600 if (outlierPad) outlierROC = outlierPad->GetCalROC(isector);
602 median[ivalue] = calROC->GetMedian();
603 mean[ivalue] = calROC->GetMean();
604 rms[ivalue] = calROC->GetRMS();
605 Double_t ltmrmsValue = 0;
606 ltm[ivalue] = calROC->GetLTM(<mrmsValue, ltmFraction);
607 ltmrms[ivalue] = ltmrmsValue;
609 medianWithOut[ivalue] = calROC->GetMedian(outlierROC);
610 meanWithOut[ivalue] = calROC->GetMean(outlierROC);
611 rmsWithOut[ivalue] = calROC->GetRMS(outlierROC);
613 ltmWithOut[ivalue] = calROC->GetLTM(<mrmsValue, ltmFraction, outlierROC);
614 ltmrmsWithOut[ivalue] = ltmrmsValue;
623 medianWithOut[ivalue] = 0.;
624 meanWithOut[ivalue] = 0.;
625 rmsWithOut[ivalue] = 0.;
626 ltmWithOut[ivalue] = 0.;
627 ltmrmsWithOut[ivalue] = 0.;
632 // fill vectors of variable per pad
634 TVectorF *posArray = new TVectorF[8];
635 for (Int_t ivalue = 0; ivalue < 8; ivalue++)
636 posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
638 Float_t posG[3] = {0};
639 Float_t posL[3] = {0};
641 for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
642 for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
643 tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
644 tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
645 posArray[0][ichannel] = irow;
646 posArray[1][ichannel] = ipad;
647 posArray[2][ichannel] = posL[0];
648 posArray[3][ichannel] = posL[1];
649 posArray[4][ichannel] = posG[0];
650 posArray[5][ichannel] = posG[1];
651 posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
652 posArray[7][ichannel] = ichannel;
654 // loop over array containing AliTPCCalPads
655 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
656 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
657 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
659 (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
661 (vectorArray[ivalue])[ichannel] = 0;
667 cstream << "calPads" <<
668 "sector=" << isector;
670 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
671 cstream << "calPads" <<
672 (Char_t*)((names[ivalue] + "_Median=").Data()) << median[ivalue] <<
673 (Char_t*)((names[ivalue] + "_Mean=").Data()) << mean[ivalue] <<
674 (Char_t*)((names[ivalue] + "_RMS=").Data()) << rms[ivalue] <<
675 (Char_t*)((names[ivalue] + "_LTM=").Data()) << ltm[ivalue] <<
676 (Char_t*)((names[ivalue] + "_RMS_LTM=").Data()) << ltmrms[ivalue];
678 cstream << "calPads" <<
679 (Char_t*)((names[ivalue] + "_Median_OutlierCutted=").Data()) << medianWithOut[ivalue] <<
680 (Char_t*)((names[ivalue] + "_Mean_OutlierCutted=").Data()) << meanWithOut[ivalue] <<
681 (Char_t*)((names[ivalue] + "_RMS_OutlierCutted=").Data()) << rmsWithOut[ivalue] <<
682 (Char_t*)((names[ivalue] + "_LTM_OutlierCutted=").Data()) << ltmWithOut[ivalue] <<
683 (Char_t*)((names[ivalue] + "_RMS_LTM_OutlierCutted=").Data()) << ltmrmsWithOut[ivalue];
687 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
688 cstream << "calPads" <<
689 (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
693 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
695 cstream << "calPads" <<
696 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
698 cstream << "calPads" <<
699 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
703 cstream << "calPads" <<
704 "row.=" << &posArray[0] <<
705 "pad.=" << &posArray[1] <<
706 "lx.=" << &posArray[2] <<
707 "ly.=" << &posArray[3] <<
708 "gx.=" << &posArray[4] <<
709 "gy.=" << &posArray[5] <<
710 "rpad.=" << &posArray[6] <<
711 "channel.=" << &posArray[7];
713 cstream << "calPads" <<
717 delete[] vectorArray;
725 delete[] mapIROCArray;
726 delete[] mapOROCArray;
733 void AliTPCcalibDB::RegisterExB(Int_t index, Float_t bz, Bool_t bdelete){
735 // Register static ExB correction map
736 // index - registration index - used for visualization
737 // bz - bz field in kGaus
739 Float_t factor = bz/(-5.); // default b filed in Cheb with minus sign
741 AliMagF* bmap = new AliMagF("MapsExB","MapsExB", 2,factor,1., 10.,AliMagF::k5kG,"$(ALICE_ROOT)/data/maps/mfchebKGI_sym.root");
743 AliTPCExBFirst *exb = new AliTPCExBFirst(bmap,0.88*2.6400e+04,50,50,50);
744 AliTPCExB::SetInstance(exb);
749 AliTPCExB::RegisterField(index,bmap);
751 if (index>=fgExBArray.GetEntries()) fgExBArray.Expand((index+1)*2+11);
752 fgExBArray.AddAt(exb,index);
756 AliTPCExB* AliTPCcalibDB::GetExB(Float_t bz, Bool_t deleteB) {
758 // bz filed in KGaus not in tesla
759 // Get ExB correction map
760 // if doesn't exist - create it
762 Int_t index = TMath::Nint(5+bz);
763 if (index>fgExBArray.GetEntries()) fgExBArray.Expand((index+1)*2+11);
764 if (!fgExBArray.At(index)) AliTPCcalibDB::RegisterExB(index,bz,deleteB);
765 return (AliTPCExB*)fgExBArray.At(index);
769 void AliTPCcalibDB::SetExBField(Float_t bz){
771 // Set magnetic filed for ExB correction
773 fExB = GetExB(bz,kFALSE);
780 void AliTPCcalibDB::UpdateRunInformations( Int_t run, Bool_t force){
782 // - > Don't use it for reconstruction - Only for Calibration studies
784 AliCDBEntry * entry = 0;
785 if (run>= fRunList.GetSize()){
786 fRunList.Set(run*2+1);
787 fGRPArray.Expand(run*2+1);
788 fGRPMaps.Expand(run*2+1);
789 fGoofieArray.Expand(run*2+1);
790 fVoltageArray.Expand(run*2+1);
791 fTemperatureArray.Expand(run*2+1);
792 fVdriftArray.Expand(run*2+1);
793 fDriftCorrectionArray.Expand(run*2+1);
794 fTimeGainSplinesArray.Expand(run*2+1);
796 if (fRunList[run]>0 &&force==kFALSE) return;
798 entry = AliCDBManager::Instance()->Get("GRP/GRP/Data",run);
800 AliGRPObject * grpRun = dynamic_cast<AliGRPObject*>(entry->GetObject());
802 TMap* map = dynamic_cast<TMap*>(entry->GetObject());
804 //grpRun = new AliGRPObject;
805 //grpRun->ReadValuesFromMap(map);
806 grpRun = MakeGRPObjectFromMap(map);
808 fGRPMaps.AddAt(map,run);
811 fGRPArray.AddAt(grpRun,run);
813 entry = AliCDBManager::Instance()->Get("TPC/Calib/Goofie",run);
815 fGoofieArray.AddAt(entry->GetObject(),run);
818 entry = AliCDBManager::Instance()->Get("TPC/Calib/HighVoltage",run);
820 fVoltageArray.AddAt(entry->GetObject(),run);
823 entry = AliCDBManager::Instance()->Get("TPC/Calib/TimeGain",run);
825 fTimeGainSplinesArray.AddAt(entry->GetObject(),run);
828 entry = AliCDBManager::Instance()->Get("TPC/Calib/TimeDrift",run);
830 fDriftCorrectionArray.AddAt(entry->GetObject(),run);
833 entry = AliCDBManager::Instance()->Get("TPC/Calib/Temperature",run);
835 fTemperatureArray.AddAt(entry->GetObject(),run);
837 fRunList[run]=1; // sign as used
839 AliDCSSensor * press = GetPressureSensor(run,0);
840 AliTPCSensorTempArray * temp = GetTemperatureSensor(run);
842 AliTPCCalibVdrift * vdrift = new AliTPCCalibVdrift(temp, press,0);
843 fVdriftArray.AddAt(vdrift,run);
848 Float_t AliTPCcalibDB::GetGain(Int_t sector, Int_t row, Int_t pad){
851 AliTPCCalPad *calPad = Instance()->fDedxGainFactor;;
852 if (!calPad) return 0;
853 return calPad->GetCalROC(sector)->GetValue(row,pad);
857 AliGRPObject *AliTPCcalibDB::GetGRP(Int_t run){
859 // Get GRP object for given run
861 AliGRPObject * grpRun = dynamic_cast<AliGRPObject *>((Instance()->fGRPArray).At(run));
863 Instance()->UpdateRunInformations(run);
864 grpRun = dynamic_cast<AliGRPObject *>(Instance()->fGRPArray.At(run));
865 if (!grpRun) return 0;
870 TMap * AliTPCcalibDB::GetGRPMap(Int_t run){
874 TMap * grpRun = dynamic_cast<TMap *>((Instance()->fGRPMaps).At(run));
876 Instance()->UpdateRunInformations(run);
877 grpRun = dynamic_cast<TMap *>(Instance()->fGRPMaps.At(run));
878 if (!grpRun) return 0;
884 AliDCSSensor * AliTPCcalibDB::GetPressureSensor(Int_t run, Int_t type){
886 // Get Pressure sensor
888 // type = 0 - Cavern pressure
889 // 1 - Suface pressure
890 // First try to get if trom map - if existing (Old format of data storing)
894 TMap *map = GetGRPMap(run);
896 AliDCSSensor * sensor = 0;
898 if (type==0) osensor = ((*map)("fCavernPressure"));
899 if (type==1) osensor = ((*map)("fP2Pressure"));
900 sensor =dynamic_cast<AliDCSSensor *>(osensor);
901 if (sensor) return sensor;
904 // If not map try to get it from the GRPObject
906 AliGRPObject * grpRun = dynamic_cast<AliGRPObject *>(fGRPArray.At(run));
908 UpdateRunInformations(run);
909 grpRun = dynamic_cast<AliGRPObject *>(fGRPArray.At(run));
910 if (!grpRun) return 0;
912 AliDCSSensor * sensor = grpRun->GetCavernAtmosPressure();
913 if (type==1) sensor = grpRun->GetSurfaceAtmosPressure();
917 AliTPCSensorTempArray * AliTPCcalibDB::GetTemperatureSensor(Int_t run){
919 // Get temperature sensor array
921 AliTPCSensorTempArray * tempArray = (AliTPCSensorTempArray *)fTemperatureArray.At(run);
923 UpdateRunInformations(run);
924 tempArray = (AliTPCSensorTempArray *)fTemperatureArray.At(run);
930 TObjArray * AliTPCcalibDB::GetTimeGainSplinesRun(Int_t run){
932 // Get temperature sensor array
934 TObjArray * gainSplines = (TObjArray *)fTimeGainSplinesArray.At(run);
936 UpdateRunInformations(run);
937 gainSplines = (TObjArray *)fTimeGainSplinesArray.At(run);
942 AliDCSSensorArray * AliTPCcalibDB::GetVoltageSensors(Int_t run){
944 // Get temperature sensor array
946 AliDCSSensorArray * voltageArray = (AliDCSSensorArray *)fVoltageArray.At(run);
948 UpdateRunInformations(run);
949 voltageArray = (AliDCSSensorArray *)fVoltageArray.At(run);
954 AliDCSSensorArray * AliTPCcalibDB::GetGoofieSensors(Int_t run){
956 // Get temperature sensor array
958 AliDCSSensorArray * goofieArray = (AliDCSSensorArray *)fGoofieArray.At(run);
960 UpdateRunInformations(run);
961 goofieArray = (AliDCSSensorArray *)fGoofieArray.At(run);
968 AliTPCCalibVdrift * AliTPCcalibDB::GetVdrift(Int_t run){
970 // Get the interface to the the vdrift
972 AliTPCCalibVdrift * vdrift = (AliTPCCalibVdrift*)fVdriftArray.At(run);
974 UpdateRunInformations(run);
975 vdrift= (AliTPCCalibVdrift*)fVdriftArray.At(run);
980 Float_t AliTPCcalibDB::GetCEdriftTime(Int_t run, Int_t sector, Double_t timeStamp, Int_t *entries)
983 // GetCE drift time information for 'sector'
984 // sector 72 is the mean drift time of the A-Side
985 // sector 73 is the mean drift time of the C-Side
986 // it timestamp==-1 return mean value
988 AliTPCcalibDB::Instance()->SetRun(run);
989 TGraph *gr=AliTPCcalibDB::Instance()->GetCErocTgraph(sector);
990 if (!gr||sector<0||sector>73) {
991 if (entries) *entries=0;
998 for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1000 gr->GetPoint(ipoint,x,y);
1001 if (x<timeStamp) continue;
1009 Float_t AliTPCcalibDB::GetCEchargeTime(Int_t run, Int_t sector, Double_t timeStamp, Int_t *entries)
1012 // GetCE mean charge for 'sector'
1013 // it timestamp==-1 return mean value
1015 AliTPCcalibDB::Instance()->SetRun(run);
1016 TGraph *gr=AliTPCcalibDB::Instance()->GetCErocQgraph(sector);
1017 if (!gr||sector<0||sector>71) {
1018 if (entries) *entries=0;
1022 if (timeStamp==-1.){
1025 for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1027 gr->GetPoint(ipoint,x,y);
1028 if (x<timeStamp) continue;
1036 Float_t AliTPCcalibDB::GetDCSSensorValue(AliDCSSensorArray *arr, Int_t timeStamp, const char * sensorName, Int_t sigDigits)
1039 // Get Value for a DCS sensor 'sensorName', run 'run' at time 'timeStamp'
1042 const TString sensorNameString(sensorName);
1043 AliDCSSensor *sensor = arr->GetSensor(sensorNameString);
1044 if (!sensor) return val;
1045 //use the dcs graph if possible
1046 TGraph *gr=sensor->GetGraph();
1048 for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1050 gr->GetPoint(ipoint,x,y);
1051 Int_t time=sensor->GetStartTime()+x*3600; //time in graph is hours
1052 if (time<timeStamp) continue;
1056 //if val is still 0, test if if the requested time if within 5min of the first/last
1057 //data point. If this is the case return the firs/last entry
1058 //the timestamps might not be syncronised for all calibration types, sometimes a 'pre'
1059 //and 'pos' period is requested. Especially to the HV this is not the case!
1063 gr->GetPoint(0,x,y);
1064 Int_t time=sensor->GetStartTime()+x*3600; //time in graph is hours
1065 if ((time-timeStamp)<5*60) val=y;
1070 gr->GetPoint(gr->GetN()-1,x,y);
1071 Int_t time=sensor->GetStartTime()+x*3600; //time in graph is hours
1072 if ((timeStamp-time)<5*60) val=y;
1075 val=sensor->GetValue(timeStamp);
1078 val=(Float_t)TMath::Floor(val * TMath::Power(10., sigDigits) + .5) / TMath::Power(10., sigDigits);
1083 Float_t AliTPCcalibDB::GetDCSSensorMeanValue(AliDCSSensorArray *arr, const char * sensorName, Int_t sigDigits)
1086 // Get mean Value for a DCS sensor 'sensorName' during run 'run'
1089 const TString sensorNameString(sensorName);
1090 AliDCSSensor *sensor = arr->GetSensor(sensorNameString);
1091 if (!sensor) return val;
1093 //use dcs graph if it exists
1094 TGraph *gr=sensor->GetGraph();
1098 //if we don't have the dcs graph, try to get some meaningful information
1099 if (!sensor->GetFit()) return val;
1100 Int_t nKnots=sensor->GetFit()->GetKnots();
1101 Double_t tMid=(sensor->GetEndTime()-sensor->GetStartTime())/2.;
1102 for (Int_t iKnot=0;iKnot<nKnots;++iKnot){
1103 if (sensor->GetFit()->GetX()[iKnot]>tMid/3600.) break;
1104 val=(Float_t)sensor->GetFit()->GetY0()[iKnot];
1109 val=(Float_t)TMath::Floor(val * TMath::Power(10., sigDigits) + .5) / TMath::Power(10., sigDigits);
1115 Float_t AliTPCcalibDB::GetChamberHighVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits) {
1117 // return the chamber HV for given run and time: 0-35 IROC, 36-72 OROC
1118 // if timeStamp==-1 return mean value
1121 TString sensorName="";
1122 TTimeStamp stamp(timeStamp);
1123 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1124 if (!voltageArray || (sector<0) || (sector>71)) return val;
1125 Char_t sideName='A';
1126 if ((sector/18)%2==1) sideName='C';
1129 sensorName=Form("TPC_ANODE_I_%c%02d_VMEAS",sideName,sector%18);
1132 sensorName=Form("TPC_ANODE_O_%c%02d_0_VMEAS",sideName,sector%18);
1135 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1137 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1141 Float_t AliTPCcalibDB::GetSkirtVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1144 // Get the skirt voltage for 'run' at 'timeStamp' and 'sector': 0-35 IROC, 36-72 OROC
1145 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1146 // if timeStamp==-1 return the mean value for the run
1149 TString sensorName="";
1150 TTimeStamp stamp(timeStamp);
1151 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1152 if (!voltageArray || (sector<0) || (sector>71)) return val;
1153 Char_t sideName='A';
1154 if ((sector/18)%2==1) sideName='C';
1155 sensorName=Form("TPC_SKIRT_%c_VMEAS",sideName);
1157 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1159 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1164 Float_t AliTPCcalibDB::GetCoverVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1167 // Get the cover voltage for run 'run' at time 'timeStamp'
1168 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1169 // if timeStamp==-1 return the mean value for the run
1172 TString sensorName="";
1173 TTimeStamp stamp(timeStamp);
1174 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1175 if (!voltageArray || (sector<0) || (sector>71)) return val;
1176 Char_t sideName='A';
1177 if ((sector/18)%2==1) sideName='C';
1180 sensorName=Form("TPC_COVER_I_%c_VMEAS",sideName);
1183 sensorName=Form("TPC_COVER_O_%c_VMEAS",sideName);
1186 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1188 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1193 Float_t AliTPCcalibDB::GetGGoffsetVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1196 // Get the GG offset voltage for run 'run' at time 'timeStamp'
1197 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1198 // if timeStamp==-1 return the mean value for the run
1201 TString sensorName="";
1202 TTimeStamp stamp(timeStamp);
1203 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1204 if (!voltageArray || (sector<0) || (sector>71)) return val;
1205 Char_t sideName='A';
1206 if ((sector/18)%2==1) sideName='C';
1209 sensorName=Form("TPC_GATE_I_%c_OFF_VMEAS",sideName);
1212 sensorName=Form("TPC_GATE_O_%c_OFF_VMEAS",sideName);
1215 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1217 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1222 Float_t AliTPCcalibDB::GetGGnegVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1225 // Get the GG offset voltage for run 'run' at time 'timeStamp'
1226 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1227 // if timeStamp==-1 return the mean value for the run
1230 TString sensorName="";
1231 TTimeStamp stamp(timeStamp);
1232 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1233 if (!voltageArray || (sector<0) || (sector>71)) return val;
1234 Char_t sideName='A';
1235 if ((sector/18)%2==1) sideName='C';
1238 sensorName=Form("TPC_GATE_I_%c_NEG_VMEAS",sideName);
1241 sensorName=Form("TPC_GATE_O_%c_NEG_VMEAS",sideName);
1244 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1246 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1251 Float_t AliTPCcalibDB::GetGGposVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1254 // Get the GG offset voltage for run 'run' at time 'timeStamp'
1255 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1256 // if timeStamp==-1 return the mean value for the run
1259 TString sensorName="";
1260 TTimeStamp stamp(timeStamp);
1261 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1262 if (!voltageArray || (sector<0) || (sector>71)) return val;
1263 Char_t sideName='A';
1264 if ((sector/18)%2==1) sideName='C';
1267 sensorName=Form("TPC_GATE_I_%c_POS_VMEAS",sideName);
1270 sensorName=Form("TPC_GATE_O_%c_POS_VMEAS",sideName);
1273 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1275 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1280 Float_t AliTPCcalibDB::GetPressure(Int_t timeStamp, Int_t run, Int_t type){
1282 // GetPressure for given time stamp and runt
1284 TTimeStamp stamp(timeStamp);
1285 AliDCSSensor * sensor = Instance()->GetPressureSensor(run,type);
1286 if (!sensor) return 0;
1287 return sensor->GetValue(stamp);
1290 Float_t AliTPCcalibDB::GetL3Current(Int_t run, Int_t statType){
1292 // return L3 current
1293 // stat type is: AliGRPObject::Stats: kMean = 0, kTruncMean = 1, kMedian = 2, kSDMean = 3, kSDMedian = 4
1296 AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1297 if (grp) current=grp->GetL3Current((AliGRPObject::Stats)statType);
1301 Float_t AliTPCcalibDB::GetBz(Int_t run){
1303 // calculate BZ in T from L3 current
1306 Float_t current=AliTPCcalibDB::GetL3Current(run);
1307 if (current>-1) bz=5*current/30000.*.1;
1311 Char_t AliTPCcalibDB::GetL3Polarity(Int_t run) {
1313 // get l3 polarity from GRP
1315 return AliTPCcalibDB::GetGRP(run)->GetL3Polarity();
1318 TString AliTPCcalibDB::GetRunType(Int_t run){
1320 // return run type from grp
1322 return AliTPCcalibDB::GetGRP(run)->GetRunType();
1325 Float_t AliTPCcalibDB::GetValueGoofie(Int_t timeStamp, Int_t run, Int_t type){
1327 // GetPressure for given time stamp and runt
1329 TTimeStamp stamp(timeStamp);
1330 AliDCSSensorArray* goofieArray = AliTPCcalibDB::Instance()->GetGoofieSensors(run);
1331 if (!goofieArray) return 0;
1332 AliDCSSensor *sensor = goofieArray->GetSensor(type);
1333 return sensor->GetValue(stamp);
1341 Bool_t AliTPCcalibDB::GetTemperatureFit(Int_t timeStamp, Int_t run, Int_t side,TVectorD& fit){
1345 TTimeStamp tstamp(timeStamp);
1346 AliTPCSensorTempArray* tempArray = Instance()->GetTemperatureSensor(run);
1347 if (! tempArray) return kFALSE;
1348 AliTPCTempMap * tempMap = new AliTPCTempMap(tempArray);
1349 TLinearFitter * fitter = tempMap->GetLinearFitter(3,side,tstamp);
1352 fitter->GetParameters(fit);
1356 if (!fitter) return kFALSE;
1360 Float_t AliTPCcalibDB::GetTemperature(Int_t timeStamp, Int_t run, Int_t side){
1366 GetTemperatureFit(timeStamp,run,0,vec);
1370 GetTemperatureFit(timeStamp,run,0,vec);
1377 Double_t AliTPCcalibDB::GetPTRelative(UInt_t timeSec, Int_t run, Int_t side){
1380 // time - absolute time
1382 // side - 0 - A side 1-C side
1383 AliTPCCalibVdrift * vdrift = Instance()->GetVdrift(run);
1384 if (!vdrift) return 0;
1385 return vdrift->GetPTRelative(timeSec,side);
1388 AliGRPObject * AliTPCcalibDB::MakeGRPObjectFromMap(TMap *map){
1390 // Function to covert old GRP run information from TMap to GRPObject
1392 // TMap * map = AliTPCcalibDB::GetGRPMap(52406);
1394 AliDCSSensor * sensor = 0;
1396 osensor = ((*map)("fP2Pressure"));
1397 sensor =dynamic_cast<AliDCSSensor *>(osensor);
1399 if (!sensor) return 0;
1401 AliDCSSensor * sensor2 = new AliDCSSensor(*sensor);
1402 osensor = ((*map)("fCavernPressure"));
1403 TGraph * gr = new TGraph(2);
1404 gr->GetX()[0]= -100000.;
1405 gr->GetX()[1]= 1000000.;
1406 gr->GetY()[0]= atof(osensor->GetName());
1407 gr->GetY()[1]= atof(osensor->GetName());
1408 sensor2->SetGraph(gr);
1412 AliGRPObject *grpRun = new AliGRPObject;
1413 grpRun->ReadValuesFromMap(map);
1414 grpRun->SetCavernAtmosPressure(sensor2);
1415 grpRun->SetSurfaceAtmosPressure(sensor);
1419 Bool_t AliTPCcalibDB::CreateGUITree(Int_t run, const char* filename)
1422 // Create a gui tree for run number 'run'
1425 if (!AliCDBManager::Instance()->GetDefaultStorage()){
1426 AliLog::Message(AliLog::kError, "Default Storage not set. Cannot create Calibration Tree!",
1427 MODULENAME(), "AliTPCcalibDB", FUNCTIONNAME(), __FILE__, __LINE__);
1431 AliTPCcalibDB *db=AliTPCcalibDB::Instance();
1432 // retrieve cal pad objects
1434 AliTPCPreprocessorOnline prep;
1435 //noise and pedestals
1436 prep.AddComponent(db->GetPedestals());
1437 prep.AddComponent(db->GetPadNoise());
1439 prep.AddComponent(db->GetPulserTmean());
1440 prep.AddComponent(db->GetPulserTrms());
1441 prep.AddComponent(db->GetPulserQmean());
1443 prep.AddComponent(db->GetCETmean());
1444 prep.AddComponent(db->GetCETrms());
1445 prep.AddComponent(db->GetCEQmean());
1447 prep.AddComponent(db->GetALTROAcqStart() );
1448 prep.AddComponent(db->GetALTROZsThr() );
1449 prep.AddComponent(db->GetALTROFPED() );
1450 prep.AddComponent(db->GetALTROAcqStop() );
1451 prep.AddComponent(db->GetALTROMasked() );
1453 TString file(filename);
1454 if (file.IsNull()) file=Form("guiTreeRun_%d.root",run);
1455 prep.DumpToFile(file.Data());