CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Types | Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | Private Attributes
GlobalTrajectoryBuilderBase Class Referenceabstract

#include <GlobalTrajectoryBuilderBase.h>

Inheritance diagram for GlobalTrajectoryBuilderBase:
MuonTrajectoryBuilder GlobalMuonTrajectoryBuilder L3MuonTrajectoryBuilder

Classes

struct  ComparatorInOut
 

Public Types

typedef
MuonTransientTrackingRecHit::ConstMuonRecHitContainer 
ConstMuonRecHitContainer
 
typedef
MuonTransientTrackingRecHit::ConstMuonRecHitPointer 
ConstMuonRecHitPointer
 
typedef
TransientTrackingRecHit::ConstRecHitContainer 
ConstRecHitContainer
 
typedef
TransientTrackingRecHit::ConstRecHitPointer 
ConstRecHitPointer
 
typedef
MuonTransientTrackingRecHit::MuonRecHitContainer 
MuonRecHitContainer
 
typedef
MuonTransientTrackingRecHit::MuonRecHitPointer 
MuonRecHitPointer
 
typedef
TransientTrackingRecHit::RecHitContainer 
RecHitContainer
 
typedef
TransientTrackingRecHit::RecHitPointer 
RecHitPointer
 
typedef std::vector< TrajectoryTC
 
typedef TC::const_iterator TI
 
- Public Types inherited from MuonTrajectoryBuilder
typedef
MuonCandidate::CandidateContainer 
CandidateContainer
 
typedef std::pair< const
Trajectory *, reco::TrackRef
TrackCand
 
typedef
MuonCandidate::TrajectoryContainer 
TrajectoryContainer
 

Public Member Functions

 GlobalTrajectoryBuilderBase (const edm::ParameterSet &, const MuonServiceProxy *, edm::ConsumesCollector &)
 constructor with Parameter Set and MuonServiceProxy More...
 
virtual void setEvent (const edm::Event &)
 pass the Event to the algo at each event More...
 
MuonTrajectoryBuilder::TrajectoryContainer trajectories (const TrajectorySeed &)
 dummy implementation, unused in this class More...
 
virtual ~GlobalTrajectoryBuilderBase ()
 destructor More...
 
- Public Member Functions inherited from MuonTrajectoryBuilder
 MuonTrajectoryBuilder ()
 constructor More...
 
virtual CandidateContainer trajectories (const TrackCand &)=0
 return a container reconstructed muons starting from a given track More...
 
virtual ~MuonTrajectoryBuilder ()
 destructor More...
 

Protected Types

enum  RefitDirection { inToOut, outToIn, undetermined }
 

Protected Member Functions

void addTraj (TrackCand &)
 if TrackCand has only a TrackRef, attempt to add Trajectory* More...
 
MuonTrajectoryBuilder::CandidateContainer build (const TrackCand &, MuonTrajectoryBuilder::CandidateContainer &) const
 build combined trajectory from sta Track and tracker RecHits More...
 
void checkMuonHits (const reco::Track &, ConstRecHitContainer &, ConstRecHitContainer &, std::vector< int > &) const
 check muon RecHits, calculate chamber occupancy and select hits to be used in the final fit More...
 
RefitDirection checkRecHitsOrdering (const ConstRecHitContainer &) const
 This does nothing now. More...
 
std::vector< TrackCandchooseRegionalTrackerTracks (const TrackCand &, const std::vector< TrackCand > &)
 choose tracker tracks within region of interest More...
 
const TrajectorychooseTrajectory (const std::vector< Trajectory * > &, int) const
 choose final trajectory More...
 
RectangularEtaPhiTrackingRegion defineRegionOfInterest (const reco::TrackRef &) const
 define region of interest with tracker More...
 
void fixTEC (ConstRecHitContainer &all, double scl_x, double scl_y) const
 rescale errors of outermost TEC RecHit More...
 
TransientTrackingRecHit::ConstRecHitContainer getTransientRecHits (const reco::Track &) const
 get transient RecHits of a Track More...
 
virtual std::vector< TrackCandmakeTkCandCollection (const TrackCand &)=0
 make a TrackCand collection using tracker Track, Trajectory information More...
 
void printHits (const ConstRecHitContainer &) const
 print all RecHits of a trajectory More...
 
ConstRecHitContainer selectMuonHits (const Trajectory &, const std::vector< int > &) const
 select muon hits compatible with trajectory; check hits in chambers with showers More...
 
ConstRecHitContainer selectTrackerHits (const ConstRecHitContainer &) const
 select tracker hits; exclude some tracker hits in the global trajectory More...
 
const MuonServiceProxyservice () const
 
GlobalMuonTrackMatchertrackMatcher () const
 
double trackProbability (const Trajectory &) const
 calculate chi2 probability (-ln(P)) More...
 

Protected Attributes

std::string theCategory
 
float thePCut
 
float thePtCut
 

Private Attributes

unsigned long long theCacheId_TRH
 
const edm::EventtheEvent
 
GlobalMuonRefittertheGlbRefitter
 
MuonDetLayerMeasurementstheLayerMeasurements
 
int theMuonHitsOption
 
edm::ESHandle
< TransientTrackingRecHitBuilder
theMuonRecHitBuilder
 
std::string theMuonRecHitBuilderName
 
MuonTrackingRegionBuildertheRegionBuilder
 
bool theRPCInTheFit
 
const MuonServiceProxytheService
 
float theTECxScale
 
float theTECyScale
 
std::string theTrackerPropagatorName
 
