CMS 3D CMS Logo

Public Types | Public Member Functions | Private Member Functions | Private Attributes

RectangularEtaPhiTrackingRegion Class Reference

#include <RectangularEtaPhiTrackingRegion.h>

Inheritance diagram for RectangularEtaPhiTrackingRegion:
TrackingRegionBase TrackingRegion

List of all members.

Public Types

typedef
TkTrackingRegionsMargin< float > 
Margin

Public Member Functions

virtual HitRZCompatibilitycheckRZ (const DetLayer *layer, const Hit &outerHit, const edm::EventSetup &iSetup) const
virtual
RectangularEtaPhiTrackingRegion
clone () const
 clone region
const RangeetaRange () const
 allowed eta range [eta_min, eta_max] interval
virtual TrackingRegion::Hits hits (const edm::Event &ev, const edm::EventSetup &es, const ctfseeding::SeedingLayer *layer) const
 get hits from layer compatible with region constraints
bool isPrecise () const
 is precise error calculation switched on
virtual std::string name () const
const MarginphiMargin () const
virtual std::string print () const
 RectangularEtaPhiTrackingRegion (const GlobalVector &dir, const GlobalPoint &vertexPos, float ptMin, float rVertex, float zVertex, float deltaEta, float deltaPhi, float whereToUseMeasurementTracker=0., bool precise=true, const std::string &measurementTrackerName="", bool etaPhiRegion=false)
 RectangularEtaPhiTrackingRegion (const GlobalVector &dir, const GlobalPoint &vertexPos, float ptMin, float rVertex, float zVertex, Margin etaMargin, Margin phiMargin, float whereToUseMeasurementTracker=0., bool precise=true, const std::string &measurementTrackerName="", bool etaPhiRegion=false)
 RectangularEtaPhiTrackingRegion ()
 dummy constructor
 RectangularEtaPhiTrackingRegion (const GlobalVector &dir, const GlobalPoint &vertexPos, Range invPtRange, float rVertex, float zVertex, Margin etaMargin, Margin phiMargin, float whereToUseMeasurementTracker=0., bool precise=true, const std::string &measurementTrackerName="", bool etaPhiRegion=false)

Private Member Functions

HitRZCompatibilitycheckRZOld (const DetLayer *layer, const TrackingRecHit *outerHit, const edm::EventSetup &iSetup) const
OuterEstimatorestimator (const BarrelDetLayer *layer, const edm::EventSetup &iSetup) const
OuterEstimatorestimator (const ForwardDetLayer *layer, const edm::EventSetup &iSetup) const
void initEtaRange (const GlobalVector &dir, const Margin &margin)
OuterHitPhiPrediction phiWindow (const edm::EventSetup &iSetup) const
HitRZConstraint rzConstraint () const

Private Attributes

Range theEtaRange
std::string theMeasurementTrackerName
double theMeasurementTrackerUsage
Margin thePhiMargin
bool thePrecise
bool theUseEtaPhi

Detailed Description

A concrete implementation of TrackingRegion. Apart of vertex constraint from TrackingRegion in this implementation the region of interest is further constrainted in phi and eta around the direction of the region

Definition at line 21 of file RectangularEtaPhiTrackingRegion.h.


Member Typedef Documentation

Definition at line 24 of file RectangularEtaPhiTrackingRegion.h.


Constructor & Destructor Documentation

RectangularEtaPhiTrackingRegion::RectangularEtaPhiTrackingRegion ( ) [inline]

dummy constructor

Definition at line 27 of file RectangularEtaPhiTrackingRegion.h.

Referenced by clone().

{ }
RectangularEtaPhiTrackingRegion::RectangularEtaPhiTrackingRegion ( const GlobalVector dir,
const GlobalPoint vertexPos,
float  ptMin,
float  rVertex,
float  zVertex,
float  deltaEta,
float  deltaPhi,
float  whereToUseMeasurementTracker = 0.,
bool  precise = true,
const std::string &  measurementTrackerName = "",
bool  etaPhiRegion = false 
) [inline]

constructor (symmetric eta and phi margins).
dir - the direction around which region is constructed
the initial direction of the momentum of the particle should be in the range
phi of dir +- deltaPhi
eta of dir +- deltaEta

vertexPos - the position of the vertex (origin) of the of the region.
It is a centre of cylinder constraind with rVertex, zVertex. The track of the particle should cross the cylinder
WARNING: in the current implementaion the vertexPos is supposed to be placed on the beam line, i.e. to be of the form (0,0,float)

