CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTPixelClusterShapeFilter.cc
Go to the documentation of this file.
2 
3 //
4 // class declaration
5 //
6 
8 public:
11 
12 private:
13 
14  edm::InputTag inputTag_; // input tag identifying product containing pixel clusters
15  double minZ_; // beginning z-vertex position
16  double maxZ_; // end z-vertex position
17  double zStep_; // size of steps in z-vertex test
18 
19  std::vector<double> clusterPars_; //pixel cluster polynomial pars for vertex compatibility cut
20  int nhitsTrunc_; //maximum pixel clusters to apply compatibility check
21  double clusterTrunc_; //maximum vertex compatibility value for event rejection
22 
23  struct VertexHit
24  {
25  float z;
26  float r;
27  float w;
28  };
29 
30  virtual bool hltFilter(edm::Event&, const edm::EventSetup&, trigger::TriggerFilterObjectWithRefs & filterproduct);
31  int getContainedHits(const std::vector<VertexHit> &hits, double z0, double &chi);
32 
33 };
34 
42 
53 
54 //
55 // constructors and destructor
56 //
57 
59  inputTag_ (config.getParameter<edm::InputTag>("inputTag")),
60  minZ_ (config.getParameter<double>("minZ")),
61  maxZ_ (config.getParameter<double>("maxZ")),
62  zStep_ (config.getParameter<double>("zStep")),
63  clusterPars_ (config.getParameter< std::vector<double> >("clusterPars")),
64  nhitsTrunc_ (config.getParameter<int>("nhitsTrunc")),
65  clusterTrunc_ (config.getParameter<double>("clusterTrunc"))
66 {
67  LogDebug("") << "Using the " << inputTag_ << " input collection";
68 }
69 
71 {
72 }
73 
74 //
75 // member functions
76 //
77 
78 // ------------ method called to produce the data ------------
80 {
81  // All HLT filters must create and fill an HLT filter object,
82  // recording any reconstructed physics objects satisfying (or not)
83  // this HLT filter, and place it in the Event.
84 
85  // The filter object
86  if (saveTags()) filterproduct.addCollectionTag(inputTag_);
87  bool accept = true;
88 
89  // get hold of products from Event
91  event.getByLabel(inputTag_, hRecHits);
92 
93  // get tracker geometry
94  if (hRecHits.isValid()) {
95  edm::ESHandle<TrackerGeometry> trackerHandle;
96  iSetup.get<TrackerDigiGeometryRecord>().get(trackerHandle);
97  const TrackerGeometry *tgeo = trackerHandle.product();
98  const SiPixelRecHitCollection *hits = hRecHits.product();
99 
100  // loop over pixel rechits
101  int nPxlHits=0;
102  std::vector<VertexHit> vhits;
104  end = hits->data().end(); hit != end; ++hit) {
105  if (!hit->isValid())
106  continue;
107  ++nPxlHits;
108  DetId id(hit->geographicalId());
109  if(id.subdetId() != int(PixelSubdetector::PixelBarrel))
110  continue;
111  const PixelGeomDetUnit *pgdu = static_cast<const PixelGeomDetUnit*>(tgeo->idToDet(id));
112  if (1) {
113  const PixelTopology *pixTopo = &(pgdu->specificTopology());
114  std::vector<SiPixelCluster::Pixel> pixels(hit->cluster()->pixels());
115  bool pixelOnEdge = false;
116  for(std::vector<SiPixelCluster::Pixel>::const_iterator pixel = pixels.begin();
117  pixel != pixels.end(); ++pixel) {
118  int pixelX = pixel->x;
119  int pixelY = pixel->y;
120  if(pixTopo->isItEdgePixelInX(pixelX) || pixTopo->isItEdgePixelInY(pixelY)) {
121  pixelOnEdge = true;
122  break;
123  }
124  }
125  if (pixelOnEdge)
126  continue;
127  }
128 
129  LocalPoint lpos = LocalPoint(hit->localPosition().x(),
130  hit->localPosition().y(),
131  hit->localPosition().z());
132  GlobalPoint gpos = pgdu->toGlobal(lpos);
133  VertexHit vh;
134  vh.z = gpos.z();
135  vh.r = gpos.perp();
136  vh.w = hit->cluster()->sizeY();
137  vhits.push_back(vh);
138  }
139 
140  // estimate z-position from cluster lengths
141  double zest = 0.0;
142  int nhits = 0, nhits_max = 0;
143  double chi = 0, chi_max = 1e+9;
144  for(double z0 = minZ_; z0 <= maxZ_; z0 += zStep_) {
145  nhits = getContainedHits(vhits, z0, chi);
146  if(nhits == 0)
147  continue;
148  if(nhits > nhits_max) {
149  chi_max = 1e+9;
150  nhits_max = nhits;
151  }
152  if(nhits >= nhits_max && chi < chi_max) {
153  chi_max = chi;
154  zest = z0;
155  }
156  }
157 
158  chi = 0;
159  int nbest=0, nminus=0, nplus=0;
160  nbest = getContainedHits(vhits,zest,chi);
161  nminus = getContainedHits(vhits,zest-10.,chi);
162  nplus = getContainedHits(vhits,zest+10.,chi);
163 
164  double clusVtxQual=0.0;
165  if ((nminus+nplus)> 0)
166  clusVtxQual = (2.0*nbest)/(nminus+nplus); // A/B
167  else if (nbest>0)
168  clusVtxQual = 1000.0; // A/0 (set to arbitrarily large number)
169  else
170  clusVtxQual = 0; // 0/0 (already the default)
171 
172  // construct polynomial cut on cluster vertex quality vs. npixelhits
173  double polyCut=0;
174  for(unsigned int i=0; i < clusterPars_.size(); i++) {
175  polyCut += clusterPars_[i]*std::pow((double)nPxlHits,(int)i);
176  }
177  if(nPxlHits < nhitsTrunc_)
178  polyCut=0; // don't use cut below nhitsTrunc_ pixel hits
179  if(polyCut > clusterTrunc_ && clusterTrunc_ > 0)
180  polyCut=clusterTrunc_; // no cut above clusterTrunc_
181 
182  if (clusVtxQual < polyCut) accept = false;
183  }
184 
185  // return with final filter decision
186  return accept;
187 }
188 
189 int HLTPixelClusterShapeFilter::getContainedHits(const std::vector<VertexHit> &hits, double z0, double &chi)
190 {
191  // Calculate number of hits contained in v-shaped window in cluster y-width vs. z-position.
192  int n = 0;
193  chi = 0.;
194 
195  for(std::vector<VertexHit>::const_iterator hit = hits.begin(); hit!= hits.end(); hit++) {
196  double p = 2 * fabs(hit->z - z0)/hit->r + 0.5; // FIXME
197  if(fabs(p - hit->w) <= 1.) {
198  chi += fabs(p - hit->w);
199  n++;
200  }
201  }
202  return n;
203 }
204 
205 
206 // define as a framework module
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:47
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:22
virtual bool isItEdgePixelInX(int ixbin) const =0
data_type const * data(size_t cell) const
#define end
Definition: vmac.h:38
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
bool isValid() const
Definition: HandleBase.h:76
virtual const GeomDet * idToDet(DetId) const
Definition: DetId.h:20
int getContainedHits(const std::vector< VertexHit > &hits, double z0, double &chi)
const T & get() const
Definition: EventSetup.h:55
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
T const * product() const
Definition: ESHandle.h:62
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
T const * product() const
Definition: Handle.h:74
virtual bool isItEdgePixelInY(int iybin) const =0
bool saveTags() const
Definition: HLTFilter.h:45
virtual bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct)
Local3DPoint LocalPoint
Definition: LocalPoint.h:11
HLTPixelClusterShapeFilter(const edm::ParameterSet &)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40