edm::ESHandle
< TransientTrackingRecHitBuilder
theTrackerRecHitBuilder
 
std::string theTrackerRecHitBuilderName
 
GlobalMuonTrackMatchertheTrackMatcher
 
TrackTransformertheTrackTransformer
 
const TrackerTopologytTopo_
 

Detailed Description

Base class for GlobalMuonTrajectoryBuilder and L3MuonTrajectoryBuilder Provide common tools and interface to reconstruct muons starting from a muon track reconstructed in the standalone muon system (with DT, CSC and RPC information). It tries to reconstruct the corresponding track in the tracker and performs matching between the reconstructed tracks in the muon system and the tracker.

Author
N. Neumeister Purdue University
C. Liu Purdue University
A. Everett Purdue University

Definition at line 49 of file GlobalTrajectoryBuilderBase.h.

Member Typedef Documentation

Definition at line 61 of file GlobalTrajectoryBuilderBase.h.

Definition at line 59 of file GlobalTrajectoryBuilderBase.h.

Definition at line 54 of file GlobalTrajectoryBuilderBase.h.

Definition at line 56 of file GlobalTrajectoryBuilderBase.h.

Definition at line 60 of file GlobalTrajectoryBuilderBase.h.

Definition at line 58 of file GlobalTrajectoryBuilderBase.h.

Definition at line 53 of file GlobalTrajectoryBuilderBase.h.

Definition at line 55 of file GlobalTrajectoryBuilderBase.h.

Definition at line 63 of file GlobalTrajectoryBuilderBase.h.

typedef TC::const_iterator GlobalTrajectoryBuilderBase::TI

Definition at line 64 of file GlobalTrajectoryBuilderBase.h.

Member Enumeration Documentation

Constructor & Destructor Documentation

GlobalTrajectoryBuilderBase::GlobalTrajectoryBuilderBase ( const edm::ParameterSet par,
const MuonServiceProxy service,
edm::ConsumesCollector iC 
)

constructor with Parameter Set and MuonServiceProxy

Definition at line 86 of file GlobalTrajectoryBuilderBase.cc.

References edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), GlobalMuonRefitter_cff::GlobalMuonRefitter, GlobalMuonTrackMatcher_cff::GlobalMuonTrackMatcher, theCacheId_TRH, theCategory, theGlbRefitter, theMuonHitsOption, theMuonRecHitBuilderName, thePCut, thePtCut, theRegionBuilder, theRPCInTheFit, theService, theTECxScale, theTECyScale, theTrackerPropagatorName, theTrackerRecHitBuilderName, theTrackMatcher, and theTrackTransformer.

88  :
90 
91  theCategory = par.getUntrackedParameter<string>("Category", "Muon|RecoMuon|GlobalMuon|GlobalTrajectoryBuilderBase");
92 
93 
94  ParameterSet trackMatcherPSet = par.getParameter<ParameterSet>("GlobalMuonTrackMatcher");
95  theTrackMatcher = new GlobalMuonTrackMatcher(trackMatcherPSet,theService);
96 
97  theTrackerPropagatorName = par.getParameter<string>("TrackerPropagator");
98 
99  ParameterSet trackTransformerPSet = par.getParameter<ParameterSet>("TrackTransformer");
100  theTrackTransformer = new TrackTransformer(trackTransformerPSet);
101 
102  ParameterSet regionBuilderPSet = par.getParameter<ParameterSet>("MuonTrackingRegionBuilder");
103 
104  theRegionBuilder = new MuonTrackingRegionBuilder(regionBuilderPSet,theService,iC);
105 
106  // TrackRefitter parameters
107  ParameterSet refitterParameters = par.getParameter<ParameterSet>("GlbRefitterParameters");
108  theGlbRefitter = new GlobalMuonRefitter(refitterParameters, theService, iC);
109 
110  theMuonHitsOption = refitterParameters.getParameter<int>("MuonHitsOption");
111 
112  theTrackerRecHitBuilderName = par.getParameter<string>("TrackerRecHitBuilder");
113  theMuonRecHitBuilderName = par.getParameter<string>("MuonRecHitBuilder");
114 
115  theRPCInTheFit = par.getParameter<bool>("RefitRPCHits");
116 
117  theTECxScale = par.getParameter<double>("ScaleTECxFactor");
118  theTECyScale = par.getParameter<double>("ScaleTECyFactor");
119  thePtCut = par.getParameter<double>("PtCut");
120  thePCut = par.getParameter<double>("PCut");
121 
122  theCacheId_TRH = 0;
123 
124 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
MuonDetLayerMeasurements * theLayerMeasurements
const MuonServiceProxy * service() const
GlobalMuonTrackMatcher * theTrackMatcher
MuonTrackingRegionBuilder * theRegionBuilder
GlobalTrajectoryBuilderBase::~GlobalTrajectoryBuilderBase ( )
virtual

destructor

Definition at line 130 of file GlobalTrajectoryBuilderBase.cc.

References theGlbRefitter, theRegionBuilder, theTrackMatcher, and theTrackTransformer.

130  {
131 
132  if (theTrackMatcher) delete theTrackMatcher;
135  if (theGlbRefitter) delete theGlbRefitter;
136 }
GlobalMuonTrackMatcher * theTrackMatcher
MuonTrackingRegionBuilder * theRegionBuilder

Member Function Documentation

void GlobalTrajectoryBuilderBase::addTraj ( TrackCand )
inlineprotected