ptMin - minimal pt of interest
rVertex - radius of the cylinder around beam line where the tracks of interest should point to.
zVertex - half height of the cylinder around the beam line where the tracks of interest should point to.
deltaEta - allowed deviation of the initial direction of particle in eta in respect to direction of the region
deltaPhi - allowed deviation of the initial direction of particle in phi in respect to direction of the region whereToUseMeasurementTracker: 1=everywhere, 0=outside pixles, -1=nowhere

Definition at line 54 of file RectangularEtaPhiTrackingRegion.h.

References initEtaRange().

    : TrackingRegionBase( dir, vertexPos, Range( -1/ptMin, 1/ptMin), 
                          rVertex, zVertex),
    thePhiMargin( Margin( fabs(deltaPhi),fabs(deltaPhi))),
    theMeasurementTrackerUsage(whereToUseMeasurementTracker), thePrecise(precise), theMeasurementTrackerName(measurementTrackerName),
    theUseEtaPhi(etaPhiRegion)
   { initEtaRange(dir, Margin( fabs(deltaEta),fabs(deltaEta))); }
RectangularEtaPhiTrackingRegion::RectangularEtaPhiTrackingRegion ( const GlobalVector dir,
const GlobalPoint vertexPos,
float  ptMin,
float  rVertex,
float  zVertex,
Margin  etaMargin,
Margin  phiMargin,
float  whereToUseMeasurementTracker = 0.,
bool  precise = true,
const std::string &  measurementTrackerName = "",
bool  etaPhiRegion = false 
) [inline]

constructor (asymmetrinc eta and phi margins).
non equal left-right eta and phi bounds around direction are possible. The ranges are defined using Margin. the meaning of other arguments is the same as in the case of symmetring bounds to direction of the region.

Definition at line 75 of file RectangularEtaPhiTrackingRegion.h.

References initEtaRange().

    : TrackingRegionBase( dir, vertexPos, Range( -1/ptMin, 1/ptMin), 
      rVertex, zVertex), thePhiMargin( phiMargin), theMeasurementTrackerUsage(whereToUseMeasurementTracker), thePrecise(precise),
      theMeasurementTrackerName(measurementTrackerName),
      theUseEtaPhi(etaPhiRegion)
    { initEtaRange(dir, etaMargin); }
RectangularEtaPhiTrackingRegion::RectangularEtaPhiTrackingRegion ( const GlobalVector dir,
const GlobalPoint vertexPos,
Range  invPtRange,
float  rVertex,
float  zVertex,
Margin  etaMargin,
Margin  phiMargin,
float  whereToUseMeasurementTracker = 0.,
bool  precise = true,
const std::string &  measurementTrackerName = "",
bool  etaPhiRegion = false 
) [inline]

constructor (explicit pt range, asymmetrinc eta and phi margins).
the meaning of other arguments is the same as in the case of symmetring bounds to direction of the region.

Definition at line 94 of file RectangularEtaPhiTrackingRegion.h.

References initEtaRange().

    : TrackingRegionBase( dir, vertexPos, invPtRange, rVertex, zVertex),
      thePhiMargin( phiMargin), theMeasurementTrackerUsage(whereToUseMeasurementTracker), thePrecise(precise),
      theMeasurementTrackerName(measurementTrackerName),
      theUseEtaPhi(etaPhiRegion)
    { initEtaRange(dir, etaMargin); }

Member Function Documentation

virtual HitRZCompatibility* RectangularEtaPhiTrackingRegion::checkRZ ( const DetLayer layer,
const Hit outerHit,
const edm::EventSetup iSetup 
) const [inline, virtual]

utility to check eta/theta hit compatibility with region constraints and outer hit constraint

Implements TrackingRegionBase.

Definition at line 126 of file RectangularEtaPhiTrackingRegion.h.

References checkRZOld().

                                         { return checkRZOld(layer,outerHit->hit(),iSetup); }
HitRZCompatibility * RectangularEtaPhiTrackingRegion::checkRZOld ( const DetLayer layer,
const TrackingRecHit outerHit,
const edm::EventSetup iSetup 
) const [private]

Definition at line 49 of file RectangularEtaPhiTrackingRegion.cc.

