]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/global/AliHLTGlobalVertexerComponent.cxx
Fix leaks in AliIsolationCut and AlianaOmegaToPi0Gamma
[u/mrichter/AliRoot.git] / HLT / global / AliHLTGlobalVertexerComponent.cxx
1 //**************************************************************************
2 //* This file is property of and copyright by the ALICE HLT Project        * 
3 //* ALICE Experiment at CERN, All rights reserved.                         *
4 //*                                                                        *
5 //* Primary Authors: S.Gorbunov <sergey.gorbunov@kip.uni-heidelberg.de>    *
6 //*                  for The ALICE HLT Project.                            *
7 //*                                                                        *
8 //* Permission to use, copy, modify and distribute this software and its   *
9 //* documentation strictly for non-commercial purposes is hereby granted   *
10 //* without fee, provided that the above copyright notice appears in all   *
11 //* copies and that both the copyright notice and this permission notice   *
12 //* appear in the supporting documentation. The authors make no claims     *
13 //* about the suitability of this software for any purpose. It is          *
14 //* provided "as is" without express or implied warranty.                  *
15 //**************************************************************************
16
17 /** @file   AliHLTGlobalVertexerComponent.cxx
18     @author Sergey Gorbunov
19     @brief  Component for reconstruct primary vertex and V0's
20 */
21
22 #if __GNUC__>= 3
23 using namespace std;
24 #endif
25
26 #include "AliHLTGlobalVertexerComponent.h"
27 #include "AliHLTDataTypes.h"
28 #include "AliHLTExternalTrackParam.h"
29 #include "AliHLTGlobalBarrelTrack.h"
30 #include "AliCDBEntry.h"
31 #include "AliCDBManager.h"
32 #include <TFile.h>
33 #include <TString.h>
34 #include "TObjString.h"
35 #include "TObjArray.h"
36 #include "AliESDEvent.h"
37 #include "AliESDtrack.h"
38 #include "AliESDVertex.h"
39 #include "AliESDv0.h"
40 #include "AliHLTMessage.h"
41 #include "TMath.h"
42 #include "AliKFParticle.h"
43 #include "AliKFVertex.h"
44 #include "TStopwatch.h"
45
46 /** ROOT macro for the implementation of ROOT specific class methods */
47 ClassImp(AliHLTGlobalVertexerComponent)
48
49 AliHLTGlobalVertexerComponent::AliHLTGlobalVertexerComponent()
50 :
51   fNTracks(0),
52   fTrackInfos(0),
53   fPrimaryVtx(),   
54   fFitTracksToVertex(1),
55   fConstrainedTrackDeviation(4.),
56   fV0DaughterPrimDeviation( 2.5 ),
57   fV0PrimDeviation( 3.5 ),
58   fV0Chi(3.5),
59   fV0DecayLengthInSigmas(3.),
60   fV0TimeLimit(10.e-3), 
61   fBenchmark("GlobalVertexer")
62 {
63   // see header file for class documentation
64   // or
65   // refer to README to build package
66   // or
67   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
68
69 }
70
71 AliHLTGlobalVertexerComponent::~AliHLTGlobalVertexerComponent()
72 {
73   // see header file for class documentation
74
75   if( fTrackInfos ) delete[] fTrackInfos;
76 }
77
78 // Public functions to implement AliHLTComponent's interface.
79 // These functions are required for the registration process
80
81 const char* AliHLTGlobalVertexerComponent::GetComponentID()
82 {
83   // see header file for class documentation
84   
85   return "GlobalVertexer";
86 }
87
88 void AliHLTGlobalVertexerComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list)
89 {
90   // see header file for class documentation
91   list.clear();
92   list.push_back( kAliHLTDataTypeESDObject );
93   list.push_back( kAliHLTDataTypeTrack|kAliHLTDataOriginITS );
94   list.push_back( kAliHLTDataTypeTrack|kAliHLTDataOriginTPC );
95 }
96
97 AliHLTComponentDataType AliHLTGlobalVertexerComponent::GetOutputDataType()
98 {
99   // see header file for class documentation
100   return kAliHLTMultipleDataType;
101 }
102
103 int AliHLTGlobalVertexerComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList)
104
105 {
106   // see header file for class documentation
107   tgtList.clear();
108   tgtList.push_back( kAliHLTDataTypeGlobalVertexer|kAliHLTDataOriginOut);
109   tgtList.push_back( kAliHLTDataTypeESDVertex|kAliHLTDataOriginOut);
110   tgtList.push_back( kAliHLTDataTypeESDObject|kAliHLTDataOriginOut);
111   return tgtList.size();
112 }
113
114 void AliHLTGlobalVertexerComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
115 {
116   // see header file for class documentation
117   // XXX TODO: Find more realistic values.
118   constBase = 80000;
119   inputMultiplier = 2.;
120 }
121
122 AliHLTComponent* AliHLTGlobalVertexerComponent::Spawn()
123 {
124   // see header file for class documentation
125   return new AliHLTGlobalVertexerComponent;
126 }
127
128 int AliHLTGlobalVertexerComponent::DoInit( int argc, const char** argv )
129 {
130   // init
131
132   fBenchmark.Reset();
133   fBenchmark.SetTimer(0,"total");
134   fBenchmark.SetTimer(1,"convert");
135   fBenchmark.SetTimer(2,"vprim");
136   fBenchmark.SetTimer(3,"v0");
137   fV0TimeLimit = 10.e-3;
138
139   int iResult=0;
140   TString configuration="";
141   TString argument="";
142   for (int i=0; i<argc && iResult>=0; i++) {
143     argument=argv[i];
144     if (!configuration.IsNull()) configuration+=" ";
145     configuration+=argument;
146   }
147   
148   if (!configuration.IsNull()) {
149     iResult=Configure(configuration.Data());
150   }  
151
152   return iResult; 
153 }
154   
155 int AliHLTGlobalVertexerComponent::DoDeinit()
156 {
157   // see header file for class documentation
158
159   if( fTrackInfos ) delete[] fTrackInfos;
160
161   fTrackInfos = 0;
162   fFitTracksToVertex = 1;
163   fConstrainedTrackDeviation = 4.;
164   fV0DaughterPrimDeviation = 2.5 ;
165   fV0PrimDeviation =3.5;
166   fV0Chi = 3.5;
167   fV0DecayLengthInSigmas = 3.;
168   fV0TimeLimit = 10.e-3;
169
170   return 0;
171 }
172
173 int AliHLTGlobalVertexerComponent::DoEvent( const AliHLTComponentEventData& /*evtData*/, 
174                                             AliHLTComponentTriggerData& /*trigData*/ )
175 {
176   //cout<<"AliHLTGlobalVertexerComponent::DoEvent called"<<endl;
177  
178   if ( GetFirstInputBlock( kAliHLTDataTypeSOR ) || GetFirstInputBlock( kAliHLTDataTypeEOR ) )
179     return 0;
180
181   int iResult = 0;
182   //cout<<"SG: GlobalVertexer:: DoEvent called"<<endl;
183   fBenchmark.StartNewEvent();
184   fBenchmark.Start(0);
185   
186   vector< AliExternalTrackParam > tracks;
187   vector< int > trackId;
188   vector< pair<int,int> > v0s;
189
190   AliESDEvent *event = 0; 
191
192   for( const AliHLTComponentBlockData *i= GetFirstInputBlock(kAliHLTDataTypeESDObject); i!=NULL; i=GetNextInputBlock() ){
193     fBenchmark.AddInput(i->fSize);
194   }
195
196
197   for ( const TObject *iter = GetFirstInputObject(kAliHLTDataTypeESDObject); iter != NULL; iter = GetNextInputObject() ) {
198     event = dynamic_cast<AliESDEvent*>(const_cast<TObject*>( iter ) );
199     if( !event ) continue;
200     
201     event->GetStdContent();
202     Int_t nESDTracks=event->GetNumberOfTracks(); 
203
204     for (Int_t iTr=0; iTr<nESDTracks; iTr++){   
205       AliESDtrack *pTrack = event->GetTrack(iTr);    
206       if( !pTrack  ) continue;
207       if( pTrack->GetKinkIndex(0)>0) continue;
208       if( !( pTrack->GetStatus()&AliESDtrack::kTPCin ) ) continue;
209       tracks.push_back(*pTrack);
210       trackId.push_back(iTr);
211     }      
212     break;
213   }
214
215   if( tracks.size()==0 ){
216     
217     for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginITS);
218          pBlock!=NULL; pBlock=GetNextInputBlock()) {
219
220       fBenchmark.AddInput(pBlock->fSize);
221
222       AliHLTTracksData* dataPtr = reinterpret_cast<AliHLTTracksData*>( pBlock->fPtr );
223       int nTracks = dataPtr->fCount;
224
225       AliHLTExternalTrackParam* currOutTrack = dataPtr->fTracklets;
226       for( int itr=0; itr<nTracks; itr++ ){
227         AliHLTGlobalBarrelTrack t(*currOutTrack);
228         tracks.push_back( t );
229         trackId.push_back( currOutTrack->fTrackID );
230         unsigned int dSize = sizeof( AliHLTExternalTrackParam ) + currOutTrack->fNPoints * sizeof( unsigned int );
231         currOutTrack = ( AliHLTExternalTrackParam* )( (( Byte_t * )currOutTrack) + dSize );
232       }
233     }
234   }
235
236   if( tracks.size()==0 ){
237     for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC);
238          pBlock!=NULL; pBlock=GetNextInputBlock()) {
239
240       fBenchmark.AddInput(pBlock->fSize);
241
242       AliHLTTracksData* dataPtr = reinterpret_cast<AliHLTTracksData*>( pBlock->fPtr );
243       int nTracks = dataPtr->fCount;      
244       AliHLTExternalTrackParam* currOutTrack = dataPtr->fTracklets;
245       for( int itr=0; itr<nTracks; itr++ ){
246         AliHLTGlobalBarrelTrack t(*currOutTrack);
247         tracks.push_back( t );
248         trackId.push_back( currOutTrack->fTrackID );
249         unsigned int dSize = sizeof( AliHLTExternalTrackParam ) + currOutTrack->fNPoints * sizeof( unsigned int );
250         currOutTrack = ( AliHLTExternalTrackParam* )( (( Byte_t * )currOutTrack) + dSize );
251       }
252     }
253   }
254
255   
256   // primary vertex & V0's 
257           
258   AliKFParticle::SetField( GetBz() );
259
260   fBenchmark.Start(1);
261
262   {  //* Fill fTrackInfo array
263
264     if( fTrackInfos ) delete[] fTrackInfos;
265     fTrackInfos = 0;
266     fNTracks=tracks.size(); 
267     fTrackInfos = new AliESDTrackInfo[ fNTracks ];
268     for (Int_t iTr=0; iTr<fNTracks; iTr++){       
269       AliESDTrackInfo &info = fTrackInfos[iTr];
270       info.fOK = 1;
271       info.fPrimUsedFlag = 0;
272       info.fParticle = AliKFParticle( tracks[iTr], 211 );
273       for( int i=0; i<8; i++ ) if( !finite(info.fParticle.GetParameter(i)) ) info.fOK = 0;
274       for( int i=0; i<36; i++ ) if( !finite(info.fParticle.GetCovariance(i)) ) info.fOK = 0;
275     }
276   }
277
278   fBenchmark.Stop(1);
279   fBenchmark.Start(2);
280   FindPrimaryVertex();
281   fBenchmark.Stop(2);
282   fBenchmark.Start(3);
283   FindV0s( v0s );
284   fBenchmark.Stop(3);
285
286   int *buf = new int[sizeof(AliHLTGlobalVertexerData)/sizeof(int)+1 + fNTracks + 2*v0s.size()];
287   AliHLTGlobalVertexerData *data = reinterpret_cast<AliHLTGlobalVertexerData*>(buf);
288
289   if( data) {  // fill the output structure
290         
291     data->fFitTracksToVertex = fFitTracksToVertex;
292     for( int i=0; i<3; i++ ) data->fPrimP[i] = fPrimaryVtx.Parameters()[i];
293     for( int i=0; i<6; i++ ) data->fPrimC[i] = fPrimaryVtx.CovarianceMatrix()[i];
294     data->fPrimChi2 = fPrimaryVtx.GetChi2();
295     data->fPrimNContributors = fPrimaryVtx.GetNContributors();
296     data->fNPrimTracks = 0;
297     for( Int_t i = 0; i<fNTracks; i++ ){
298       if( !fTrackInfos[i].fPrimUsedFlag ) continue;       
299       if( fTrackInfos[i].fPrimDeviation > fConstrainedTrackDeviation ) continue;
300       data->fTrackIndices[ (data->fNPrimTracks)++ ] =  trackId[i];
301     }
302     int *listV0 = data->fTrackIndices + data->fNPrimTracks;
303     data->fNV0s = v0s.size();
304     for( int i=0; i<data->fNV0s; i++ ){
305       listV0[2*i] = trackId[v0s[i].first];
306       listV0[2*i+1] = trackId[v0s[i].second];
307     }
308
309     unsigned int blockSize = sizeof(AliHLTGlobalVertexerData) + (data->fNPrimTracks + 2*data->fNV0s)*sizeof(int);
310
311     iResult = PushBack( reinterpret_cast<void*>(data), blockSize, kAliHLTDataTypeGlobalVertexer|kAliHLTDataOriginOut );  
312     fBenchmark.AddOutput(blockSize);
313   }  
314   
315   
316   // output the vertex if found
317   {
318     if( iResult==0 && data && data->fPrimNContributors >=3 ){
319       AliESDVertex vESD( data->fPrimP, data->fPrimC, data->fPrimChi2, data->fPrimNContributors );
320       iResult = PushBack( (TObject*) &vESD, kAliHLTDataTypeESDVertex|kAliHLTDataOriginOut,0 );
321       fBenchmark.AddOutput(GetLastObjectSize());
322     }
323   }
324
325   // output the ESD event
326   if( iResult==0 && event && data ){  
327     FillESD( event, data ); 
328     iResult = PushBack( event, kAliHLTDataTypeESDObject|kAliHLTDataOriginOut, 0);
329     fBenchmark.AddOutput(GetLastObjectSize());
330   }
331
332   delete[] buf;
333
334   fBenchmark.Stop(0);
335   HLTInfo(fBenchmark.GetStatistics());
336   return iResult;
337 }
338
339
340 int AliHLTGlobalVertexerComponent::Configure(const char* arguments)
341 {
342   
343   int iResult=0;
344   if (!arguments) return iResult;
345   
346   TString allArgs=arguments;
347   TString argument;
348   
349   TObjArray* pTokens=allArgs.Tokenize(" ");
350   int bMissingParam=0;
351
352   if (pTokens) {
353     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
354       argument=((TObjString*)pTokens->At(i))->GetString();
355       if (argument.IsNull()) continue;
356       
357       if (argument.CompareTo("-fitTracksToVertex")==0) {
358         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
359         HLTInfo("fitTracksToVertex is set set to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
360         fFitTracksToVertex=((TObjString*)pTokens->At(i))->GetString().Atoi();
361         continue;
362       }
363       else if (argument.CompareTo("-constrainedTrackDeviation")==0) {
364         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
365         HLTInfo("constrainedTrackDeviation is set set to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
366         fConstrainedTrackDeviation=((TObjString*)pTokens->At(i))->GetString().Atof();
367         continue;
368       }
369       else if (argument.CompareTo("-v0DaughterPrimDeviation")==0) {
370         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
371         HLTInfo("v0DaughterPrimDeviation is set set to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
372         fV0DaughterPrimDeviation=((TObjString*)pTokens->At(i))->GetString().Atof();
373         continue;
374       }
375       else if (argument.CompareTo("-v0PrimDeviation")==0) {
376         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
377         HLTInfo("v0PrimDeviation is set set to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
378         fV0PrimDeviation=((TObjString*)pTokens->At(i))->GetString().Atof();
379         continue;
380       }
381       else if (argument.CompareTo("-v0Chi")==0) {
382         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
383         HLTInfo("v0Chi is set set to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
384         fV0Chi=((TObjString*)pTokens->At(i))->GetString().Atof();
385         continue;
386       }
387       else if (argument.CompareTo("-v0DecayLengthInSigmas")==0) {
388         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
389         HLTInfo("v0DecayLengthInSigmas is set set to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
390         fV0DecayLengthInSigmas=((TObjString*)pTokens->At(i))->GetString().Atof();
391         continue;
392       }
393       else if (argument.CompareTo("-v0MinEventRate")==0) {
394         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
395         HLTInfo("Minimum event rate for V0 finder is set set to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
396         Double_t rate = ((TObjString*)pTokens->At(i))->GetString().Atof();
397         fV0TimeLimit = (rate >0 ) ?1./rate :60; // 1 minute maximum time
398         continue;
399       }
400       else {
401         HLTError("unknown argument %s", argument.Data());
402         iResult=-EINVAL;
403         break;
404       }
405     }
406     delete pTokens;
407   }
408   if (bMissingParam) {
409     HLTError("missing parameter for argument %s", argument.Data());
410     iResult=-EINVAL;
411   }  
412   
413   return iResult;
414 }
415
416 int AliHLTGlobalVertexerComponent::Reconfigure(const char* cdbEntry, const char* chainId)
417 {
418   // see header file for class documentation
419
420   return 0; // no CDB path is set so far
421
422   int iResult=0;  
423   const char* path="HLT/ConfigTPC/KryptonHistoComponent";
424   const char* defaultNotify="";
425   if (cdbEntry) {
426     path=cdbEntry;
427     defaultNotify=" (default)";
428   }
429   if (path) {
430     HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
431     AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
432     if (pEntry) {
433       TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
434       if (pString) {
435         HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
436         iResult=Configure(pString->GetString().Data());
437       } else {
438         HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
439       }
440     } else {
441       HLTError("can not fetch object \"%s\" from CDB", path);
442     }
443   }
444
445   return iResult;
446 }
447
448
449 struct AliHLTGlobalVertexerDeviation
450 {
451   int fI; // track index
452   float fD; // deviation from vertex
453
454   bool operator<(const AliHLTGlobalVertexerDeviation &a) const { return fD<a.fD; }
455 };
456
457
458 void AliHLTGlobalVertexerComponent::FindPrimaryVertex()
459 {
460   //* Find event primary vertex
461
462   fPrimaryVtx.Initialize();
463   //fPrimaryVtx.SetBeamConstraint(fESD->GetDiamondX(),fESD->GetDiamondY(),0,
464   //TMath::Sqrt(fESD->GetSigma2DiamondX()),TMath::Sqrt(fESD->GetSigma2DiamondY()),5.3);
465
466   // select rough region (in sigmas) in which the vertex could be found, all tracks outside these limits are rejected
467   // from the primary vertex finding
468   fPrimaryVtx.SetBeamConstraint( 0, 0, 0, 3., 3., 5.3 );
469
470   const AliKFParticle **vSelected = new const AliKFParticle*[fNTracks]; //* Selected particles for the vertex fit
471   AliHLTGlobalVertexerDeviation *dev = new AliHLTGlobalVertexerDeviation[fNTracks];  
472
473   Int_t nSelected = 0;
474   
475   for( Int_t i = 0; i<fNTracks; i++){ 
476     if(!fTrackInfos[i].fOK ) continue;
477     //if( fESD->GetTrack(i)->GetTPCNcls()<60  ) continue;
478     const AliKFParticle &p = fTrackInfos[i].fParticle;
479     Double_t chi = p.GetDeviationFromVertex( fPrimaryVtx );      
480     if( chi > fConstrainedTrackDeviation ) continue;
481     dev[nSelected].fI = i; 
482     dev[nSelected].fD = chi; 
483     vSelected[nSelected] = &(fTrackInfos[i].fParticle);
484     nSelected++;  
485   }    
486
487   // fit
488
489   while( nSelected>2 ){ 
490
491     //* Primary vertex finder with rejection of outliers
492
493     for( Int_t i = 0; i<nSelected; i++){ 
494       vSelected[i] = &(fTrackInfos[dev[i].fI].fParticle);      
495     }
496
497     double xv = fPrimaryVtx.GetX();
498     double yv = fPrimaryVtx.GetY();
499     double zv = fPrimaryVtx.GetZ(); // values from previous iteration of calculations          
500     fPrimaryVtx.Initialize();
501     fPrimaryVtx.SetBeamConstraint( 0, 0, 0, 3., 3., 5.3 );
502     fPrimaryVtx.SetVtxGuess( xv, yv, zv );    
503
504     fPrimaryVtx.Construct( vSelected, nSelected, 0, -1, 1 ); // refilled for every iteration
505     
506     for( Int_t it=0; it<nSelected; it++ ){ 
507       const AliKFParticle &p = fTrackInfos[dev[it].fI].fParticle;
508       if( nSelected <= 20 ){
509         AliKFVertex tmp =  fPrimaryVtx - p; // exclude the current track from the sample and recalculate the vertex
510         dev[it].fD = p.GetDeviationFromVertex( tmp );
511       } else {
512         dev[it].fD = p.GetDeviationFromVertex( fPrimaryVtx );   
513       }
514     }
515     sort(dev,dev+nSelected); // sort tracks with increasing chi2 (used for rejection)
516        
517     int nRemove = (int) ( 0.3*nSelected );  //remove 30% of the tracks (done for performance, only if there are more than 20 tracks)  
518     if( nSelected - nRemove <=20 ) nRemove = 1;  // removal based on the chi2 of every track   
519     int firstRemove = nSelected - nRemove;
520     while( firstRemove<nSelected ){
521       if( dev[firstRemove].fD >= fConstrainedTrackDeviation ) break;
522       firstRemove++;
523     }
524     if( firstRemove>=nSelected ) break;
525     nSelected = firstRemove;
526   }
527
528   for( Int_t i = 0; i<fNTracks; i++){ 
529     fTrackInfos[i].fPrimUsedFlag = 0;
530   }
531
532   if( nSelected < 3 ){  // no vertex for fewer than 3 contributors  
533     fPrimaryVtx.NDF() = -3;
534     fPrimaryVtx.Chi2() = 0;
535     nSelected = 0;
536   }
537   
538   for( Int_t i = 0; i<nSelected; i++){ 
539     AliESDTrackInfo &info = fTrackInfos[dev[i].fI];
540     info.fPrimUsedFlag = 1;
541     info.fPrimDeviation = dev[i].fD;
542   }
543
544   for( Int_t i = 0; i<fNTracks; i++ ){
545     AliESDTrackInfo &info = fTrackInfos[i];
546     if( info.fPrimUsedFlag ) continue;
547     info.fPrimDeviation = info.fParticle.GetDeviationFromVertex( fPrimaryVtx );   
548   }
549
550   delete[] vSelected;
551   delete[] dev;
552 }
553
554
555
556 void AliHLTGlobalVertexerComponent::FindV0s( vector<pair<int,int> > &v0s  )
557 {
558   //* V0 finder
559
560   AliKFVertex &primVtx = fPrimaryVtx;
561   if( primVtx.GetNContributors()<3 ) return;
562
563   TStopwatch timer;
564   Int_t statN = 0;
565   Bool_t run = 1;
566
567   for( Int_t iTr = 0; iTr<fNTracks && run; iTr++ ){ //* first daughter
568
569     AliESDTrackInfo &info = fTrackInfos[iTr];
570     if( !info.fOK ) continue;    
571     if( info.fParticle.GetQ() >0 ) continue;    
572     if( info.fPrimDeviation < fV0DaughterPrimDeviation ) continue;
573
574     for( Int_t jTr = 0; jTr<fNTracks; jTr++ ){  //* second daughter
575       
576       
577       AliESDTrackInfo &jnfo = fTrackInfos[jTr];
578       if( !jnfo.fOK ) continue;
579       if( jnfo.fParticle.GetQ() < 0 ) continue;
580       if( jnfo.fPrimDeviation < fV0DaughterPrimDeviation ) continue;
581
582       // check the time once a while...
583
584       if( (++statN)%100 ==0 ){ 
585         if( timer.RealTime()>= fV0TimeLimit ){  run = 0; break; }
586         timer.Continue();
587       }
588
589       //* check if the particles fit
590
591       if( info.fParticle.GetDeviationFromParticle(jnfo.fParticle) > fV0Chi ) continue;
592
593       //* construct V0 mother
594
595       AliKFParticle v0( info.fParticle, jnfo.fParticle );     
596
597       //* check V0 Chi^2
598       
599       if( v0.GetChi2()<0 || v0.GetChi2() > fV0Chi*fV0Chi*v0.GetNDF() ) continue;
600
601       //* subtruct daughters from primary vertex 
602
603       AliKFVertex primVtxCopy = primVtx;    
604        
605       if( info.fPrimUsedFlag ){ 
606         if( primVtxCopy.GetNContributors()<=2 ) continue;
607         primVtxCopy -= info.fParticle;
608       }
609       if( jnfo.fPrimUsedFlag ){
610         if( primVtxCopy.GetNContributors()<=2 ) continue;
611         primVtxCopy -= jnfo.fParticle;
612       }
613       //* Check v0 Chi^2 deviation from primary vertex 
614
615       if( v0.GetDeviationFromVertex( primVtxCopy ) > fV0PrimDeviation ) continue;
616
617       //* Add V0 to primary vertex to improve the primary vertex resolution
618
619       primVtxCopy += v0;      
620
621       //* Set production vertex for V0
622
623       v0.SetProductionVertex( primVtxCopy );
624
625       //* Get V0 decay length with estimated error
626
627       Double_t length, sigmaLength;
628       if( v0.GetDecayLength( length, sigmaLength ) ) continue;
629
630       //* Reject V0 if it decays too close[sigma] to the primary vertex
631
632       if( length  < fV0DecayLengthInSigmas*sigmaLength ) continue;
633       
634       //* keep v0 
635
636       v0s.push_back(pair<int,int>(iTr,jTr));
637     }
638   }
639 }
640
641
642
643
644 void AliHLTGlobalVertexerComponent::FillESD( AliESDEvent *event, AliHLTGlobalVertexerData *data
645 )
646 {
647   //* put output of a vertexer to the esd event
648
649   Int_t nESDTracks = event->GetNumberOfTracks();
650
651   const int *listPrim = data->fTrackIndices;
652   const int *listV0 = data->fTrackIndices + data->fNPrimTracks;
653
654   std::map<int,int> mapId;
655   bool *constrainedToVtx   = new bool[nESDTracks];
656
657   for( int i=0; i<nESDTracks; i++ ){
658     constrainedToVtx[i] = 0;
659     if( !event->GetTrack(i) ) continue;
660     mapId[ event->GetTrack(i)->GetID() ] = i;
661   }
662
663   if( data->fPrimNContributors >=3 ){
664
665     AliESDVertex vESD( data->fPrimP, data->fPrimC, data->fPrimChi2, data->fPrimNContributors );
666     event->SetPrimaryVertexTPC( &vESD );
667     event->SetPrimaryVertexTracks( &vESD );
668
669     // relate tracks to the primary vertex
670
671     if( data->fFitTracksToVertex ){
672       for( Int_t i = 0; i<data->fNPrimTracks; i++ ){
673         Int_t id = listPrim[ i ];
674         map<int,int>::iterator it = mapId.find(id);
675         if( it==mapId.end() ) continue;
676         Int_t itr = it->second;
677         event->GetTrack(itr)->RelateToVertex( &vESD, event->GetMagneticField(),100. );
678         constrainedToVtx[ itr ] = 1;
679       }
680     }
681   }
682
683   //* add ESD v0s and relate tracks to v0s
684  
685
686   for( int i=0; i<data->fNV0s; i++ ){
687
688     Int_t id1 = listV0[ 2*i ];
689     Int_t id2 = listV0[ 2*i + 1];
690     map<int,int>::iterator it = mapId.find(id1);
691     if( it==mapId.end() ) continue;
692     Int_t iTr = it->second;
693     it = mapId.find(id2);
694     if( it==mapId.end() ) continue;
695     Int_t jTr = it->second;
696    
697     AliESDv0 v0( *event->GetTrack( iTr ), iTr, *event->GetTrack( jTr ), jTr );  
698     event->AddV0( &v0 );
699
700     // relate the tracks to the vertex  
701
702     if( data->fFitTracksToVertex ){
703       if( constrainedToVtx[iTr] || constrainedToVtx[jTr] ) continue;
704       double pos[3];
705       double sigma[3] = {.1,.1,.1};
706       v0.XvYvZv(pos);
707       AliESDVertex vESD(pos, sigma);
708       event->GetTrack(iTr)->RelateToVertex( &vESD, event->GetMagneticField(),100. );
709       event->GetTrack(jTr)->RelateToVertex( &vESD, event->GetMagneticField(),100. );
710       constrainedToVtx[iTr] = 1;
711       constrainedToVtx[jTr] = 1;    
712     }
713   }
714
715   delete[] constrainedToVtx;
716 }