]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/ITS/AliHLTITSVertexerSPDComponent.cxx
compilation error in debug mode fixed
[u/mrichter/AliRoot.git] / HLT / ITS / AliHLTITSVertexerSPDComponent.cxx
1 // $Id: AliHLTITSVertexerSPDComponent.cxx 32659 2009-06-02 16:08:40Z sgorbuno $
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: Sergey Gorbunov <sergey.gorbunov@kip.uni-heidelberg.de> *
7 //                  Ivan Kisel <kisel@kip.uni-heidelberg.de>                *
8 //                  for The ALICE HLT Project.                              *
9 //                                                                          *
10 // Permission to use, copy, modify and distribute this software and its     *
11 // documentation strictly for non-commercial purposes is hereby granted     *
12 // without fee, provided that the above copyright notice appears in all     *
13 // copies and that both the copyright notice and this permission notice     *
14 // appear in the supporting documentation. The authors make no claims       *
15 // about the suitability of this software for any purpose. It is            *
16 // provided "as is" without express or implied warranty.                    *
17 //                                                                          *
18 //***************************************************************************
19
20 ///  @file   AliHLTITSVertexerSPDComponent.cxx
21 ///  @author Sergey Gorbunov <sergey.gorbunov@kip.uni-heidelberg.de>
22 ///  @date   June 2009
23 ///  @brief  An ITS tracker processing component for the HLT
24
25
26 /////////////////////////////////////////////////////
27 //                                                 //
28 // a ITS tracker processing component for the HLT  //
29 //                                                 //
30 /////////////////////////////////////////////////////
31
32 #if __GNUC__>= 3
33 using namespace std;
34 #endif
35
36 #include "AliHLTITSVertexerSPDComponent.h"
37 #include "AliHLTArray.h"
38 #include "AliExternalTrackParam.h"
39 #include "TStopwatch.h"
40 #include "TMath.h"
41 #include "AliCDBEntry.h"
42 #include "AliCDBManager.h"
43 #include "AliGeomManager.h"
44 #include "TObjString.h"
45 #include "TObjArray.h"
46 #include "AliITStrackerHLT.h"
47 #include "AliHLTITSSpacePointData.h"
48 #include "AliHLTITSClusterDataFormat.h"
49 #include "AliHLTDataTypes.h"
50 #include "AliHLTExternalTrackParam.h"
51 #include "AliHLTGlobalBarrelTrack.h"
52 #include "AliGeomManager.h"
53 #include "TH1F.h"
54 #include "TH2F.h"
55 #include "AliESDVertex.h"
56
57 /** ROOT macro for the implementation of ROOT specific class methods */
58 ClassImp( AliHLTITSVertexerSPDComponent )
59 AliHLTITSVertexerSPDComponent::AliHLTITSVertexerSPDComponent()
60     :
61     fSolenoidBz( 0 ),
62     fProduceHistos(1),
63     fAutoCalibration(1000),
64     fFullTime( 0 ),
65     fRecoTime( 0 ),
66     fNEvents( 0 ),    
67     fHistoVertexXY(0),
68     fHistoVertexX(0),
69     fHistoVertexY(0),
70   fHistoVertexZ(0)
71
72 {
73   // see header file for class documentation
74   // or
75   // refer to README to build package
76   // or
77   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
78
79   fDefRunVtx[0] = 0;
80   fDefRunVtx[1] = 0;
81   fDefRunVtx[2] = 0;
82 }
83
84 AliHLTITSVertexerSPDComponent::AliHLTITSVertexerSPDComponent( const AliHLTITSVertexerSPDComponent& )
85     :
86     AliHLTProcessor(),
87     fSolenoidBz( 0 ),
88     fProduceHistos(1),
89     fAutoCalibration(1000),
90     fFullTime( 0 ),
91     fRecoTime( 0 ),
92     fNEvents( 0 ),    
93     fHistoVertexXY(0),
94     fHistoVertexX(0),
95     fHistoVertexY(0),
96   fHistoVertexZ(0)
97 {
98   // see header file for class documentation
99   HLTFatal( "copy constructor untested" );
100 }
101
102 AliHLTITSVertexerSPDComponent& AliHLTITSVertexerSPDComponent::operator=( const AliHLTITSVertexerSPDComponent& )
103 {
104   // see header file for class documentation
105   HLTFatal( "assignment operator untested" );
106   return *this;
107 }
108
109 AliHLTITSVertexerSPDComponent::~AliHLTITSVertexerSPDComponent()
110 {
111   // see header file for class documentation
112   
113   delete fHistoVertexXY;
114   delete fHistoVertexX;
115   delete fHistoVertexY;
116   delete fHistoVertexZ;
117 }
118
119 //
120 // Public functions to implement AliHLTComponent's interface.
121 // These functions are required for the registration process
122 //
123
124 const char* AliHLTITSVertexerSPDComponent::GetComponentID()
125 {
126   // see header file for class documentation
127   return "ITSVertexerSPD";
128 }
129
130 void AliHLTITSVertexerSPDComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list )
131 {
132   // see header file for class documentation
133   list.clear();
134   list.push_back( kAliHLTDataTypeClusters|kAliHLTDataOriginITSSPD );
135 }
136
137 AliHLTComponentDataType AliHLTITSVertexerSPDComponent::GetOutputDataType()
138 {
139   // see header file for class documentation  
140   return kAliHLTMultipleDataType;
141 }
142
143 int AliHLTITSVertexerSPDComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList)
144
145 {
146   // see header file for class documentation
147   tgtList.clear();
148   tgtList.push_back( kAliHLTDataTypeESDVertex|kAliHLTDataOriginITS);
149   tgtList.push_back(kAliHLTDataTypeHistogram);
150   return tgtList.size();
151 }
152
153 void AliHLTITSVertexerSPDComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
154 {
155   // define guess for the output data size
156   constBase = 10000;       // minimum size
157   inputMultiplier = 0.5; // size relative to input
158 }
159
160 AliHLTComponent* AliHLTITSVertexerSPDComponent::Spawn()
161 {
162   // see header file for class documentation
163   return new AliHLTITSVertexerSPDComponent;
164 }
165
166 void AliHLTITSVertexerSPDComponent::SetDefaultConfiguration()
167 {
168   // Set default configuration for the CA tracker component
169   // Some parameters can be later overwritten from the OCDB
170
171   fSolenoidBz = -5.00668;
172   fProduceHistos=1;
173   fAutoCalibration = 1000;
174   fDefRunVtx[0] = 0;
175   fDefRunVtx[1] = 0;
176   fDefRunVtx[2] = 0;
177   fFullTime = 0;
178   fRecoTime = 0;
179   fNEvents = 0;
180 }
181
182 int AliHLTITSVertexerSPDComponent::ReadConfigurationString(  const char* arguments )
183 {
184   // Set configuration parameters for the CA tracker component from the string
185
186   int iResult = 0;
187   if ( !arguments ) return iResult;
188
189   TString allArgs = arguments;
190   TString argument;
191   int bMissingParam = 0;
192
193   TObjArray* pTokens = allArgs.Tokenize( " " );
194
195   int nArgs =  pTokens ? pTokens->GetEntries() : 0;
196
197   for ( int i = 0; i < nArgs; i++ ) {
198     argument = ( ( TObjString* )pTokens->At( i ) )->GetString();
199     if ( argument.IsNull() ) continue;
200
201     if ( argument.CompareTo( "-solenoidBz" ) == 0 ) {
202       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
203       HLTWarning("argument -solenoidBz is deprecated, magnetic field set up globally (%f)", GetBz());
204       continue;
205     }
206
207     if ( argument.CompareTo( "-produceHistos" ) == 0 ) {
208       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
209       fProduceHistos = ( ( TObjString* )pTokens->At( i ) )->GetString().Atoi();
210       HLTInfo( "fProduceHistos set to: %d", fProduceHistos );
211       continue;
212     }
213
214     if ( argument.CompareTo( "-runVertex" ) == 0 ) {
215       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
216       fDefRunVtx[0] = ( ( TObjString* )pTokens->At( i ) )->GetString().Atof();
217       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
218       fDefRunVtx[1] = ( ( TObjString* )pTokens->At( i ) )->GetString().Atof();
219       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
220       fDefRunVtx[2] = ( ( TObjString* )pTokens->At( i ) )->GetString().Atof();
221       HLTInfo( "Default run vertex is set to (%f,%f,%f)",fDefRunVtx[0],
222                fDefRunVtx[1],fDefRunVtx[2] );
223       continue;
224     }
225
226     if ( argument.CompareTo( "-beamDiamondCalibration" ) == 0 ) {
227       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
228       fAutoCalibration = ( ( TObjString* )pTokens->At( i ) )->GetString().Atoi();
229       HLTInfo( "N events for recalibration of the run vertex is set to: %d", fAutoCalibration );
230       continue;
231     }
232     
233
234     HLTError( "Unknown option \"%s\"", argument.Data() );
235     iResult = -EINVAL;
236   }
237   delete pTokens;
238
239   if ( bMissingParam ) {
240     HLTError( "Specifier missed for parameter \"%s\"", argument.Data() );
241     iResult = -EINVAL;
242   }
243
244   return iResult;
245 }
246
247
248 int AliHLTITSVertexerSPDComponent::ReadCDBEntry( const char* cdbEntry, const char* chainId )
249 {
250   // see header file for class documentation
251
252   const char* defaultNotify = "";
253
254   if ( !cdbEntry ) {
255     return 0;// need to add the HLT/ConfigITS/ITSTracker directory to cdb SG!!!
256     cdbEntry = "HLT/ConfigITS/ITSTracker";
257     defaultNotify = " (default)";
258     chainId = 0;
259   }
260
261   HLTInfo( "configure from entry \"%s\"%s, chain id %s", cdbEntry, defaultNotify, ( chainId != NULL && chainId[0] != 0 ) ? chainId : "<none>" );
262   AliCDBEntry *pEntry = AliCDBManager::Instance()->Get( cdbEntry );//,GetRunNo());
263
264   if ( !pEntry ) {
265     HLTError( "cannot fetch object \"%s\" from CDB", cdbEntry );
266     return -EINVAL;
267   }
268
269   TObjString* pString = dynamic_cast<TObjString*>( pEntry->GetObject() );
270
271   if ( !pString ) {
272     HLTError( "configuration object \"%s\" has wrong type, required TObjString", cdbEntry );
273     return -EINVAL;
274   }
275
276   HLTInfo( "received configuration object string: \"%s\"", pString->GetString().Data() );
277
278   return  ReadConfigurationString( pString->GetString().Data() );
279 }
280
281
282 int AliHLTITSVertexerSPDComponent::Configure( const char* cdbEntry, const char* chainId, const char *commandLine )
283 {
284   // Configure the component
285   // There are few levels of configuration,
286   // parameters which are set on one step can be overwritten on the next step
287
288   //* read hard-coded values
289
290   SetDefaultConfiguration();
291
292   //* read the default CDB entry
293
294   int iResult1 = ReadCDBEntry( NULL, chainId );
295
296   //* read magnetic field
297
298   int iResult2 = 0; //ReadCDBEntry( kAliHLTCDBSolenoidBz, chainId );
299   fSolenoidBz = GetBz();
300
301   //* read the actual CDB entry if required
302
303   int iResult3 = ( cdbEntry ) ? ReadCDBEntry( cdbEntry, chainId ) : 0;
304
305   //* read extra parameters from input (if they are)
306
307   int iResult4 = 0;
308
309   if ( commandLine && commandLine[0] != '\0' ) {
310     HLTInfo( "received configuration string from HLT framework: \"%s\"", commandLine );
311     iResult4 = ReadConfigurationString( commandLine );
312   }
313
314   // Initialise the tracker here
315
316   return iResult1 ? iResult1 : ( iResult2 ? iResult2 : ( iResult3 ? iResult3 : iResult4 ) );
317 }
318
319
320
321 int AliHLTITSVertexerSPDComponent::DoInit( int argc, const char** argv )
322 {
323   // Configure the ITS tracker component
324   
325   fProduceHistos = 1;
326   fAutoCalibration = 1000;
327
328   // Check field
329   if (!TGeoGlobalMagField::Instance()) {
330     HLTError("magnetic field not initialized, please set up TGeoGlobalMagField and AliMagF");
331     return -ENODEV;
332   }
333   fSolenoidBz=GetBz();
334
335   if(AliGeomManager::GetGeometry()==NULL){
336     AliGeomManager::LoadGeometry("");
337   }
338   AliGeomManager::ApplyAlignObjsFromCDB("ITS");
339   
340   TString arguments = "";
341   for ( int i = 0; i < argc; i++ ) {
342     if ( !arguments.IsNull() ) arguments += " ";
343     arguments += argv[i];
344   }
345
346   int ret = Configure( NULL, NULL, arguments.Data() );
347
348   if( fProduceHistos ){
349     fHistoVertexXY = new TH2F("hITSvertexXY", "ITSSPD vertex in XY", 100,-2,2,100,-2,2);
350     fHistoVertexX = new TH1F("hITSvertexX", "ITSSPD vertex X", 100,-2,2);
351     fHistoVertexY = new TH1F("hITSvertexY", "ITSSPD vertex Y", 100,-2,2);
352     fHistoVertexZ = new TH1F("hITSvertexZ", "ITSSPD vertex Z", 100,-15,15);
353   }
354
355   for( int i=0; i<3; i++){
356     fRunVtx[i] = fDefRunVtx[i];
357     fRunVtxNew[i] = fDefRunVtx[i];
358   }
359   fRunVtx[3] = 0.;
360   fRunVtxNew[3] = 0.;
361
362   return ret;
363 }
364
365
366 int AliHLTITSVertexerSPDComponent::DoDeinit()
367 {
368   // see header file for class documentation
369
370   delete fHistoVertexXY;
371   delete fHistoVertexX;
372   delete fHistoVertexY;
373   delete fHistoVertexZ;
374
375   fHistoVertexXY = 0;
376   fHistoVertexX = 0;
377   fHistoVertexY = 0;
378   fHistoVertexZ = 0;
379   fAutoCalibration = 1000;
380   fProduceHistos = 1;
381   return 0;
382 }
383
384
385
386 int AliHLTITSVertexerSPDComponent::Reconfigure( const char* cdbEntry, const char* chainId )
387 {
388   // Reconfigure the component from OCDB .
389
390   return Configure( cdbEntry, chainId, NULL );
391 }
392
393
394 int AliHLTITSVertexerSPDComponent::DoEvent
395 (
396   const AliHLTComponentEventData& evtData,
397   const AliHLTComponentBlockData* blocks,
398   AliHLTComponentTriggerData& /*trigData*/,
399   AliHLTUInt8_t* /*outputPtr*/,
400   AliHLTUInt32_t& size,
401   vector<AliHLTComponentBlockData>& /*outputBlocks*/ )
402 {
403   //* process event
404
405   //AliHLTUInt32_t maxBufferSize = size;
406   size = 0; // output size
407
408   if (!IsDataEvent()) return 0;
409
410   if ( evtData.fBlockCnt <= 0 ) {
411     HLTWarning( "no blocks in event" );
412     return 0;
413   }
414
415
416   TStopwatch timer;
417
418   // Event reconstruction in ITS
419
420   int iResult=0;
421
422   int nBlocks = evtData.fBlockCnt;
423   int nClustersTotal = 0;
424  
425
426   const int kNPhiBins = 20;
427   const double kDPhi = TMath::TwoPi() / kNPhiBins;
428   double vtxX = fRunVtx[0], vtxY = fRunVtx[1], vtxZ = fRunVtx[2];
429
430
431   std::vector<AliHLTITSVertexerSPDComponent::AliHLTITSVZCluster> clusters[2][ kNPhiBins ];
432
433   for (int ndx=0; ndx<nBlocks && iResult>=0; ndx++) {
434
435     const AliHLTComponentBlockData* iter = blocks+ndx;
436  
437     // Read ITS SPD clusters
438
439     if ( iter->fDataType == (kAliHLTDataTypeClusters|kAliHLTDataOriginITSSPD) ){
440
441       AliHLTITSClusterData *inPtr=reinterpret_cast<AliHLTITSClusterData*>( iter->fPtr );
442       int nClusters = inPtr->fSpacePointCnt;
443       nClustersTotal+=nClusters;
444       for( int icl=0; icl<nClusters; icl++ ){   
445         AliHLTITSSpacePointData &d = inPtr->fSpacePoints[icl];
446         if( d.fLayer>1 ) continue;// SPD only;
447         if( d.fLayer<0 ) continue;// SPD only;
448         Int_t lab[4] = { d.fTracks[0], d.fTracks[1], d.fTracks[2], d.fIndex };
449         Int_t info[3] = { d.fNy, d.fNz, d.fLayer };
450         Float_t hit[6] = { d.fY, d.fZ, d.fSigmaY2, d.fSigmaZ2, d.fQ, d.fSigmaYZ };
451         if( d.fLayer==4 ) hit[5] = -hit[5];
452
453
454         AliITSRecPoint p( lab, hit, info );
455
456         if (!p.Misalign()){
457           HLTWarning("Can not misalign an SPD cluster");
458         }
459         Float_t xyz[3];      
460         p.GetGlobalXYZ(xyz);      
461
462         // get phi bin
463         double phi = TMath::ATan2( xyz[1]-vtxY, xyz[0]-vtxX );
464         if( phi<0 ) phi+=TMath::TwoPi();
465         int iphi = (int ) phi/kDPhi;
466         if( iphi<0 ) iphi = 0;
467         if( iphi>=kNPhiBins ) iphi = kNPhiBins-1;
468         AliHLTITSVZCluster c;
469         c.fX = xyz[0];
470         c.fY = xyz[1];
471         c.fZ = xyz[2];
472         clusters[d.fLayer][iphi].push_back( c );        
473      }   
474     }    
475   }// end read input blocks
476   
477
478   // Reconstruct the vertex
479
480   TStopwatch timerReco;
481     
482
483   const double kZBinMin = -10;
484   const int kNZBins = 40;
485   const double kZStep = 1./kNZBins;
486
487   double histX[kNZBins];
488   double histY[kNZBins];
489   double histZ[kNZBins];
490   double histW[kNZBins];
491   int histN[kNZBins];
492   
493   // fit few times with decreasing of the radius cut
494
495   double vtxDR2[2] = { 2.*2., .2*.2};
496   double  maxW = 0;  
497   int bestBin = -1;
498   
499   for( int iter = 0; iter<2; iter++ ){
500     
501     for( int i=0; i<kNZBins; i++ ){
502       histX[i] = 0;
503       histY[i] = 0;
504       histZ[i] = 0;
505       histW[i] = 0;
506       histN[i] = 0;
507     }
508  
509     for( int iPhi=0; iPhi<kNPhiBins; iPhi++ ){
510       int nCluUp = clusters[1][iPhi].size();    
511       int nCluDn = clusters[0][iPhi].size();    
512       for( int icUp=0; icUp<nCluUp; icUp++ ){
513         AliHLTITSVZCluster &cu = clusters[1][iPhi][icUp];
514         double x0 = cu.fX - vtxX;
515         double y0 = cu.fY - vtxY;
516         double z0 = cu.fZ - vtxZ;
517         double bestR2 = 1.e10;
518         int bestDn=-1, bestBinDn =0;
519         double bestV[3]={0,0,0};
520         for( int icDn=0; icDn<nCluDn; icDn++ ){
521           AliHLTITSVZCluster &cd = clusters[0][iPhi][icDn];
522           double x1 = cd.fX - vtxX;
523           double y1 = cd.fY - vtxY;
524           double z1 = cd.fZ - vtxZ;
525           double dx = x1 - x0;
526           double dy = y1 - y0;
527           double l2 = 1./(dx*dx + dy*dy);
528           double a = x1*y0 - x0*y1;
529           double r2 = a*a*l2; 
530           if( r2>vtxDR2[iter] ) continue;
531           double xv = -dy*a*l2;
532           double yv =  dx*a*l2;
533           double zv = ( (x1*z0-x0*z1)*dx + (y1*z0-y0*z1)*dy )*l2;
534           int zbin = (int)((zv - kZBinMin)*kZStep);
535           if( zbin<0 || zbin>=kNZBins ) continue;
536           if( r2<bestR2 ){
537             bestR2 = r2;
538             bestDn = icDn;
539             bestV[0] = xv;
540             bestV[1] = yv;
541             bestV[2] = zv;
542             bestBinDn = zbin;
543           }
544         }
545         if( bestDn>=0 ){
546           double w = 1./(1.+bestR2);
547           histX[bestBinDn]+=bestV[0]*w;
548           histY[bestBinDn]+=bestV[1]*w;
549           histZ[bestBinDn]+=bestV[2]*w;
550           histW[bestBinDn]+=w;
551           histN[bestBinDn]+=1;
552         }
553       }
554     }
555     
556     maxW = 0;  
557     bestBin = -1;
558     for( int i=0; i<kNZBins; i++ ){
559       if( histW[i]>maxW ){
560         bestBin = i;
561         maxW = histW[i];
562       }
563     }
564     if( bestBin<0 || histN[bestBin] <3 ){
565       bestBin = -1;
566       break;
567     }
568     vtxX +=histX[bestBin]/maxW;
569     vtxY +=histY[bestBin]/maxW;
570     vtxZ +=histZ[bestBin]/maxW;    
571   }
572
573   if( bestBin>=0 ){
574     double nv = 1;
575     double nNew = fRunVtxNew[3] + nv;
576     double v[3] = {vtxX, vtxY, vtxZ};
577     for( int i=0; i<3; i++){
578       fRunVtxNew[i] = ( fRunVtxNew[i]*fRunVtxNew[3] + v[i]*nv )/nNew;
579     }
580     fRunVtxNew[3] = nNew;
581   }  
582   
583   if( fAutoCalibration>0 && fRunVtxNew[3] >= fAutoCalibration ){
584     for( int i=0; i<4; i++ ){
585       fRunVtx[i] = fRunVtxNew[i];
586       fRunVtxNew[i] = 0;
587     }
588     //cout<<"ITSVertexerSPD: set run vtx to "<<fRunVtx[0]<<" "<<fRunVtx[1]<<" "<<fRunVtx[2]<<endl;  
589     if( fRunVtx[0]>3. ) fRunVtx[0] = 3;
590     if( fRunVtx[0]<-3. ) fRunVtx[0] = -3;
591     if( fRunVtx[1]>3. ) fRunVtx[1] = 3;
592     if( fRunVtx[1]<-3. ) fRunVtx[1] = -3;
593     if( fRunVtx[2]>30. ) fRunVtx[2] = 30;
594     if( fRunVtx[2]<30. ) fRunVtx[2] = -30;
595   }
596
597   timerReco.Stop();
598   
599   // Fill output 
600
601   if( bestBin>=0 ){
602     double pos[3] = {vtxX, vtxY, vtxZ};
603     double s = 400.E-4;
604     double cov[6] = {s*s,0,s*s,0,0,s*s};
605     AliESDVertex v(pos, cov, 0, histN[bestBin]);
606     PushBack( (TObject*) &v, kAliHLTDataTypeESDVertex|kAliHLTDataOriginITS,0 );
607
608     //cout<<"ITSVertexerSPD: vtx found "<<vtxX<<" "<<vtxY<<" "<<vtxZ<<endl;
609
610     if( fHistoVertexXY ) fHistoVertexXY->Fill( vtxX, vtxY );
611     if( fHistoVertexX ) fHistoVertexX->Fill( vtxX );
612     if( fHistoVertexY ) fHistoVertexY->Fill( vtxY );
613     if( fHistoVertexZ ) fHistoVertexZ->Fill( vtxZ );
614   }
615
616   if( fHistoVertexXY ) PushBack( (TObject*) fHistoVertexXY, kAliHLTDataTypeHistogram,0);
617   if( fHistoVertexX ) PushBack( (TObject*) fHistoVertexX, kAliHLTDataTypeHistogram,0);
618   if( fHistoVertexY ) PushBack( (TObject*) fHistoVertexY, kAliHLTDataTypeHistogram,0);
619   if( fHistoVertexZ ) PushBack( (TObject*) fHistoVertexZ, kAliHLTDataTypeHistogram,0);
620
621
622   timer.Stop();
623   fFullTime += timer.RealTime();
624   fRecoTime += timerReco.RealTime();
625   fNEvents++;
626
627   // Set log level to "Warning" for on-line system monitoring
628   int hz = ( int ) ( fFullTime > 1.e-10 ? fNEvents / fFullTime : 100000 );
629   int hz1 = ( int ) ( fRecoTime > 1.e-10 ? fNEvents / fRecoTime : 100000 );
630   HLTInfo( "ITS Z Vertexer: input %d clusters; time: full %d / reco %d Hz",
631            nClustersTotal, hz, hz1 );
632
633   return iResult;
634 }