]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/ITS/AliHLTITSVertexerSPDComponent.cxx
Compilation warnings 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       fSolenoidBz = ( ( TObjString* )pTokens->At( i ) )->GetString().Atof();
204       HLTInfo( "Magnetic Field set to: %f", fSolenoidBz );
205       continue;
206     }
207
208     if ( argument.CompareTo( "-produceHistos" ) == 0 ) {
209       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
210       fProduceHistos = ( ( TObjString* )pTokens->At( i ) )->GetString().Atoi();
211       HLTInfo( "fProduceHistos set to: %d", fProduceHistos );
212       continue;
213     }
214
215     if ( argument.CompareTo( "-runVertex" ) == 0 ) {
216       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
217       fDefRunVtx[0] = ( ( TObjString* )pTokens->At( i ) )->GetString().Atof();
218       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
219       fDefRunVtx[1] = ( ( TObjString* )pTokens->At( i ) )->GetString().Atof();
220       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
221       fDefRunVtx[2] = ( ( TObjString* )pTokens->At( i ) )->GetString().Atof();
222       HLTInfo( "Default run vertex is set to (%f,%f,%f)",fDefRunVtx[0],
223                fDefRunVtx[1],fDefRunVtx[2] );
224       continue;
225     }
226
227     if ( argument.CompareTo( "-beamDiamondCalibration" ) == 0 ) {
228       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
229       fAutoCalibration = ( ( TObjString* )pTokens->At( i ) )->GetString().Atoi();
230       HLTInfo( "N events for recalibration of the run vertex is set to: %d", fAutoCalibration );
231       continue;
232     }
233     
234
235     HLTError( "Unknown option \"%s\"", argument.Data() );
236     iResult = -EINVAL;
237   }
238   delete pTokens;
239
240   if ( bMissingParam ) {
241     HLTError( "Specifier missed for parameter \"%s\"", argument.Data() );
242     iResult = -EINVAL;
243   }
244
245   return iResult;
246 }
247
248
249 int AliHLTITSVertexerSPDComponent::ReadCDBEntry( const char* cdbEntry, const char* chainId )
250 {
251   // see header file for class documentation
252
253   const char* defaultNotify = "";
254
255   if ( !cdbEntry ) {
256     return 0;// need to add the HLT/ConfigITS/ITSTracker directory to cdb SG!!!
257     cdbEntry = "HLT/ConfigITS/ITSTracker";
258     defaultNotify = " (default)";
259     chainId = 0;
260   }
261
262   HLTInfo( "configure from entry \"%s\"%s, chain id %s", cdbEntry, defaultNotify, ( chainId != NULL && chainId[0] != 0 ) ? chainId : "<none>" );
263   AliCDBEntry *pEntry = AliCDBManager::Instance()->Get( cdbEntry );//,GetRunNo());
264
265   if ( !pEntry ) {
266     HLTError( "cannot fetch object \"%s\" from CDB", cdbEntry );
267     return -EINVAL;
268   }
269
270   TObjString* pString = dynamic_cast<TObjString*>( pEntry->GetObject() );
271
272   if ( !pString ) {
273     HLTError( "configuration object \"%s\" has wrong type, required TObjString", cdbEntry );
274     return -EINVAL;
275   }
276
277   HLTInfo( "received configuration object string: \"%s\"", pString->GetString().Data() );
278
279   return  ReadConfigurationString( pString->GetString().Data() );
280 }
281
282
283 int AliHLTITSVertexerSPDComponent::Configure( const char* cdbEntry, const char* chainId, const char *commandLine )
284 {
285   // Configure the component
286   // There are few levels of configuration,
287   // parameters which are set on one step can be overwritten on the next step
288
289   //* read hard-coded values
290
291   SetDefaultConfiguration();
292
293   //* read the default CDB entry
294
295   int iResult1 = ReadCDBEntry( NULL, chainId );
296
297   //* read magnetic field
298
299   int iResult2 = ReadCDBEntry( kAliHLTCDBSolenoidBz, chainId );
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   if(AliGeomManager::GetGeometry()==NULL){
329     AliGeomManager::LoadGeometry("");
330   }
331   AliGeomManager::ApplyAlignObjsFromCDB("ITS");
332   
333   TString arguments = "";
334   for ( int i = 0; i < argc; i++ ) {
335     if ( !arguments.IsNull() ) arguments += " ";
336     arguments += argv[i];
337   }
338
339   int ret = Configure( NULL, NULL, arguments.Data() );
340
341   if( fProduceHistos ){
342     fHistoVertexXY = new TH2F("hITSvertexXY", "ITSSPD vertex in XY", 100,-2,2,100,-2,2);
343     fHistoVertexX = new TH1F("hITSvertexX", "ITSSPD vertex X", 100,-2,2);
344     fHistoVertexY = new TH1F("hITSvertexY", "ITSSPD vertex Y", 100,-2,2);
345     fHistoVertexZ = new TH1F("hITSvertexZ", "ITSSPD vertex Z", 100,-15,15);
346   }
347
348   for( int i=0; i<3; i++){
349     fRunVtx[i] = fDefRunVtx[i];
350     fRunVtxNew[i] = fDefRunVtx[i];
351   }
352   fRunVtx[3] = 0.;
353   fRunVtxNew[3] = 0.;
354
355   return ret;
356 }
357
358
359 int AliHLTITSVertexerSPDComponent::DoDeinit()
360 {
361   // see header file for class documentation
362
363   delete fHistoVertexXY;
364   delete fHistoVertexX;
365   delete fHistoVertexY;
366   delete fHistoVertexZ;
367
368   fHistoVertexXY = 0;
369   fHistoVertexX = 0;
370   fHistoVertexY = 0;
371   fHistoVertexZ = 0;
372   fAutoCalibration = 1000;
373   fProduceHistos = 1;
374   return 0;
375 }
376
377
378
379 int AliHLTITSVertexerSPDComponent::Reconfigure( const char* cdbEntry, const char* chainId )
380 {
381   // Reconfigure the component from OCDB .
382
383   return Configure( cdbEntry, chainId, NULL );
384 }
385
386
387 int AliHLTITSVertexerSPDComponent::DoEvent
388 (
389   const AliHLTComponentEventData& evtData,
390   const AliHLTComponentBlockData* blocks,
391   AliHLTComponentTriggerData& /*trigData*/,
392   AliHLTUInt8_t* /*outputPtr*/,
393   AliHLTUInt32_t& size,
394   vector<AliHLTComponentBlockData>& /*outputBlocks*/ )
395 {
396   //* process event
397
398   //AliHLTUInt32_t maxBufferSize = size;
399   size = 0; // output size
400
401   if (!IsDataEvent()) return 0;
402
403   if ( evtData.fBlockCnt <= 0 ) {
404     HLTWarning( "no blocks in event" );
405     return 0;
406   }
407
408
409   TStopwatch timer;
410
411   // Event reconstruction in ITS
412
413   int iResult=0;
414
415   int nBlocks = evtData.fBlockCnt;
416   int nClustersTotal = 0;
417  
418
419   const int kNPhiBins = 20;
420   const double kDPhi = TMath::TwoPi() / kNPhiBins;
421   double vtxX = fRunVtx[0], vtxY = fRunVtx[1], vtxZ = fRunVtx[2];
422
423
424   std::vector<AliHLTITSVertexerSPDComponent::AliHLTITSVZCluster> clusters[2][ kNPhiBins ];
425
426   for (int ndx=0; ndx<nBlocks && iResult>=0; ndx++) {
427
428     const AliHLTComponentBlockData* iter = blocks+ndx;
429  
430     // Read ITS SPD clusters
431
432     if ( iter->fDataType == (kAliHLTDataTypeClusters|kAliHLTDataOriginITSSPD) ){
433
434       AliHLTITSClusterData *inPtr=reinterpret_cast<AliHLTITSClusterData*>( iter->fPtr );
435       int nClusters = inPtr->fSpacePointCnt;
436       nClustersTotal+=nClusters;
437       for( int icl=0; icl<nClusters; icl++ ){   
438         AliHLTITSSpacePointData &d = inPtr->fSpacePoints[icl];
439         if( d.fLayer>1 ) continue;// SPD only;
440         if( d.fLayer<0 ) continue;// SPD only;
441         Int_t lab[4] = { d.fTracks[0], d.fTracks[1], d.fTracks[2], d.fIndex };
442         Int_t info[3] = { d.fNy, d.fNz, d.fLayer };
443         Float_t hit[6] = { d.fY, d.fZ, d.fSigmaY2, d.fSigmaZ2, d.fQ, d.fSigmaYZ };
444         if( d.fLayer==4 ) hit[5] = -hit[5];
445
446
447         AliITSRecPoint p( lab, hit, info );
448
449         if (!p.Misalign()){
450           HLTWarning("Can not misalign an SPD cluster");
451         }
452         Float_t xyz[3];      
453         p.GetGlobalXYZ(xyz);      
454
455         // get phi bin
456         double phi = TMath::ATan2( xyz[1]-vtxY, xyz[0]-vtxX );
457         if( phi<0 ) phi+=TMath::TwoPi();
458         int iphi = (int ) phi/kDPhi;
459         if( iphi<0 ) iphi = 0;
460         if( iphi>=kNPhiBins ) iphi = kNPhiBins-1;
461         AliHLTITSVZCluster c;
462         c.fX = xyz[0];
463         c.fY = xyz[1];
464         c.fZ = xyz[2];
465         clusters[d.fLayer][iphi].push_back( c );        
466      }   
467     }    
468   }// end read input blocks
469   
470
471   // Reconstruct the vertex
472
473   TStopwatch timerReco;
474     
475
476   const double kZBinMin = -10;
477   const int kNZBins = 40;
478   const double kZStep = 1./kNZBins;
479
480   double histX[kNZBins];
481   double histY[kNZBins];
482   double histZ[kNZBins];
483   double histW[kNZBins];
484   int histN[kNZBins];
485   
486   // fit few times with decreasing of the radius cut
487
488   double vtxDR2[2] = { 2.*2., .2*.2};
489   double  maxW = 0;  
490   int bestBin = -1;
491   
492   for( int iter = 0; iter<2; iter++ ){
493     
494     for( int i=0; i<kNZBins; i++ ){
495       histX[i] = 0;
496       histY[i] = 0;
497       histZ[i] = 0;
498       histW[i] = 0;
499       histN[i] = 0;
500     }
501  
502     for( int iPhi=0; iPhi<kNPhiBins; iPhi++ ){
503       int nCluUp = clusters[1][iPhi].size();    
504       int nCluDn = clusters[0][iPhi].size();    
505       for( int icUp=0; icUp<nCluUp; icUp++ ){
506         AliHLTITSVZCluster &cu = clusters[1][iPhi][icUp];
507         double x0 = cu.fX - vtxX;
508         double y0 = cu.fY - vtxY;
509         double z0 = cu.fZ - vtxZ;
510         double bestR2 = 1.e10;
511         int bestDn=-1, bestBinDn =0;
512         double bestV[3]={0,0,0};
513         for( int icDn=0; icDn<nCluDn; icDn++ ){
514           AliHLTITSVZCluster &cd = clusters[0][iPhi][icDn];
515           double x1 = cd.fX - vtxX;
516           double y1 = cd.fY - vtxY;
517           double z1 = cd.fZ - vtxZ;
518           double dx = x1 - x0;
519           double dy = y1 - y0;
520           double l2 = 1./(dx*dx + dy*dy);
521           double a = x1*y0 - x0*y1;
522           double r2 = a*a*l2; 
523           if( r2>vtxDR2[iter] ) continue;
524           double xv = -dy*a*l2;
525           double yv =  dx*a*l2;
526           double zv = ( (x1*z0-x0*z1)*dx + (y1*z0-y0*z1)*dy )*l2;
527           int zbin = (int)((zv - kZBinMin)*kZStep);
528           if( zbin<0 || zbin>=kNZBins ) continue;
529           if( r2<bestR2 ){
530             bestR2 = r2;
531             bestDn = icDn;
532             bestV[0] = xv;
533             bestV[1] = yv;
534             bestV[2] = zv;
535             bestBinDn = zbin;
536           }
537         }
538         if( bestDn>=0 ){
539           double w = 1./(1.+bestR2);
540           histX[bestBinDn]+=bestV[0]*w;
541           histY[bestBinDn]+=bestV[1]*w;
542           histZ[bestBinDn]+=bestV[2]*w;
543           histW[bestBinDn]+=w;
544           histN[bestBinDn]+=1;
545         }
546       }
547     }
548     
549     maxW = 0;  
550     bestBin = -1;
551     for( int i=0; i<kNZBins; i++ ){
552       if( histW[i]>maxW ){
553         bestBin = i;
554         maxW = histW[i];
555       }
556     }
557     if( bestBin<0 || histN[bestBin] <3 ){
558       bestBin = -1;
559       break;
560     }
561     vtxX +=histX[bestBin]/maxW;
562     vtxY +=histY[bestBin]/maxW;
563     vtxZ +=histZ[bestBin]/maxW;    
564   }
565
566   if( bestBin>=0 ){
567     double nv = 1;
568     double nNew = fRunVtxNew[3] + nv;
569     double v[3] = {vtxX, vtxY, vtxZ};
570     for( int i=0; i<3; i++){
571       fRunVtxNew[i] = ( fRunVtxNew[i]*fRunVtxNew[3] + v[i]*nv )/nNew;
572     }
573     fRunVtxNew[3] = nNew;
574   }  
575   
576   if( fAutoCalibration>0 && fRunVtxNew[3] >= fAutoCalibration ){
577     for( int i=0; i<4; i++ ){
578       fRunVtx[i] = fRunVtxNew[i];
579       fRunVtxNew[i] = 0;
580     }
581     //cout<<"ITSVertexerSPD: set run vtx to "<<fRunVtx[0]<<" "<<fRunVtx[1]<<" "<<fRunVtx[2]<<endl;  
582     if( fRunVtx[0]>3. ) fRunVtx[0] = 3;
583     if( fRunVtx[0]<-3. ) fRunVtx[0] = -3;
584     if( fRunVtx[1]>3. ) fRunVtx[1] = 3;
585     if( fRunVtx[1]<-3. ) fRunVtx[1] = -3;
586     if( fRunVtx[2]>30. ) fRunVtx[2] = 30;
587     if( fRunVtx[2]<30. ) fRunVtx[2] = -30;
588   }
589
590   timerReco.Stop();
591   
592   // Fill output 
593
594   if( bestBin>=0 ){
595     double pos[3] = {vtxX, vtxY, vtxZ};
596     double s = 400.E-4;
597     double cov[6] = {s*s,0,s*s,0,0,s*s};
598     AliESDVertex v(pos, cov, 0, histN[bestBin]);
599     PushBack( (TObject*) &v, kAliHLTDataTypeESDVertex|kAliHLTDataOriginITS,0 );
600
601     //cout<<"ITSVertexerSPD: vtx found "<<vtxX<<" "<<vtxY<<" "<<vtxZ<<endl;
602
603     if( fHistoVertexXY ) fHistoVertexXY->Fill( vtxX, vtxY );
604     if( fHistoVertexX ) fHistoVertexX->Fill( vtxX );
605     if( fHistoVertexY ) fHistoVertexY->Fill( vtxY );
606     if( fHistoVertexZ ) fHistoVertexZ->Fill( vtxZ );
607   }
608
609   if( fHistoVertexXY ) PushBack( (TObject*) fHistoVertexXY, kAliHLTDataTypeHistogram,0);
610   if( fHistoVertexX ) PushBack( (TObject*) fHistoVertexX, kAliHLTDataTypeHistogram,0);
611   if( fHistoVertexY ) PushBack( (TObject*) fHistoVertexY, kAliHLTDataTypeHistogram,0);
612   if( fHistoVertexZ ) PushBack( (TObject*) fHistoVertexZ, kAliHLTDataTypeHistogram,0);
613
614
615   timer.Stop();
616   fFullTime += timer.RealTime();
617   fRecoTime += timerReco.RealTime();
618   fNEvents++;
619
620   // Set log level to "Warning" for on-line system monitoring
621   int hz = ( int ) ( fFullTime > 1.e-10 ? fNEvents / fFullTime : 100000 );
622   int hz1 = ( int ) ( fRecoTime > 1.e-10 ? fNEvents / fRecoTime : 100000 );
623   HLTInfo( "ITS Z Vertexer: input %d clusters; time: full %d / reco %d Hz",
624            nClustersTotal, hz, hz1 );
625
626   return iResult;
627 }