CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Protected Member Functions | Protected Attributes
StripSubClusterShapeFilterBase Class Reference

#include <StripSubClusterShapeTrajectoryFilter.h>

Inheritance diagram for StripSubClusterShapeFilterBase:
StripSubClusterShapeSeedFilter StripSubClusterShapeTrajectoryFilter

Public Member Functions

 StripSubClusterShapeFilterBase (const edm::ParameterSet &iConfig, edm::ConsumesCollector &iC)
 
virtual ~StripSubClusterShapeFilterBase ()
 

Protected Member Functions

void setEventBase (const edm::Event &, const edm::EventSetup &)
 
bool testLastHit (const TrackingRecHit *hit, const TrajectoryStateOnSurface &tsos, bool mustProject=false) const
 
bool testLastHit (const TrackingRecHit *hit, const GlobalPoint &gpos, const GlobalVector &gdir, bool mustProject=false) const
 

Protected Attributes

std::string label_
 
std::array< std::array
< uint8_t, 10 >, 7 > 
layerMask_
 
uint32_t maxNSat_
 
float maxTrimmedSizeDiffNeg_
 
float maxTrimmedSizeDiffPos_
 
float seedCutMIPs_
 
float seedCutSN_
 
float subclusterCutMIPs_
 
float subclusterCutSN_
 
float subclusterWindow_
 
edm::ESHandle
< ClusterShapeHitFilter
theFilter
 
edm::ESHandle< SiStripNoisestheNoise
 
edm::ESHandle< TrackerTopologytheTopology
 
edm::ESHandle< TrackerGeometrytheTracker
 
uint8_t trimMaxADC_
 
float trimMaxFracNeigh_
 
float trimMaxFracTotal_
 

Detailed Description

Definition at line 26 of file StripSubClusterShapeTrajectoryFilter.h.

Constructor & Destructor Documentation

StripSubClusterShapeFilterBase::StripSubClusterShapeFilterBase ( const edm::ParameterSet iConfig,
edm::ConsumesCollector iC 
)

Definition at line 125 of file StripSubClusterShapeTrajectoryFilter.cc.

References begin, end, edm::ParameterSet::exists(), edm::ParameterSet::existsAs(), lumiContext::fill, edm::ParameterSet::getParameter(), and i.

125  :
126  label_(iCfg.getUntrackedParameter<std::string>("label","")),
127  maxNSat_(iCfg.getParameter<uint32_t>("maxNSat")),
128  trimMaxADC_(iCfg.getParameter<double>("trimMaxADC")),
129  trimMaxFracTotal_(iCfg.getParameter<double>("trimMaxFracTotal")),
130  trimMaxFracNeigh_(iCfg.getParameter<double>("trimMaxFracNeigh")),
131  maxTrimmedSizeDiffPos_(iCfg.getParameter<double>("maxTrimmedSizeDiffPos")),
132  maxTrimmedSizeDiffNeg_(iCfg.getParameter<double>("maxTrimmedSizeDiffNeg")),
133  subclusterWindow_(iCfg.getParameter<double>("subclusterWindow")),
134  seedCutMIPs_(iCfg.getParameter<double>("seedCutMIPs")),
135  seedCutSN_(iCfg.getParameter<double>("seedCutSN")),
136  subclusterCutMIPs_(iCfg.getParameter<double>("subclusterCutMIPs")),
137  subclusterCutSN_(iCfg.getParameter<double>("subclusterCutSN"))
138 #ifdef StripSubClusterShapeFilterBase_COUNTERS
139  ,called_(0),saturated_(0),test_(0),passTrim_(0),failTooLarge_(0),passSC_(0),failTooNarrow_(0)
140 #endif
141 {
142  if (iCfg.exists("layerMask")) {
143  const edm::ParameterSet &iLM = iCfg.getParameter<edm::ParameterSet>("layerMask");
144  const char *ndets[4] = { "TIB", "TID", "TOB", "TEC" };
145  const int idets[4] = { 3, 4, 5, 6 };
146  for (unsigned int i = 0; i < 4; ++i) {
147  if (iLM.existsAs<bool>(ndets[i])) {
148  std::fill(layerMask_[idets[i]].begin(), layerMask_[idets[i]].end(), iLM.getParameter<bool>(ndets[i]));
149  } else {
150  layerMask_[idets[i]][0] = 2;
151  std::fill(layerMask_[idets[i]].begin()+1, layerMask_[idets[i]].end(), 0);
152  for (uint32_t lay : iLM.getParameter<std::vector<uint32_t>>(ndets[i])) {
153  layerMask_[idets[i]][lay] = 1;
154  }
155  }
156  }
157  } else {
158  for (auto & arr : layerMask_) {
159  std::fill(arr.begin(), arr.end(), 1);
160  }
161  }
162 }
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
string fill
Definition: lumiContext.py:319
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:186
std::array< std::array< uint8_t, 10 >, 7 > layerMask_
#define end
Definition: vmac.h:37
#define begin
Definition: vmac.h:30
StripSubClusterShapeFilterBase::~StripSubClusterShapeFilterBase ( )
virtual

