CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

DTSegmentSelector Class Reference

#include <DTSegmentSelector.h>

List of all members.

Public Member Functions

 DTSegmentSelector (edm::ParameterSet const &pset)
bool operator() (DTRecSegment4D const &, edm::Event const &, edm::EventSetup const &)
 ~DTSegmentSelector ()

Private Member Functions

bool checkNoisySegment (edm::ESHandle< DTStatusFlag > const &, std::vector< DTRecHit1D > const &)

Private Attributes

bool checkNoisyChannels_
double maxAnglePhi_
double maxAngleZ_
double maxChi2_
int minHitsPhi_
int minHitsZ_

Detailed Description

Definition at line 21 of file DTSegmentSelector.h.


Constructor & Destructor Documentation

DTSegmentSelector::DTSegmentSelector ( edm::ParameterSet const &  pset) [inline]

Definition at line 23 of file DTSegmentSelector.h.

                                                    :
         checkNoisyChannels_(pset.getParameter<bool>("checkNoisyChannels")),
         minHitsPhi_(pset.getParameter<int>("minHitsPhi")),
         minHitsZ_(pset.getParameter<int>("minHitsZ")),
         maxChi2_(pset.getParameter<double>("maxChi2")),
         maxAnglePhi_(pset.getParameter<double>("maxAnglePhi")),
         maxAngleZ_(pset.getParameter<double>("maxAngleZ")) {
      }
DTSegmentSelector::~DTSegmentSelector ( ) [inline]

Definition at line 31 of file DTSegmentSelector.h.

{}

Member Function Documentation

bool DTSegmentSelector::checkNoisySegment ( edm::ESHandle< DTStatusFlag > const &  statusMap,
std::vector< DTRecHit1D > const &  dtHits 
) [private]

Definition at line 84 of file DTSegmentSelector.cc.

References LogTrace.

Referenced by operator()().

                                                                                                                        {

  bool segmentNoisy = false;

  std::vector<DTRecHit1D>::const_iterator dtHit = dtHits.begin();
  std::vector<DTRecHit1D>::const_iterator dtHits_end = dtHits.end();
  for(; dtHit != dtHits_end; ++dtHit){
     //DTRecHit1D const* dtHit1D = dynamic_cast<DTRecHit1D const*>(*recHit);
     DTWireId wireId = dtHit->wireId();
     // Check for noisy channels to skip them
     bool isNoisy = false, isFEMasked = false, isTDCMasked = false, isTrigMask = false,
          isDead = false, isNohv = false;
     statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
     if(isNoisy) {
        LogTrace("Calibration") << "Wire: " << wireId << " is noisy, skipping!";
        segmentNoisy = true; break;
     }
  }
  return segmentNoisy;
}
bool DTSegmentSelector::operator() ( DTRecSegment4D const &  segment,
edm::Event const &  event,
edm::EventSetup const &  setup 
)

Definition at line 16 of file DTSegmentSelector.cc.

References checkNoisyChannels_, checkNoisySegment(), DTRecSegment4D::chi2(), DTRecSegment4D::degreesOfFreedom(), edm::EventSetup::get(), DTRecSegment4D::hasPhi(), DTRecSegment4D::hasZed(), DTRecSegment4D::localDirection(), DTRecSegment4D::localPosition(), maxAnglePhi_, maxAngleZ_, maxChi2_, minHitsPhi_, minHitsZ_, DTRecSegment4D::phiSegment(), Geom::pi(), query::result, DTRecSegment2D::specificRecHits(), crabStatusFromReport::statusMap, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), PV3DBase< T, PVType, FrameType >::z(), and DTRecSegment4D::zSegment().

                                                                                                                   {

  bool result = true;

  /*
  // Get the DT Geometry
  ESHandle<DTGeometry> dtGeom;
  eventSetup.get<MuonGeometryRecord>().get(dtGeom);
  */
 
  edm::ESHandle<DTStatusFlag> statusMap;
  if(checkNoisyChannels_) setup.get<DTStatusFlagRcd>().get(statusMap);

  // Get the Phi 2D segment
  int nPhiHits = -1;
  bool segmentNoisyPhi = false;
  if( segment.hasPhi() ){
     const DTChamberRecSegment2D* phiSeg = segment.phiSegment();  // phiSeg lives in the chamber RF
     //LocalPoint phiSeg2DPosInCham = phiSeg->localPosition();  
     //LocalVector phiSeg2DDirInCham = phiSeg->localDirection();
     std::vector<DTRecHit1D> phiRecHits = phiSeg->specificRecHits();
     nPhiHits = phiRecHits.size(); 
     if(checkNoisyChannels_) segmentNoisyPhi = checkNoisySegment(statusMap,phiRecHits);
  }
  // Get the Theta 2D segment
  int nZHits = -1;
  bool segmentNoisyZ = false;
  if( segment.hasZed() ){
     const DTSLRecSegment2D* zSeg = segment.zSegment();  // zSeg lives in the SL RF
     //const DTSuperLayer* sl = chamber->superLayer(zSeg->superLayerId());
     //LocalPoint zSeg2DPosInCham = chamber->toLocal(sl->toGlobal((*zSeg).localPosition())); 
     //LocalVector zSeg2DDirInCham = chamber->toLocal(sl->toGlobal((*zSeg).localDirection()));
     std::vector<DTRecHit1D> zRecHits = zSeg->specificRecHits();
     nZHits = zRecHits.size();
     if(checkNoisyChannels_) segmentNoisyZ = checkNoisySegment(statusMap,zRecHits);
  } 

  // Segment selection 
  // Discard segment if it has a noisy cell
  if(segmentNoisyPhi || segmentNoisyZ)
     result = false;

  // 2D-segment number of hits
  if(segment.hasPhi() && nPhiHits < minHitsPhi_)
     result = false;

  if(segment.hasZed() && nZHits < minHitsZ_)
     result = false;

  // Segment chi2
  double chiSquare = segment.chi2()/segment.degreesOfFreedom();
  if(chiSquare > maxChi2_)
     result = false;

  // Segment angle
  LocalPoint segment4DLocalPos = segment.localPosition();
  LocalVector segment4DLocalDir = segment.localDirection();
  double angleZ = fabs( atan(segment4DLocalDir.y()/segment4DLocalDir.z())*180./Geom::pi() ); 
  if( angleZ > maxAngleZ_)
     result = false;

  double anglePhi = fabs( atan(segment4DLocalDir.x()/segment4DLocalDir.z())*180./Geom::pi() );
  if( anglePhi > maxAnglePhi_)
     result = false;

  return result;
}

Member Data Documentation

Definition at line 37 of file DTSegmentSelector.h.

Referenced by operator()().

Definition at line 41 of file DTSegmentSelector.h.

Referenced by operator()().

Definition at line 42 of file DTSegmentSelector.h.

Referenced by operator()().

double DTSegmentSelector::maxChi2_ [private]

Definition at line 40 of file DTSegmentSelector.h.

Referenced by operator()().

Definition at line 38 of file DTSegmentSelector.h.

Referenced by operator()().

Definition at line 39 of file DTSegmentSelector.h.

Referenced by operator()().