References Reference_intrackfit_cff::barrel, PixelRecoLineRZ::cotLine(), TrackingRecHit::geographicalId(), edm::EventSetup::get(), TrackingRecHit::localPosition(), DetLayer::location(), max(), min, SurfaceOrientation::outer, PtMinSelector_cfg::ptMin, funct::sqr(), mathSSE::sqrt(), patCandidatesForDimuonsSequences_cff::tracker, vtxMean(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), PixelRecoPointRZ::z(), PV3DBase< T, PVType, FrameType >::z(), and PixelRecoLineRZ::zAtR().

Referenced by checkRZ().

{
  edm::ESHandle<TrackerGeometry> tracker;
  iSetup.get<TrackerDigiGeometryRecord>().get(tracker);

  bool isBarrel = (layer->location() == GeomDetEnumerators::barrel);
  GlobalPoint ohit =  tracker->idToDet(outerHit->geographicalId())->surface().toGlobal(outerHit->localPosition());
  float outerred_r = sqrt( sqr(ohit.x()-origin().x())+sqr(ohit.y()-origin().y()) );
  //PixelRecoPointRZ outer(ohit.perp(), ohit.z());
  PixelRecoPointRZ outer(outerred_r, ohit.z());
  
  float zMinOrigin = origin().z() - originZBound();
  float zMaxOrigin = origin().z() + originZBound();

  if (!thePrecise) {
    double vcotMin = (outer.z() > zMaxOrigin) ?
        (outer.z()-zMaxOrigin)/(outer.r()+originRBound())
      : (outer.z()-zMaxOrigin)/(outer.r()-originRBound());
    double vcotMax = (outer.z() > zMinOrigin) ?
        (outer.z()-zMinOrigin)/(outer.r()-originRBound())
      : (outer.z()-zMinOrigin)/(outer.r()+originRBound());
    float cotRight  = max(vcotMin,(double) sinh(theEtaRange.min()));
    float cotLeft = min(vcotMax, (double) sinh(theEtaRange.max()));
    return new HitEtaCheck( isBarrel, outer, cotLeft, cotRight);
  }
  float hitZErr = 0.;
  float hitRErr = 0.;

  PixelRecoPointRZ  outerL, outerR;
  if (layer->location() == GeomDetEnumerators::barrel) {
    outerL = PixelRecoPointRZ(outer.r(), outer.z()-hitZErr);
    outerR = PixelRecoPointRZ(outer.r(), outer.z()+hitZErr);
  } else if (outer.z() > 0) {
    outerL = PixelRecoPointRZ(outer.r()+hitRErr, outer.z());
    outerR = PixelRecoPointRZ(outer.r()-hitRErr, outer.z());
  } else {
    outerL = PixelRecoPointRZ(outer.r()-hitRErr, outer.z());
    outerR = PixelRecoPointRZ(outer.r()+hitRErr, outer.z());
  }
  //CHECK
  MultipleScatteringParametrisation oSigma(layer,iSetup);
  float cotThetaOuter = sinh(theEtaRange.mean());
  float sinThetaOuter = 1/sqrt(1+sqr(cotThetaOuter)); 
  float outerZscatt = 3*oSigma(ptMin(),cotThetaOuter) / sinThetaOuter;

  PixelRecoLineRZ boundL(outerL, sinh(theEtaRange.max()));
  PixelRecoLineRZ boundR(outerR, sinh(theEtaRange.min()));
  float zMinLine = boundL.zAtR(0.)-outerZscatt;
  float zMaxLine = boundR.zAtR(0.)+outerZscatt;
  PixelRecoPointRZ vtxL(0.,max(zMinLine, zMinOrigin));
  PixelRecoPointRZ vtxR(0.,min(zMaxLine, zMaxOrigin)); 
  PixelRecoPointRZ vtxMean(0.,(vtxL.z()+vtxR.z())/2.);
  //CHECK
  MultipleScatteringParametrisation iSigma(layer,iSetup);
  float innerScatt = 3 * iSigma(ptMin(),vtxMean, outer);
  
  PixelRecoLineRZ leftLine( vtxL, outerL);
  PixelRecoLineRZ rightLine( vtxR, outerR);

  HitRZConstraint rzConstraint(leftLine, rightLine);
  float cotTheta = fabs(leftLine.cotLine()+rightLine.cotLine())/2;

//  float bendR = longitudinalBendingCorrection(outer.r(),ptMin());

  if (isBarrel) {
    float sinTheta = 1/sqrt(1+sqr(cotTheta));
    float corrZ = innerScatt/sinTheta + hitZErr;
    return new HitZCheck(rzConstraint, HitZCheck::Margin(corrZ,corrZ));
  } else {
    float cosTheta = 1/sqrt(1+sqr(1/cotTheta));
    float corrR = innerScatt/cosTheta + hitRErr;
    return new HitRCheck( rzConstraint, HitRCheck::Margin(corrR,corrR));
  }
}
virtual RectangularEtaPhiTrackingRegion* RectangularEtaPhiTrackingRegion::clone ( ) const [inline, virtual]

