CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

OuterDetCompatibility Class Reference

#include <OuterDetCompatibility.h>

List of all members.

Public Member Functions

GlobalPoint center () const
MeasurementEstimator::Local2DVector maximalLocalDisplacement (const GlobalPoint &ts, const BoundPlane &plane) const
MeasurementEstimator::Local2DVector maximalLocalDisplacement (const TrajectoryStateOnSurface &ts, const BoundPlane &plane) const
bool operator() (const BoundPlane &plane) const
 OuterDetCompatibility (const ForwardDetLayer *layer, const OuterHitPhiPrediction::Range &phiRange, const HitRZConstraint::Range &rRange, const HitRZConstraint::Range &zRange)
 OuterDetCompatibility (const BarrelDetLayer *layer, const OuterHitPhiPrediction::Range &phiRange, const HitRZConstraint::Range &rRange, const HitRZConstraint::Range &zRange)
const
OuterHitPhiPrediction::Range
phiRange () const
const HitRZConstraint::RangerRange () const
const HitRZConstraint::RangezRange () const

Private Member Functions

bool checkPhi (const OuterHitPhiPrediction::Range &detPhiRange) const
bool checkR (const HitRZConstraint::Range &detRRange) const
bool checkZ (const HitRZConstraint::Range &detZRange) const
double loc_dist (double radius, double ts_phi, double range_phi, double cosGamma) const

Private Attributes

bool barrel
OuterHitPhiPrediction::Range hitDetPhiRange
HitRZConstraint::Range hitDetRRange
HitRZConstraint::Range hitDetZRange
const DetLayertheLayer

Detailed Description

check det compatibility by comparistion of det BoundPlane ranges with phi,r,z ranges (given at construction).

Definition at line 13 of file OuterDetCompatibility.h.


Constructor & Destructor Documentation

OuterDetCompatibility::OuterDetCompatibility ( const BarrelDetLayer layer,
const OuterHitPhiPrediction::Range phiRange,
const HitRZConstraint::Range rRange,
const HitRZConstraint::Range zRange 
) [inline]

Definition at line 16 of file OuterDetCompatibility.h.

    : theLayer(layer), barrel(true), 
      hitDetPhiRange(phiRange), hitDetRRange(rRange), hitDetZRange(zRange) { } 
OuterDetCompatibility::OuterDetCompatibility ( const ForwardDetLayer layer,
const OuterHitPhiPrediction::Range phiRange,
const HitRZConstraint::Range rRange,
const HitRZConstraint::Range zRange 
) [inline]

Definition at line 23 of file OuterDetCompatibility.h.

    : theLayer(layer), barrel(false), 
      hitDetPhiRange(phiRange), hitDetRRange(rRange), hitDetZRange(zRange) { } 

Member Function Documentation

GlobalPoint OuterDetCompatibility::center ( ) const

Definition at line 37 of file OuterDetCompatibility.cc.

References funct::cos(), phi, alignCSCRings::r, and funct::sin().

Referenced by OuterEstimator::center().

{
  float phi = hitDetPhiRange.mean();
  float r   = hitDetRRange.mean();
  return GlobalPoint( r*cos(phi), r*sin(phi), hitDetZRange.mean() );
}
bool OuterDetCompatibility::checkPhi ( const OuterHitPhiPrediction::Range detPhiRange) const [private]

Definition at line 24 of file OuterDetCompatibility.cc.

References rangesIntersect().

{ return rangesIntersect(detPhiRange, hitDetPhiRange, PhiLess()); }
bool OuterDetCompatibility::checkR ( const HitRZConstraint::Range detRRange) const [private]

Definition at line 28 of file OuterDetCompatibility.cc.

References rangesIntersect().

{ return rangesIntersect(detRRange, hitDetRRange); }
bool OuterDetCompatibility::checkZ ( const HitRZConstraint::Range detZRange) const [private]

Definition at line 32 of file OuterDetCompatibility.cc.

References rangesIntersect().

{ return rangesIntersect(detZRange, hitDetZRange); }
double OuterDetCompatibility::loc_dist ( double  radius,
double  ts_phi,
double  range_phi,
double  cosGamma 
) const [private]

Definition at line 117 of file OuterDetCompatibility.cc.

References funct::sin(), and mathSSE::sqrt().

{
    double sinDphi = sin(ts_phi - range_phi);
    double cosDphi = sqrt(1-sinDphi*sinDphi);
    double sinGamma = sqrt(1-cosGamma*cosGamma);
    double sinBeta = fabs(cosDphi*cosGamma -  sinDphi* sinGamma); 
    return radius * fabs(sinDphi) / sinBeta; 
}
MeasurementEstimator::Local2DVector OuterDetCompatibility::maximalLocalDisplacement ( const TrajectoryStateOnSurface ts,
const BoundPlane plane 
) const
MeasurementEstimator::Local2DVector OuterDetCompatibility::maximalLocalDisplacement ( const GlobalPoint ts,
const BoundPlane plane 
) const

Definition at line 45 of file OuterDetCompatibility.cc.

References Reference_intrackfit_cff::barrel, funct::cos(), Vector3DBase< T, FrameTag >::dot(), M_PI, max(), Plane::normalVector(), PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), CosmicsPD_Skims::radius, funct::sin(), GloballyPositioned< T >::toLocal(), csvLumiCalc::unit, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