Definition at line 165 of file StripSubClusterShapeTrajectoryFilter.cc.

References gather_cfg::cout, and label_.

166 {
167 #if 0
168  std::cout << "StripSubClusterShapeFilterBase " << label_ <<": called " << called_ << std::endl;
169  std::cout << "StripSubClusterShapeFilterBase " << label_ <<": saturated " << saturated_ << std::endl;
170  std::cout << "StripSubClusterShapeFilterBase " << label_ <<": test " << test_ << std::endl;
171  std::cout << "StripSubClusterShapeFilterBase " << label_ <<": failTooNarrow " << failTooNarrow_ << std::endl;
172  std::cout << "StripSubClusterShapeFilterBase " << label_ <<": passTrim " << passTrim_ << std::endl;
173  std::cout << "StripSubClusterShapeFilterBase " << label_ <<": passSC " << passSC_ << std::endl;
174  std::cout << "StripSubClusterShapeFilterBase " << label_ <<": failTooLarge " << failTooLarge_ << std::endl;
175 #endif
176 }
tuple cout
Definition: gather_cfg.py:121

Member Function Documentation

void StripSubClusterShapeFilterBase::setEventBase ( const edm::Event event,
const edm::EventSetup es 
)
protected

Definition at line 288 of file StripSubClusterShapeTrajectoryFilter.cc.

References edm::EventSetup::get().

Referenced by StripSubClusterShapeSeedFilter::init(), and StripSubClusterShapeTrajectoryFilter::setEvent().

289 {
290  // Get tracker geometry
292 
293  es.get<CkfComponentsRecord>().get("ClusterShapeHitFilter",theFilter);
294 
295  //Retrieve tracker topology from geometry
297 
298  es.get<SiStripNoisesRcd>().get(theNoise);
299 
300 }
edm::ESHandle< TrackerTopology > theTopology
edm::ESHandle< ClusterShapeHitFilter > theFilter
const T & get() const
Definition: EventSetup.h:56
bool StripSubClusterShapeFilterBase::testLastHit ( const TrackingRecHit hit,
const TrajectoryStateOnSurface tsos,
bool  mustProject = false 
) const
protected

Definition at line 180 of file StripSubClusterShapeTrajectoryFilter.cc.

References TrajectoryStateOnSurface::globalDirection(), and TrajectoryStateOnSurface::globalPosition().

181 {
182  return testLastHit(hit, tsos.globalPosition(), tsos.globalDirection(), mustProject);
183 
184 }
GlobalPoint globalPosition() const
bool testLastHit(const TrackingRecHit *hit, const TrajectoryStateOnSurface &tsos, bool mustProject=false) const
GlobalVector globalDirection() const
bool StripSubClusterShapeFilterBase::testLastHit ( const TrackingRecHit hit,
const GlobalPoint gpos,
const GlobalVector gdir,
bool  mustProject = false 
) const
protected

Definition at line 186 of file StripSubClusterShapeTrajectoryFilter.cc.

References funct::abs(), SiStripCluster::amplitudes(), begin, Surface::bounds(), f, SiStripCluster::firstStrip(), TrackingRecHit::geographicalId(), i, INC_COUNTER, prof2calltree::last, HLT_25ns14e33_v1_cff::MeVperADCStrip, SiStripMatchedRecHit2D::monoHit(), gen::n, ProjectedSiStripRecHit2D::originalHit(), DetId::rawId(), SiStripMatchedRecHit2D::stereoHit(), TrackerSingleRecHit::stripCluster(), DetId::subdetId(), GeomDet::surface(), run_regression::test, Bounds::thickness(), GeomDet::toLocal(), and PV3DBase< T, PVType, FrameType >::z().