if TrackCand has only a TrackRef, attempt to add Trajectory*

Definition at line 126 of file GlobalTrajectoryBuilderBase.h.

126 {}
MuonCandidate::CandidateContainer GlobalTrajectoryBuilderBase::build ( const TrackCand staCand,
MuonTrajectoryBuilder::CandidateContainer tkTrajs 
) const
protected

build combined trajectory from sta Track and tracker RecHits

Definition at line 172 of file GlobalTrajectoryBuilderBase.cc.

References checkRecHitsOrdering(), fixTEC(), getTransientRecHits(), getDQMSummary::iter, LogDebug, LogTrace, MuonCandidate::muonTrack(), outToIn, mix_2012_Summer_inTimeOnly_cff::prob, GlobalMuonRefitter::refit(), selectTrackerHits(), Trajectory::setSeedRef(), theCategory, theGlbRefitter, theMuonHitsOption, thePCut, thePtCut, theService, theTECxScale, theTECyScale, theTrackTransformer, MuonCandidate::trackerTrack(), MuonCandidate::trackerTrajectory(), trackProbability(), MuonCandidate::trajectory(), TrackTransformer::transform(), and tTopo_.

Referenced by GlobalMuonTrajectoryBuilder::trajectories(), and L3MuonTrajectoryBuilder::trajectories().

173  {
174 
175  LogTrace(theCategory) << " Begin Build" << endl;
176 
177  // tracker trajectory should be built and refit before this point
178  if ( tkTrajs.empty() ) return CandidateContainer();
179 
180  // add muon hits and refit/smooth trajectories
181  CandidateContainer refittedResult;
182  ConstRecHitContainer muonRecHits = getTransientRecHits(*(staCand.second));
183 
184  // check order of muon measurements
185  if ( (muonRecHits.size() > 1) &&
186  ( muonRecHits.front()->globalPosition().mag() >
187  muonRecHits.back()->globalPosition().mag() ) ) {
188  LogTrace(theCategory)<< " reverse order: ";
189  }
190 
191  for ( CandidateContainer::const_iterator it = tkTrajs.begin(); it != tkTrajs.end(); it++ ) {
192 
193  // cut on tracks with low momenta
194  LogTrace(theCategory)<< " Track p and pT " << (*it)->trackerTrack()->p() << " " << (*it)->trackerTrack()->pt();
195  if( (*it)->trackerTrack()->p() < thePCut || (*it)->trackerTrack()->pt() < thePtCut ) continue;
196 
197  ConstRecHitContainer trackerRecHits;
198  if ((*it)->trackerTrack().isNonnull()) {
199  trackerRecHits = getTransientRecHits(*(*it)->trackerTrack());
200  } else {
201  LogDebug(theCategory)<<" NEED HITS FROM TRAJ";
202  //trackerRecHits = (*it)->trackerTrajectory()->recHits();
203  }
204 
205  // check for single TEC RecHits in trajectories in the overalp region
206  if ( fabs((*it)->trackerTrack()->eta()) > 0.95 && fabs((*it)->trackerTrack()->eta()) < 1.15 && (*it)->trackerTrack()->pt() < 60 ) {
207  if ( theTECxScale < 0 || theTECyScale < 0 )
208  trackerRecHits = selectTrackerHits(trackerRecHits);
209  else
210  fixTEC(trackerRecHits,theTECxScale,theTECyScale);
211  }
212 
213  RefitDirection recHitDir = checkRecHitsOrdering(trackerRecHits);
214  if ( recHitDir == outToIn ) reverse(trackerRecHits.begin(),trackerRecHits.end());
215 
216  reco::TransientTrack tTT((*it)->trackerTrack(),&*theService->magneticField(),theService->trackingGeometry());
217  TrajectoryStateOnSurface innerTsos = tTT.innermostMeasurementState();
218 
220  if((*it)->trackerTrack()->seedRef().isAvailable()) tmpSeed = (*it)->trackerTrack()->seedRef();
221 
222  if ( !innerTsos.isValid() ) {
223  LogTrace(theCategory) << " inner Trajectory State is invalid. ";
224  continue;
225  }
226 
227  innerTsos.rescaleError(100.);
228 
229  TC refitted0,refitted1;
230  MuonCandidate* finalTrajectory = 0;
231  Trajectory *tkTrajectory = 0;
232 
233  // tracker only track
234  if ( ! ((*it)->trackerTrajectory() && (*it)->trackerTrajectory()->isValid()) ) {
235  refitted0 = theTrackTransformer->transform((*it)->trackerTrack()) ;
236  if (!refitted0.empty()) tkTrajectory = new Trajectory(*(refitted0.begin()));
237  else LogWarning(theCategory)<< " Failed to load tracker track trajectory";
238  } else tkTrajectory = (*it)->trackerTrajectory();
239  if (tkTrajectory) tkTrajectory->setSeedRef(tmpSeed);
240 
241  // full track with all muon hits using theGlbRefitter
242  ConstRecHitContainer allRecHits = trackerRecHits;
243  allRecHits.insert(allRecHits.end(), muonRecHits.begin(),muonRecHits.end());
244  refitted1 = theGlbRefitter->refit( *(*it)->trackerTrack(), tTT, allRecHits,theMuonHitsOption, tTopo_);
245  LogTrace(theCategory)<<" This track-sta refitted to " << refitted1.size() << " trajectories";
246 
247  Trajectory *glbTrajectory1 = 0;
248  if (!refitted1.empty()) glbTrajectory1 = new Trajectory(*(refitted1.begin()));
249  else LogDebug(theCategory)<< " Failed to load global track trajectory 1";
250  if (glbTrajectory1) glbTrajectory1->setSeedRef(tmpSeed);
251 
252  finalTrajectory = 0;
253  if(glbTrajectory1 && tkTrajectory) finalTrajectory = new MuonCandidate(glbTrajectory1, (*it)->muonTrack(), (*it)->trackerTrack(),
254  tkTrajectory? new Trajectory(*tkTrajectory) : 0);
255 
256  if ( finalTrajectory )
257  refittedResult.push_back(finalTrajectory);
258 
259  if(tkTrajectory) delete tkTrajectory;
260  }
261 
262  // choose the best global fit for this Standalone Muon based on the track probability
263  CandidateContainer selectedResult;
264  MuonCandidate* tmpCand = 0;
265  if ( refittedResult.size() > 0 ) tmpCand = *(refittedResult.begin());
266  double minProb = 9999;
267 
268  for (CandidateContainer::const_iterator iter=refittedResult.begin(); iter != refittedResult.end(); iter++) {
269  double prob = trackProbability(*(*iter)->trajectory());
270  LogTrace(theCategory)<<" refitted-track-sta with pT " << (*iter)->trackerTrack()->pt() << " has probability " << prob;
271 
272  if (prob < minProb) {
273  minProb = prob;
274  tmpCand = (*iter);
275  }
276  }
277 
278  if ( tmpCand ) selectedResult.push_back(new MuonCandidate(new Trajectory(*(tmpCand->trajectory())), tmpCand->muonTrack(), tmpCand->trackerTrack(),
279  (tmpCand->trackerTrajectory())? new Trajectory( *(tmpCand->trackerTrajectory()) ):0 ) );
280 
281  for (CandidateContainer::const_iterator it = refittedResult.begin(); it != refittedResult.end(); ++it) {
282  if ( (*it)->trajectory() ) delete (*it)->trajectory();
283  if ( (*it)->trackerTrajectory() ) delete (*it)->trackerTrajectory();
284  if ( *it ) delete (*it);
285  }
286  refittedResult.clear();
287 
288  return selectedResult;
289 
290 }
#define LogDebug(id)
const reco::TrackRef muonTrack() const
return muon track
Definition: MuonCandidate.h:43
MuonCandidate::CandidateContainer CandidateContainer
ConstRecHitContainer selectTrackerHits(const ConstRecHitContainer &) const
select tracker hits; exclude some tracker hits in the global trajectory
const reco::TrackRef trackerTrack() const
return tracker track
Definition: MuonCandidate.h:46
void fixTEC(ConstRecHitContainer &all, double scl_x, double scl_y) const
rescale errors of outermost TEC RecHit
double trackProbability(const Trajectory &) const
calculate chi2 probability (-ln(P))
Trajectory * trackerTrajectory() const
return tracker trajectory
Definition: MuonCandidate.h:49
TransientTrackingRecHit::ConstRecHitContainer getTransientRecHits(const reco::Track &) const
get transient RecHits of a Track
#define LogTrace(id)
void setSeedRef(const edm::RefToBase< TrajectorySeed > &seedRef)
Definition: Trajectory.h:308
virtual std::vector< Trajectory > transform(const reco::Track &) const
Convert a reco::Track into Trajectory.
RefitDirection checkRecHitsOrdering(const ConstRecHitContainer &) const
This does nothing now.
TransientTrackingRecHit::ConstRecHitContainer ConstRecHitContainer
Trajectory * trajectory() const
return trajectory
Definition: MuonCandidate.h:40
std::vector< Trajectory > refit(const reco::Track &globalTrack, const int theMuonHitsOption, const TrackerTopology *tTopo) const
build combined trajectory from sta Track and tracker RecHits
void GlobalTrajectoryBuilderBase::checkMuonHits ( const reco::Track ,
ConstRecHitContainer ,
ConstRecHitContainer ,
std::vector< int > &   
) const
protected

