CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

L3MuonTrajectoryBuilder Class Reference

#include <L3MuonTrajectoryBuilder.h>

Inheritance diagram for L3MuonTrajectoryBuilder:
GlobalTrajectoryBuilderBase MuonTrajectoryBuilder

List of all members.

Public Member Functions

 L3MuonTrajectoryBuilder (const edm::ParameterSet &, const MuonServiceProxy *)
 constructor with Parameter Set and MuonServiceProxy
virtual void setEvent (const edm::Event &)
 pass the Event to the algo at each event
MuonTrajectoryBuilder::CandidateContainer trajectories (const TrackCand &)
 reconstruct trajectories from standalone and tracker only Tracks
 ~L3MuonTrajectoryBuilder ()
 destructor

Private Member Functions

std::vector< TrackCandmakeTkCandCollection (const TrackCand &)
 make a TrackCand collection using tracker Track, Trajectory information

Private Attributes

edm::Handle
< reco::TrackCollection
allTrackerTracks
edm::InputTag theTkCollName
TrajectoryCleanertheTrajectoryCleaner

Detailed Description

class to build muon trajectory

Date:
2009/02/24 07:07:17
Revision:
1.9
Author:
N. Neumeister Purdue University
C. Liu Purdue University
A. Everett Purdue University

Definition at line 30 of file L3MuonTrajectoryBuilder.h.


Constructor & Destructor Documentation

L3MuonTrajectoryBuilder::L3MuonTrajectoryBuilder ( const edm::ParameterSet par,
const MuonServiceProxy *  service 
)
L3MuonTrajectoryBuilder::~L3MuonTrajectoryBuilder ( )

destructor

Definition at line 85 of file L3MuonTrajectoryBuilder.cc.

References theTrajectoryCleaner.


Member Function Documentation

vector< L3MuonTrajectoryBuilder::TrackCand > L3MuonTrajectoryBuilder::makeTkCandCollection ( const TrackCand staCand) [private, virtual]

make a TrackCand collection using tracker Track, Trajectory information

Implements GlobalTrajectoryBuilderBase.

Definition at line 178 of file L3MuonTrajectoryBuilder.cc.

References allTrackerTracks, python::rootplot::argparse::category, and position.

Referenced by trajectories().

                                                                                                               {

  const std::string category = "Muon|RecoMuon|L3MuonTrajectoryBuilder|makeTkCandCollection";

  vector<TrackCand> tkCandColl;
  
  vector<TrackCand> tkTrackCands;
  
  for ( unsigned int position = 0; position != allTrackerTracks->size(); ++position ) {
    reco::TrackRef tkTrackRef(allTrackerTracks,position);
    TrackCand tkCand = TrackCand((Trajectory*)(0),tkTrackRef);
    tkCandColl.push_back(tkCand);
  }

  for(vector<TrackCand>::const_iterator tk = tkCandColl.begin(); tk != tkCandColl.end() ; ++tk) { 
    edm::Ref<L3MuonTrajectorySeedCollection> l3seedRef = (*tk).second->seedRef().castTo<edm::Ref<L3MuonTrajectorySeedCollection> >() ;
    reco::TrackRef staTrack = l3seedRef->l2Track();
    if(staTrack == (staCand.second) ) tkTrackCands.push_back(*tk);
  }

  return tkTrackCands;
  
}
void L3MuonTrajectoryBuilder::setEvent ( const edm::Event event) [virtual]

pass the Event to the algo at each event

Reimplemented from GlobalTrajectoryBuilderBase.

Definition at line 92 of file L3MuonTrajectoryBuilder.cc.

References allTrackerTracks, python::rootplot::argparse::category, LogDebug, and theTkCollName.

                                                            {
  
  const std::string category = "Muon|RecoMuon|L3MuonTrajectoryBuilder|setEvent";
  
  GlobalTrajectoryBuilderBase::setEvent(event);
      
  // get tracker TrackCollection from Event
  event.getByLabel(theTkCollName,allTrackerTracks);
  LogDebug(category) 
      << "Found " << allTrackerTracks->size() 
      << " tracker Tracks with label "<< theTkCollName;  
  
}
MuonCandidate::CandidateContainer L3MuonTrajectoryBuilder::trajectories ( const TrackCand staCandIn) [virtual]