187 {
188  const TrackerSingleRecHit *stripHit = 0;
189  if (typeid(*hit) == typeid(SiStripMatchedRecHit2D)) {
190  const SiStripMatchedRecHit2D & mhit = static_cast<const SiStripMatchedRecHit2D &>(*hit);
191  SiStripRecHit2D mono = mhit.monoHit();
192  SiStripRecHit2D stereo = mhit.stereoHit();
193  return testLastHit(&mono, gpos, gdir, true) && testLastHit(&stereo, gpos, gdir, true);
194  } else if (typeid(*hit) == typeid(ProjectedSiStripRecHit2D)) {
195  const ProjectedSiStripRecHit2D & mhit = static_cast<const ProjectedSiStripRecHit2D &>(*hit);
196  const SiStripRecHit2D &orig = mhit.originalHit();
197  return testLastHit(&orig, gpos, gdir, true);
198  } else if ((stripHit = dynamic_cast<const TrackerSingleRecHit *>(hit)) != 0) {
199  DetId detId = hit->geographicalId();
200 
201  if (layerMask_[detId.subdetId()][0] == 0) {
202  return true; // no filtering here
203  } else if (layerMask_[detId.subdetId()][0] == 2) {
204  unsigned int ilayer = theTopology->layer(detId);
205  if (layerMask_[detId.subdetId()][ilayer] == 0) {
206  return true; // no filtering here
207  }
208  }
209 
210  const GeomDet *det = theTracker->idToDet(detId);
211  LocalVector ldir = det->toLocal(gdir);
212  LocalPoint lpos = det->toLocal(gpos);
213  if (mustProject) {
214  lpos -= ldir * lpos.z()/ldir.z();
215  }
216  int hitStrips; float hitPredPos;
217  const SiStripCluster &cluster = stripHit->stripCluster();
218  bool usable = theFilter->getSizes(detId, cluster, lpos, ldir, hitStrips, hitPredPos);
219  if (!usable) return true;
220 
221  INC_COUNTER(called_)
222  const auto &ampls = cluster.amplitudes();
223 
224  // pass-through of trivial case
225  if (std::abs(hitPredPos) < 1.5f && hitStrips <= 2) {
226  return true;
227  }
228 
229 
230  // compute number of consecutive saturated strips
231  // (i.e. with adc count >= 254, see SiStripCluster class for documentation)
232  unsigned int thisSat = (ampls[0] >= 254), maxSat = thisSat;
233  for (unsigned int i = 1, n = ampls.size(); i < n; ++i) {
234  if (ampls[i] >= 254) {
235  thisSat++;
236  } else if (thisSat > 0) {
237  maxSat = std::max<int>(maxSat, thisSat);
238  thisSat = 0;
239  }
240  }
241  if (thisSat > 0) {
242  maxSat = std::max<int>(maxSat, thisSat);
243  }
244  if (maxSat >= maxNSat_) {
245  INC_COUNTER(saturated_)
246  return true;
247  }
248 
249  // trimming
250  INC_COUNTER(test_)
251  unsigned int hitStripsTrim = ampls.size();
252  int sum = std::accumulate(ampls.begin(), ampls.end(), 0);
253  uint8_t trimCut = std::min<uint8_t>(trimMaxADC_, std::floor(trimMaxFracTotal_ * sum));
254  auto begin = ampls.begin();
255  auto last = ampls.end()-1;
256  while (hitStripsTrim > 1 && (*begin < std::max<uint8_t>(trimCut, trimMaxFracNeigh_*(*(begin+1)))) ) { hitStripsTrim--; ++begin; }
257  while (hitStripsTrim > 1 && (*last < std::max<uint8_t>(trimCut, trimMaxFracNeigh_*(*(last -1)))) ) { hitStripsTrim--; --last; }
258 
259  if (hitStripsTrim < std::floor(std::abs(hitPredPos)-maxTrimmedSizeDiffNeg_)) {
260  INC_COUNTER(failTooNarrow_)
261  return false;
262  } else if (hitStripsTrim <= std::ceil(std::abs(hitPredPos)+maxTrimmedSizeDiffPos_)) {
263  INC_COUNTER(passTrim_)
264  return true;
265  }
266 
267  const StripGeomDetUnit* stripDetUnit = dynamic_cast<const StripGeomDetUnit *>(det);
268  if (det == 0) { edm::LogError("Strip not a StripGeomDetUnit?") << " on " << detId.rawId() << "\n"; return true; }
269 
270  float MeVperADCStrip = 9.5665E-4; // conversion constant from ADC counts to MeV for the SiStrip detector
271  float mip = 3.9 / ( MeVperADCStrip/stripDetUnit->surface().bounds().thickness() ); // 3.9 MeV/cm = ionization in silicon
272  float mipnorm = mip/std::abs(ldir.z());
273  ::SlidingPeakFinder pf(std::max<int>(2,std::ceil(std::abs(hitPredPos)+subclusterWindow_)));
274  ::PeakFinderTest test(mipnorm, detId(), cluster.firstStrip(), &*theNoise, seedCutMIPs_, seedCutSN_, subclusterCutMIPs_, subclusterCutSN_);
275  if (pf.apply(cluster.amplitudes(), test)) {
276  INC_COUNTER(passSC_)
277  return true;
278  } else {
279  INC_COUNTER(failTooLarge_)
280  return false;
281  }
282 
283  }
284  return true;
285 }
int i
Definition: DBlmapReader.cc:9
edm::ESHandle< TrackerTopology > theTopology
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:67
uint16_t firstStrip() const
std::array< std::array< uint8_t, 10 >, 7 > layerMask_
return((rh^lh)&mask)
else
Definition: XrdSource.cc:104
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
edm::ESHandle< ClusterShapeHitFilter > theFilter
T z() const
Definition: PV3DBase.h:64
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
if(c.getParameter< edm::InputTag >("puppiValueMap").label().size()!=0)
double f[11][100]
bool testLastHit(const TrackingRecHit *hit, const TrajectoryStateOnSurface &tsos, bool mustProject=false) const
#define end
Definition: vmac.h:37
T min(T a, T b)
Definition: MathUtil.h:58
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
SiStripRecHit2D originalHit() const
Definition: DetId.h:18
SiStripRecHit2D stereoHit() const
SiStripCluster const & stripCluster() const
string const
Definition: compareJSON.py:14
SiStripRecHit2D monoHit() const
#define begin
Definition: vmac.h:30
DetId geographicalId() const
volatile std::atomic< bool > shutdown_flag false
tuple size
Write out results.
const std::vector< uint8_t > & amplitudes() const