check muon RecHits, calculate chamber occupancy and select hits to be used in the final fit

GlobalTrajectoryBuilderBase::RefitDirection GlobalTrajectoryBuilderBase::checkRecHitsOrdering ( const ConstRecHitContainer recHits) const
protected

This does nothing now.

check order of RechIts on a trajectory

Definition at line 402 of file GlobalTrajectoryBuilderBase.cc.

References inToOut, outToIn, theCategory, and undetermined.

Referenced by build().

402  {
403 
404  if ( !recHits.empty() ) {
405  ConstRecHitContainer::const_iterator frontHit = recHits.begin();
406  ConstRecHitContainer::const_iterator backHit = recHits.end() - 1;
407  while ( !(*frontHit)->isValid() && frontHit != backHit ) {frontHit++;}
408  while ( !(*backHit)->isValid() && backHit != frontHit ) {backHit--;}
409 
410  double rFirst = (*frontHit)->globalPosition().mag();
411  double rLast = (*backHit) ->globalPosition().mag();
412 
413  if ( rFirst < rLast ) return inToOut;
414  else if (rFirst > rLast) return outToIn;
415  else {
416  LogError(theCategory) << "Impossible to determine the rechits order" << endl;
417  return undetermined;
418  }
419  }
420  else {
421  LogError(theCategory) << "Impossible to determine the rechits order" << endl;
422  return undetermined;
423  }
424 }
vector< GlobalTrajectoryBuilderBase::TrackCand > GlobalTrajectoryBuilderBase::chooseRegionalTrackerTracks ( const TrackCand staCand,
const std::vector< TrackCand > &  tkTs 
)
protected

