2 /**************************************************************************
3 * This file is property of and copyright by the ALICE HLT Project *
4 * ALICE Experiment at CERN, All rights reserved. *
6 * Primary Authors: Jochen Thaeder <thaeder@kip.uni-heidelberg.de> *
7 * for The ALICE HLT Project. *
9 * Permission to use, copy, modify and distribute this software and its *
10 * documentation strictly for non-commercial purposes is hereby granted *
11 * without fee, provided that the above copyright notice appears in all *
12 * copies and that both the copyright notice and this permission notice *
13 * appear in the supporting documentation. The authors make no claims *
14 * about the suitability of this software for any purpose. It is *
15 * provided "as is" without express or implied warranty. *
16 **************************************************************************/
18 /** @file AliHLTTPCCalibTimeComponent.cxx
19 @author Kalliopi Kanaki
21 @brief A calibration component for interfacing the offline calculation of TPC drift velocity correction
24 // see header file for class documentation
26 // refer to README to build package
28 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
30 #include "AliHLTTPCCalibTimeComponent.h"
31 #include "AliHLTTPCDefinitions.h"
32 #include "AliHLTMisc.h"
34 #include "AliESDEvent.h"
35 #include "AliESDtrack.h"
36 #include "AliESDfriend.h"
38 #include "AliTPCcalibTime.h"
39 #include "AliTPCcalibCalib.h"
40 #include "AliTPCseed.h"
41 #include "AliTPCcalibDB.h"
42 #include "AliTPCClusterParam.h"
44 #include "AliHLTTPCOfflineCluster.h"
45 #include "AliHLTTPCSpacePointData.h"
46 #include "AliHLTTPCTrackletDataFormat.h"
47 #include "AliHLTExternalTrackParam.h"
48 #include "AliHLTGlobalBarrelTrack.h"
49 #include "AliHLTTPCTransform.h"
51 #include "TObjArray.h"
55 #include "THnSparse.h"
56 #include "TGraphErrors.h"
61 #include "AliHLTReadoutList.h"
65 ClassImp(AliHLTTPCCalibTimeComponent) // ROOT macro for the implementation of ROOT specific class methods
67 AliHLTTPCCalibTimeComponent::AliHLTTPCCalibTimeComponent()
77 // see header file for class documentation
79 // refer to README to build package
81 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
83 for(int i=0; i<fkNPartition; i++){
84 fPartitionClusters[i] = 0;
85 fNPartitionClusters[i] = 0;
89 const char* AliHLTTPCCalibTimeComponent::fgkOCDBEntry="HLT/ConfigTPC/TPCCalibTime";
91 AliHLTTPCCalibTimeComponent::~AliHLTTPCCalibTimeComponent(){
92 // see header file for class documentation
94 for(int i=0; i<fkNPartition; i++){
95 delete[] fPartitionClusters[i];
99 const char* AliHLTTPCCalibTimeComponent::GetComponentID() {
100 // see header file for class documentation
101 return "TPCCalibTime";
104 void AliHLTTPCCalibTimeComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list){
105 // see header file for class documentation
108 list.push_back( kAliHLTDataTypeTObjArray|kAliHLTDataOriginTPC ); // output of TPCCalibSeedMaker
109 list.push_back( AliHLTTPCDefinitions::fgkClustersDataType ); // output of the TPC CF
110 list.push_back( kAliHLTDataTypeTrack|kAliHLTDataOriginTPC ); // output of the global merger
111 list.push_back( kAliHLTDataTypeESDObject|kAliHLTDataOriginOut ); // output of global esd converter
114 AliHLTComponentDataType AliHLTTPCCalibTimeComponent::GetOutputDataType() {
115 // see header file for class documentation
117 return AliHLTTPCDefinitions::CalibCEDataType()|kAliHLTDataOriginOut;
120 void AliHLTTPCCalibTimeComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ){
121 // see header file for class documentation
123 constBase = fOutputSize;
124 inputMultiplier = 0; // to be estimated
127 AliHLTComponent* AliHLTTPCCalibTimeComponent::Spawn(){
128 // see header file for class documentation
130 return new AliHLTTPCCalibTimeComponent();
134 Int_t AliHLTTPCCalibTimeComponent::ScanConfigurationArgument( Int_t argc, const char** argv ){
135 // see header file for class documentation
137 if (argc<=0) return 0;
139 TString argument=argv[i];
142 if (argument.CompareTo("-output-size")==0) {
143 if (++i>=argc) return -EPROTO;
145 fOutputSize=argument.Atoi();
151 Int_t AliHLTTPCCalibTimeComponent::InitCalibration(){
152 // see header file for class documentation
154 AliTPCcalibDB::Instance()->SetRun(AliHLTMisc::Instance().GetCDBRunNo());
155 AliTPCcalibDB::Instance()->GetClusterParam()->SetInstance(AliTPCcalibDB::Instance()->GetClusterParam());
158 // AliTPCcalibDB *calib = AliTPCcalibDB::Instance();
161 // HLTError("AliTPCcalibDB does not exist");
165 // AliTPCClusterParam *clusPar = calib->GetClusterParam();
167 // HLTError("OCDB entry TPC/Calib/ClusterParam (AliTPCcalibDB::GetClusterParam()) is not available.");
171 // first configure the default
173 if (iResult>=0) iResult=ConfigureFromCDBTObjString(fgkOCDBEntry);
175 if(fCalibTime) return EINPROGRESS;
176 fCal = new AliTPCcalibCalib();
178 fSeedArray = new TObjArray();
184 Int_t AliHLTTPCCalibTimeComponent::DeinitCalibration() {
185 // see header file for class documentation
187 if(fCalibTime) delete fCalibTime; fCalibTime = NULL;
188 if(fCal) delete fCal; fCal = NULL;
189 if(fSeedArray) delete fSeedArray; fSeedArray = NULL;
191 //if(fESDfriend) delete fESDfriend; fESDfriend = NULL;
193 //if(arr) delete arr; arr = NULL;
198 int AliHLTTPCCalibTimeComponent::Reconfigure(const char* cdbEntry, const char* /*chainId*/){
199 // see header file for class documentation
201 // configure from the specified antry or the default one
202 const char* entry=cdbEntry;
203 if(!entry || entry[0]==0){
207 return ConfigureFromCDBTObjString(entry);
210 Int_t AliHLTTPCCalibTimeComponent::ProcessCalibration( const AliHLTComponentEventData& /*evtData*/, AliHLTComponentTriggerData& /*trigData*/ ){
211 // see header file for class documentation
213 if(GetFirstInputBlock( kAliHLTDataTypeSOR ) || GetFirstInputBlock( kAliHLTDataTypeEOR )) return 0;
215 const AliHLTComponentBlockData *iter = NULL;
217 //--------------- output over TObjArray of AliTPCseed objects (output of TPCSeedMaker) -------------------//
219 // A previous component in the chain (TPCSeedMaker) has processed the TPC clusters and tracks and created a TClonesArray of AliTPCseed objects
220 // In this loop the iterator accesses this array stored in memory, in order to use it in the next loop over the AliESDevent of the HLT
222 // for(iter = (TObject*)GetFirstInputObject(kAliHLTDataTypeTObjArray|kAliHLTDataOriginTPC); iter != NULL; iter = (TObject*)GetNextInputObject()){
224 // if(GetDataType(iter) != (kAliHLTDataTypeTObjArray | kAliHLTDataOriginTPC)) continue;
225 // fSeedArray = dynamic_cast<TClonesArray*>(iter);
228 // int nInputClusters = 0;
229 // int nInputTracks = 0;
231 //TObjArray *arr = new TObjArray(1000);
232 //arr->SetOwner(kTRUE);
236 for( int i=0; i<fkNPartition; i++ ){
237 delete[] fPartitionClusters[i];
238 fPartitionClusters[i] = 0;
239 fNPartitionClusters[i] = 0;
244 //------------------- loop over clusters -------------//
246 for(iter = GetFirstInputBlock(AliHLTTPCDefinitions::fgkClustersDataType); iter != NULL; iter = GetNextInputBlock()){
248 if ( iter->fDataType != AliHLTTPCDefinitions::fgkClustersDataType ) continue;
250 Int_t slice = AliHLTTPCDefinitions::GetMinSliceNr(iter->fSpecification);
251 Int_t partition = AliHLTTPCDefinitions::GetMinPatchNr(iter->fSpecification);
253 Int_t slicepartition = slice*6+partition;
255 if(slicepartition > fkNPartition){
256 HLTWarning("Wrong header of TPC cluster data, slice %d, partition %d", slice, partition );
260 AliHLTTPCClusterData* inPtrSP = ( AliHLTTPCClusterData* )( iter->fPtr );
261 // nInputClusters += inPtrSP->fSpacePointCnt;
263 delete[] fPartitionClusters[slicepartition];
264 fPartitionClusters[slicepartition] = new AliTPCclusterMI[inPtrSP->fSpacePointCnt];
265 fNPartitionClusters[slicepartition] = inPtrSP->fSpacePointCnt;
267 // create offline clusters out of the HLT clusters
268 // todo: check which cluster information is really needed for the dEdx
269 for(unsigned int i = 0; i < inPtrSP->fSpacePointCnt; i++){
270 AliHLTTPCSpacePointData *chlt = &( inPtrSP->fSpacePoints[i] );
271 AliTPCclusterMI *c = fPartitionClusters[slicepartition]+i;
275 c->SetSigmaY2(chlt->fSigmaY2);
277 c->SetSigmaZ2(chlt->fSigmaZ2);
278 c->SetQ( chlt->fCharge );
279 c->SetMax( chlt->fQMax );
281 Float_t padtime[3] = {0,chlt->fY,chlt->fZ};
282 AliHLTTPCTransform::Slice2Sector(slice,chlt->fPadRow, sector, row);
283 AliHLTTPCTransform::Local2Raw( padtime, sector, row);
284 c->SetDetector( sector );
286 c->SetPad( (Int_t) padtime[1] );
287 c->SetTimeBin( (Int_t) padtime[2] );
289 } // end of loop over blocks of clusters
294 //---------- loop over merged tracks ------------------ //
296 for(const AliHLTComponentBlockData *pBlock = GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC); pBlock != NULL; pBlock = GetNextInputBlock()){
298 AliHLTTracksData *dataPtr = (AliHLTTracksData*) pBlock->fPtr;
299 nTracks = dataPtr->fCount;
301 AliHLTExternalTrackParam *currTrack = dataPtr->fTracklets;
303 //nInputTracks += nTracks;
305 for(int itr=0; itr<nTracks && ( (AliHLTUInt8_t *)currTrack < ((AliHLTUInt8_t *) pBlock->fPtr)+pBlock->fSize); itr++){
307 // create an offline track
308 AliHLTGlobalBarrelTrack gb(*currTrack);
310 tTPC.Set(gb.GetX(), gb.GetAlpha(), gb.GetParameter(), gb.GetCovariance());
312 // set the cluster pointers
313 for(UInt_t ic=0; ic<currTrack->fNPoints; ic++){
315 tTPC.SetNumberOfClusters(currTrack->fNPoints);
317 UInt_t id = currTrack->fPointIDs[ic];
318 int iSlice = AliHLTTPCSpacePointData::GetSlice(id);
319 int iPartition = AliHLTTPCSpacePointData::GetPatch(id);
320 int iCluster = AliHLTTPCSpacePointData::GetNumber(id);
322 if(iSlice<0 || iSlice>36 || iPartition<0 || iPartition>5){
323 HLTError("Corrupted TPC cluster Id: slice %d, partition %d, cluster %d", iSlice, iPartition, iCluster);
327 AliTPCclusterMI *partitionClusters = fPartitionClusters[iSlice*6 + iPartition];
329 if(!partitionClusters){
330 HLTError("Clusters are missed for slice %d, partition %d", iSlice, iPartition );
334 if(iCluster >= fNPartitionClusters[iSlice*6 + iPartition]){
335 HLTError("TPC slice %d, partition %d: ClusterID==%d >= N Clusters==%d ", iSlice, iPartition,iCluster, fNPartitionClusters[iSlice*6 + iPartition]);
339 AliTPCclusterMI *c = &(partitionClusters[iCluster]);
340 int sec = c->GetDetector();
341 int row = c->GetRow();
342 if(sec >= 36) row = row + AliHLTTPCTransform::GetNRowLow();
344 tTPC.SetClusterPointer(row, c);
346 AliTPCTrackerPoint &point = *( tTPC.GetTrackPoint( row ) );
347 //tTPC.Propagate( TMath::DegToRad()*(sec%18*20.+10.), c->GetX(), fSolenoidBz );
348 Double_t angle2 = tTPC.GetSnp()*tTPC.GetSnp();
349 angle2 = (angle2<1) ?TMath::Sqrt(angle2/(1-angle2)) :10.;
350 point.SetAngleY( angle2 );
351 point.SetAngleZ( tTPC.GetTgl() );
352 } // end of associated cluster loop
354 AliTPCseed *seed = &(tTPC);
355 fSeedArray->Add(seed);
357 unsigned int step = sizeof( AliHLTExternalTrackParam ) + currTrack->fNPoints * sizeof( unsigned int );
358 currTrack = ( AliHLTExternalTrackParam* )( (( Byte_t * )currTrack) + step );
360 }// end of vector track loop
361 } // end of loop over blocks of merged tracks
363 HLTInfo("Number of reconstructed tracks %d, number of produced seeds %d\n", nTracks, fSeedArray->GetEntries());
366 //----------- loop over output of global esd converter ----------------//
368 // In this loop we access the AliESDevent that was produced by the HLT and is stored in memory. There should exist 1 object
369 // of type kAliHLTDataTypeESDObject per event.
371 TObject *iterOb = NULL;
372 for(iterOb = (TObject*)GetFirstInputObject(kAliHLTDataTypeESDObject | kAliHLTDataOriginOut); iterOb != NULL; iterOb = (TObject*)GetNextInputObject()){
374 if(GetDataType(iterOb) != (kAliHLTDataTypeESDObject | kAliHLTDataOriginOut)) continue;
376 fESDevent = dynamic_cast<AliESDEvent*>(iterOb);
377 fESDevent->GetStdContent();
379 HLTInfo("Number of seeds: %i\n", fSeedArray->GetEntriesFast()); // access of the info from the previous loop over the AliTPCseed array
381 fCal->UpdateEventInfo(fESDevent);
382 for(Int_t i=0; i<fSeedArray->GetEntriesFast(); i++){ // loop over TObjArray with seeds
384 AliTPCseed *seed = (AliTPCseed*)fSeedArray->UncheckedAt(i);
385 fESDtrack = fESDevent->GetTrack(i);
386 if(!fESDtrack || !seed) continue;
388 //if(fESDtrack->GetID() != seed->GetLabel()) {
389 // HLTWarning("Mismatch of track id between seed and ESD track: %i, %i\n", fESDtrack->GetID(), seed->GetLabel());
393 if(seed->GetNumberOfClusters()==0) continue;
394 fCal->RefitTrack(fESDtrack, seed, GetBz()); // update AliESDtrack and AliTPCseed info, acccording to Marian's request
396 AliTPCseed *seedCopy = new AliTPCseed(*seed, kTRUE);
397 fESDtrack->AddCalibObject(seedCopy); // add the AliTPCseed as a friend track to the AliESDtrack (to be accessed in TPC/AliTPCcalibTime.cxx)
399 //fESDfriendTrack = const_cast<AliESDfriendTrack*>(fESDtrack->GetFriendTrack());
403 Int_t startTime = fESDevent->GetTimeStamp()-60*60*1; //Start time one hour before first event, will make precise cuts later.
404 Int_t endTime = fESDevent->GetTimeStamp()+60*60*23; //End time 23 hours after first event.
405 // create the calibration object that will call the offline functions
406 fCalibTime = new AliTPCcalibTime("calibTime","time dependent Vdrift calibration", startTime, endTime, 20*60);
408 printf("fCalibTime pointer is NULL\n");
412 fCalibTime->SetStreamLevel(20);
413 fCalibTime->SetDebugLevel(20);
414 printf("fCalibTime = %i, startTime = %i, endTime = %i \n", fCalibTime!=0, startTime, endTime);
415 fESDfriend = new AliESDfriend();
416 fESDevent->GetESDfriend(fESDfriend);
417 fESDevent->SetESDfriend(fESDfriend);
418 fESDevent->AddObject(fESDfriend);
419 // create the AliESDfriend and add it to the event, now both the friend tracks and the friends are available for the offline functions to be called
421 fCalibTime->UpdateEventInfo(fESDevent); // needed for getting the run number and time stamp information correct on the offline side
422 fCalibTime->Process(fESDevent); // first offline function called
424 // delete fESDfriend;
426 //PushBack( (TObject*)fCalibTime, AliHLTTPCDefinitions::fgkCalibCEDataType | kAliHLTDataOriginOut, 0x0);
431 Int_t AliHLTTPCCalibTimeComponent::ShipDataToFXS( const AliHLTComponentEventData& /*evtData*/, AliHLTComponentTriggerData& /*trigData*/ ){
432 // see header file for class documentation
434 HLTInfo("Shipping data to FXS...\n");
436 fCalibTime->Analyze(); // called at the end of the run or event modulo
438 // the rest of the histogram and graph declarations were copied by Dag as a first attempt to get the start/end time bin "automatically". Perhaps we need some more
439 // thinking here to avoid copying all these lines that might chain in offline without HLT realizing.
441 THnSparse* addHist = fCalibTime->GetHistoDrift("all");
442 if(!addHist) return -1;
444 //Identifying used range of histogram
446 Int_t startTimeBin = 0;
447 Int_t endTimeBin = 0;
449 TH1D *histoTime = addHist->Projection(0);
451 startTimeBin = histoTime->FindFirstBinAbove(0);
452 endTimeBin = histoTime->FindLastBinAbove(0);
453 printf("startTimeBin = %i endTimeBin = %i\n", startTimeBin, endTimeBin);
454 printf("startTimeBinCentre = %f endTimeBinCentre = %f\n", histoTime->GetBinCenter(startTimeBin), histoTime->GetBinCenter(endTimeBin));
455 printf("startTimeBinWidth = %f endTimeBinWidth = %f\n", histoTime->GetBinWidth(startTimeBin), histoTime->GetBinWidth(endTimeBin));
456 delete histoTime; histoTime = 0;
459 Int_t startPtBin = 0;
461 TH1D *histoPt = addHist->Projection(1);
463 startPtBin = histoPt->FindFirstBinAbove(0);
464 endPtBin = histoPt->FindLastBinAbove(0);
465 printf("startPtBin = %i endPtBin = %i\n", startPtBin, endPtBin);
466 printf("startPtBinCentre = %f endPtBinCentre = %f\n", histoPt->GetBinCenter(startPtBin), histoPt->GetBinCenter(endPtBin));
467 printf("startPtinWidth = %f endPtBinWidth = %f\n", histoPt->GetBinWidth(startPtBin), histoPt->GetBinWidth(endPtBin));
468 delete histoPt; histoPt = 0;
471 Int_t startVdBin = 0;
473 TH1D *histoVd = addHist->Projection(2);
475 startVdBin = histoVd->FindFirstBinAbove(0);
476 endVdBin = histoVd->FindLastBinAbove(0);
477 printf("startVdBin = %i endVdBin = %i\n", startVdBin, endVdBin);
478 printf("startVdBinCentre = %f endVdBinCentre = %f\n", histoVd->GetBinCenter(startVdBin), histoVd->GetBinCenter(endVdBin));
479 printf("startVdBinWidth = %f endVdBinWidth = %f\n", histoVd->GetBinWidth(startVdBin), histoVd->GetBinWidth(endVdBin));
480 delete histoVd; histoVd = 0;
483 Int_t startRunBin = 0;
485 TH1D *histoRun = addHist->Projection(3);
487 startRunBin = histoRun->FindFirstBinAbove(0);
488 endRunBin = histoRun->FindLastBinAbove(0);
489 printf("startRunBin = %i endRunBin = %i\n", startRunBin, endRunBin);
490 printf("startRunBinCentre = %f endRunBinCentre = %f\n", histoRun->GetBinCenter(startRunBin), histoRun->GetBinCenter(endRunBin));
491 printf("startRunBinWidth = %f endRunBinWidth = %f\n", histoRun->GetBinWidth(startRunBin), histoRun->GetBinWidth(endRunBin));
492 delete histoRun; histoRun = 0;
495 TObjArray *vdriftArray = new TObjArray();
496 if(!vdriftArray) return -2;
498 TObjArray *array = fCalibTime->GetHistoDrift();
499 if(!array) return -3;
501 TIterator *iterator = array->MakeIterator();
502 if(!iterator) return -4;
505 THnSparse *hist = NULL;
506 while((hist = (THnSparseF*)iterator->Next())){
508 //if(!hist) continue;
510 hist->GetAxis(0)->SetRange(startTimeBin, endTimeBin);
511 hist->GetAxis(1)->SetRange(startPtBin, endPtBin);
512 hist->GetAxis(0)->SetRange(startVdBin, endVdBin);
513 hist->GetAxis(3)->SetRange(startRunBin, endRunBin);
515 TString name = hist->GetName();
516 Int_t dim[4] = {0,1,2,3};
517 THnSparse *newHist = hist->Projection(4,dim);
518 newHist->SetName(name);
519 vdriftArray->Add(newHist);
521 TGraphErrors *graph = AliTPCcalibBase::FitSlices(newHist,2,0,400,100,0.05,0.95, kTRUE);
522 printf("name = %s graph = %i\n", name.Data(), graph==0);
523 if(!graph || !graph->GetN()) continue;
524 printf("name = %s graph = %i, N = %i\n", name.Data(), graph==0, graph->GetN());
525 Int_t pos = name.Index("_");
526 name = name(pos,name.Capacity()-pos);
527 TString graphName = graph->ClassName();
530 graph->SetName(graphName);
531 printf("name = %s\n", graphName.Data());
532 vdriftArray->Add(graph);
534 //Currently, AliSplineFits can not be given names...
535 //AliSplineFit* fit=new AliSplineFit();
536 //fit->SetGraph(graph);
537 //fit->SetMinPoints(graph->GetN()+1);
538 //fit->InitKnots(graph,2,0,0.001);
540 //TString fiName=fit->ClassName();
544 //fit->SetName(fiName.Data());
545 //printf("name=%s\n", fiName.Data());
546 //vdriftArray->Add(fit);
549 THnSparse *laserHist = NULL;
550 TGraphErrors *laserGraph = NULL;
551 TString laserName = "";
553 //Histograms and graphs for A side lasers
554 laserHist = fCalibTime->GetHistVdriftLaserA(1);
557 laserName=laserHist->ClassName();
558 laserName+="_MEAN_DRIFT_LASER_ALL_A";
560 laserHist->SetName(laserName);
561 vdriftArray->Add(laserHist);
562 laserGraph=AliTPCcalibBase::FitSlices(laserHist,2,0,400,100,0.05,0.95, kTRUE);
563 if(laserGraph && laserGraph->GetN()){
564 laserName=laserGraph->GetName();
565 laserName+="_MEAN_DRIFT_LASER_ALL_A";
567 laserGraph->SetName(laserName);
568 vdriftArray->Add(laserGraph);
572 //Histograms and graphs for C side lasers
573 laserHist=fCalibTime->GetHistVdriftLaserC(1);
575 laserName=laserHist->ClassName();
576 laserName+="_MEAN_DRIFT_LASER_ALL_C";
578 laserHist->SetName(laserName);
579 vdriftArray->Add(laserHist);
580 laserGraph=AliTPCcalibBase::FitSlices(laserHist,2,0,400,100,0.05,0.95, kTRUE);
581 if(laserGraph && laserGraph->GetN()){
582 laserName=laserGraph->GetName();
583 laserName+="_MEAN_DRIFT_LASER_ALL_C";
585 laserGraph->SetName(laserName);
586 vdriftArray->Add(laserGraph);
590 //Meatdata set in off-line...
591 //AliCDBMetaData *metaData= new AliCDBMetaData();
592 //metaData->SetObjectClassName("TObjArray");
593 //metaData->SetResponsible("Dag Toppe Larsen");
594 //metaData->SetBeamPeriod(1);
595 //metaData->SetAliRootVersion("05-25-01"); //root version
596 //metaData->SetComment("Calibration of the time dependence of the drift velocity due to pressure and temperature changes");
597 //AliCDBId* id1=NULL;
598 //if(end) id1=new AliCDBId("TPC/Calib/TimeDrift", runNumber, end);
599 //else id1=new AliCDBId("TPC/Calib/TimeDrift", runNumber, runNumber);
600 //AliCDBStorage* gStorage = AliCDBManager::Instance()->GetStorage("local://$ALICE_ROOT/OCDB");
601 //gStorage->Put(vdriftArray, (*id1), metaData);
602 //printf("done runNumber=%i, end=%i\n", runNumber, end);
604 static AliHLTReadoutList rdList(AliHLTReadoutList::kTPC);
607 TFile *file = TFile::Open("vdrift.root", "RECREATE");
608 vdriftArray->Write();
612 file = TFile::Open("calibTime.root", "RECREATE");
617 // the vdriftArray is pushed to the HLT-FXSsubscriber
618 PushToFXS( (TObject*)vdriftArray, "TPC", "TIMEDRIFT", &rdList );
620 //Should array be deleted now?
622 // vdriftArray.Clear();
623 // delete vdriftArray;