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