CMS 3D CMS Logo

GsfTrackProducerBase Class Reference

Produce Tracks from TrackCandidates. More...

#include <RecoTracker/TrackProducer/interface/GsfTrackProducerBase.h>

Inheritance diagram for GsfTrackProducerBase:

TrackProducerBase< reco::GsfTrack > GsfTrackProducer GsfTrackRefitter

List of all members.

Public Member Functions

 GsfTrackProducerBase (bool trajectoryInEvent, bool split)
 Constructor.
virtual void putInEvt (edm::Event &, std::auto_ptr< TrackingRecHitCollection > &, std::auto_ptr< reco::GsfTrackCollection > &, std::auto_ptr< reco::TrackExtraCollection > &, std::auto_ptr< reco::GsfTrackExtraCollection > &, std::auto_ptr< std::vector< Trajectory > > &, AlgoProductCollection &, const reco::BeamSpot &)
 Put produced collections in the event.

Protected Member Functions

void fillMode (reco::GsfTrack &track, const TrajectoryStateOnSurface innertsos, const Propagator &gsfProp, const TransverseImpactPointExtrapolator &tipExtrapolator, TrajectoryStateClosestToBeamLineBuilder &tscblBuilder, const reco::BeamSpot &bs) const
void fillStates (TrajectoryStateOnSurface tsos, std::vector< reco::GsfComponent5D > &states) const

Private Attributes

bool useSplitting


Detailed Description

Produce Tracks from TrackCandidates.

Date
2008/03/03 15:06:40
Revision
1.8
Author:
cerati

Definition at line 27 of file GsfTrackProducerBase.h.


Constructor & Destructor Documentation

GsfTrackProducerBase::GsfTrackProducerBase ( bool  trajectoryInEvent,
bool  split 
) [inline, explicit]

Constructor.

Definition at line 31 of file GsfTrackProducerBase.h.

00031                                                                     :
00032     TrackProducerBase<reco::GsfTrack>(trajectoryInEvent),
00033     useSplitting(split){}


Member Function Documentation

void GsfTrackProducerBase::fillMode ( reco::GsfTrack track,
const TrajectoryStateOnSurface  innertsos,
const Propagator gsfProp,
const TransverseImpactPointExtrapolator tipExtrapolator,
TrajectoryStateClosestToBeamLineBuilder tscblBuilder,
const reco::BeamSpot bs 
) const [protected]

Definition at line 227 of file GsfTrackProducerBase.cc.