choose tracker tracks within region of interest

Definition at line 297 of file GlobalTrajectoryBuilderBase.cc.

References defineRegionOfInterest(), deltaR(), TrackingRegion::direction(), PV3DBase< T, PVType, FrameType >::eta(), PV3DBase< T, PVType, FrameType >::phi(), and query::result.

Referenced by GlobalMuonTrajectoryBuilder::makeTkCandCollection().

298  {
299 
300  // define eta-phi region
301  RectangularEtaPhiTrackingRegion regionOfInterest = defineRegionOfInterest(staCand.second);
302 
303  // get region's etaRange and phiMargin
304  //UNUSED: PixelRecoRange<float> etaRange = regionOfInterest.etaRange();
305  //UNUSED: TkTrackingRegionsMargin<float> phiMargin = regionOfInterest.phiMargin();
306 
307  vector<TrackCand> result;
308 
309  double deltaR_max = 1.0;
310 
311  for ( vector<TrackCand>::const_iterator is = tkTs.begin(); is != tkTs.end(); ++is ) {
312  // check if each trackCand is in region of interest
313 // bool inEtaRange = etaRange.inside(is->second->eta());
314 // bool inPhiRange = (fabs(Geom::Phi<float>(is->second->phi()) - Geom::Phi<float>(regionOfInterest.direction().phi())) < phiMargin.right() ) ? true : false ;
315 
316  double deltaR_tmp = deltaR(static_cast<double>(regionOfInterest.direction().eta()),
317  static_cast<double>(regionOfInterest.direction().phi()),
318  is->second->eta(), is->second->phi());
319 
320  // for each trackCand in region, add trajectory and add to result
321  //if ( inEtaRange && inPhiRange ) {
322  if (deltaR_tmp < deltaR_max) {
323  TrackCand tmpCand = TrackCand(*is);
324  result.push_back(tmpCand);
325  }
326  }
327 
328  return result;
329 
330 }
std::pair< const Trajectory *, reco::TrackRef > TrackCand
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
GlobalVector const & direction() const
the direction around which region is constructed
tuple result
Definition: query.py:137
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
RectangularEtaPhiTrackingRegion defineRegionOfInterest(const reco::TrackRef &) const
define region of interest with tracker
T eta() const
Definition: PV3DBase.h:76
const Trajectory* GlobalTrajectoryBuilderBase::chooseTrajectory ( const std::vector< Trajectory * > &  ,
int   
) const
protected

choose final trajectory

RectangularEtaPhiTrackingRegion GlobalTrajectoryBuilderBase::defineRegionOfInterest ( const reco::TrackRef staTrack) const
protected

define region of interest with tracker

Definition at line 337 of file GlobalTrajectoryBuilderBase.cc.

References TrackingRegion::direction(), RectangularEtaPhiTrackingRegion::etaRange(), PixelRecoRange< T >::max(), PixelRecoRange< T >::mean(), PixelRecoRange< T >::min(), TrackingRegion::origin(), TrackingRegion::originRBound(), TrackingRegion::originZBound(), RectangularEtaPhiTrackingRegion::phiMargin(), TrackingRegion::ptMin(), MuonTrackingRegionBuilder::region(), and theRegionBuilder.

Referenced by chooseRegionalTrackerTracks().

337  {
338 
340 
341  TkTrackingRegionsMargin<float> etaMargin(fabs(region1->etaRange().min() - region1->etaRange().mean()),
342  fabs(region1->etaRange().max() - region1->etaRange().mean()));
343 
344  RectangularEtaPhiTrackingRegion region2(region1->direction(),
345  region1->origin(),
346  region1->ptMin(),
347  region1->originRBound(),
348  region1->originZBound(),
349  etaMargin,
350  region1->phiMargin());
351 
352  delete region1;
353  return region2;
354 
355 }
RectangularEtaPhiTrackingRegion * region(const reco::TrackRef &) const
define tracking region
float originRBound() const
bounds the particle vertex in the transverse plane
GlobalPoint const & origin() const
T max() const
T min() const
GlobalVector const & direction() const
the direction around which region is constructed
float originZBound() const
bounds the particle vertex in the longitudinal plane
T mean() const
float ptMin() const
minimal pt of interest
MuonTrackingRegionBuilder * theRegionBuilder
const Range & etaRange() const
allowed eta range [eta_min, eta_max] interval
void GlobalTrajectoryBuilderBase::fixTEC ( ConstRecHitContainer all,
double  scl_x,
double  scl_y 
) const
protected

rescale errors of outermost TEC RecHit

 freeze this hit, make sure it will not be recomputed during fitting

/ the implemetantion below works with cloning / to get a RecHitPointer to SiStripRecHit2D, the only method that works is / RecHitPointer MuonTransientTrackingRecHit::build(const GeomDet*,const TrackingRecHit*)

Definition at line 455 of file GlobalTrajectoryBuilderBase.cc.

References angle(), MuonTransientTrackingRecHit::build(), SiStripRecHit2D::cluster(), TrackingRecHit::det(), relativeConstraints::error, i, BaseTrackerRecHit::localPosition(), BaseTrackerRecHit::localPositionError(), LocalError::rotate(), StripTopology::strip(), StripTopology::stripAngle(), StripSubdetector::TEC, ecaldqm::topology(), DetId::Tracker, LocalError::xx(), and LocalError::yy().

