]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/calibration/AliHLTTPCCalibTimeComponent.cxx
Fixing coverity defects 24797 and 24795
[u/mrichter/AliRoot.git] / HLT / TPCLib / calibration / AliHLTTPCCalibTimeComponent.cxx
1 // $Id$
2 /**************************************************************************
3  * This file is property of and copyright by the ALICE HLT Project        * 
4  * ALICE Experiment at CERN, All rights reserved.                         *
5  *                                                                        *
6  * Primary Authors: Jochen Thaeder <thaeder@kip.uni-heidelberg.de>        *
7  *                  for The ALICE HLT Project.                            *
8  *                                                                        *
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  **************************************************************************/
17
18 /** @file   AliHLTTPCCalibTimeComponent.cxx
19     @author Kalliopi Kanaki
20     @date   2009-07-08
21     @brief  A calibration component for interfacing the offline calculation of TPC drift velocity correction 
22 */
23
24 // see header file for class documentation
25 // or
26 // refer to README to build package
27 // or
28 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt   
29
30 #include "AliHLTTPCCalibTimeComponent.h"
31 #include "AliHLTTPCDefinitions.h"
32 #include "AliHLTMisc.h"
33
34 #include "AliESDEvent.h"
35 #include "AliESDtrack.h"
36 #include "AliESDfriend.h"
37
38 #include "AliTPCcalibTime.h"
39 #include "AliTPCcalibCalib.h"
40 #include "AliTPCseed.h"
41 #include "AliTPCcalibDB.h"
42 #include "AliTPCClusterParam.h"
43
44 #include "AliHLTTPCOfflineCluster.h"
45 #include "AliHLTTPCSpacePointData.h"
46 #include "AliHLTTPCTrackletDataFormat.h"
47 #include "AliHLTExternalTrackParam.h"
48 #include "AliHLTGlobalBarrelTrack.h"
49 #include "AliHLTTPCTransform.h"
50
51 #include "TObjArray.h"
52 #include "TString.h"
53 #include "TFile.h"
54
55 #include "THnSparse.h"
56 #include "TGraphErrors.h"
57
58 #include <cstdlib>
59 #include <cerrno>
60
61 #include "AliHLTReadoutList.h"
62
63 using namespace std;
64
65 ClassImp(AliHLTTPCCalibTimeComponent) // ROOT macro for the implementation of ROOT specific class methods
66
67 AliHLTTPCCalibTimeComponent::AliHLTTPCCalibTimeComponent()
68   :
69    fCalibTime(NULL)
70   ,fCal(NULL)
71   ,fESDevent(NULL)
72   ,fESDtrack(NULL)
73   ,fESDfriend(NULL)
74   ,fSeedArray(NULL)
75   ,fOutputSize(50000)
76 {
77   // see header file for class documentation
78   // or
79   // refer to README to build package
80   // or
81   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt  
82   
83   for(int i=0; i<fkNPartition; i++){
84       fPartitionClusters[i]  = 0;    
85       fNPartitionClusters[i] = 0;    
86   }
87 }
88
89 const char* AliHLTTPCCalibTimeComponent::fgkOCDBEntry="HLT/ConfigTPC/TPCCalibTime";
90
91 AliHLTTPCCalibTimeComponent::~AliHLTTPCCalibTimeComponent(){
92 // see header file for class documentation 
93  
94   for(int i=0; i<fkNPartition; i++){
95       delete[] fPartitionClusters[i];
96   }
97 }
98
99 const char* AliHLTTPCCalibTimeComponent::GetComponentID() {
100 // see header file for class documentation
101   return "TPCCalibTime";
102 }
103
104 void AliHLTTPCCalibTimeComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list){
105 // see header file for class documentation
106
107   list.clear(); 
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
112 }
113
114 AliHLTComponentDataType AliHLTTPCCalibTimeComponent::GetOutputDataType() {
115 // see header file for class documentation
116
117   return AliHLTTPCDefinitions::CalibCEDataType()|kAliHLTDataOriginOut;
118 }
119
120 void AliHLTTPCCalibTimeComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ){
121 // see header file for class documentation
122
123   constBase = fOutputSize;
124   inputMultiplier = 0; // to be estimated
125 }
126
127 AliHLTComponent* AliHLTTPCCalibTimeComponent::Spawn(){
128 // see header file for class documentation
129
130   return new AliHLTTPCCalibTimeComponent();
131 }  
132
133
134 Int_t AliHLTTPCCalibTimeComponent::ScanConfigurationArgument( Int_t argc, const char** argv ){
135 // see header file for class documentation
136  
137   if (argc<=0) return 0;
138   int i=0;
139   TString argument=argv[i];
140
141   // -output-size
142   if (argument.CompareTo("-output-size")==0) {
143     if (++i>=argc) return -EPROTO;
144     argument=argv[i];
145     fOutputSize=argument.Atoi();
146     return 2;
147   }
148   return -EINVAL;
149 }
150
151 Int_t AliHLTTPCCalibTimeComponent::InitCalibration(){
152 // see header file for class documentation
153   
154   AliTPCcalibDB::Instance()->SetRun(AliHLTMisc::Instance().GetCDBRunNo());
155   AliTPCcalibDB::Instance()->GetClusterParam()->SetInstance(AliTPCcalibDB::Instance()->GetClusterParam());
156   
157
158 //   AliTPCcalibDB *calib = AliTPCcalibDB::Instance();
159 // 
160 //   if(!calib){
161 //     HLTError("AliTPCcalibDB does not exist");
162 //     return -ENOENT;
163 //   }
164 //   
165 //   AliTPCClusterParam *clusPar = calib->GetClusterParam();
166 //   if(!clusPar){
167 //     HLTError("OCDB entry TPC/Calib/ClusterParam (AliTPCcalibDB::GetClusterParam()) is not available.");
168 //     return -ENOENT;
169 //   }
170
171   // first configure the default
172   int iResult=0;
173   if (iResult>=0) iResult=ConfigureFromCDBTObjString(fgkOCDBEntry);
174     
175   if(fCalibTime) return EINPROGRESS;
176   fCal = new AliTPCcalibCalib();
177   
178   fSeedArray = new TObjArray();
179   
180   return iResult;
181  
182 }
183
184 Int_t AliHLTTPCCalibTimeComponent::DeinitCalibration() {
185 // see header file for class documentation
186
187   if(fCalibTime) delete fCalibTime; fCalibTime = NULL;
188   if(fCal)       delete fCal;       fCal       = NULL;
189   if(fSeedArray) delete fSeedArray; fSeedArray = NULL;
190
191   //if(fESDfriend) delete fESDfriend; fESDfriend = NULL;
192   
193   //if(arr) delete arr; arr = NULL;
194   
195   return 0;
196 }
197
198 int AliHLTTPCCalibTimeComponent::Reconfigure(const char* cdbEntry, const char* /*chainId*/){
199 // see header file for class documentation
200
201   // configure from the specified antry or the default one
202   const char* entry=cdbEntry;
203   if(!entry || entry[0]==0){
204      entry=fgkOCDBEntry;
205   }
206
207   return ConfigureFromCDBTObjString(entry);
208 }
209
210 Int_t AliHLTTPCCalibTimeComponent::ProcessCalibration( const AliHLTComponentEventData& /*evtData*/,  AliHLTComponentTriggerData& /*trigData*/ ){
211 // see header file for class documentation
212
213   if(GetFirstInputBlock( kAliHLTDataTypeSOR ) || GetFirstInputBlock( kAliHLTDataTypeEOR )) return 0;
214
215   const AliHLTComponentBlockData *iter = NULL;      
216
217   //--------------- output over TObjArray of AliTPCseed objects (output of TPCSeedMaker) -------------------//
218   
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
221   
222 //   for(iter = (TObject*)GetFirstInputObject(kAliHLTDataTypeTObjArray|kAliHLTDataOriginTPC); iter != NULL; iter = (TObject*)GetNextInputObject()){  
223 //               
224 //       if(GetDataType(iter) != (kAliHLTDataTypeTObjArray | kAliHLTDataOriginTPC)) continue;      
225 //       fSeedArray = dynamic_cast<TClonesArray*>(iter);      
226 //    }
227  
228  // int nInputClusters = 0;
229  // int nInputTracks = 0;
230
231   //TObjArray *arr = new TObjArray(1000);
232   //arr->SetOwner(kTRUE);
233   fSeedArray->Clear();
234
235   
236   for( int i=0; i<fkNPartition; i++ ){
237     delete[] fPartitionClusters[i];    
238     fPartitionClusters[i] = 0;
239     fNPartitionClusters[i] = 0;    
240   }
241
242   
243   
244   //------------------- loop over clusters -------------//
245   
246   for(iter = GetFirstInputBlock(AliHLTTPCDefinitions::fgkClustersDataType); iter != NULL; iter = GetNextInputBlock()){
247       
248       if ( iter->fDataType != AliHLTTPCDefinitions::fgkClustersDataType ) continue;
249     
250       Int_t slice     = AliHLTTPCDefinitions::GetMinSliceNr(iter->fSpecification);
251       Int_t partition = AliHLTTPCDefinitions::GetMinPatchNr(iter->fSpecification);
252        
253       Int_t slicepartition = slice*6+partition;
254       
255       if(slicepartition > fkNPartition){
256          HLTWarning("Wrong header of TPC cluster data, slice %d, partition %d", slice, partition );
257          continue;
258       }
259       
260       AliHLTTPCClusterData* inPtrSP = ( AliHLTTPCClusterData* )( iter->fPtr );
261       // nInputClusters += inPtrSP->fSpacePointCnt;
262
263       delete[] fPartitionClusters[slicepartition];
264       fPartitionClusters[slicepartition]  = new AliTPCclusterMI[inPtrSP->fSpacePointCnt];
265       fNPartitionClusters[slicepartition] = inPtrSP->fSpacePointCnt;
266     
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;
272           c->SetX(chlt->fX);
273           c->SetY(chlt->fY);
274           c->SetZ(chlt->fZ);
275           c->SetSigmaY2(chlt->fSigmaY2);
276           c->SetSigmaYZ( 0 );
277           c->SetSigmaZ2(chlt->fSigmaZ2);
278           c->SetQ( chlt->fCharge );
279           c->SetMax( chlt->fQMax );
280           Int_t sector, row;
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 );
285           c->SetRow( row );
286           c->SetPad( (Int_t) padtime[1] );
287           c->SetTimeBin( (Int_t) padtime[2] );
288       }      
289   } // end of loop over blocks of clusters
290   
291   
292   
293   
294   //---------- loop over merged tracks ------------------ //
295   int nTracks = 0;
296   for(const AliHLTComponentBlockData *pBlock = GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC); pBlock != NULL; pBlock = GetNextInputBlock()){
297  
298       AliHLTTracksData *dataPtr = (AliHLTTracksData*) pBlock->fPtr;
299       nTracks = dataPtr->fCount;
300     
301       AliHLTExternalTrackParam *currTrack = dataPtr->fTracklets;
302       
303       //nInputTracks += nTracks;
304     
305       for(int itr=0;  itr<nTracks && ( (AliHLTUInt8_t *)currTrack < ((AliHLTUInt8_t *) pBlock->fPtr)+pBlock->fSize); itr++){    
306           
307           // create an offline track
308           AliHLTGlobalBarrelTrack gb(*currTrack);
309           AliTPCseed tTPC;
310           tTPC.Set(gb.GetX(), gb.GetAlpha(), gb.GetParameter(), gb.GetCovariance());
311             
312           // set the cluster pointers     
313           for(UInt_t ic=0; ic<currTrack->fNPoints; ic++){       
314       
315               tTPC.SetNumberOfClusters(currTrack->fNPoints);
316           
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);    
321         
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);
324                  continue;
325               }
326         
327               AliTPCclusterMI *partitionClusters = fPartitionClusters[iSlice*6 + iPartition];
328               
329               if(!partitionClusters){
330                   HLTError("Clusters are missed for slice %d, partition %d", iSlice, iPartition );
331                   continue;
332               }
333         
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]);
336                  continue;
337               }
338         
339               AliTPCclusterMI *c = &(partitionClusters[iCluster]);              
340               int sec = c->GetDetector();
341               int row = c->GetRow();
342               if(sec >= 36) row = row + AliHLTTPCTransform::GetNRowLow();
343         
344               tTPC.SetClusterPointer(row, c);   
345         
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
353        
354       AliTPCseed *seed = &(tTPC);
355       fSeedArray->Add(seed);
356      
357       unsigned int step = sizeof( AliHLTExternalTrackParam ) + currTrack->fNPoints * sizeof( unsigned int );
358       currTrack = ( AliHLTExternalTrackParam* )( (( Byte_t * )currTrack) + step );  
359
360       }// end of vector track loop           
361   } // end of loop over blocks of merged tracks  
362   
363   HLTInfo("Number of reconstructed tracks %d, number of produced seeds %d\n", nTracks, fSeedArray->GetEntries());
364
365  
366   //----------- loop over output of global esd converter ----------------//
367   
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.
370  
371   TObject *iterOb = NULL; 
372   for(iterOb = (TObject*)GetFirstInputObject(kAliHLTDataTypeESDObject | kAliHLTDataOriginOut); iterOb != NULL; iterOb = (TObject*)GetNextInputObject()){   
373       
374       if(GetDataType(iterOb) != (kAliHLTDataTypeESDObject | kAliHLTDataOriginOut)) continue;
375             
376       fESDevent = dynamic_cast<AliESDEvent*>(iterOb);
377       fESDevent->GetStdContent();    
378                  
379       HLTInfo("Number of seeds: %i\n", fSeedArray->GetEntriesFast()); // access of the info from the previous loop over the AliTPCseed array        
380    
381       fCal->UpdateEventInfo(fESDevent);     
382       for(Int_t i=0; i<fSeedArray->GetEntriesFast(); i++){  // loop over TObjArray with seeds
383          
384           AliTPCseed *seed = (AliTPCseed*)fSeedArray->UncheckedAt(i);                    
385           fESDtrack = fESDevent->GetTrack(i);
386           if(!fESDtrack || !seed) continue; 
387
388           //if(fESDtrack->GetID() != seed->GetLabel()) { 
389           //   HLTWarning("Mismatch of track id between seed and ESD track: %i, %i\n", fESDtrack->GetID(), seed->GetLabel());
390           //   continue;            
391           //}
392           
393           if(seed->GetNumberOfClusters()==0) continue;      
394           fCal->RefitTrack(fESDtrack, seed, GetBz()); // update AliESDtrack and AliTPCseed info, acccording to Marian's request
395                
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)              
398         
399           //fESDfriendTrack = const_cast<AliESDfriendTrack*>(fESDtrack->GetFriendTrack());        
400       }
401   } 
402   
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);
407   if(!fCalibTime){
408      printf("fCalibTime pointer is NULL\n");
409      return 0;
410   }
411   
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
420
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
423   
424   // delete fESDfriend;
425   
426   //PushBack( (TObject*)fCalibTime, AliHLTTPCDefinitions::fgkCalibCEDataType | kAliHLTDataOriginOut, 0x0);
427   
428   return 0;
429 }
430
431 Int_t AliHLTTPCCalibTimeComponent::ShipDataToFXS( const AliHLTComponentEventData& /*evtData*/,  AliHLTComponentTriggerData& /*trigData*/ ){
432 // see header file for class documentation
433
434   HLTInfo("Shipping data to FXS...\n");
435  
436   fCalibTime->Analyze(); // called at the end of the run or event modulo
437   
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.
440   
441   THnSparse* addHist = fCalibTime->GetHistoDrift("all");
442   if(!addHist) return -1;
443
444   //Identifying used range of histogram
445  
446   Int_t startTimeBin = 0;
447   Int_t endTimeBin   = 0;
448
449   TH1D *histoTime = addHist->Projection(0);
450   if(histoTime){
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;
457   }
458
459   Int_t startPtBin = 0;
460   Int_t endPtBin   = 0;
461   TH1D *histoPt = addHist->Projection(1);
462   if(histoPt){
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;
469   }
470
471   Int_t startVdBin = 0;
472   Int_t endVdBin   = 0;
473   TH1D *histoVd = addHist->Projection(2);
474   if(histoVd){
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;
481   }
482
483   Int_t startRunBin = 0;
484   Int_t endRunBin   = 0;
485   TH1D *histoRun = addHist->Projection(3);
486   if(histoRun){
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;
493   }
494
495   TObjArray *vdriftArray = new TObjArray();
496   if(!vdriftArray) return -2;
497
498   TObjArray *array = fCalibTime->GetHistoDrift();
499   if(!array) return -3;
500
501   TIterator *iterator = array->MakeIterator();
502   if(!iterator) return -4;
503
504   iterator->Reset();
505   THnSparse *hist = NULL;
506   while((hist = (THnSparseF*)iterator->Next())){
507        
508          //if(!hist) continue;
509          hist->Print();
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);
514          
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);
520          
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();
528          graphName+=name;
529          graphName.ToUpper();
530          graph->SetName(graphName);
531          printf("name = %s\n", graphName.Data());
532          vdriftArray->Add(graph);
533
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);
539          //fit->SplineFit(0);
540          //TString fiName=fit->ClassName();
541          //fiName+=type;
542          //fiName+=trigger;
543          //fiName.ToUpper();
544          //fit->SetName(fiName.Data());
545          //printf("name=%s\n", fiName.Data());
546          //vdriftArray->Add(fit);
547   }
548   
549   THnSparse    *laserHist  = NULL;
550   TGraphErrors *laserGraph = NULL;
551   TString laserName = "";
552
553   //Histograms and graphs for A side lasers
554   laserHist = fCalibTime->GetHistVdriftLaserA(1);
555   if(laserHist){
556     
557      laserName=laserHist->ClassName();
558      laserName+="_MEAN_DRIFT_LASER_ALL_A";
559      laserName.ToUpper();
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";
566        laserName.ToUpper();
567        laserGraph->SetName(laserName);
568        vdriftArray->Add(laserGraph);
569      }
570   }
571
572   //Histograms and graphs for C side lasers
573   laserHist=fCalibTime->GetHistVdriftLaserC(1);
574   if(laserHist){
575      laserName=laserHist->ClassName();
576      laserName+="_MEAN_DRIFT_LASER_ALL_C";
577      laserName.ToUpper();
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";
584        laserName.ToUpper();
585        laserGraph->SetName(laserName);
586        vdriftArray->Add(laserGraph);
587      }
588   }
589
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);
603
604   static AliHLTReadoutList rdList(AliHLTReadoutList::kTPC);
605   
606
607   TFile *file = TFile::Open("vdrift.root", "RECREATE");
608   vdriftArray->Write();
609   file->Close();
610   delete file;
611
612   file = TFile::Open("calibTime.root", "RECREATE");
613   fCalibTime->Write();
614   file->Close();
615   delete file;
616  
617   // the vdriftArray is pushed to the HLT-FXSsubscriber 
618   PushToFXS( (TObject*)vdriftArray, "TPC", "TIMEDRIFT", &rdList );
619  
620   //Should array be deleted now?
621   //  if(vdriftArray){
622   //      vdriftArray.Clear();
623   //    delete vdriftArray;
624   //    vdriftArray=0;
625   //  }
626   
627   return 0;
628
629