CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
DTSegmentSelector Class Reference

#include <DTSegmentSelector.h>

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 19 of file DTSegmentSelector.h.

Constructor & Destructor Documentation

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

Definition at line 21 of file DTSegmentSelector.h.

21  :
22  checkNoisyChannels_(pset.getParameter<bool>("checkNoisyChannels")),
23  minHitsPhi_(pset.getParameter<int>("minHitsPhi")),
24  minHitsZ_(pset.getParameter<int>("minHitsZ")),
25  maxChi2_(pset.getParameter<double>("maxChi2")),
26  maxAnglePhi_(pset.getParameter<double>("maxAnglePhi")),
27  maxAngleZ_(pset.getParameter<double>("maxAngleZ")) {
28  }
DTSegmentSelector::~DTSegmentSelector ( )
inline

Definition at line 29 of file DTSegmentSelector.h.

References checkNoisySegment(), and operator()().

29 {}

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 DTStatusFlag::cellStatus(), and LogTrace.

Referenced by operator()(), and ~DTSegmentSelector().

84  {
85 
86  bool segmentNoisy = false;
87 
88  std::vector<DTRecHit1D>::const_iterator dtHit = dtHits.begin();
89  std::vector<DTRecHit1D>::const_iterator dtHits_end = dtHits.end();
90  for(; dtHit != dtHits_end; ++dtHit){
91  //DTRecHit1D const* dtHit1D = dynamic_cast<DTRecHit1D const*>(*recHit);
92  DTWireId wireId = dtHit->wireId();
93  // Check for noisy channels to skip them
94  bool isNoisy = false, isFEMasked = false, isTDCMasked = false, isTrigMask = false,
95  isDead = false, isNohv = false;
96  statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
97  if(isNoisy) {
98  LogTrace("Calibration") << "Wire: " << wireId << " is noisy, skipping!";
99  segmentNoisy = true; break;
100  }
101  }
102  return segmentNoisy;
103 }
#define LogTrace(id)
int cellStatus(int wheelId, int stationId, int sectorId, int slId, int layerId, int cellId, bool &noiseFlag, bool &feMask, bool &tdcMask, bool &trigMask, bool &deadFlag, bool &nohvFlag) const
get content
Definition: DTStatusFlag.h:100
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(), maxAnglePhi_, maxAngleZ_, maxChi2_, minHitsPhi_, minHitsZ_, DTRecSegment4D::phiSegment(), Geom::pi(), mps_fire::result, DTRecSegment2D::specificRecHits(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), PV3DBase< T, PVType, FrameType >::z(), and DTRecSegment4D::zSegment().

Referenced by ~DTSegmentSelector().

16  {
17 
18  bool result = true;
19 
20  /*
21  // Get the DT Geometry
22  ESHandle<DTGeometry> dtGeom;
23  eventSetup.get<MuonGeometryRecord>().get(dtGeom);
24  */
25 
27  if(checkNoisyChannels_) setup.get<DTStatusFlagRcd>().get(statusMap);
28 
29  // Get the Phi 2D segment
30  int nPhiHits = -1;
31  bool segmentNoisyPhi = false;
32  if( segment.hasPhi() ){
33  const DTChamberRecSegment2D* phiSeg = segment.phiSegment(); // phiSeg lives in the chamber RF
34  //LocalPoint phiSeg2DPosInCham = phiSeg->localPosition();
35  //LocalVector phiSeg2DDirInCham = phiSeg->localDirection();
36  std::vector<DTRecHit1D> phiRecHits = phiSeg->specificRecHits();
37  nPhiHits = phiRecHits.size();
38  if(checkNoisyChannels_) segmentNoisyPhi = checkNoisySegment(statusMap,phiRecHits);
39  //} else result = false;
40  }
41  // Get the Theta 2D segment
42  int nZHits = -1;
43  bool segmentNoisyZ = false;
44  if( segment.hasZed() ){
45  const DTSLRecSegment2D* zSeg = segment.zSegment(); // zSeg lives in the SL RF
46  //const DTSuperLayer* sl = chamber->superLayer(zSeg->superLayerId());
47  //LocalPoint zSeg2DPosInCham = chamber->toLocal(sl->toGlobal((*zSeg).localPosition()));
48  //LocalVector zSeg2DDirInCham = chamber->toLocal(sl->toGlobal((*zSeg).localDirection()));
49  std::vector<DTRecHit1D> zRecHits = zSeg->specificRecHits();
50  nZHits = zRecHits.size();
51  if(checkNoisyChannels_) segmentNoisyZ = checkNoisySegment(statusMap,zRecHits);
52  }
53 
54  // Segment selection
55  // Discard segment if it has a noisy cell
56  if(segmentNoisyPhi || segmentNoisyZ)
57  result = false;
58 
59  // 2D-segment number of hits
60  if(segment.hasPhi() && nPhiHits < minHitsPhi_)
61  result = false;
62 
63  if(segment.hasZed() && nZHits < minHitsZ_)
64  result = false;
65 
66  // Segment chi2
67  double chiSquare = segment.chi2()/segment.degreesOfFreedom();
68  if(chiSquare > maxChi2_)
69  result = false;
70 
71  // Segment angle
72  LocalVector segment4DLocalDir = segment.localDirection();
73  double angleZ = fabs( atan(segment4DLocalDir.y()/segment4DLocalDir.z())*180./Geom::pi() );
74  if( angleZ > maxAngleZ_)
75  result = false;
76 
77  double anglePhi = fabs( atan(segment4DLocalDir.x()/segment4DLocalDir.z())*180./Geom::pi() );
78  if( anglePhi > maxAnglePhi_)
79  result = false;
80 
81  return result;
82 }
T y() const
Definition: PV3DBase.h:63
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:2
T z() const
Definition: PV3DBase.h:64
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
bool checkNoisySegment(edm::ESHandle< DTStatusFlag > const &, std::vector< DTRecHit1D > const &)
constexpr double pi()
Definition: Pi.h:31
T x() const
Definition: PV3DBase.h:62

Member Data Documentation

bool DTSegmentSelector::checkNoisyChannels_
private

Definition at line 35 of file DTSegmentSelector.h.

Referenced by operator()().

double DTSegmentSelector::maxAnglePhi_
private

Definition at line 39 of file DTSegmentSelector.h.

Referenced by operator()().

double DTSegmentSelector::maxAngleZ_
private

Definition at line 40 of file DTSegmentSelector.h.

Referenced by operator()().

double DTSegmentSelector::maxChi2_
private

Definition at line 38 of file DTSegmentSelector.h.

Referenced by operator()().

int DTSegmentSelector::minHitsPhi_
private

Definition at line 36 of file DTSegmentSelector.h.

Referenced by operator()().

int DTSegmentSelector::minHitsZ_
private

Definition at line 37 of file DTSegmentSelector.h.

Referenced by operator()().