Referenced by build().

457  {
458 
459  int nTEC(0);
460  ConstRecHitContainer::iterator lone_tec;
461 
462  for ( ConstRecHitContainer::iterator i = all.begin(); i != all.end(); i++) {
463  if ( !(*i)->isValid() ) continue;
464 
465  if ( (*i)->det()->geographicalId().det() == DetId::Tracker &&
466  (*i)->det()->geographicalId().subdetId() == StripSubdetector::TEC) {
467  lone_tec = i;
468  nTEC++;
469 
470  if ( (i+1) != all.end() && (*(i+1))->isValid() &&
471  (*(i+1))->det()->geographicalId().det() == DetId::Tracker &&
472  (*(i+1))->det()->geographicalId().subdetId() == StripSubdetector::TEC) {
473  nTEC++;
474  break;
475  }
476  }
477 
478  if (nTEC > 1) break;
479  }
480 
481  int hitDet = (*lone_tec)->hit()->geographicalId().det();
482  int hitSubDet = (*lone_tec)->hit()->geographicalId().subdetId();
483  if ( nTEC == 1 && (*lone_tec)->hit()->isValid() &&
484  hitDet == DetId::Tracker && hitSubDet == StripSubdetector::TEC) {
485 
486  // rescale the TEC rechit error matrix in its rotated frame
487  const SiStripRecHit2D* strip = dynamic_cast<const SiStripRecHit2D*>((*lone_tec)->hit());
488  if (strip && strip->det() ) {
489  LocalPoint pos = strip->localPosition();
490  if ((*lone_tec)->detUnit()) {
491  const StripTopology* topology = dynamic_cast<const StripTopology*>(&(*lone_tec)->detUnit()->topology());
492  if (topology) {
493  // rescale the local error along/perp the strip by a factor
494  float angle = topology->stripAngle(topology->strip((*lone_tec)->hit()->localPosition()));
496  LocalError rotError = error.rotate(angle);
497  LocalError scaledError(rotError.xx() * scl_x * scl_x, 0, rotError.yy() * scl_y * scl_y);
498  error = scaledError.rotate(-angle);
503  SiStripRecHit2D* st = new SiStripRecHit2D(pos,error,
504  *strip->det(),
505  strip->cluster());
506  *lone_tec = MuonTransientTrackingRecHit::build((*lone_tec)->det(),st);
507 
508  }
509  }
510  }
511  }
512 
513 }
int i
Definition: DBlmapReader.cc:9
float xx() const
Definition: LocalError.h:24
virtual float stripAngle(float strip) const =0
CaloTopology const * topology(0)
virtual float strip(const LocalPoint &) const =0
virtual LocalError localPositionError() const
float yy() const
Definition: LocalError.h:26
const GeomDet * det() const
ClusterRef cluster() const
LocalError rotate(float x, float y) const
Return a new LocalError, rotated by an angle defined by the direction (x,y)
Definition: LocalError.h:39
virtual LocalPoint localPosition() const
static RecHitPointer build(const GeomDet *geom, const TrackingRecHit *rh)
FIXME virtual ConstMuonRecHitContainer specificTransientHits() const;.
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11
TransientTrackingRecHit::ConstRecHitContainer GlobalTrajectoryBuilderBase::getTransientRecHits ( const reco::Track track) const
protected

get transient RecHits of a Track

Definition at line 521 of file GlobalTrajectoryBuilderBase.cc.

References DetId::det(), trajectoryStateTransform::innerStateOnSurface(), TrajectoryStateOnSurface::isValid(), LogDebug, DetId::Muon, edm::ESHandle< class >::product(), reco::Track::recHitsBegin(), reco::Track::recHitsEnd(), query::result, theCategory, theMuonRecHitBuilder, theRPCInTheFit, theService, theTrackerPropagatorName, theTrackerRecHitBuilder, and DetId::Tracker.

Referenced by build().

521  {
522 
524 
525 
526 
527  TrajectoryStateOnSurface currTsos = trajectoryStateTransform::innerStateOnSurface(track, *theService->trackingGeometry(), &*theService->magneticField());
528 
529  auto hitCloner = static_cast<TkTransientTrackingRecHitBuilder const *>(theTrackerRecHitBuilder.product())->cloner();
530  for (trackingRecHit_iterator hit = track.recHitsBegin(); hit != track.recHitsEnd(); ++hit) {
531  if((*hit)->isValid()) {
532  DetId recoid = (*hit)->geographicalId();
533  if ( recoid.det() == DetId::Tracker ) {
534  if (!(*hit)->hasPositionAndError()){
535  TrajectoryStateOnSurface predTsos = theService->propagator(theTrackerPropagatorName)->propagate(currTsos, theService->trackingGeometry()->idToDet(recoid)->surface());
536 
537  if ( !predTsos.isValid() ) {
538  edm::LogError("MissingTransientHit")
539  <<"Could not get a tsos on the hit surface. We will miss a tracking hit.";
540  continue;
541  }
542  currTsos = predTsos;
543  result.emplace_back(hitCloner(**hit,predTsos));
544  }else{
545  result.push_back((*hit)->cloneSH());
546  }
547  } else if ( recoid.det() == DetId::Muon ) {
548  if ( (*hit)->geographicalId().subdetId() == 3 && !theRPCInTheFit) {
549  LogDebug(theCategory) << "RPC Rec Hit discarded";
550  continue;
551  }
552  result.push_back(theMuonRecHitBuilder->build(&**hit));
553  }
554  }
555  }
556 
557  return result;
558 }
#define LogDebug(id)
TrajectoryStateOnSurface innerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field)
edm::ESHandle< TransientTrackingRecHitBuilder > theMuonRecHitBuilder
tuple result
Definition: query.py:137
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:63
std::vector< ConstRecHitPointer > ConstRecHitContainer
Definition: DetId.h:18
T const * product() const
Definition: ESHandle.h:62
Detector det() const
get the detector field from this detid
Definition: DetId.h:35
edm::ESHandle< TransientTrackingRecHitBuilder > theTrackerRecHitBuilder
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:65
virtual std::vector<TrackCand> GlobalTrajectoryBuilderBase::makeTkCandCollection ( const TrackCand )
protectedpure virtual