clone region

Implements TrackingRegionBase.

Definition at line 131 of file RectangularEtaPhiTrackingRegion.h.

References RectangularEtaPhiTrackingRegion().

                                                         { 
    return new RectangularEtaPhiTrackingRegion(*this);
  }
OuterEstimator * RectangularEtaPhiTrackingRegion::estimator ( const BarrelDetLayer layer,
const edm::EventSetup iSetup 
) const [private]

Definition at line 125 of file RectangularEtaPhiTrackingRegion.cc.

References BoundSurface::bounds(), PixelRecoRange< T >::empty(), reco::helper::VirtualJetProducerHelper::intersection(), PixelRecoRange< T >::intersection(), Bounds::length(), PixelRecoUtilities::longitudinalBendingCorrection(), max(), PixelRecoRange< T >::max(), PixelRecoRange< T >::mean(), min, PixelRecoRange< T >::min(), GeometricSearchDet::position(), PtMinSelector_cfg::ptMin, CosmicsPD_Skims::radius, HitZCheck::range(), HitZCheck::setTolerance(), OuterHitPhiPrediction::setTolerance(), BarrelDetLayer::specificSurface(), funct::sqr(), mathSSE::sqrt(), BarrelDetLayer::surface(), Bounds::thickness(), w2, and PV3DBase< T, PVType, FrameType >::z().

{

  // det dimensions 
  float halfLength = layer->surface().bounds().length()/2;
  float halfThickness  = layer->surface().bounds().thickness()/2;
  float z0 = layer->position().z();
  float radius = layer->specificSurface().radius();

  // det ranges
  Range detRWindow (radius-halfThickness, radius+halfThickness);
  Range detZWindow(z0-halfLength,z0+halfLength);

  // z prediction, skip if not intersection
  HitZCheck zPrediction(rzConstraint());
  Range hitZWindow = zPrediction.range(detRWindow.min()).
                                               intersection(detZWindow);
  if (hitZWindow.empty()) return 0;

  // phi prediction
  OuterHitPhiPrediction phiPrediction = phiWindow(iSetup);

  //
  // optional corrections for tolerance (mult.scatt, error, bending)
  //
  OuterHitPhiPrediction::Range phiRange;
  if (thePrecise) {
    float cotTheta = (hitZWindow.mean()-origin().z()) / radius;
    float sinTheta = 1/sqrt(1+sqr(cotTheta));
    MultipleScatteringParametrisation msSigma(layer,iSetup);
    float scatt = 3 * msSigma(ptMin(), cotTheta);
    float bendR = longitudinalBendingCorrection(radius,ptMin(),iSetup);

    float hitErrRPhi = 0.;
    float hitErrZ = 0.;
    float corrPhi = (scatt+ hitErrRPhi)/radius;
    float corrZ = scatt/sinTheta + bendR*fabs(cotTheta) + hitErrZ;

    phiPrediction.setTolerance(OuterHitPhiPrediction::Margin(corrPhi,corrPhi));
    zPrediction.setTolerance(HitZCheck::Margin(corrZ,corrZ));

    //
    // hit ranges in det
    //
    OuterHitPhiPrediction::Range phi1 = phiPrediction(detRWindow.min());
    OuterHitPhiPrediction::Range phi2 = phiPrediction(detRWindow.max());
    phiRange = Range( min(phi1.min(),phi2.min()), max(phi1.max(),phi2.max()));
    Range w1 = zPrediction.range(detRWindow.min());
    Range w2 = zPrediction.range(detRWindow.max());
    hitZWindow = Range(
      min(w1.min(),w2.min()), max(w1.max(),w2.max())).intersection(detZWindow);
  }
  else {
    phiRange = phiPrediction(detRWindow.mean()); 
  }

  return new OuterEstimator(
                            OuterDetCompatibility( layer, phiRange, detRWindow, hitZWindow),
                            OuterHitCompatibility( phiPrediction, zPrediction ),
                            iSetup);
}
OuterEstimator * RectangularEtaPhiTrackingRegion::estimator ( const ForwardDetLayer layer,
const edm::EventSetup iSetup 
) const [private]

