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. *
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. *
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. *
18 //***************************************************************************
20 /// @file AliHLTITSVertexerSPDComponent.cxx
21 /// @author Sergey Gorbunov <sergey.gorbunov@kip.uni-heidelberg.de>
23 /// @brief An ITS tracker processing component for the HLT
26 /////////////////////////////////////////////////////
28 // a ITS tracker processing component for the HLT //
30 /////////////////////////////////////////////////////
32 #include "AliHLTITSVertexerSPDComponent.h"
33 #include "AliHLTArray.h"
34 #include "AliExternalTrackParam.h"
35 #include "TStopwatch.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"
53 /** ROOT macro for the implementation of ROOT specific class methods */
54 ClassImp( AliHLTITSVertexerSPDComponent )
55 AliHLTITSVertexerSPDComponent::AliHLTITSVertexerSPDComponent()
67 // see header file for class documentation
69 // refer to README to build package
71 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
73 for( int i=0; i<9; i++ ) fSum[i] = 0;
80 AliHLTITSVertexerSPDComponent::AliHLTITSVertexerSPDComponent( const AliHLTITSVertexerSPDComponent& )
93 // see header file for class documentation
95 for( int i=0; i<9; i++ ) fSum[i] = 0;
101 HLTFatal( "copy constructor untested" );
104 AliHLTITSVertexerSPDComponent& AliHLTITSVertexerSPDComponent::operator=( const AliHLTITSVertexerSPDComponent& )
106 // see header file for class documentation
107 HLTFatal( "assignment operator untested" );
111 AliHLTITSVertexerSPDComponent::~AliHLTITSVertexerSPDComponent()
113 // see header file for class documentation
114 for( int i=0; i<9; i++ ){
123 // Public functions to implement AliHLTComponent's interface.
124 // These functions are required for the registration process
127 const char* AliHLTITSVertexerSPDComponent::GetComponentID()
129 // see header file for class documentation
130 return "ITSVertexerSPD";
133 void AliHLTITSVertexerSPDComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list )
135 // see header file for class documentation
137 list.push_back( kAliHLTDataTypeClusters|kAliHLTDataOriginITSSPD );
140 AliHLTComponentDataType AliHLTITSVertexerSPDComponent::GetOutputDataType()
142 // see header file for class documentation
143 return kAliHLTMultipleDataType;
146 int AliHLTITSVertexerSPDComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList)
149 // see header file for class documentation
151 tgtList.push_back( kAliHLTDataTypeESDVertex|kAliHLTDataOriginITS);
152 return tgtList.size();
155 void AliHLTITSVertexerSPDComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
157 // define guess for the output data size
158 constBase = 10000; // minimum size
159 inputMultiplier = 0.5; // size relative to input
162 AliHLTComponent* AliHLTITSVertexerSPDComponent::Spawn()
164 // see header file for class documentation
165 return new AliHLTITSVertexerSPDComponent;
168 void AliHLTITSVertexerSPDComponent::SetDefaultConfiguration()
170 // Set default configuration for the CA tracker component
171 // Some parameters can be later overwritten from the OCDB
183 int AliHLTITSVertexerSPDComponent::ReadConfigurationString( const char* arguments )
185 // Set configuration parameters for the CA tracker component from the string
188 if ( !arguments ) return iResult;
190 TString allArgs = arguments;
192 int bMissingParam = 0;
194 TObjArray* pTokens = allArgs.Tokenize( " " );
196 int nArgs = pTokens ? pTokens->GetEntries() : 0;
198 for ( int i = 0; i < nArgs; i++ ) {
199 argument = ( ( TObjString* )pTokens->At( i ) )->GetString();
200 if ( argument.IsNull() ) continue;
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] );
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 );
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 );
228 HLTError( "Unknown option \"%s\"", argument.Data() );
233 if ( bMissingParam ) {
234 HLTError( "Specifier missed for parameter \"%s\"", argument.Data() );
242 int AliHLTITSVertexerSPDComponent::ReadCDBEntry( const char* cdbEntry, const char* chainId )
244 // see header file for class documentation
246 const char* defaultNotify = "";
249 return 0;// need to add the HLT/ConfigITS/ITSTracker directory to cdb SG!!!
250 //cdbEntry = "HLT/ConfigITS/ITSTracker";
251 //defaultNotify = " (default)";
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());
259 HLTError( "cannot fetch object \"%s\" from CDB", cdbEntry );
263 TObjString* pString = dynamic_cast<TObjString*>( pEntry->GetObject() );
266 HLTError( "configuration object \"%s\" has wrong type, required TObjString", cdbEntry );
270 HLTInfo( "received configuration object string: \"%s\"", pString->GetString().Data() );
272 return ReadConfigurationString( pString->GetString().Data() );
276 int AliHLTITSVertexerSPDComponent::Configure( const char* cdbEntry, const char* chainId, const char *commandLine )
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
282 //* read hard-coded values
284 SetDefaultConfiguration();
286 //* read the default CDB entry
288 int iResult1 = ReadCDBEntry( NULL, chainId );
290 //* read the actual CDB entry if required
292 int iResult2 = ( cdbEntry ) ? ReadCDBEntry( cdbEntry, chainId ) : 0;
294 //* read extra parameters from input (if they are)
298 if ( commandLine && commandLine[0] != '\0' ) {
299 HLTInfo( "received configuration string from HLT framework: \"%s\"", commandLine );
300 iResult3 = ReadConfigurationString( commandLine );
303 // Initialise the tracker here
305 return iResult1 ? iResult1 : ( iResult2 ? iResult2 : iResult3 );
310 int AliHLTITSVertexerSPDComponent::DoInit( int argc, const char** argv )
312 // Configure the ITS tracker component
314 SetDefaultConfiguration();
316 if(AliGeomManager::GetGeometry()==NULL){
317 AliGeomManager::LoadGeometry("");
319 AliGeomManager::ApplyAlignObjsFromCDB("ITS");
321 TString arguments = "";
322 for ( int i = 0; i < argc; i++ ) {
323 if ( !arguments.IsNull() ) arguments += " ";
324 arguments += argv[i];
327 int ret = Configure( NULL, NULL, arguments.Data() );
329 for( int i=0; i<9; i++ ) delete[] fSum[i];
334 fNZBins = ( (int) (2*fZRange/fZBinSize ))+1;
336 for( int i=0; i<9; i++ ) fSum[i] = new double [fNZBins];
338 fSumW = new double [fNZBins];
339 fSumN = new int [fNZBins];
345 int AliHLTITSVertexerSPDComponent::DoDeinit()
347 // see header file for class documentation
349 for( int i=0; i<9; i++ ){
358 SetDefaultConfiguration();
365 int AliHLTITSVertexerSPDComponent::Reconfigure( const char* cdbEntry, const char* chainId )
367 // Reconfigure the component from OCDB .
369 return Configure( cdbEntry, chainId, NULL );
373 int AliHLTITSVertexerSPDComponent::DoEvent
375 const AliHLTComponentEventData& evtData,
376 const AliHLTComponentBlockData* blocks,
377 AliHLTComponentTriggerData& /*trigData*/,
378 AliHLTUInt8_t* /*outputPtr*/,
379 AliHLTUInt32_t& size,
380 vector<AliHLTComponentBlockData>& /*outputBlocks*/ )
384 //AliHLTUInt32_t maxBufferSize = size;
385 size = 0; // output size
387 if (!IsDataEvent()) return 0;
389 if ( evtData.fBlockCnt <= 0 ) {
390 HLTWarning( "no blocks in event" );
396 // Event reconstruction in ITS
400 int nBlocks = evtData.fBlockCnt;
401 int nClustersTotal = 0;
404 const int kNPhiBins = 20;
405 const double kDPhi = TMath::TwoPi() / kNPhiBins;
406 double vtxX = fRunVtx[0], vtxY = fRunVtx[1], vtxZ = fRunVtx[2];
408 std::vector<AliHLTITSVertexerSPDComponent::AliHLTITSVZCluster> clusters[2][ kNPhiBins ];
410 for (int ndx=0; ndx<nBlocks && iResult>=0; ndx++) {
412 const AliHLTComponentBlockData* iter = blocks+ndx;
414 // Read ITS SPD clusters
416 if ( iter->fDataType == (kAliHLTDataTypeClusters|kAliHLTDataOriginITSSPD) ){
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];
431 AliITSRecPoint p( lab, hit, info );
434 HLTWarning("Can not misalign an SPD cluster");
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;
449 clusters[d.fLayer][iphi].push_back( c );
452 }// end read input blocks
455 // Reconstruct the vertex
457 TStopwatch timerReco;
459 double zScale = 1./fZBinSize;
461 // fit few times decreasing the radius cut
463 double vtxDR2[5] = { 1.*1., .5*.5, .4*.4, .4*.4, .2*.2};
464 double vtxDZ [5] = { 1., .5, .3, .3, .3};
468 for( int iter = 0; iter<3; iter++ ){
470 bool doZBins = (iter<2);
471 int nSearchBins = doZBins ?fNZBins :1;
473 for( int i=0; i<nSearchBins; i++ ){
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;
506 double l2 = 1./(dx*dx + dy*dy);
507 double a = x1*y0 - x0*y1;
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;
516 zbin = (int)((zv - fZMin)*zScale);
517 if( zbin<0 || zbin>=fNZBins ) continue;
520 if( TMath::Abs( zv - vtxZ ) > vtxDZ[ iter ] ) continue;
529 if( bestDn < 0 ) continue;
531 double dx = bestP[0] - x0;
532 double dy = bestP[1] - y0;
533 double dz = bestP[2] - z0;
534 double w = 1./(1.+bestR2);
538 // A_i*x + B_i*y + C_i*z = D_i
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
550 double n = w / TMath::Sqrt(dx*dx + dy*dy + dz*dz);
555 a = dy; b = -dx; c = 0; d = dy*x0 - dx*y0;
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;
567 a = dz; b = 0; c = -dx; d = dz*x0 - dx*z0;
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;
580 a = 0; b = dz; c = -dy; d = dz*y0 - dy*z0;
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;
593 a = dz; b = 0; c = -dx; d = dz*x0 - dx*z0;
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;
613 for( int i=0; i<nSearchBins; i++ ){
619 if( bestBin<0 || fSumN[bestBin] < 3 ){
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] };
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];
637 Double_t s = ( w[0]*wI[0] + w[1]*wI[1] + w[3]*wI[3] );
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;
663 double pos[3] = {vtxX, vtxY, vtxZ};
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 );
669 //cout<<"ITSVertexerSPD: vtx found "<<vtxX<<" "<<vtxY<<" "<<vtxZ<<endl;
673 fFullTime += timer.RealTime();
674 fRecoTime += timerReco.RealTime();
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 );