make a TrackCand collection using tracker Track, Trajectory information

Implemented in L3MuonTrajectoryBuilder, and GlobalMuonTrajectoryBuilder.

void GlobalTrajectoryBuilderBase::printHits ( const ConstRecHitContainer hits) const
protected

print all RecHits of a trajectory

Definition at line 376 of file GlobalTrajectoryBuilderBase.cc.

References LogTrace, mathSSE::sqrt(), theCategory, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

376  {
377 
378  LogTrace(theCategory) << "Used RecHits: " << hits.size();
379  for (ConstRecHitContainer::const_iterator ir = hits.begin(); ir != hits.end(); ir++ ) {
380  if ( !(*ir)->isValid() ) {
381  LogTrace(theCategory) << "invalid RecHit";
382  continue;
383  }
384 
385  const GlobalPoint& pos = (*ir)->globalPosition();
386 
388  << "r = " << sqrt(pos.x() * pos.x() + pos.y() * pos.y())
389  << " z = " << pos.z()
390  << " dimension = " << (*ir)->dimension()
391  << " " << (*ir)->det()->geographicalId().det()
392  << " " << (*ir)->det()->subDetector();
393 
394  }
395 
396 }
T y() const
Definition: PV3DBase.h:63
T sqrt(T t)
Definition: SSEVec.h:48
T z() const
Definition: PV3DBase.h:64
#define LogTrace(id)
T x() const
Definition: PV3DBase.h:62
ConstRecHitContainer GlobalTrajectoryBuilderBase::selectMuonHits ( const Trajectory ,
const std::vector< int > &   
) const
protected

select muon hits compatible with trajectory; check hits in chambers with showers

GlobalTrajectoryBuilderBase::ConstRecHitContainer GlobalTrajectoryBuilderBase::selectTrackerHits ( const ConstRecHitContainer all) const
protected

select tracker hits; exclude some tracker hits in the global trajectory

Definition at line 431 of file GlobalTrajectoryBuilderBase.cc.

References cond::ecalcond::all, i, StripSubdetector::TEC, and DetId::Tracker.

Referenced by build().

431  {
432 
433  int nTEC(0);
434 
436  for (ConstRecHitContainer::const_iterator i = all.begin(); i != all.end(); i++) {
437  if ( !(*i)->isValid() ) continue;
438  if ( (*i)->det()->geographicalId().det() == DetId::Tracker &&
439  (*i)->det()->geographicalId().subdetId() == StripSubdetector::TEC) {
440  nTEC++;
441  } else {
442  hits.push_back(*i);
443  }
444  if ( nTEC > 1 ) return all;
445  }
446 
447  return hits;
448 
449 }
int i
Definition: DBlmapReader.cc:9
TransientTrackingRecHit::ConstRecHitContainer ConstRecHitContainer
const MuonServiceProxy* GlobalTrajectoryBuilderBase::service ( ) const
inlineprotected

Definition at line 139 of file GlobalTrajectoryBuilderBase.h.

References theService.

139 { return theService; }
void GlobalTrajectoryBuilderBase::setEvent ( const edm::Event event)
virtual

pass the Event to the algo at each event

Implements MuonTrajectoryBuilder.

Reimplemented in L3MuonTrajectoryBuilder, and GlobalMuonTrajectoryBuilder.

Definition at line 142 of file GlobalTrajectoryBuilderBase.cc.

References event(), LogDebug, edm::ESHandle< class >::product(), MuonTrackingRegionBuilder::setEvent(), GlobalMuonRefitter::setEvent(), TrackTransformer::setServices(), GlobalMuonRefitter::setServices(), theCacheId_TRH, theCategory, theEvent, theGlbRefitter, theMuonRecHitBuilder, theMuonRecHitBuilderName, theRegionBuilder, theService, theTrackerRecHitBuilder, theTrackerRecHitBuilderName, theTrackTransformer, and tTopo_.

Referenced by GlobalMuonTrajectoryBuilder::setEvent(), and L3MuonTrajectoryBuilder::setEvent().