References FreeTrajectoryState::charge(), FreeTrajectoryState::curvilinearError(), reco::GsfTrack::dimensionMode, reco::BeamSpot::dxdz(), reco::BeamSpot::dydz(), TransverseImpactPointExtrapolator::extrapolate(), CurvilinearTrajectoryError::matrix(), SingleGaussianState1D::mean(), GaussianSumUtilities1D::mean(), GaussianSumUtilities1D::mode(), GaussianSumUtilities1D::modeIsValid(), FreeTrajectoryState::momentum(), MultiGaussianStateTransform::multiState1D(), reco::BeamSpot::position(), Propagator::propagate(), reco::GsfTrack::setMode(), funct::sqrt(), TrajectoryStateOnSurface::surface(), GaussianSumUtilities1D::variance(), SingleGaussianState1D::variance(), reco::TrackBase::vz(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by putInEvt().

00232 {
00233   // Get transverse impact parameter plane (from mean). This is a first approximation;
00234   // the mode is then extrapolated to the
00235   // final position closest to the beamline.
00236   GlobalPoint bsPos(bs.position().x()+(track.vz()-bs.position().z())*bs.dxdz(),
00237                     bs.position().y()+(track.vz()-bs.position().z())*bs.dydz(),
00238                     track.vz());
00239   TrajectoryStateOnSurface vtxTsos = tipExtrapolator.extrapolate(innertsos,bsPos);
00240   if ( !vtxTsos.isValid() )  vtxTsos = innertsos;
00241  // extrapolate mixture
00242   vtxTsos = gsfProp.propagate(innertsos,vtxTsos.surface());
00243   if ( vtxTsos.isValid() ) {
00244     // extract mode
00245     // build perigee parameters (for covariance to be stored)
00246     AlgebraicVector5 modeParameters;
00247     AlgebraicSymMatrix55 modeCovariance;
00248     // set parameters and variances for "mode" state (local parameters)
00249     for ( unsigned int iv=0; iv<5; ++iv ) {
00250       MultiGaussianState1D state1D = MultiGaussianStateTransform::multiState1D(vtxTsos,iv);
00251       GaussianSumUtilities1D utils(state1D);
00252       modeParameters(iv) = utils.mode().mean();
00253       modeCovariance(iv,iv) = utils.mode().variance();
00254       if ( !utils.modeIsValid() ) {
00255         // if mode calculation fails: use mean
00256         modeParameters(iv) = utils.mean();
00257         modeCovariance(iv,iv) = utils.variance();
00258       }
00259     }
00260     // complete covariance matrix
00261     // approximation: use correlations from mean
00262     const AlgebraicSymMatrix55& meanCovariance(vtxTsos.localError().matrix());
00263     for ( unsigned int iv1=0; iv1<5; ++iv1 ) {
00264       for ( unsigned int iv2=0; iv2<iv1; ++iv2 ) {
00265         double cov12 = meanCovariance(iv1,iv2) * 
00266           sqrt(modeCovariance(iv1,iv1)/meanCovariance(iv1,iv1)*
00267                modeCovariance(iv2,iv2)/meanCovariance(iv2,iv2));
00268         modeCovariance(iv1,iv2) = modeCovariance(iv2,iv1) = cov12;
00269       }
00270     }
00271     TrajectoryStateOnSurface modeTsos(LocalTrajectoryParameters(modeParameters,
00272                                                                 vtxTsos.localParameters().pzSign()),
00273                                       LocalTrajectoryError(modeCovariance),
00274                                       vtxTsos.surface(),
00275                                       vtxTsos.magneticField(),
00276                                       vtxTsos.surfaceSide());
00277     TrajectoryStateClosestToBeamLine tscbl = tscblBuilder(*modeTsos.freeState(),bs);
00278     if ( tscbl.isValid() ) {
00279       FreeTrajectoryState fts = tscbl.trackStateAtPCA();
00280       GlobalVector tscblMom = fts.momentum();
00281       reco::GsfTrack::Vector mom(tscblMom.x(),tscblMom.y(),tscblMom.z());
00282       reco::GsfTrack::CovarianceMatrixMode cov;
00283       const AlgebraicSymMatrix55& tscblCov = fts.curvilinearError().matrix();
00284       for ( unsigned int iv1=0; iv1<reco::GsfTrack::dimensionMode; ++iv1 ) {
00285         for ( unsigned int iv2=0; iv2<reco::GsfTrack::dimensionMode; ++iv2 ) {
00286           cov(iv1,iv2) = tscblCov(iv1,iv2);
00287         }
00288       } 
00289       track.setMode(fts.charge(),mom,cov);
00290     }
00291   }
00292 }

void GsfTrackProducerBase::fillStates ( TrajectoryStateOnSurface  tsos,
std::vector< reco::GsfComponent5D > &  states 
) const [protected]

Definition at line 196 of file GsfTrackProducerBase.cc.

References TrajectoryStateOnSurface::components(), and i.

Referenced by putInEvt().

00198 {
00199   //   std::cout << "in fill states" << std::endl;
00200   //   if ( !tsos.isValid() ) {
00201   //     std::cout << std::endl << std::endl << "invalid tsos" << std::endl;
00202   //     return;
00203   //   }
00204   reco::GsfComponent5D::ParameterVector pLocS;
00205   reco::GsfComponent5D::CovarianceMatrix cLocS;
00206   std::vector<TrajectoryStateOnSurface> components(tsos.components());
00207   for ( std::vector<TrajectoryStateOnSurface>::const_iterator i=components.begin();
00208         i!=components.end(); ++i ) {
00209     //     if ( !(*i).isValid() ) {
00210     //       std::cout << std::endl << "invalid component" << std::endl;
00211     //       continue;
00212     //     }
00213     // Unneeded hack ... now we have SMatrix in tracking too
00214     // const AlgebraicVector& pLoc = i->localParameters().vector();
00215     // for ( int j=0; j<reco::GsfTrackExtra::dimension; ++j )  pLocS(j) = pLoc[j];
00216     // const AlgebraicSymMatrix& cLoc = i->localError().matrix();
00217     // for ( int j1=0; j1<reco::GsfTrack::dimension; ++j1 )
00218     // for ( int j2=0; j2<=j1; ++j2 )  cLocS(j1,j2) = cLoc[j1][j2];
00219     // states.push_back(reco::GsfComponent5D(i->weight(),pLocS,cLocS));
00220     
00221     states.push_back(reco::GsfComponent5D(i->weight(),i->localParameters().vector(),i->localError().matrix()));
00222   }
00223   //   std::cout << "end fill states" << std::endl;
00224 }

void GsfTrackProducerBase::putInEvt ( edm::Event evt,
std::auto_ptr< TrackingRecHitCollection > &  selHits,
std::auto_ptr< reco::GsfTrackCollection > &  selTracks,
std::auto_ptr< reco::TrackExtraCollection > &  selTrackExtras,
std::auto_ptr< reco::GsfTrackExtraCollection > &  selGsfTrackExtras,
std::auto_ptr< std::vector< Trajectory > > &  selTrajectories,
AlgoProductCollection algoResults,
const reco::BeamSpot bs 
) [virtual]

Put produced collections in the event.

Definition at line 25 of file GsfTrackProducerBase.cc.

References reco::TrackExtraBase::add(), alongMomentum, AnalyticalPropagator_cfi::AnalyticalPropagator, anyDirection, TrajectoryStateOnSurface::components(), TrajectoryStateOnSurface::curvilinearError(), Trajectory::direction(), fillMode(), fillStates(), Trajectory::firstMeasurement(), edm::Event::getRefBeforePut(), TrajectoryStateOnSurface::globalParameters(), i, TrajectoryStateOnSurface::isValid(), it, j, Trajectory::lastMeasurement(), TrajectoryStateOnSurface::localParameters(), LogDebug, LogTrace, TrajectoryStateOnSurface::magneticField(), GlobalTrajectoryParameters::momentum(), p, GlobalTrajectoryParameters::position(), edm::Event::put(), LocalTrajectoryParameters::pzSign(), TrajectoryMeasurement::recHit(), Trajectory::recHits(), TrackProducerBase< reco::GsfTrack >::rTracks_, Trajectory::seedRef(), reco::Track::setExtra(), reco::GsfTrack::setGsfExtra(), reco::TrackBase::setHitPattern(), t, track, TrackProducerBase< reco::GsfTrack >::trajectoryInEvent_, TrajectoryMeasurement::updatedState(), useSplitting, v, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by GsfTrackProducer::produce(), and GsfTrackRefitter::produce().

00033 {
00034 
00035   TrackingRecHitRefProd rHits = evt.getRefBeforePut<TrackingRecHitCollection>();
00036   reco::TrackExtraRefProd rTrackExtras = evt.getRefBeforePut<reco::TrackExtraCollection>();
00037   reco::GsfTrackExtraRefProd rGsfTrackExtras = evt.getRefBeforePut<reco::GsfTrackExtraCollection>();
00038   reco::GsfTrackRefProd rTracks = evt.getRefBeforePut<reco::GsfTrackCollection>();
00039 
00040   edm::Ref<reco::TrackExtraCollection>::key_type idx = 0;
00041   edm::Ref<reco::TrackExtraCollection>::key_type hidx = 0;
00042   edm::Ref<reco::GsfTrackExtraCollection>::key_type idxGsf = 0;
00043 //   edm::Ref<reco::TrackCollection>::key_type iTkRef = 0;
00044 //   edm::Ref< std::vector<Trajectory> >::key_type iTjRef = 0;
00045 //   std::map<unsigned int, unsigned int> tjTkMap;
00046 
00047   TrajectoryStateClosestToBeamLineBuilder tscblBuilder;
00048 
00049   for(AlgoProductCollection::iterator i=algoResults.begin(); i!=algoResults.end();i++){
00050     Trajectory * theTraj = (*i).first;
00051     if(trajectoryInEvent_) {
00052       selTrajectories->push_back(*theTraj);
00053 //       iTjRef++;
00054     }
00055     // const TrajectoryFitter::RecHitContainer& transHits = theTraj->recHits(useSplitting);  // NO: the return type in Trajectory is by VALUE
00056     TrajectoryFitter::RecHitContainer transHits = theTraj->recHits(useSplitting);
00057     reco::GsfTrack * theTrack = (*i).second.first;
00058     PropagationDirection seedDir = (*i).second.second;  
00059     
00060     LogDebug("TrackProducer") << "In GsfTrackProducerBase::putInEvt - seedDir=" << seedDir;
00061 
00062     reco::GsfTrack t = * theTrack;
00063     selTracks->push_back( t );
00064 //     iTkRef++;
00065 
00066 //     // Store indices in local map (starts at 0)
00067 //     if(trajectoryInEvent_) tjTkMap[iTjRef-1] = iTkRef-1;
00068 
00069     //sets the outermost and innermost TSOSs
00070     TrajectoryStateOnSurface outertsos;
00071     TrajectoryStateOnSurface innertsos;
00072     unsigned int innerId, outerId;
00073 
00074     // ---  NOTA BENE: the convention is to sort hits and measurements "along the momentum".
00075     // This is consistent with innermost and outermost labels only for tracks from LHC collision
00076     if (theTraj->direction() == alongMomentum) {
00077       outertsos = theTraj->lastMeasurement().updatedState();
00078       innertsos = theTraj->firstMeasurement().updatedState();
00079       outerId = theTraj->lastMeasurement().recHit()->geographicalId().rawId();
00080       innerId = theTraj->firstMeasurement().recHit()->geographicalId().rawId();
00081     } else { 
00082       outertsos = theTraj->firstMeasurement().updatedState();
00083       innertsos = theTraj->lastMeasurement().updatedState();
00084       outerId = theTraj->firstMeasurement().recHit()->geographicalId().rawId();
00085       innerId = theTraj->lastMeasurement().recHit()->geographicalId().rawId();
00086     }
00087     //build the TrackExtra
00088     GlobalPoint v = outertsos.globalParameters().position();
00089     GlobalVector p = outertsos.globalParameters().momentum();
00090     math::XYZVector outmom( p.x(), p.y(), p.z() );
00091     math::XYZPoint  outpos( v.x(), v.y(), v.z() );
00092     v = innertsos.globalParameters().position();
00093     p = innertsos.globalParameters().momentum();
00094     math::XYZVector inmom( p.x(), p.y(), p.z() );
00095     math::XYZPoint  inpos( v.x(), v.y(), v.z() );
00096 
00097     reco::TrackExtraRef teref= reco::TrackExtraRef ( rTrackExtras, idx ++ );
00098     reco::GsfTrack & track = selTracks->back();
00099     track.setExtra( teref );
00100 
00101 
00102     selTrackExtras->push_back( reco::TrackExtra (outpos, outmom, true, inpos, inmom, true,
00103                                                  outertsos.curvilinearError(), outerId,
00104                                                  innertsos.curvilinearError(), innerId,
00105                                                  seedDir,theTraj->seedRef()));
00106 
00107 
00108     reco::TrackExtra & tx = selTrackExtras->back();
00109 
00110 
00111     size_t i = 0;
00112     // ---  NOTA BENE: the convention is to sort hits and measurements "along the momentum".
00113     // This is consistent with innermost and outermost labels only for tracks from LHC collisions
00114     if (theTraj->direction() == alongMomentum) {
00115       for( TrajectoryFitter::RecHitContainer::const_iterator j = transHits.begin();
00116            j != transHits.end(); j ++ ) {
00117         if ((**j).hit()!=0){
00118           TrackingRecHit * hit = (**j).hit()->clone();
00119           track.setHitPattern( * hit, i ++ );
00120           selHits->push_back( hit );
00121           tx.add( TrackingRecHitRef( rHits, hidx ++ ) );
00122         }
00123       }
00124     }else{
00125       for( TrajectoryFitter::RecHitContainer::const_iterator j = transHits.end()-1;
00126            j != transHits.begin()-1; --j ) {
00127         if ((**j).hit()!=0){
00128           TrackingRecHit * hit = (**j).hit()->clone();
00129           track.setHitPattern( * hit, i ++ );
00130           selHits->push_back( hit );
00131         tx.add( TrackingRecHitRef( rHits, hidx ++ ) );
00132         }
00133       }
00134     }
00135     // ----
00136 
00137 
00138     //build the GsfTrackExtra
00139     std::vector<reco::GsfComponent5D> outerStates;
00140     outerStates.reserve(outertsos.components().size());
00141     fillStates(outertsos,outerStates);
00142     std::vector<reco::GsfComponent5D> innerStates;
00143     innerStates.reserve(innertsos.components().size());
00144     fillStates(innertsos,innerStates);
00145 
00146     reco::GsfTrackExtraRef terefGsf = reco::GsfTrackExtraRef ( rGsfTrackExtras, idxGsf ++ );
00147     track.setGsfExtra( terefGsf );
00148     selGsfTrackExtras->push_back( reco::GsfTrackExtra (outerStates, outertsos.localParameters().pzSign(),
00149                                                        innerStates, innertsos.localParameters().pzSign()));
00150 
00151     if ( innertsos.isValid() ) {
00152       GsfPropagatorAdapter gsfProp(AnalyticalPropagator(innertsos.magneticField(),anyDirection));
00153       TransverseImpactPointExtrapolator tipExtrapolator(gsfProp);
00154       fillMode(track,innertsos,gsfProp,tipExtrapolator,tscblBuilder,bs);
00155     }
00156 
00157     delete theTrack;
00158     delete theTraj;
00159   }
00160 
00161   LogTrace("TrackingRegressionTest") << "========== TrackProducer Info ===================";
00162   LogTrace("TrackingRegressionTest") << "number of finalGsfTracks: " << selTracks->size();
00163   for (reco::GsfTrackCollection::const_iterator it = selTracks->begin(); it != selTracks->end(); it++) {
00164     LogTrace("TrackingRegressionTest") << "track's n valid and invalid hit, chi2, pt : " 
00165                                        << it->found() << " , " 
00166                                        << it->lost()  <<" , " 
00167                                        << it->normalizedChi2() << " , "
00168                                        << it->pt() << " , "
00169                                        << it->eta() ;
00170   }
00171   LogTrace("TrackingRegressionTest") << "=================================================";
00172   
00173 
00174   rTracks_ = evt.put( selTracks );
00175   evt.put( selTrackExtras );
00176   evt.put( selGsfTrackExtras );
00177   evt.put( selHits );
00178 
00179   if(trajectoryInEvent_) {
00180     edm::OrphanHandle<std::vector<Trajectory> > rTrajs = evt.put(selTrajectories);
00181 
00182 //     // Now Create traj<->tracks association map
00183 //     std::auto_ptr<TrajTrackAssociationCollection> trajTrackMap( new TrajTrackAssociationCollection() );
00184 //     for ( std::map<unsigned int, unsigned int>::iterator i = tjTkMap.begin(); 
00185 //           i != tjTkMap.end(); i++ ) {
00186 //       edm::Ref<std::vector<Trajectory> > trajRef( rTrajs, (*i).first );
00187 //       edm::Ref<reco::GsfTrackCollection>    tkRef( rTracks_, (*i).second );
00188 //       trajTrackMap->insert( edm::Ref<std::vector<Trajectory> >( rTrajs, (*i).first ),
00189 //                             edm::Ref<reco::GsfTrackCollection>( rTracks_, (*i).second ) );
00190 //     }
00191 //     evt.put( trajTrackMap );
00192   }
00193 }


Member Data Documentation

bool GsfTrackProducerBase::useSplitting [private]

Definition at line 55 of file GsfTrackProducerBase.h.

Referenced by putInEvt().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:23:18 2009 for CMSSW by  doxygen 1.5.4