reconstruct trajectories from standalone and tracker only Tracks

Implements MuonTrajectoryBuilder.

Definition at line 109 of file L3MuonTrajectoryBuilder.cc.

References GlobalTrajectoryBuilderBase::build(), python::rootplot::argparse::category, LogDebug, makeTkCandCollection(), GlobalMuonTrackMatcher::match(), query::result, edm::second(), GlobalTrajectoryBuilderBase::thePtCut, and GlobalTrajectoryBuilderBase::trackMatcher().

                                                                                                {

  const std::string category = "Muon|RecoMuon|L3MuonTrajectoryBuilder|trajectories";
  
  // cut on muons with low momenta
  if ( (staCandIn).second->pt() < thePtCut || (staCandIn).second->innerMomentum().Rho() < thePtCut || (staCandIn).second->innerMomentum().R() < 2.5 ) return CandidateContainer();
  
  // convert the STA track into a Trajectory if Trajectory not already present
  TrackCand staCand(staCandIn);
  
  vector<TrackCand> trackerTracks;
  
  vector<TrackCand> regionalTkTracks = makeTkCandCollection(staCand);
  LogDebug(category) << "Found " << regionalTkTracks.size() << " tracks within region of interest";  
  
  // match tracker tracks to muon track
  trackerTracks = trackMatcher()->match(staCand, regionalTkTracks);
  
  LogDebug(category) << "Found " << trackerTracks.size() << " matching tracker tracks within region of interest";
  if ( trackerTracks.empty() ) return CandidateContainer();
  
  // build a combined tracker-muon MuonCandidate
  //
  // turn tkMatchedTracks into MuonCandidates
  //
  LogDebug(category) << "turn tkMatchedTracks into MuonCandidates";
  CandidateContainer tkTrajs;
  for (vector<TrackCand>::const_iterator tkt = trackerTracks.begin(); tkt != trackerTracks.end(); tkt++) {
    if ((*tkt).first != 0 && (*tkt).first->isValid()) {
      //
      MuonCandidate* muonCand = new MuonCandidate( 0 ,staCand.second,(*tkt).second, new Trajectory(*(*tkt).first));
      tkTrajs.push_back(muonCand);
      //      LogTrace(category) << "tpush";
      //
    } else {
      MuonCandidate* muonCand = new MuonCandidate( 0 ,staCand.second,(*tkt).second, 0);
      tkTrajs.push_back(muonCand);
    }
  }
    
  if ( tkTrajs.empty() )  {
    LogDebug(category) << "tkTrajs empty";
    return CandidateContainer();
  }
  
  CandidateContainer result = build(staCand, tkTrajs);  
  LogDebug(category) << "Found "<< result.size() << " L3Muons from one L2Cand";

  // free memory
  if ( staCandIn.first == 0) delete staCand.first;

  for( CandidateContainer::const_iterator it = tkTrajs.begin(); it != tkTrajs.end(); ++it) {
    if ( (*it)->trajectory() ) delete (*it)->trajectory();
    if ( (*it)->trackerTrajectory() ) delete (*it)->trackerTrajectory();
    if ( *it ) delete (*it);
  }
  tkTrajs.clear();  

  for ( vector<TrackCand>::const_iterator is = regionalTkTracks.begin(); is != regionalTkTracks.end(); ++is) {
    delete (*is).first;   
  }
  
  return result;
  
}

Member Data Documentation

Definition at line 56 of file L3MuonTrajectoryBuilder.h.

Referenced by makeTkCandCollection(), and setEvent().

Definition at line 55 of file L3MuonTrajectoryBuilder.h.

Referenced by L3MuonTrajectoryBuilder(), and setEvent().

Definition at line 53 of file L3MuonTrajectoryBuilder.h.

Referenced by L3MuonTrajectoryBuilder(), and ~L3MuonTrajectoryBuilder().