Member Data Documentation

std::string StripSubClusterShapeFilterBase::label_
protected
std::array<std::array<uint8_t,10>, 7> StripSubClusterShapeFilterBase::layerMask_
protected

Definition at line 57 of file StripSubClusterShapeTrajectoryFilter.h.

uint32_t StripSubClusterShapeFilterBase::maxNSat_
protected

Definition at line 42 of file StripSubClusterShapeTrajectoryFilter.h.

float StripSubClusterShapeFilterBase::maxTrimmedSizeDiffNeg_
protected

Definition at line 49 of file StripSubClusterShapeTrajectoryFilter.h.

float StripSubClusterShapeFilterBase::maxTrimmedSizeDiffPos_
protected

Definition at line 49 of file StripSubClusterShapeTrajectoryFilter.h.

float StripSubClusterShapeFilterBase::seedCutMIPs_
protected

Definition at line 53 of file StripSubClusterShapeTrajectoryFilter.h.

float StripSubClusterShapeFilterBase::seedCutSN_
protected

Definition at line 53 of file StripSubClusterShapeTrajectoryFilter.h.

float StripSubClusterShapeFilterBase::subclusterCutMIPs_
protected

Definition at line 54 of file StripSubClusterShapeTrajectoryFilter.h.

float StripSubClusterShapeFilterBase::subclusterCutSN_
protected

Definition at line 54 of file StripSubClusterShapeTrajectoryFilter.h.

float StripSubClusterShapeFilterBase::subclusterWindow_
protected

Definition at line 52 of file StripSubClusterShapeTrajectoryFilter.h.

edm::ESHandle<ClusterShapeHitFilter> StripSubClusterShapeFilterBase::theFilter
protected

Definition at line 64 of file StripSubClusterShapeTrajectoryFilter.h.

edm::ESHandle<SiStripNoises> StripSubClusterShapeFilterBase::theNoise
protected

Definition at line 65 of file StripSubClusterShapeTrajectoryFilter.h.

edm::ESHandle<TrackerTopology> StripSubClusterShapeFilterBase::theTopology
protected

Definition at line 66 of file StripSubClusterShapeTrajectoryFilter.h.

edm::ESHandle<TrackerGeometry> StripSubClusterShapeFilterBase::theTracker
protected

Definition at line 63 of file StripSubClusterShapeTrajectoryFilter.h.

uint8_t StripSubClusterShapeFilterBase::trimMaxADC_
protected

Definition at line 45 of file StripSubClusterShapeTrajectoryFilter.h.

float StripSubClusterShapeFilterBase::trimMaxFracNeigh_
protected

Definition at line 46 of file StripSubClusterShapeTrajectoryFilter.h.

float StripSubClusterShapeFilterBase::trimMaxFracTotal_
protected

Definition at line 46 of file StripSubClusterShapeTrajectoryFilter.h.