142  {
143 
144  theEvent = &event;
145 
147  theRegionBuilder->setEvent(event);
148 
149  theGlbRefitter->setEvent(event);
150  theGlbRefitter->setServices(theService->eventSetup());
151 
152  unsigned long long newCacheId_TRH = theService->eventSetup().get<TransientRecHitRecord>().cacheIdentifier();
153  if ( newCacheId_TRH != theCacheId_TRH ) {
154  LogDebug(theCategory) << "TransientRecHitRecord changed!";
155  theCacheId_TRH = newCacheId_TRH;
158  }
159 
160  //Retrieve tracker topology from geometry
162  theService->eventSetup().get<IdealGeometryRecord>().get(tTopoHand);
163  tTopo_=tTopoHand.product();
164 
165 }
#define LogDebug(id)
edm::ESHandle< TransientTrackingRecHitBuilder > theMuonRecHitBuilder
virtual void setEvent(const edm::Event &)
pass the Event to the algo at each event
void setServices(const edm::EventSetup &)
set the services needed by the TrackTransformer
virtual void setServices(const edm::EventSetup &)
set the services needed by the TrackTransformer
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
virtual void setEvent(const edm::Event &)
pass the Event to the algo at each event
T const * product() const
Definition: ESHandle.h:62
MuonTrackingRegionBuilder * theRegionBuilder
edm::ESHandle< TransientTrackingRecHitBuilder > theTrackerRecHitBuilder
GlobalMuonTrackMatcher* GlobalTrajectoryBuilderBase::trackMatcher ( ) const
inlineprotected
double GlobalTrajectoryBuilderBase::trackProbability ( const Trajectory track) const
protected

calculate chi2 probability (-ln(P))

Definition at line 362 of file GlobalTrajectoryBuilderBase.cc.

References Trajectory::chiSquared(), LnChiSquaredProbability(), and Trajectory::ndof().

Referenced by build().

362  {
363 
364  if ( track.ndof() > 0 && track.chiSquared() > 0 ) {
365  return -LnChiSquaredProbability(track.chiSquared(), track.ndof());
366  } else {
367  return 0.0;
368  }
369 
370 }
float LnChiSquaredProbability(double chiSquared, double nrDOF)
int ndof(bool bon=true) const
Definition: Trajectory.cc:78
float chiSquared() const
Definition: Trajectory.h:252
MuonTrajectoryBuilder::TrajectoryContainer GlobalTrajectoryBuilderBase::trajectories ( const TrajectorySeed )
inlinevirtual

dummy implementation, unused in this class

Implements MuonTrajectoryBuilder.

Definition at line 75 of file GlobalTrajectoryBuilderBase.h.

MuonCandidate::TrajectoryContainer TrajectoryContainer

Member Data Documentation

unsigned long long GlobalTrajectoryBuilderBase::theCacheId_TRH
private

Definition at line 175 of file GlobalTrajectoryBuilderBase.h.

Referenced by GlobalTrajectoryBuilderBase(), and setEvent().

std::string GlobalTrajectoryBuilderBase::theCategory
protected
const edm::Event* GlobalTrajectoryBuilderBase::theEvent
private

Definition at line 183 of file GlobalTrajectoryBuilderBase.h.

Referenced by setEvent().

GlobalMuonRefitter* GlobalTrajectoryBuilderBase::theGlbRefitter
private
MuonDetLayerMeasurements* GlobalTrajectoryBuilderBase::theLayerMeasurements
private

Definition at line 170 of file GlobalTrajectoryBuilderBase.h.

int GlobalTrajectoryBuilderBase::theMuonHitsOption
private

Definition at line 178 of file GlobalTrajectoryBuilderBase.h.

Referenced by build(), and GlobalTrajectoryBuilderBase().

edm::ESHandle<TransientTrackingRecHitBuilder> GlobalTrajectoryBuilderBase::theMuonRecHitBuilder
private

Definition at line 189 of file GlobalTrajectoryBuilderBase.h.

Referenced by getTransientRecHits(), and setEvent().

std::string GlobalTrajectoryBuilderBase::theMuonRecHitBuilderName
private

Definition at line 188 of file GlobalTrajectoryBuilderBase.h.

Referenced by GlobalTrajectoryBuilderBase(), and setEvent().

float GlobalTrajectoryBuilderBase::thePCut
protected

Definition at line 165 of file GlobalTrajectoryBuilderBase.h.

Referenced by build(), and GlobalTrajectoryBuilderBase().

float GlobalTrajectoryBuilderBase::thePtCut
protected
MuonTrackingRegionBuilder* GlobalTrajectoryBuilderBase::theRegionBuilder
private
bool GlobalTrajectoryBuilderBase::theRPCInTheFit
private
const MuonServiceProxy* GlobalTrajectoryBuilderBase::theService
private
float GlobalTrajectoryBuilderBase::theTECxScale
private

Definition at line 179 of file GlobalTrajectoryBuilderBase.h.

Referenced by build(), and GlobalTrajectoryBuilderBase().

float GlobalTrajectoryBuilderBase::theTECyScale
private

Definition at line 180 of file GlobalTrajectoryBuilderBase.h.

Referenced by build(), and GlobalTrajectoryBuilderBase().

std::string GlobalTrajectoryBuilderBase::theTrackerPropagatorName
private
edm::ESHandle<TransientTrackingRecHitBuilder> GlobalTrajectoryBuilderBase::theTrackerRecHitBuilder
private

Definition at line 186 of file GlobalTrajectoryBuilderBase.h.

Referenced by getTransientRecHits(), and setEvent().

std::string GlobalTrajectoryBuilderBase::theTrackerRecHitBuilderName
private

Definition at line 185 of file GlobalTrajectoryBuilderBase.h.

Referenced by GlobalTrajectoryBuilderBase(), and setEvent().

GlobalMuonTrackMatcher* GlobalTrajectoryBuilderBase::theTrackMatcher
private
TrackTransformer* GlobalTrajectoryBuilderBase::theTrackTransformer
private
const TrackerTopology* GlobalTrajectoryBuilderBase::tTopo_
private

Definition at line 191 of file GlobalTrajectoryBuilderBase.h.

Referenced by build(), and setEvent().