{
  float x_loc = 0.;
  float y_loc = 0.;
  if(barrel) {
    double  radius = ts.perp();
    GlobalVector planeNorm = plane.normalVector();
    GlobalVector tsDir = GlobalVector( ts.x(), ts.y(),0. ).unit();
    double ts_phi = tsDir.phi();
    if (! hitDetPhiRange.inside(ts_phi) ) {
     while (ts_phi >= hitDetPhiRange.max() ) ts_phi -= 2*M_PI;
     while (ts_phi < hitDetPhiRange.min() ) ts_phi += 2*M_PI;
     if (!hitDetPhiRange.inside(ts_phi)) return MeasurementEstimator::Local2DVector(0.,0.);
    }
    double cosGamma = tsDir.dot(planeNorm);

    double dx1 = loc_dist( radius, ts_phi, hitDetPhiRange.min(), cosGamma);
    double dx2 = loc_dist( radius, ts_phi, hitDetPhiRange.max(), cosGamma);

    double ts_z = ts.z();
    double dy1 = ts_z - hitDetZRange.min();
    double dy2 = hitDetZRange.max() - ts_z;

    x_loc = max(dx1,dx2);
    y_loc = max(dy1,dy2);

    // debug only
/*
    double r1 = dx1 * fabs(cosGamma) / sin(ts_phi-hitDetPhiRange.min());
    double r2 = dx2 * fabs(cosGamma) / sin(hitDetPhiRange.max()-ts_phi);
    GlobalPoint p1( r1* cos(hitDetPhiRange.min()), r1 * sin(hitDetPhiRange.min()), hitDetZRange.min());
    GlobalPoint p2( r2* cos(hitDetPhiRange.max()), r2 * sin(hitDetPhiRange.max()), hitDetZRange.min());
    GlobalPoint p3( r1* cos(hitDetPhiRange.min()), r1 * sin(hitDetPhiRange.min()), hitDetZRange.max());
    GlobalPoint p4( r2* cos(hitDetPhiRange.max()), r2 * sin(hitDetPhiRange.max()), hitDetZRange.max());
    cout << " Local1: " << plane.toLocal(ts-p1) << endl;
    cout << " Local2: " << plane.toLocal(ts-p2) << endl;
    cout << " Local3: " << plane.toLocal(ts-p3) << endl;
    cout << " Local4: " << plane.toLocal(ts-p4) << endl;
*/
  }
  else {
    LocalPoint ts_loc = plane.toLocal(ts);
    GlobalVector planeNorm = plane.normalVector();

    double x_glob[4], y_glob[4], z_glob[4];
    x_glob[0] = hitDetRRange.min()*cos(hitDetPhiRange.min());
    y_glob[0] = hitDetRRange.min()*sin(hitDetPhiRange.min());
    x_glob[1] = hitDetRRange.max()*cos(hitDetPhiRange.min());
    y_glob[1] = hitDetRRange.max()*sin(hitDetPhiRange.min());
    x_glob[2] = hitDetRRange.min()*cos(hitDetPhiRange.max());
    y_glob[2] = hitDetRRange.min()*sin(hitDetPhiRange.max());
    x_glob[3] = hitDetRRange.max()*cos(hitDetPhiRange.max());
    y_glob[3] = hitDetRRange.max()*sin(hitDetPhiRange.max());

    for (int idx = 0; idx < 4; idx++) {
      double dx_glob = x_glob[idx] - ts.x();
      double dy_glob = y_glob[idx] - ts.y();
      double dz_glob = -(dx_glob * planeNorm.x() + dy_glob*planeNorm.y()) / planeNorm.z();
      z_glob[idx] = dz_glob + ts.z();
    }

    for (int idx=0; idx <4; idx++) {
      LocalPoint lp = plane.toLocal( GlobalPoint( x_glob[idx], y_glob[idx], z_glob[idx]));
      x_loc = max(x_loc, fabs(lp.x()-ts_loc.x()));
      y_loc = max(y_loc, fabs(lp.y()-ts_loc.y())); 
    }
  }
  MeasurementEstimator::Local2DVector distance(x_loc,y_loc);
  return distance;
}
bool OuterDetCompatibility::operator() ( const BoundPlane plane) const

Definition at line 9 of file OuterDetCompatibility.cc.

References Reference_intrackfit_cff::barrel, GlobalDetRangeRPhi::phiRange(), GlobalDetRangeZPhi::phiRange(), GlobalDetRangeRPhi::rRange(), and GlobalDetRangeZPhi::zRange().

{
  if (barrel) {
    GlobalDetRangeZPhi detRange(plane);
    if (!checkPhi(detRange.phiRange())) return 0;
    if (!checkZ(detRange.zRange())) return 0;
  } else {
    GlobalDetRangeRPhi detRange(plane);
    if (!checkPhi(detRange.phiRange())) return 0;
    if (!checkR(detRange.rRange())) return 0;
  }
  return 1;
}
const OuterHitPhiPrediction::Range& OuterDetCompatibility::phiRange ( ) const [inline]

Definition at line 39 of file OuterDetCompatibility.h.

References hitDetPhiRange.

{return hitDetPhiRange;}
const HitRZConstraint::Range& OuterDetCompatibility::rRange ( ) const [inline]

Definition at line 40 of file OuterDetCompatibility.h.

References hitDetRRange.

{ return hitDetRRange; }
const HitRZConstraint::Range& OuterDetCompatibility::zRange ( ) const [inline]

Definition at line 41 of file OuterDetCompatibility.h.

References hitDetZRange.

{ return hitDetZRange; }

Member Data Documentation

Definition at line 53 of file OuterDetCompatibility.h.

Definition at line 54 of file OuterDetCompatibility.h.

Referenced by phiRange().

Definition at line 55 of file OuterDetCompatibility.h.

Referenced by rRange().

Definition at line 55 of file OuterDetCompatibility.h.

Referenced by zRange().

Definition at line 52 of file OuterDetCompatibility.h.