Definition at line 188 of file RectangularEtaPhiTrackingRegion.cc.

References BoundSurface::bounds(), PixelRecoRange< T >::empty(), reco::helper::VirtualJetProducerHelper::intersection(), PixelRecoRange< T >::intersection(), PixelRecoUtilities::longitudinalBendingCorrection(), PixelRecoRange< T >::max(), PixelRecoRange< T >::mean(), PixelRecoRange< T >::min(), GeometricSearchDet::position(), PtMinSelector_cfg::ptMin, OuterHitPhiPrediction::setTolerance(), ForwardDetLayer::specificSurface(), funct::sqr(), mathSSE::sqrt(), ForwardDetLayer::surface(), Bounds::thickness(), w2, and PV3DBase< T, PVType, FrameType >::z().

{

  // det dimensions, ranges
  float halfThickness  = layer->surface().bounds().thickness()/2;
  float zLayer = layer->position().z() ;
  Range detZWindow( zLayer-halfThickness, zLayer+halfThickness);
  Range detRWindow( layer->specificSurface().innerRadius(), 
                    layer->specificSurface().outerRadius());
  
  // r prediction, skip if not intersection
  HitRCheck rPrediction(rzConstraint());
  Range hitRWindow = rPrediction.range(zLayer).intersection(detRWindow);
  if (hitRWindow.empty()) return 0;

  // phi prediction
  OuterHitPhiPrediction phiPrediction = phiWindow(iSetup);
  OuterHitPhiPrediction::Range phiRange = phiPrediction(detRWindow.max());

  //
  // optional corrections for tolerance (mult.scatt, error, bending)
  //
  if (thePrecise) {
    float cotTheta = (detZWindow.mean()-origin().z())/hitRWindow.mean();
    float cosTheta = cotTheta/sqrt(1+sqr(cotTheta)); 
    MultipleScatteringParametrisation msSigma(layer,iSetup);
    float scatt = 3 * msSigma(ptMin(),cotTheta);
    float bendR = longitudinalBendingCorrection(hitRWindow.max(),ptMin(),iSetup);
    float hitErrRPhi = 0.;
    float hitErrR = 0.;
    float corrPhi = (scatt+hitErrRPhi)/detRWindow.min();
    float corrR   = scatt/fabs(cosTheta) + bendR + hitErrR;

    phiPrediction.setTolerance(OuterHitPhiPrediction::Margin(corrPhi,corrPhi));
    rPrediction.setTolerance(HitRCheck::Margin(corrR,corrR));

    //
    // hit ranges in det
    //
    Range w1,w2;
    if (zLayer > 0) {
      w1 = rPrediction.range(detZWindow.min());
      w2 = rPrediction.range(detZWindow.max());
    } else {
      w1 = rPrediction.range(detZWindow.max());
      w2 = rPrediction.range(detZWindow.min());
    }
    hitRWindow = Range(w1.min(),w2.max()).intersection(detRWindow);
  }

  return new OuterEstimator(
    OuterDetCompatibility( layer, phiRange, hitRWindow, detZWindow),
    OuterHitCompatibility( phiPrediction, rPrediction),iSetup );
}
const Range& RectangularEtaPhiTrackingRegion::etaRange ( ) const [inline]

allowed eta range [eta_min, eta_max] interval

Definition at line 112 of file RectangularEtaPhiTrackingRegion.h.

References theEtaRange.

Referenced by TrackerSeedCleaner::clean(), FastTSGFromL2Muon::clean(), FastTSGFromIOHit::clean(), and GlobalTrajectoryBuilderBase::defineRegionOfInterest().

{ return theEtaRange; }
TrackingRegion::Hits RectangularEtaPhiTrackingRegion::hits ( const edm::Event ev,
const edm::EventSetup es,
const ctfseeding::SeedingLayer layer 
) const [virtual]

get hits from layer compatible with region constraints

Implements TrackingRegion.

Definition at line 280 of file RectangularEtaPhiTrackingRegion.cc.

