]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/ITS/AliHLTITSVertexerSPDComponent.cxx
coverity warnings 15382 15100 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 "AliESDVertex.h"
54
55 /** ROOT macro for the implementation of ROOT specific class methods */
56 ClassImp( AliHLTITSVertexerSPDComponent )
57 AliHLTITSVertexerSPDComponent::AliHLTITSVertexerSPDComponent()
58     :
59     fZRange(40),
60     fZBinSize(1),
61     fFullTime( 0 ),
62     fRecoTime( 0 ),
63     fNEvents( 0 ),
64     fSumW(0),
65     fSumN(0),
66     fZMin(0),
67     fNZBins(0)
68 {
69   // see header file for class documentation
70   // or
71   // refer to README to build package
72   // or
73   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
74
75   for( int i=0; i<9; i++ ) fSum[i] = 0;
76
77   fRunVtx[0] = 0;
78   fRunVtx[1] = 0;
79   fRunVtx[2] = 0;
80 }
81
82 AliHLTITSVertexerSPDComponent::AliHLTITSVertexerSPDComponent( const AliHLTITSVertexerSPDComponent& )
83     :
84     AliHLTProcessor(),
85     fZRange(40),
86     fZBinSize(1),
87     fFullTime( 0 ),
88     fRecoTime( 0 ),
89     fNEvents( 0 ),
90     fSumW(0),
91     fSumN(0),
92     fZMin(0),
93     fNZBins(0)
94 {
95   // see header file for class documentation
96
97   for( int i=0; i<9; i++ ) fSum[i] = 0;
98
99   fRunVtx[0] = 0;
100   fRunVtx[1] = 0;
101   fRunVtx[2] = 0;
102
103   HLTFatal( "copy constructor untested" );
104 }
105
106 AliHLTITSVertexerSPDComponent& AliHLTITSVertexerSPDComponent::operator=( const AliHLTITSVertexerSPDComponent& )
107 {
108   // see header file for class documentation
109   HLTFatal( "assignment operator untested" );
110   return *this;
111 }
112
113 AliHLTITSVertexerSPDComponent::~AliHLTITSVertexerSPDComponent()
114 {
115   // see header file for class documentation  
116   for( int i=0; i<9; i++ ){
117     delete[] fSum[i];
118     fSum[i] = 0;
119   }
120   delete[] fSumW;
121   delete[] fSumN;
122 }
123
124 //
125 // Public functions to implement AliHLTComponent's interface.
126 // These functions are required for the registration process
127 //
128
129 const char* AliHLTITSVertexerSPDComponent::GetComponentID()
130 {
131   // see header file for class documentation
132   return "ITSVertexerSPD";
133 }
134
135 void AliHLTITSVertexerSPDComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list )
136 {
137   // see header file for class documentation
138   list.clear();
139   list.push_back( kAliHLTDataTypeClusters|kAliHLTDataOriginITSSPD );
140 }
141
142 AliHLTComponentDataType AliHLTITSVertexerSPDComponent::GetOutputDataType()
143 {
144   // see header file for class documentation  
145   return kAliHLTMultipleDataType;
146 }
147
148 int AliHLTITSVertexerSPDComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList)
149
150 {
151   // see header file for class documentation
152   tgtList.clear();
153   tgtList.push_back( kAliHLTDataTypeESDVertex|kAliHLTDataOriginITS);
154   return tgtList.size();
155 }
156
157 void AliHLTITSVertexerSPDComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
158 {
159   // define guess for the output data size
160   constBase = 10000;       // minimum size
161   inputMultiplier = 0.5; // size relative to input
162 }
163
164 AliHLTComponent* AliHLTITSVertexerSPDComponent::Spawn()
165 {
166   // see header file for class documentation
167   return new AliHLTITSVertexerSPDComponent;
168 }
169
170 void AliHLTITSVertexerSPDComponent::SetDefaultConfiguration()
171 {
172   // Set default configuration for the CA tracker component
173   // Some parameters can be later overwritten from the OCDB
174
175   fRunVtx[0] = 0;
176   fRunVtx[1] = 0;
177   fRunVtx[2] = 0;
178   fZRange = 40;
179   fZBinSize = 1;
180   fFullTime = 0;
181   fRecoTime = 0;
182   fNEvents = 0;
183 }
184
185 int AliHLTITSVertexerSPDComponent::ReadConfigurationString(  const char* arguments )
186 {
187   // Set configuration parameters for the CA tracker component from the string
188
189   int iResult = 0;
190   if ( !arguments ) return iResult;
191
192   TString allArgs = arguments;
193   TString argument;
194   int bMissingParam = 0;
195
196   TObjArray* pTokens = allArgs.Tokenize( " " );
197
198   int nArgs =  pTokens ? pTokens->GetEntries() : 0;
199
200   for ( int i = 0; i < nArgs; i++ ) {
201     argument = ( ( TObjString* )pTokens->At( i ) )->GetString();
202     if ( argument.IsNull() ) continue;
203
204     if ( argument.CompareTo( "-runVertex" ) == 0 ) {
205       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
206       fRunVtx[0] = ( ( TObjString* )pTokens->At( i ) )->GetString().Atof();
207       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
208       fRunVtx[1] = ( ( TObjString* )pTokens->At( i ) )->GetString().Atof();
209       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
210       fRunVtx[2] = ( ( TObjString* )pTokens->At( i ) )->GetString().Atof();
211       HLTInfo( "Default run vertex is set to (%f,%f,%f)",fRunVtx[0],
212                fRunVtx[1],fRunVtx[2] );
213       continue;
214     }
215
216     if ( argument.CompareTo( "-zRange" ) == 0 ) {
217       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
218       fZRange = ( ( TObjString* )pTokens->At( i ) )->GetString().Atof();
219       HLTInfo( "Z range for the vertex search is set to +-%f cm", fZRange );
220       continue;
221     }
222
223     if ( argument.CompareTo( "-zBinSize" ) == 0 ) {
224       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
225       fZBinSize = ( ( TObjString* )pTokens->At( i ) )->GetString().Atof();
226       HLTInfo( "Size of the Z bin for the vertex search is set to %f cm", fZBinSize );
227       continue;
228     }
229
230     HLTError( "Unknown option \"%s\"", argument.Data() );
231     iResult = -EINVAL;
232   }
233   delete pTokens;
234
235   if ( bMissingParam ) {
236     HLTError( "Specifier missed for parameter \"%s\"", argument.Data() );
237     iResult = -EINVAL;
238   }
239
240   return iResult;
241 }
242
243
244 int AliHLTITSVertexerSPDComponent::ReadCDBEntry( const char* cdbEntry, const char* chainId )
245 {
246   // see header file for class documentation
247
248   const char* defaultNotify = "";
249
250   if ( !cdbEntry ) {
251     return 0;// need to add the HLT/ConfigITS/ITSTracker directory to cdb SG!!!
252     //cdbEntry = "HLT/ConfigITS/ITSTracker";
253     //defaultNotify = " (default)";
254     //chainId = 0;
255   }
256
257   HLTInfo( "configure from entry \"%s\"%s, chain id %s", cdbEntry, defaultNotify, ( chainId != NULL && chainId[0] != 0 ) ? chainId : "<none>" );
258   AliCDBEntry *pEntry = AliCDBManager::Instance()->Get( cdbEntry );//,GetRunNo());
259
260   if ( !pEntry ) {
261     HLTError( "cannot fetch object \"%s\" from CDB", cdbEntry );
262     return -EINVAL;
263   }
264
265   TObjString* pString = dynamic_cast<TObjString*>( pEntry->GetObject() );
266
267   if ( !pString ) {
268     HLTError( "configuration object \"%s\" has wrong type, required TObjString", cdbEntry );
269     return -EINVAL;
270   }
271
272   HLTInfo( "received configuration object string: \"%s\"", pString->GetString().Data() );
273
274   return  ReadConfigurationString( pString->GetString().Data() );
275 }
276
277
278 int AliHLTITSVertexerSPDComponent::Configure( const char* cdbEntry, const char* chainId, const char *commandLine )
279 {
280   // Configure the component
281   // There are few levels of configuration,
282   // parameters which are set on one step can be overwritten on the next step
283
284   //* read hard-coded values
285
286   SetDefaultConfiguration();
287
288   //* read the default CDB entry
289
290   int iResult1 = ReadCDBEntry( NULL, chainId );
291
292   //* read the actual CDB entry if required
293
294   int iResult2 = ( cdbEntry ) ? ReadCDBEntry( cdbEntry, chainId ) : 0;
295   
296   //* read extra parameters from input (if they are)
297
298   int iResult3 = 0;
299
300   if ( commandLine && commandLine[0] != '\0' ) {
301     HLTInfo( "received configuration string from HLT framework: \"%s\"", commandLine );
302     iResult3 = ReadConfigurationString( commandLine );
303   }
304
305   // Initialise the tracker here
306
307   return iResult1 ? iResult1 : ( iResult2 ? iResult2 : iResult3 );
308 }
309
310
311
312 int AliHLTITSVertexerSPDComponent::DoInit( int argc, const char** argv )
313 {
314   // Configure the ITS tracker component
315   
316   SetDefaultConfiguration();
317
318   if(AliGeomManager::GetGeometry()==NULL){
319     AliGeomManager::LoadGeometry("");
320   }
321   AliGeomManager::ApplyAlignObjsFromCDB("ITS");
322   
323   TString arguments = "";
324   for ( int i = 0; i < argc; i++ ) {
325     if ( !arguments.IsNull() ) arguments += " ";
326     arguments += argv[i];
327   }
328
329   int ret = Configure( NULL, NULL, arguments.Data() );
330
331   for( int i=0; i<9; i++ ) delete[] fSum[i];
332   delete[] fSumW;
333   delete[] fSumN;
334
335   fZMin = -fZRange;
336   fNZBins = ( (int) (2*fZRange/fZBinSize ))+1;
337
338   for( int i=0; i<9; i++ ) fSum[i] = new double [fNZBins];
339
340   fSumW = new double [fNZBins];
341   fSumN = new int [fNZBins];
342
343   return ret;
344 }
345
346
347 int AliHLTITSVertexerSPDComponent::DoDeinit()
348 {
349   // see header file for class documentation
350
351   for( int i=0; i<9; i++ ){    
352     delete[] fSum[i];
353     fSum[i] = 0;
354   }
355   delete[] fSumW;
356   delete[] fSumN;
357   fSumW = 0;
358   fSumN = 0;
359   
360   SetDefaultConfiguration();
361
362   return 0;
363 }
364
365
366
367 int AliHLTITSVertexerSPDComponent::Reconfigure( const char* cdbEntry, const char* chainId )
368 {
369   // Reconfigure the component from OCDB .
370
371   return Configure( cdbEntry, chainId, NULL );
372 }
373
374
375 int AliHLTITSVertexerSPDComponent::DoEvent
376 (
377   const AliHLTComponentEventData& evtData,
378   const AliHLTComponentBlockData* blocks,
379   AliHLTComponentTriggerData& /*trigData*/,
380   AliHLTUInt8_t* /*outputPtr*/,
381   AliHLTUInt32_t& size,
382   vector<AliHLTComponentBlockData>& /*outputBlocks*/ )
383 {
384   //* process event
385
386   //AliHLTUInt32_t maxBufferSize = size;
387   size = 0; // output size
388
389   if (!IsDataEvent()) return 0;
390
391   if ( evtData.fBlockCnt <= 0 ) {
392     HLTWarning( "no blocks in event" );
393     return 0;
394   }
395
396   TStopwatch timer;
397
398   // Event reconstruction in ITS
399
400   int iResult=0;
401
402   int nBlocks = evtData.fBlockCnt;
403   int nClustersTotal = 0;
404  
405
406   const int kNPhiBins = 20;
407   const double kDPhi = TMath::TwoPi() / kNPhiBins;
408   double vtxX = fRunVtx[0], vtxY = fRunVtx[1], vtxZ = fRunVtx[2];
409
410   std::vector<AliHLTITSVertexerSPDComponent::AliHLTITSVZCluster> clusters[2][ kNPhiBins ];
411
412   for (int ndx=0; ndx<nBlocks && iResult>=0; ndx++) {
413
414     const AliHLTComponentBlockData* iter = blocks+ndx;
415  
416     // Read ITS SPD clusters
417
418     if ( iter->fDataType == (kAliHLTDataTypeClusters|kAliHLTDataOriginITSSPD) ){
419
420       AliHLTITSClusterData *inPtr=reinterpret_cast<AliHLTITSClusterData*>( iter->fPtr );
421       int nClusters = inPtr->fSpacePointCnt;
422       nClustersTotal+=nClusters;
423       for( int icl=0; icl<nClusters; icl++ ){   
424         AliHLTITSSpacePointData &d = inPtr->fSpacePoints[icl];
425         if( d.fLayer>1 ) continue;// SPD only;
426         if( d.fLayer<0 ) continue;// SPD only;
427         Int_t lab[4] = { d.fTracks[0], d.fTracks[1], d.fTracks[2], d.fIndex };
428         Int_t info[3] = { d.fNy, d.fNz, d.fLayer };
429         Float_t hit[6] = { d.fY, d.fZ, d.fSigmaY2, d.fSigmaZ2, d.fQ, d.fSigmaYZ };
430         if( d.fLayer==4 ) hit[5] = -hit[5];
431
432
433         AliITSRecPoint p( lab, hit, info );
434
435         if (!p.Misalign()){
436           HLTWarning("Can not misalign an SPD cluster");
437         }
438         Float_t xyz[3];      
439         p.GetGlobalXYZ(xyz);      
440
441         // get phi bin
442         double phi = TMath::ATan2( xyz[1]-vtxY, xyz[0]-vtxX );
443         if( phi<0 ) phi+=TMath::TwoPi();
444         int iphi = (int ) (phi/kDPhi);
445         if( iphi<0 ) iphi = 0;
446         if( iphi>=kNPhiBins ) iphi = kNPhiBins-1;
447         AliHLTITSVZCluster c;
448         c.fX = xyz[0];
449         c.fY = xyz[1];
450         c.fZ = xyz[2];
451         clusters[d.fLayer][iphi].push_back( c );        
452      }   
453     }    
454   }// end read input blocks
455   
456
457   // Reconstruct the vertex
458
459   TStopwatch timerReco;
460     
461   double zScale = 1./fZBinSize;
462
463   // fit few times decreasing the radius cut
464
465   double vtxDR2[5] = { 1.*1., .5*.5, .4*.4, .4*.4, .2*.2};
466   double vtxDZ [5] = { 1., .5, .3, .3, .3};
467   double  maxW = 0;
468   int bestBin = -1;
469   
470   for( int iter = 0; iter<3; iter++ ){
471
472     bool doZBins = (iter<2);
473     int nSearchBins = doZBins ?fNZBins :1;
474     {
475       for( int i=0; i<nSearchBins; i++ ){
476         fSum[0][i] = 0;
477         fSum[1][i] = 0;
478         fSum[2][i] = 0;
479         fSum[3][i] = 0;
480         fSum[4][i] = 0;
481         fSum[5][i] = 0;
482         fSum[6][i] = 0;
483         fSum[7][i] = 0;
484         fSum[8][i] = 0;
485         fSumW[i] = 0;
486         fSumN[i] = 0;
487       }
488     }
489
490     for( int iPhi=0; iPhi<kNPhiBins; iPhi++ ){
491       int nCluUp = clusters[1][iPhi].size();    
492       int nCluDn = clusters[0][iPhi].size();    
493       for( int icUp=0; icUp<nCluUp; icUp++ ){
494         AliHLTITSVZCluster &cu = clusters[1][iPhi][icUp];
495         double x0 = cu.fX - vtxX;
496         double y0 = cu.fY - vtxY;
497         double z0 = cu.fZ;// - vtxZ;
498         double bestR2 = 1.e10;
499         int bestDn=-1, bestBinDn =0;
500         double bestP[3]={0,0,0};
501         for( int icDn=0; icDn<nCluDn; icDn++ ){
502           AliHLTITSVZCluster &cd = clusters[0][iPhi][icDn];
503           double x1 = cd.fX - vtxX;
504           double y1 = cd.fY - vtxY;
505           double z1 = cd.fZ;// - vtxZ;
506           double dx = x1 - x0;
507           double dy = y1 - y0;
508           double l2 = 1./(dx*dx + dy*dy);
509           double a = x1*y0 - x0*y1;
510           double r2 = a*a*l2; 
511           if( r2 > vtxDR2[iter] ) continue;
512           if( r2 > bestR2 ) continue;
513           //double xv = -dy*a*l2;
514           //double yv =  dx*a*l2;
515           double zv = ( (x1*z0-x0*z1)*dx + (y1*z0-y0*z1)*dy )*l2;
516           int zbin;
517           if( doZBins ){
518             zbin = (int)((zv - fZMin)*zScale);
519             if( zbin<0 || zbin>=fNZBins ) continue;       
520           } else {
521             zbin = 0;
522             if( TMath::Abs( zv - vtxZ ) > vtxDZ[ iter ] ) continue;
523           }
524           bestR2 = r2;
525           bestDn = icDn;
526           bestP[0] = x1;
527           bestP[1] = y1;
528           bestP[2] = z1;
529           bestBinDn = zbin;       
530         }
531         if( bestDn < 0 ) continue;
532
533         double dx = bestP[0] - x0;
534         double dy = bestP[1] - y0;
535         double dz = bestP[2] - z0;        
536         double w = 1./(1.+bestR2);
537           
538         // Equations:
539         //
540         // A_i*x + B_i*y + C_i*z = D_i
541         //
542         // Sum[0] = sum A_i*A_i
543         // Sum[1] = sum B_i*A_i
544         // Sum[2] = sum B_i*B_i
545         // Sum[3] = sum C_i*A_i
546         // Sum[4] = sum C_i*B_i
547         // Sum[5] = sum C_i*C_i   
548         // Sum[6] = sum A_i*D_i
549         // Sum[7] = sum B_i*D_i
550         // Sum[8] = sum C_i*D_i
551         
552         double n = w / TMath::Sqrt(dx*dx + dy*dy + dz*dz);        
553         dy*=n;
554         dx*=n;
555         dz*=n;
556         double a, b, c, d;
557         a = dy; b = -dx; c = 0; d = dy*x0 - dx*y0;
558         {
559           fSum[0][bestBinDn]+= a*a;
560           fSum[1][bestBinDn]+= b*a;
561           fSum[2][bestBinDn]+= b*b;
562           //fSum[3][bestBinDn]+= c*a;
563           //fSum[4][bestBinDn]+= c*b;
564           //fSum[5][bestBinDn]+= c*c;
565           fSum[6][bestBinDn]+= a*d;
566           fSum[7][bestBinDn]+= b*d;
567           //fSum[8][bestBinDn]+= c*d;
568         }
569         a = dz; b = 0; c = -dx; d = dz*x0 - dx*z0;
570         {
571           fSum[0][bestBinDn]+= a*a;
572           //fSum[1][bestBinDn]+= b*a;
573           //fSum[2][bestBinDn]+= b*b;
574           fSum[3][bestBinDn]+= c*a;
575           //fSum[4][bestBinDn]+= c*b;
576           fSum[5][bestBinDn]+= c*c;
577           fSum[6][bestBinDn]+= a*d;
578           //fSum[7][bestBinDn]+= b*d;
579           fSum[8][bestBinDn]+= c*d;
580         }
581         
582         a = 0; b = dz; c = -dy; d = dz*y0 - dy*z0;
583         {
584           //fSum[0][bestBinDn]+= a*a;
585           //fSum[1][bestBinDn]+= b*a;
586           fSum[2][bestBinDn]+= b*b;
587           //fSum[3][bestBinDn]+= c*a;
588           fSum[4][bestBinDn]+= c*b;
589           fSum[5][bestBinDn]+= c*c;
590           //fSum[6][bestBinDn]+= a*d;
591           fSum[7][bestBinDn]+= b*d;
592           fSum[8][bestBinDn]+= c*d;
593         }
594
595         a = dz; b = 0; c = -dx; d = dz*x0 - dx*z0;
596         {
597           fSum[0][bestBinDn]+= a*a;
598           //fSum[1][bestBinDn]+= b*a;
599           //fSum[2][bestBinDn]+= b*b;
600           fSum[3][bestBinDn]+= c*a;
601           //fSum[4][bestBinDn]+= c*b;
602           fSum[5][bestBinDn]+= c*c;
603           fSum[6][bestBinDn]+= a*d;
604           //fSum[7][bestBinDn]+= b*d;
605           fSum[8][bestBinDn]+= c*d;
606         }
607         
608         fSumW[bestBinDn]+=w;
609         fSumN[bestBinDn]+=1;
610       }    
611     }
612     
613     maxW = 0;  
614     bestBin = -1;
615     for( int i=0; i<nSearchBins; i++ ){
616       if( fSumN[i]>maxW ){
617         bestBin = i;
618         maxW = fSumN[i];
619       }
620     }
621     if( bestBin<0 || fSumN[bestBin] < 3 ){
622       bestBin = -1;
623       break;
624     }
625     
626     // calculate the vertex position in the best bin
627     Double_t w[6] = { fSum[0][bestBin],
628                       fSum[1][bestBin], fSum[2][bestBin],
629                       fSum[3][bestBin], fSum[4][bestBin], fSum[5][bestBin] };
630     Double_t wI[6];
631     {
632       wI[0] = w[2]*w[5] - w[4]*w[4];
633       wI[1] = w[3]*w[4] - w[1]*w[5];
634       wI[2] = w[0]*w[5] - w[3]*w[3];
635       wI[3] = w[1]*w[4] - w[2]*w[3];
636       wI[4] = w[1]*w[3] - w[0]*w[4];
637       wI[5] = w[0]*w[2] - w[1]*w[1];     
638       
639       Double_t s = ( w[0]*wI[0] + w[1]*wI[1] + w[3]*wI[3] );
640       
641       if( s<1.e-10 ){
642         bestBin = -1;
643         break;
644       }
645       s = 1./s;   
646       wI[0]*=s;
647       wI[1]*=s;
648       wI[2]*=s;
649       wI[3]*=s;
650       wI[4]*=s;
651       wI[5]*=s;
652     }    
653
654     vtxX += wI[0]*fSum[6][bestBin] + wI[1]*fSum[7][bestBin] + wI[3]*fSum[8][bestBin];
655     vtxY += wI[1]*fSum[6][bestBin] + wI[2]*fSum[7][bestBin] + wI[4]*fSum[8][bestBin];
656     vtxZ = wI[3]*fSum[6][bestBin] + wI[4]*fSum[7][bestBin] + wI[5]*fSum[8][bestBin];
657     //cout<<"SG: "<<iter<<": "<<vtxX<<" "<<vtxY<<" "<<vtxZ<<endl;
658   }
659
660   timerReco.Stop();
661   
662   // Fill output 
663
664   if( bestBin>=0 ){
665     double pos[3] = {vtxX, vtxY, vtxZ};
666     double s = 400.E-4;
667     double cov[6] = {s*s,0,s*s,0,0,s*s};
668     AliESDVertex v(pos, cov, 0, fSumN[bestBin]);
669     PushBack( (TObject*) &v, kAliHLTDataTypeESDVertex|kAliHLTDataOriginITS,0 );
670
671     //cout<<"ITSVertexerSPD: vtx found "<<vtxX<<" "<<vtxY<<" "<<vtxZ<<endl;
672   }
673
674   timer.Stop();
675   fFullTime += timer.RealTime();
676   fRecoTime += timerReco.RealTime();
677   fNEvents++;
678
679   // Set log level to "Warning" for on-line system monitoring
680   int hz = ( int ) ( fFullTime > 1.e-10 ? fNEvents / fFullTime : 100000 );
681   int hz1 = ( int ) ( fRecoTime > 1.e-10 ? fNEvents / fRecoTime : 100000 );
682   HLTInfo( "ITS Z Vertexer: input %d clusters; time: full %d / reco %d Hz",
683            nClustersTotal, hz, hz1 );
684
685   return iResult;
686 }