CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTTrackClusterRemover.cc
Go to the documentation of this file.
6 
23 
26 
30 //
31 // class decleration
32 //
33 
35  public:
38  void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) ;
39  private:
40  struct ParamBlock {
42  ParamBlock(const edm::ParameterSet& iConfig) :
43  isSet_(true),
44  usesCharge_(iConfig.exists("maxCharge")),
45  usesSize_(iConfig.exists("maxSize")),
46  maxChi2_(iConfig.getParameter<double>("maxChi2")),
47  maxCharge_(usesCharge_ ? iConfig.getParameter<double>("maxCharge") : 0),
48  maxSize_(usesSize_ ? iConfig.getParameter<uint32_t>("maxSize") : 0) { }
51  size_t maxSize_;
52  };
53  static const unsigned int NumberOfParamBlocks = 6;
54 
58  bool mergeOld_;
60 
62  void readPSet(const edm::ParameterSet& iConfig, const std::string &name,
63  int id1=-1, int id2=-1, int id3=-1, int id4=-1, int id5=-1, int id6=-1) ;
64 
65  std::vector<uint8_t> pixels, strips; // avoid unneed alloc/dealloc of this
66  edm::ProductID pixelSourceProdID, stripSourceProdID; // ProdIDs refs must point to (for consistency tests)
67 
68  inline void process(const TrackingRecHit *hit, float chi2);
69  inline void process(const OmniClusterRef & cluRef, uint32_t subdet);
70 
71  template<typename T>
72  std::auto_ptr<edmNew::DetSetVector<T> >
73  cleanup(const edmNew::DetSetVector<T> &oldClusters, const std::vector<uint8_t> &isGood,
75 
76  // Carries in full removal info about a given det from oldRefs
78 
80  std::vector<bool> collectedRegStrips_;
81  std::vector<bool> collectedPixels_;
82 };
83 
84 
85 using namespace std;
86 using namespace edm;
87 using namespace reco;
88 
89 void
90 HLTTrackClusterRemover::readPSet(const edm::ParameterSet& iConfig, const std::string &name,
91  int id1, int id2, int id3, int id4, int id5, int id6)
92 {
93  if (iConfig.exists(name)) {
94  ParamBlock pblock(iConfig.getParameter<ParameterSet>(name));
95  if (id1 == -1) {
96  fill(pblocks_, pblocks_+NumberOfParamBlocks, pblock);
97  } else {
98  pblocks_[id1] = pblock;
99  if (id2 != -1) pblocks_[id2] = pblock;
100  if (id3 != -1) pblocks_[id3] = pblock;
101  if (id4 != -1) pblocks_[id4] = pblock;
102  if (id5 != -1) pblocks_[id5] = pblock;
103  if (id6 != -1) pblocks_[id6] = pblock;
104  }
105  }
106 }
107 
109  trajectories_(iConfig.getParameter<InputTag>("trajectories")),
110  doStrip_(iConfig.existsAs<bool>("doStrip") ? iConfig.getParameter<bool>("doStrip") : true),
111  doPixel_(iConfig.existsAs<bool>("doPixel") ? iConfig.getParameter<bool>("doPixel") : true),
112  stripClusters_(doStrip_ ? iConfig.getParameter<InputTag>("stripClusters") : InputTag("NONE")),
113  pixelClusters_(doPixel_ ? iConfig.getParameter<InputTag>("pixelClusters") : InputTag("NONE")),
114  mergeOld_(false),
115  oldRemovalInfo_(edm::InputTag()),
116  makeProducts_(true)
117 {
118 
119  if (iConfig.exists("oldClusterRemovalInfo"))
120  {
121  oldRemovalInfo_=iConfig.getParameter<InputTag>("oldClusterRemovalInfo");
122  if (not (oldRemovalInfo_== edm::InputTag())) mergeOld_=true;
123  }
124 
126  readPSet(iConfig, "Common",-1);
127  if (doPixel_) {
128  readPSet(iConfig, "Pixel" ,0,1);
129  readPSet(iConfig, "PXB" ,0);
130  readPSet(iConfig, "PXE" ,1);
131  }
132  if (doStrip_) {
133  readPSet(iConfig, "Strip" ,2,3,4,5);
134  readPSet(iConfig, "StripInner" ,2,3);
135  readPSet(iConfig, "StripOuter" ,4,5);
136  readPSet(iConfig, "TIB" ,2);
137  readPSet(iConfig, "TID" ,3);
138  readPSet(iConfig, "TOB" ,4);
139  readPSet(iConfig, "TEC" ,5);
140  }
141 
142  bool usingCharge = false;
143  for (size_t i = 0; i < NumberOfParamBlocks; ++i) {
144  if (!pblocks_[i].isSet_) throw cms::Exception("Configuration Error") << "HLTTrackClusterRemover: Missing configuration for detector with subDetID = " << (i+1);
145  if (pblocks_[i].usesCharge_ && !usingCharge) {
146  throw cms::Exception("Configuration Error") << "HLTTrackClusterRemover: Configuration for subDetID = " << (i+1) << " uses cluster charge, which is not enabled.";
147  }
148  }
149 
150  //produces<edmNew::DetSetVector<SiPixelClusterRefNew> >();
151  //produces<edmNew::DetSetVector<SiStripRecHit1D::ClusterRegionalRef> >();
152 
153  produces< edm::ContainerMask<edmNew::DetSetVector<SiPixelCluster> > >();
154  produces< edm::ContainerMask<edm::LazyGetter<SiStripCluster> > >();
155 
156 }
157 
158 
160 {
161 }
162 
164  const ClusterRemovalInfo::Indices &oldRefs)
165 {
166  for (size_t i = 0, n = refs.size(); i < n; ++i) {
167  refs[i] = oldRefs[refs[i]];
168  }
169 }
170 
171 
172 template<typename T>
173 auto_ptr<edmNew::DetSetVector<T> >
174 HLTTrackClusterRemover::cleanup(const edmNew::DetSetVector<T> &oldClusters, const std::vector<uint8_t> &isGood,
176  typedef typename edmNew::DetSetVector<T> DSV;
177  typedef typename edmNew::DetSetVector<T>::FastFiller DSF;
178  typedef typename edmNew::DetSet<T> DS;
179  auto_ptr<DSV> output(new DSV());
180  output->reserve(oldClusters.size(), oldClusters.dataSize());
181 
182  unsigned int countOld=0;
183  unsigned int countNew=0;
184 
185  // cluster removal loop
186  const T * firstOffset = & oldClusters.data().front();
187  for (typename DSV::const_iterator itdet = oldClusters.begin(), enddet = oldClusters.end(); itdet != enddet; ++itdet) {
188  DS oldDS = *itdet;
189 
190  if (oldDS.empty()) continue; // skip empty detsets
191 
192  uint32_t id = oldDS.detId();
193  DSF outds(*output, id);
194 
195  for (typename DS::const_iterator it = oldDS.begin(), ed = oldDS.end(); it != ed; ++it) {
196  uint32_t index = ((&*it) - firstOffset);
197  countOld++;
198  if (isGood[index]) {
199  outds.push_back(*it);
200  countNew++;
201  refs.push_back(index);
202  //std::cout << "HLTTrackClusterRemover::cleanup " << typeid(T).name() << " reference " << index << " to " << (refs.size() - 1) << std::endl;
203  }
204 
205  }
206  if (outds.empty()) outds.abort(); // not write in an empty DSV
207  }
208  // double fraction = countNew / (double) countOld;
209  // std::cout<<"fraction: "<<fraction<<std::endl;
210  if (oldRefs != 0) mergeOld(refs, *oldRefs);
211  return output;
212 }
213 
214 void HLTTrackClusterRemover::process(OmniClusterRef const & clusterReg, uint32_t subdet) {
215  if (clusterReg.id() != stripSourceProdID) throw cms::Exception("Inconsistent Data") <<
216  "HLTTrackClusterRemover: strip cluster ref from Product ID = " << clusterReg.id() <<
217  " does not match with source cluster collection (ID = " << stripSourceProdID << ")\n.";
218  if (collectedRegStrips_.size()<=clusterReg.key()){
219  edm::LogError("BadCollectionSize")<<collectedRegStrips_.size()<<" is smaller than "<<clusterReg.key();
220  assert(collectedRegStrips_.size()>clusterReg.key());
221  }
222  collectedRegStrips_[clusterReg.key()]=true;
223 }
224 
225 
226 
228  DetId detid = hit->geographicalId();
229  uint32_t subdet = detid.subdetId();
230 
231  assert ((subdet > 0) && (subdet <= NumberOfParamBlocks));
232 
233  // chi2 cut
234  if (chi2 > pblocks_[subdet-1].maxChi2_) return;
235 
236  if ((subdet == PixelSubdetector::PixelBarrel) || (subdet == PixelSubdetector::PixelEndcap)) {
237  // std::cout<<"process pxl hit"<<std::endl;
238  if (!doPixel_) return;
239  // this is a pixel, and i *know* it is
240  const SiPixelRecHit *pixelHit = static_cast<const SiPixelRecHit *>(hit);
241 
242  SiPixelRecHit::ClusterRef cluster = pixelHit->cluster();
243 
244  if (cluster.id() != pixelSourceProdID) throw cms::Exception("Inconsistent Data") <<
245  "HLTTrackClusterRemover: pixel cluster ref from Product ID = " << cluster.id() <<
246  " does not match with source cluster collection (ID = " << pixelSourceProdID << ")\n.";
247 
248  assert(cluster.id() == pixelSourceProdID);
249 //DBG// cout << "HIT NEW PIXEL DETID = " << detid.rawId() << ", Cluster [ " << cluster.key().first << " / " << cluster.key().second << " ] " << endl;
250 
251  // if requested, cut on cluster size
252  if (pblocks_[subdet-1].usesSize_ && (cluster->pixels().size() > pblocks_[subdet-1].maxSize_)) return;
253 
254  // mark as used
255  //pixels[cluster.key()] = false;
256  assert(collectedPixels_.size() > cluster.key());
257  collectedPixels_[cluster.key()]=true;
258  } else { // aka Strip
259  if (!doStrip_) return;
260  const type_info &hitType = typeid(*hit);
261  if (hitType == typeid(SiStripRecHit2D)) {
262  const SiStripRecHit2D *stripHit = static_cast<const SiStripRecHit2D *>(hit);
263 //DBG// cout << "Plain RecHit 2D: " << endl;
264  process(stripHit->omniClusterRef(),subdet);}
265  else if (hitType == typeid(SiStripRecHit1D)) {
266  const SiStripRecHit1D *hit1D = static_cast<const SiStripRecHit1D *>(hit);
267  process(hit1D->omniClusterRef(),subdet);
268  } else if (hitType == typeid(SiStripMatchedRecHit2D)) {
269  const SiStripMatchedRecHit2D *matchHit = static_cast<const SiStripMatchedRecHit2D *>(hit);
270 //DBG// cout << "Matched RecHit 2D: " << endl;
271  process(matchHit->monoClusterRef(),subdet);
272  process(matchHit->stereoClusterRef(),subdet);
273  } else if (hitType == typeid(ProjectedSiStripRecHit2D)) {
274  const ProjectedSiStripRecHit2D *projHit = static_cast<const ProjectedSiStripRecHit2D *>(hit);
275 //DBG// cout << "Projected RecHit 2D: " << endl;
276  process(projHit->originalHit().omniClusterRef(),subdet);
277  } else throw cms::Exception("NOT IMPLEMENTED") << "Don't know how to handle " << hitType.name() << " on detid " << detid.rawId() << "\n";
278  }
279 }
280 
281 /* Schematic picture of n-th step Iterative removal
282  * (that os removing clusters after n-th step tracking)
283  * clusters: [ C1 ] -> [ C2 ] -> ... -> [ Cn ] -> [ Cn + 1 ]
284  * ^ ^ ^--- OUTPUT "new" ID
285  * |-- before any removal |----- Source clusters
286  * |-- OUTPUT "old" ID |----- Hits in Traj. point here
287  * | \----- Old ClusterRemovalInfo "new" ID
288  * \-- Old ClusterRemovalInfo "old" ID
289  */
290 
291 
292 void
294 {
295  ProductID pixelOldProdID, stripOldProdID;
296 
298  if (doPixel_) {
299  iEvent.getByLabel(pixelClusters_, pixelClusters);
300  pixelSourceProdID = pixelClusters.id();
301  }
302 
304  if (doStrip_) {
305  iEvent.getByLabel(stripClusters_, stripClusters);
306  stripSourceProdID = stripClusters.id();
307  }
308 
309  //Handle<TrajTrackAssociationCollection> trajectories;
310  Handle<vector<Trajectory> > trajectories;
311  iEvent.getByLabel(trajectories_, trajectories);
312 
313  typedef edm::ContainerMask<edmNew::DetSetVector<SiPixelCluster> > PixelMaskContainer;
314  typedef edm::ContainerMask<edm::LazyGetter<SiStripCluster> > StripMaskContainer;
315  if(mergeOld_) {
318  iEvent.getByLabel(oldRemovalInfo_,oldPxlMask);
319  iEvent.getByLabel(oldRemovalInfo_,oldStrMask);
320  LogDebug("TrackClusterRemover")<<"to merge in, "<<oldStrMask->size()<<" strp and "<<oldPxlMask->size()<<" pxl";
321  oldStrMask->copyMaskTo(collectedRegStrips_);
322  oldPxlMask->copyMaskTo(collectedPixels_);
323  collectedRegStrips_.resize(stripClusters->size());
324  }else {
325  collectedRegStrips_.resize(stripClusters->size()); fill(collectedRegStrips_.begin(), collectedRegStrips_.end(), false);
326  collectedPixels_.resize(pixelClusters->dataSize()); fill(collectedPixels_.begin(), collectedPixels_.end(), false);
327  }
328 
329 
330  //for (TrajTrackAssociationCollection::const_iterator it = trajectories->begin(), ed = trajectories->end(); it != ed; ++it) {
331  // const Trajectory &tj = * it->key;
332  for (vector<Trajectory>::const_iterator it = trajectories->begin(), ed = trajectories->end(); it != ed; ++it) {
333  const Trajectory &tj = * it;
334  const vector<TrajectoryMeasurement> &tms = tj.measurements();
335  vector<TrajectoryMeasurement>::const_iterator itm, endtm;
336  for (itm = tms.begin(), endtm = tms.end(); itm != endtm; ++itm) {
337  const TrackingRecHit *hit = itm->recHit()->hit();
338  if (!hit->isValid()) continue;
339  // std::cout<<"process hit"<<std::endl;
340  process( hit, itm->estimate() );
341  }
342  }
343 
344  std::auto_ptr<StripMaskContainer> removedStripClusterMask(
345  new StripMaskContainer(edm::RefProd<edm::LazyGetter<SiStripCluster> >(stripClusters),collectedRegStrips_));
346  LogDebug("TrackClusterRemover")<<"total strip to skip: "<<std::count(collectedRegStrips_.begin(),collectedRegStrips_.end(),true);
347  iEvent.put( removedStripClusterMask );
348 
349  std::auto_ptr<PixelMaskContainer> removedPixelClusterMask(
350  new PixelMaskContainer(edm::RefProd<edmNew::DetSetVector<SiPixelCluster> >(pixelClusters),collectedPixels_));
351  LogDebug("TrackClusterRemover")<<"total pxl to skip: "<<std::count(collectedPixels_.begin(),collectedPixels_.end(),true);
352  iEvent.put( removedPixelClusterMask );
353 
354 }
355 
#define LogDebug(id)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
const_iterator begin() const
size_type dataSize() const
string fill
Definition: lumiContext.py:319
ProductID id() const
Definition: HandleBase.cc:15
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void readPSet(const edm::ParameterSet &iConfig, const std::string &name, int id1=-1, int id2=-1, int id3=-1, int id4=-1, int id5=-1, int id6=-1)
std::vector< uint8_t > strips
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::vector< uint8_t > pixels
std::vector< bool > collectedRegStrips_
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
DataContainer const & measurements() const
Definition: Trajectory.h:203
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
static const unsigned int NumberOfParamBlocks
data_type const * data(size_t cell) const
const_iterator end() const
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
std::vector< uint32_t > Indices
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
ParamBlock pblocks_[NumberOfParamBlocks]
void mergeOld(reco::ClusterRemovalInfo::Indices &refs, const reco::ClusterRemovalInfo::Indices &oldRefs)
Definition: DetId.h:20
ParamBlock(const edm::ParameterSet &iConfig)
bool isValid() const
std::vector< bool > collectedPixels_
size_type size() const
HLTTrackClusterRemover(const edm::ParameterSet &iConfig)
edm::ProductID id() const
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)
std::auto_ptr< edmNew::DetSetVector< T > > cleanup(const edmNew::DetSetVector< T > &oldClusters, const std::vector< uint8_t > &isGood, reco::ClusterRemovalInfo::Indices &refs, const reco::ClusterRemovalInfo::Indices *oldRefs)
void process(const TrackingRecHit *hit, float chi2)
DetId geographicalId() const
long double T
unsigned int key() const
const SiStripRecHit2D & originalHit() const
Pixel Reconstructed Hit.