References alongMomentum, Reference_intrackfit_cff::barrel, newFWLiteAna::build, funct::cos(), ctfseeding::SeedingLayer::detLayer(), dir, Reference_intrackfit_cff::endcap, edm::EventSetup::get(), OuterEstimator::hitCompatibility(), ctfseeding::SeedingLayer::hits(), DetLayer::location(), LogDebug, LayerMeasurements::measurements(), MeasurementTrackerESProducer_cfi::MeasurementTracker, PV3DBase< T, PVType, FrameType >::phi(), phi, GeomDetEnumerators::PixelBarrel, GeomDetEnumerators::PixelEndcap, edm::ESHandle< T >::product(), query::result, makeMuonMisalignmentScenario::rot, funct::sin(), StraightLinePropagator_cfi::StraightLinePropagator, and DetLayer::subDetector().

{


  //ESTIMATOR
  TrackingRegion::Hits result;

  const DetLayer * detLayer = layer->detLayer();
  OuterEstimator * est = 0;

  bool measurementMethod = false;
  if ( theMeasurementTrackerUsage > 0.5) measurementMethod = true;
  if ( theMeasurementTrackerUsage > -0.5 &&
       !(detLayer->subDetector() == GeomDetEnumerators::PixelBarrel ||
         detLayer->subDetector() == GeomDetEnumerators::PixelEndcap) ) measurementMethod = true;

  if(measurementMethod) {
  edm::ESHandle<MagneticField> field;
  es.get<IdealMagneticFieldRecord>().get(field);
  const MagneticField * magField = field.product();

  const GlobalPoint vtx = origin();
  GlobalVector dir = direction();
  
  if (detLayer->subDetector() == GeomDetEnumerators::PixelBarrel || (!theUseEtaPhi  && detLayer->location() == GeomDetEnumerators::barrel)){
    const BarrelDetLayer& bl = dynamic_cast<const BarrelDetLayer&>(*detLayer);
    est = estimator(&bl,es);
  } else if (detLayer->subDetector() == GeomDetEnumerators::PixelEndcap || (!theUseEtaPhi  && detLayer->location() == GeomDetEnumerators::endcap)) {
    const ForwardDetLayer& fl = dynamic_cast<const ForwardDetLayer&>(*detLayer);
    est = estimator(&fl,es);
  }
  
  EtaPhiMeasurementEstimator etaPhiEstimator ((theEtaRange.second-theEtaRange.first)/2.,
                                              (thePhiMargin.left()+thePhiMargin.right())/2.);
  MeasurementEstimator & findDetAndHits = etaPhiEstimator;
  if (est){
    LogDebug("RectangularEtaPhiTrackingRegion")<<"use pixel specific estimator.";
    findDetAndHits =*est;
  }
  else{
    LogDebug("RectangularEtaPhiTrackingRegion")<<"use generic etat phi estimator.";
  }
   
  // TSOS
  float phi = dir.phi();
  Surface::RotationType rot( sin(phi), -cos(phi),           0,
                             0,                0,          -1,
                             cos(phi),  sin(phi),           0);

  Plane::PlanePointer surface = Plane::build(GlobalPoint(0.,0.,0.), rot);
  //TrajectoryStateOnSurface tsos(lpar, *surface, magField);

  FreeTrajectoryState fts( GlobalTrajectoryParameters(vtx, dir, 1, magField) );
  TrajectoryStateOnSurface tsos(fts, *surface);

  // propagator
  StraightLinePropagator prop( magField, alongMomentum);

  edm::ESHandle<MeasurementTracker> measurementTrackerESH;
  es.get<CkfComponentsRecord>().get(theMeasurementTrackerName, measurementTrackerESH); 
  const MeasurementTracker * measurementTracker = measurementTrackerESH.product(); 
  measurementTracker->update(ev);

  LayerMeasurements lm(measurementTracker);
   
  vector<TrajectoryMeasurement> meas = lm.measurements(*detLayer, tsos, prop, findDetAndHits);
  typedef vector<TrajectoryMeasurement>::const_iterator IM;
  for (IM im = meas.begin(); im != meas.end(); im++) {
    TrajectoryMeasurement::ConstRecHitPointer ptrHit = im->recHit();
    if (ptrHit->isValid())  result.push_back( ptrHit );
  }

  LogDebug("RectangularEtaPhiTrackingRegion")<<" found "<< meas.size()<<" minus one measurements on layer: "<<detLayer->subDetector();

  } else {
  //
  // temporary solution 
  //
  if (detLayer->location() == GeomDetEnumerators::barrel) {
    const BarrelDetLayer& bl = dynamic_cast<const BarrelDetLayer&>(*detLayer);
    est = estimator(&bl,es);
  } else {
    const ForwardDetLayer& fl = dynamic_cast<const ForwardDetLayer&>(*detLayer);
    est = estimator(&fl,es);
  }
  if (!est) return result;

  TrackingRegion::Hits layerHits = layer->hits(ev,es);
  for (TrackingRegion::Hits::const_iterator ih= layerHits.begin(); ih != layerHits.end(); ih++) {
    if ( est->hitCompatibility()(ih->get()) ) {
       result.push_back( *ih );
    }
  }
  }
  
  delete est;
  return result;
}
void RectangularEtaPhiTrackingRegion::initEtaRange ( const GlobalVector dir,
const Margin margin 
) [private]
bool RectangularEtaPhiTrackingRegion::isPrecise ( ) const [inline]

