00001 #include "Fireworks/ParticleFlow/interface/FWPFTrackUtils.h"
00002
00003 FWPFTrackSingleton *FWPFTrackSingleton::pInstance = NULL;
00004 bool FWPFTrackSingleton::instanceFlag = false;
00005
00006
00007 FWPFTrackSingleton *
00008 FWPFTrackSingleton::Instance()
00009 {
00010 if( !instanceFlag )
00011 {
00012 pInstance = new FWPFTrackSingleton();
00013 instanceFlag = true;
00014 }
00015
00016 return pInstance;
00017 }
00018
00019
00020 void
00021 FWPFTrackSingleton::initPropagator()
00022 {
00023 m_magField = new FWMagField();
00024
00025
00026 m_trackPropagator = new TEveTrackPropagator();
00027 m_trackPropagator->SetMagFieldObj( m_magField, false );
00028 m_trackPropagator->SetMaxR( FWPFGeom::caloR3() );
00029 m_trackPropagator->SetMaxZ( FWPFGeom::caloZ2() );
00030 m_trackPropagator->SetDelta( 0.01 );
00031 m_trackPropagator->SetProjTrackBreaking( TEveTrackPropagator::kPTB_UseLastPointPos );
00032 m_trackPropagator->SetRnrPTBMarkers( kTRUE );
00033 m_trackPropagator->IncDenyDestroy();
00034
00035
00036 m_trackerTrackPropagator = new TEveTrackPropagator();
00037 m_trackerTrackPropagator->SetStepper( TEveTrackPropagator::kRungeKutta );
00038 m_trackerTrackPropagator->SetMagFieldObj( m_magField, false );
00039 m_trackerTrackPropagator->SetDelta( 0.01 );
00040 m_trackerTrackPropagator->SetMaxR( FWPFGeom::caloR3() );
00041 m_trackerTrackPropagator->SetMaxZ( FWPFGeom::caloZ2() );
00042 m_trackerTrackPropagator->SetProjTrackBreaking( TEveTrackPropagator::kPTB_UseLastPointPos );
00043 m_trackerTrackPropagator->SetRnrPTBMarkers( kTRUE );
00044 m_trackerTrackPropagator->IncDenyDestroy();
00045 }
00046
00047
00048 FWPFTrackUtils::FWPFTrackUtils()
00049 {
00050 m_singleton = FWPFTrackSingleton::Instance();
00051 }
00052
00053
00054
00055 TEveTrack *
00056 FWPFTrackUtils::getTrack( const reco::Track &iData )
00057 {
00058 TEveTrackPropagator *propagator = ( !iData.extra().isAvailable() ) ?
00059 m_singleton->getTrackerTrackPropagator() : m_singleton->getTrackPropagator();
00060
00061 TEveRecTrack t;
00062 t.fBeta = 1;
00063 t.fP = TEveVector( iData.px(), iData.py(), iData.pz() );
00064 t.fV = TEveVector( iData.vertex().x(), iData.vertex().y(), iData.vertex().z() );
00065 t.fSign = iData.charge();
00066 TEveTrack *trk = new TEveTrack( &t, propagator );
00067 trk->MakeTrack();
00068
00069 return trk;
00070 }
00071
00072
00073 TEveStraightLineSet *
00074 FWPFTrackUtils::setupLegoTrack( const reco::Track &iData )
00075 {
00076 using namespace FWPFGeom;
00077
00078
00079 int wraps[3] = { -1, -1, -1 };
00080 bool ECAL = false;
00081 TEveTrack *trk = getTrack( iData );
00082 std::vector<TEveVector> trackPoints( trk->GetN() - 1 );
00083 const Float_t *points = trk->GetP();
00084 TEveStraightLineSet *legoTrack = new TEveStraightLineSet();
00085
00086 if( m_singleton->getField()->getSource() == FWMagField::kNone )
00087 {
00088 if( fabs( iData.eta() ) < 2.0 && iData.pt() > 0.5 && iData.pt() < 30 )
00089 {
00090 double estimate = fw::estimate_field( iData, true );
00091 if( estimate >= 0 ) m_singleton->getField()->guessField( estimate );
00092 }
00093 }
00094
00095
00096 for( Int_t i = 1; i < trk->GetN(); ++i )
00097 {
00098 int j = i * 3;
00099 TEveVector temp = TEveVector( points[j], points[j+1], points[j+2] );
00100 TEveVector vec = TEveVector( temp.Eta(), temp.Phi(), 0.001 );
00101
00102 trackPoints[i-1] = vec;
00103 }
00104
00105
00106 for( Int_t i = 1; i < trk->GetN(); ++i )
00107 {
00108 int j = i * 3;
00109 TEveVector v1 = TEveVector( points[j], points[j+1], points[j+2] );
00110
00111 if( !ECAL )
00112 {
00113 if( FWPFMaths::checkIntersect( v1, caloR1() ) )
00114 {
00115 TEveVector v2 = TEveVector( points[j-3], points[j-2], points[j-1] );
00116 TEveVector xyPoint = FWPFMaths::lineCircleIntersect( v1, v2, caloR1() );
00117 TEveVector zPoint;
00118 if( v1.fZ < 0 )
00119 zPoint = TEveVector( xyPoint.fX, xyPoint.fY, v1.fZ - 50.f );
00120 else
00121 zPoint = TEveVector( xyPoint.fX, xyPoint.fY, v1.fZ + 50.f );
00122
00123 TEveVector vec = FWPFMaths::lineLineIntersect( v1, v2, xyPoint, zPoint );
00124 legoTrack->AddMarker( vec.Eta(), vec.Phi(), 0.001, 0 );
00125
00126 wraps[0] = i;
00127 ECAL = true;
00128 }
00129 else if( fabs( v1.fZ ) >= caloZ1() )
00130 {
00131 TEveVector p1, p2;
00132 TEveVector vec, v2 = TEveVector( points[j-3], points[j-2], points[j-1] );
00133 float z, y = FWPFMaths::linearInterpolation( v2, v1, caloZ1() );
00134
00135 if( v2.fZ > 0 )
00136 z = caloZ1();
00137 else
00138 z = caloZ1() * -1;
00139
00140 p1 = TEveVector( v2.fX + 50, y, z );
00141 p2 = TEveVector( v2.fX - 50, y, z );
00142 vec = FWPFMaths::lineLineIntersect( v1, v2, p1, p2 );
00143
00144 legoTrack->AddMarker( vec.Eta(), vec.Phi(), 0.001, 0 );
00145 wraps[0] = i;
00146 ECAL = true;
00147 }
00148 }
00149 else if( FWPFMaths::checkIntersect( v1, caloR2() ) )
00150 {
00151 TEveVector v2 = TEveVector( points[j-3], points[j-2], points[j-1] );
00152 TEveVector xyPoint = FWPFMaths::lineCircleIntersect( v1, v2, caloR2() );
00153 TEveVector zPoint;
00154 if( v1.fZ < 0 )
00155 zPoint = TEveVector( xyPoint.fX, xyPoint.fY, v1.fZ - 50.f );
00156 else
00157 zPoint = TEveVector( xyPoint.fX, xyPoint.fY, v1.fZ + 50.f );
00158
00159 TEveVector vec = FWPFMaths::lineLineIntersect( v1, v2, xyPoint, zPoint );
00160 legoTrack->AddMarker( vec.Eta(), vec.Phi(), 0.001, 1 );
00161
00162 wraps[1] = i;
00163 break;
00164 }
00165 }
00166
00167 if( wraps[0] != -1 )
00168 {
00169 int i = ( trk->GetN() - 1 ) * 3;
00170 int j = trk->GetN() - 2;
00171 TEveVector vec = TEveVector( points[i], points[i+1], points[i+2] );
00172
00173 if( FWPFMaths::checkIntersect( vec, caloR3() - 1 ) )
00174 {
00175 legoTrack->AddMarker( vec.Eta(), vec.Phi(), 0.001, 2 );
00176 wraps[2] = j;
00177 }
00178 else if( fabs( vec.fZ ) >= caloZ2() )
00179 {
00180 legoTrack->AddMarker( vec.Eta(), vec.Phi(), 0.001, 2 );
00181 wraps[2] = j;
00182 }
00183 }
00184
00185
00186 for( int i = 0; i < static_cast<int>( trackPoints.size() - 1 ); ++i )
00187 {
00188 if( ( trackPoints[i+1].fY - trackPoints[i].fY ) > 1 )
00189 {
00190 trackPoints[i+1].fY -= TMath::TwoPi();
00191 if( i == wraps[0] )
00192 {
00193 TEveChunkManager::iterator mi( legoTrack->GetMarkerPlex() );
00194 mi.next();
00195 TEveStraightLineSet::Marker_t &m = * ( TEveStraightLineSet::Marker_t* ) mi();
00196 m.fV[0] = trackPoints[i+1].fX; m.fV[1] = trackPoints[i+1].fY; m.fV[2] = 0.001;
00197 }
00198 else if( i == wraps[1] )
00199 {
00200 TEveChunkManager::iterator mi( legoTrack->GetMarkerPlex() );
00201 mi.next(); mi.next();
00202 TEveStraightLineSet::Marker_t &m = * ( TEveStraightLineSet::Marker_t* ) mi();
00203 m.fV[0] = trackPoints[i+1].fX; m.fV[1] = trackPoints[i+1].fY; m.fV[2] = 0.001;
00204 }
00205 }
00206
00207 if( ( trackPoints[i].fY - trackPoints[i+1].fY ) > 1 )
00208 {
00209 trackPoints[i+1].fY += TMath::TwoPi();
00210 if( i == wraps[0] )
00211 {
00212 TEveChunkManager::iterator mi( legoTrack->GetMarkerPlex() );
00213 mi.next();
00214 TEveStraightLineSet::Marker_t &m = * ( TEveStraightLineSet::Marker_t* ) mi();
00215 m.fV[0] = trackPoints[i+1].fX; m.fV[1] = trackPoints[i+1].fY; m.fV[2] = 0.001;
00216 }
00217 else if( i == wraps[1] )
00218 {
00219 TEveChunkManager::iterator mi( legoTrack->GetMarkerPlex() );
00220 mi.next(); mi.next();
00221 TEveStraightLineSet::Marker_t &m = * ( TEveStraightLineSet::Marker_t* ) mi();
00222 m.fV[0] = trackPoints[i+1].fX; m.fV[1] = trackPoints[i+1].fY; m.fV[2] = 0.001;
00223 }
00224 }
00225 }
00226
00227 int end = static_cast<int>( trackPoints.size() - 1 );
00228 if( wraps[2] == end )
00229 {
00230 TEveChunkManager::iterator mi( legoTrack->GetMarkerPlex() );
00231 mi.next(); mi.next(); mi.next();
00232 TEveStraightLineSet::Marker_t &m = * ( TEveStraightLineSet::Marker_t* ) mi();
00233 m.fV[0] = trackPoints[end].fX; m.fV[1] = trackPoints[end].fY; m.fV[2] = 0.001;
00234 }
00235
00236
00237 for( unsigned int i = 0;i < trackPoints.size() - 1; ++i )
00238 legoTrack->AddLine( trackPoints[i], trackPoints[i+1] );
00239
00240 legoTrack->SetDepthTest( false );
00241 legoTrack->SetMarkerStyle( 4 );
00242 legoTrack->SetMarkerSize( 1 );
00243 legoTrack->SetRnrMarkers( true );
00244
00245 delete trk;
00246
00247 return legoTrack;
00248 }
00249
00250
00251 TEveTrack *
00252 FWPFTrackUtils::setupTrack( const reco::Track &iData )
00253 {
00254 if( m_singleton->getField()->getSource() == FWMagField::kNone )
00255 {
00256 if( fabs( iData.eta() ) < 2.0 && iData.pt() > 0.5 && iData.pt() < 30 )
00257 {
00258 double estimate = fw::estimate_field( iData, true );
00259 if( estimate >= 0 ) m_singleton->getField()->guessField( estimate );
00260 }
00261 }
00262
00263 TEveTrack *trk = getTrack( iData );
00264
00265 return trk;
00266 }
00267
00268
00269
00270 TEvePointSet *
00271 FWPFTrackUtils::getCollisionMarkers( const TEveTrack *trk )
00272 {
00273 using namespace FWPFGeom;
00274
00275 bool ECAL = false;
00276 const Float_t *points = trk->GetP();
00277 TEvePointSet *ps = new TEvePointSet();
00278
00279 for( Int_t i = 1; i < trk->GetN(); ++i )
00280 {
00281 int j = i * 3;
00282 TEveVector v1 = TEveVector( points[j], points[j+1], points[j+2] );
00283
00284 if( !ECAL )
00285 {
00286 if( FWPFMaths::checkIntersect( v1, caloR1() ) )
00287 {
00288 TEveVector v2 = TEveVector( points[j-3], points[j-2], points[j-1] );
00289 TEveVector xyPoint = FWPFMaths::lineCircleIntersect( v1, v2, caloR1() );
00290 TEveVector zPoint;
00291 if( v1.fZ < 0 )
00292 zPoint = TEveVector( xyPoint.fX, xyPoint.fY, v1.fZ - 50.f );
00293 else
00294 zPoint = TEveVector( xyPoint.fX, xyPoint.fY, v1.fZ + 50.f );
00295
00296 TEveVector vec = FWPFMaths::lineLineIntersect( v1, v2, xyPoint, zPoint );
00297 ps->SetNextPoint( vec.fX, vec.fY, vec.fZ );
00298
00299 ECAL = true;
00300 }
00301 else if( fabs( v1.fZ ) >= caloZ1() )
00302 {
00303 TEveVector p1, p2;
00304 TEveVector vec, v2 = TEveVector( points[j-3], points[j-2], points[j-1] );
00305 float z, y = FWPFMaths::linearInterpolation( v2, v1, caloZ1() );
00306
00307 if( v2.fZ > 0 )
00308 z = caloZ1();
00309 else
00310 z = caloZ1() * -1;
00311
00312 p1 = TEveVector( v2.fX + 50, y, z );
00313 p2 = TEveVector( v2.fX - 50, y, z );
00314 vec = FWPFMaths::lineLineIntersect( v1, v2, p1, p2 );
00315
00316 ps->SetNextPoint( vec.fX, vec.fY, vec.fZ );
00317 ECAL = true;
00318 }
00319 }
00320 else if ( FWPFMaths::checkIntersect( v1, caloR2() ) )
00321 {
00322 TEveVector v2 = TEveVector( points[j-3], points[j-2], points[j-1] );
00323 TEveVector xyPoint = FWPFMaths::lineCircleIntersect( v1, v2, caloR2() );
00324 TEveVector zPoint;
00325 if( v1.fZ < 0 )
00326 zPoint = TEveVector( xyPoint.fX, xyPoint.fY, v1.fZ - 50.f );
00327 else
00328 zPoint = TEveVector( xyPoint.fX, xyPoint.fY, v1.fZ + 50.f );
00329
00330 TEveVector vec = FWPFMaths::lineLineIntersect( v1, v2, xyPoint, zPoint );
00331 ps->SetNextPoint( vec.fX, vec.fY, vec.fZ );
00332 break;
00333 }
00334 }
00335
00336
00337 int i = ( trk->GetN() - 1 ) * 3;
00338 TEveVector vec = TEveVector( points[i], points[i+1], points[i+2] );
00339
00340 if( FWPFMaths::checkIntersect( vec, caloR3() - 1 ) )
00341 ps->SetNextPoint( vec.fX, vec.fY, vec.fZ );
00342 else if( fabs( vec.fZ ) >= caloZ2() )
00343 ps->SetNextPoint( vec.fX, vec.fY, vec.fZ );
00344
00345 return ps;
00346 }