is precise error calculation switched on

Definition at line 119 of file RectangularEtaPhiTrackingRegion.h.

References thePrecise.

{ return thePrecise; }
virtual std::string RectangularEtaPhiTrackingRegion::name ( void  ) const [inline, virtual]

Reimplemented from TrackingRegion.

Definition at line 135 of file RectangularEtaPhiTrackingRegion.h.

{ return "RectangularEtaPhiTrackingRegion"; }
const Margin& RectangularEtaPhiTrackingRegion::phiMargin ( ) const [inline]

defined phi range around phi0, margin is [phi_left,phi_right]. region is defined in a range: [phi0-phi_left, phi0+phi_right]

Definition at line 116 of file RectangularEtaPhiTrackingRegion.h.

References thePhiMargin.

Referenced by TrackerSeedCleaner::clean(), FastTSGFromIOHit::clean(), FastTSGFromL2Muon::clean(), and GlobalTrajectoryBuilderBase::defineRegionOfInterest().

{ return thePhiMargin; }
OuterHitPhiPrediction RectangularEtaPhiTrackingRegion::phiWindow ( const edm::EventSetup iSetup) const [private]
std::string RectangularEtaPhiTrackingRegion::print ( void  ) const [virtual]

Reimplemented from TrackingRegionBase.

Definition at line 382 of file RectangularEtaPhiTrackingRegion.cc.

References reco::print().

                                                       {
  std::ostringstream str;
  str << TrackingRegionBase::print() 
      <<" eta: "<<theEtaRange<<" phi:"<<thePhiMargin
      << "precise: "<<thePrecise;
  return str.str();
}
HitRZConstraint RectangularEtaPhiTrackingRegion::rzConstraint ( ) const [private]

Definition at line 259 of file RectangularEtaPhiTrackingRegion.cc.

References PixelRecoPointRZ::z().

{
  HitRZConstraint::LineOrigin pLeft,pRight;
  float zMin = origin().z() - originZBound();
  float zMax = origin().z() + originZBound();
  float rMin = -originRBound();
  float rMax =  originRBound();
  if(theEtaRange.max() > 0) {
    pRight = HitRZConstraint::LineOrigin(rMin,zMax);
  } else { 
    pRight = HitRZConstraint::LineOrigin(rMax,zMax);
  } 
  if (theEtaRange.min() > 0.) {
    pLeft = HitRZConstraint::LineOrigin(rMax, zMin);
  } else {
    pLeft = HitRZConstraint::LineOrigin(rMin, zMin);
  } 
  return HitRZConstraint(pLeft, sinh(theEtaRange.min()),
                              pRight, sinh(theEtaRange.max()) );
}

Member Data Documentation

Definition at line 154 of file RectangularEtaPhiTrackingRegion.h.

Referenced by etaRange().

Definition at line 158 of file RectangularEtaPhiTrackingRegion.h.

Definition at line 156 of file RectangularEtaPhiTrackingRegion.h.

Definition at line 155 of file RectangularEtaPhiTrackingRegion.h.

Referenced by phiMargin().

Definition at line 157 of file RectangularEtaPhiTrackingRegion.h.

Referenced by isPrecise().

Definition at line 159 of file RectangularEtaPhiTrackingRegion.h.