00001 #include "FWCore/Framework/interface/Frameworkfwd.h"
00002 #include "FWCore/Framework/interface/EDProducer.h"
00003 #include "FWCore/Framework/interface/Event.h"
00004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00005 #include "FWCore/Utilities/interface/InputTag.h"
00006
00007 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
00008 #include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
00009 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
00010 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit1D.h"
00011 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
00012 #include "DataFormats/TrackerRecHit2D/interface/ProjectedSiStripRecHit2D.h"
00013 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
00014 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHitCollection.h"
00015 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2DCollection.h"
00016 #include "DataFormats/Common/interface/Handle.h"
00017 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00018 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00019 #include "DataFormats/Common/interface/DetSetVector.h"
00020 #include "DataFormats/Common/interface/DetSetVectorNew.h"
00021 #include "DataFormats/Provenance/interface/ProductID.h"
00022
00023 #include "DataFormats/TrackReco/interface/Track.h"
00024 #include "DataFormats/TrackerRecHit2D/interface/ClusterRemovalInfo.h"
00025
00026 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00027 #include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
00028 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00029
00030
00031
00032
00033 class HLTTrackClusterRemover : public edm::EDProducer {
00034 public:
00035 HLTTrackClusterRemover(const edm::ParameterSet& iConfig) ;
00036 ~HLTTrackClusterRemover() ;
00037 void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) ;
00038 private:
00039 struct ParamBlock {
00040 ParamBlock() : isSet_(false), usesCharge_(false) {}
00041 ParamBlock(const edm::ParameterSet& iConfig) :
00042 isSet_(true),
00043 usesCharge_(iConfig.exists("maxCharge")),
00044 usesSize_(iConfig.exists("maxSize")),
00045 maxChi2_(iConfig.getParameter<double>("maxChi2")),
00046 maxCharge_(usesCharge_ ? iConfig.getParameter<double>("maxCharge") : 0),
00047 maxSize_(usesSize_ ? iConfig.getParameter<uint32_t>("maxSize") : 0) { }
00048 bool isSet_, usesCharge_, usesSize_;
00049 float maxChi2_, maxCharge_;
00050 size_t maxSize_;
00051 };
00052 static const unsigned int NumberOfParamBlocks = 6;
00053
00054 edm::InputTag trajectories_;
00055 bool doStrip_, doPixel_;
00056 edm::InputTag stripClusters_, pixelClusters_;
00057 bool mergeOld_;
00058 edm::InputTag oldRemovalInfo_;
00059
00060 ParamBlock pblocks_[NumberOfParamBlocks];
00061 void readPSet(const edm::ParameterSet& iConfig, const std::string &name,
00062 int id1=-1, int id2=-1, int id3=-1, int id4=-1, int id5=-1, int id6=-1) ;
00063
00064 std::vector<uint8_t> pixels, strips;
00065 edm::ProductID pixelSourceProdID, stripSourceProdID;
00066
00067 inline void process(const TrackingRecHit *hit, float chi2);
00068 inline void process(const SiStripRecHit2D *hit2D, uint32_t subdet);
00069 inline void process(const SiStripRecHit1D *hit1D, uint32_t subdet);
00070
00071 template<typename T>
00072 std::auto_ptr<edmNew::DetSetVector<T> >
00073 cleanup(const edmNew::DetSetVector<T> &oldClusters, const std::vector<uint8_t> &isGood,
00074 reco::ClusterRemovalInfo::Indices &refs, const reco::ClusterRemovalInfo::Indices *oldRefs) ;
00075
00076
00077 void mergeOld(reco::ClusterRemovalInfo::Indices &refs, const reco::ClusterRemovalInfo::Indices &oldRefs) ;
00078
00079 bool makeProducts_;
00080 std::map<uint32_t, std::set< SiStripRecHit1D::ClusterRegionalRef > > collectedRegStrip;
00081 std::map<uint32_t, std::set< SiPixelRecHit::ClusterRef > > collectedPixel;
00082 };
00083
00084
00085 using namespace std;
00086 using namespace edm;
00087 using namespace reco;
00088
00089 void
00090 HLTTrackClusterRemover::readPSet(const edm::ParameterSet& iConfig, const std::string &name,
00091 int id1, int id2, int id3, int id4, int id5, int id6)
00092 {
00093 if (iConfig.exists(name)) {
00094 ParamBlock pblock(iConfig.getParameter<ParameterSet>(name));
00095 if (id1 == -1) {
00096 fill(pblocks_, pblocks_+NumberOfParamBlocks, pblock);
00097 } else {
00098 pblocks_[id1] = pblock;
00099 if (id2 != -1) pblocks_[id2] = pblock;
00100 if (id3 != -1) pblocks_[id3] = pblock;
00101 if (id4 != -1) pblocks_[id4] = pblock;
00102 if (id5 != -1) pblocks_[id5] = pblock;
00103 if (id6 != -1) pblocks_[id6] = pblock;
00104 }
00105 }
00106 }
00107
00108 HLTTrackClusterRemover::HLTTrackClusterRemover(const ParameterSet& iConfig):
00109 trajectories_(iConfig.getParameter<InputTag>("trajectories")),
00110 doStrip_(iConfig.existsAs<bool>("doStrip") ? iConfig.getParameter<bool>("doStrip") : true),
00111 doPixel_(iConfig.existsAs<bool>("doPixel") ? iConfig.getParameter<bool>("doPixel") : true),
00112 stripClusters_(doStrip_ ? iConfig.getParameter<InputTag>("stripClusters") : InputTag("NONE")),
00113 pixelClusters_(doPixel_ ? iConfig.getParameter<InputTag>("pixelClusters") : InputTag("NONE")),
00114 mergeOld_(false),
00115 oldRemovalInfo_(edm::InputTag()),
00116 makeProducts_(true)
00117 {
00118
00119 if (iConfig.exists("oldClusterRemovalInfo"))
00120 {
00121 oldRemovalInfo_=iConfig.getParameter<InputTag>("oldClusterRemovalInfo");
00122 if (not (oldRemovalInfo_== edm::InputTag())) mergeOld_=true;
00123 }
00124
00125 fill(pblocks_, pblocks_+NumberOfParamBlocks, ParamBlock());
00126 readPSet(iConfig, "Common",-1);
00127 if (doPixel_) {
00128 readPSet(iConfig, "Pixel" ,0,1);
00129 readPSet(iConfig, "PXB" ,0);
00130 readPSet(iConfig, "PXE" ,1);
00131 }
00132 if (doStrip_) {
00133 readPSet(iConfig, "Strip" ,2,3,4,5);
00134 readPSet(iConfig, "StripInner" ,2,3);
00135 readPSet(iConfig, "StripOuter" ,4,5);
00136 readPSet(iConfig, "TIB" ,2);
00137 readPSet(iConfig, "TID" ,3);
00138 readPSet(iConfig, "TOB" ,4);
00139 readPSet(iConfig, "TEC" ,5);
00140 }
00141
00142 bool usingCharge = false;
00143 for (size_t i = 0; i < NumberOfParamBlocks; ++i) {
00144 if (!pblocks_[i].isSet_) throw cms::Exception("Configuration Error") << "HLTTrackClusterRemover: Missing configuration for detector with subDetID = " << (i+1);
00145 if (pblocks_[i].usesCharge_ && !usingCharge) {
00146 throw cms::Exception("Configuration Error") << "HLTTrackClusterRemover: Configuration for subDetID = " << (i+1) << " uses cluster charge, which is not enabled.";
00147 }
00148 }
00149
00150 produces<edmNew::DetSetVector<SiPixelClusterRefNew> >();
00151 produces<edmNew::DetSetVector<SiStripRecHit1D::ClusterRegionalRef> >();
00152
00153 }
00154
00155
00156 HLTTrackClusterRemover::~HLTTrackClusterRemover()
00157 {
00158 }
00159
00160 void HLTTrackClusterRemover::mergeOld(ClusterRemovalInfo::Indices &refs,
00161 const ClusterRemovalInfo::Indices &oldRefs)
00162 {
00163 for (size_t i = 0, n = refs.size(); i < n; ++i) {
00164 refs[i] = oldRefs[refs[i]];
00165 }
00166 }
00167
00168
00169 template<typename T>
00170 auto_ptr<edmNew::DetSetVector<T> >
00171 HLTTrackClusterRemover::cleanup(const edmNew::DetSetVector<T> &oldClusters, const std::vector<uint8_t> &isGood,
00172 reco::ClusterRemovalInfo::Indices &refs, const reco::ClusterRemovalInfo::Indices *oldRefs){
00173 typedef typename edmNew::DetSetVector<T> DSV;
00174 typedef typename edmNew::DetSetVector<T>::FastFiller DSF;
00175 typedef typename edmNew::DetSet<T> DS;
00176 auto_ptr<DSV> output(new DSV());
00177 output->reserve(oldClusters.size(), oldClusters.dataSize());
00178
00179 unsigned int countOld=0;
00180 unsigned int countNew=0;
00181
00182
00183 const T * firstOffset = & oldClusters.data().front();
00184 for (typename DSV::const_iterator itdet = oldClusters.begin(), enddet = oldClusters.end(); itdet != enddet; ++itdet) {
00185 DS oldDS = *itdet;
00186
00187 if (oldDS.empty()) continue;
00188
00189 uint32_t id = oldDS.detId();
00190 DSF outds(*output, id);
00191
00192 for (typename DS::const_iterator it = oldDS.begin(), ed = oldDS.end(); it != ed; ++it) {
00193 uint32_t index = ((&*it) - firstOffset);
00194 countOld++;
00195 if (isGood[index]) {
00196 outds.push_back(*it);
00197 countNew++;
00198 refs.push_back(index);
00199
00200 }
00201
00202 }
00203 if (outds.empty()) outds.abort();
00204 }
00205
00206
00207 if (oldRefs != 0) mergeOld(refs, *oldRefs);
00208 return output;
00209 }
00210
00211
00212 void HLTTrackClusterRemover::process(const SiStripRecHit2D *hit, uint32_t subdet) {
00213 SiStripRecHit2D::ClusterRegionalRef clusterReg = hit->cluster_regional();
00214 if (clusterReg.id() != stripSourceProdID) throw cms::Exception("Inconsistent Data") <<
00215 "HLTTrackClusterRemover: strip cluster ref from Product ID = " << clusterReg.id() <<
00216 " does not match with source cluster collection (ID = " << stripSourceProdID << ")\n.";
00217 collectedRegStrip[hit->geographicalId()].insert(clusterReg);
00218 }
00219
00220 void HLTTrackClusterRemover::process(const SiStripRecHit1D *hit, uint32_t subdet) {
00221 SiStripRecHit2D::ClusterRegionalRef clusterReg = hit->cluster_regional();
00222 if (clusterReg.id() != stripSourceProdID) throw cms::Exception("Inconsistent Data") <<
00223 "HLTTrackClusterRemover: strip cluster ref from Product ID = " << clusterReg.id() <<
00224 " does not match with source cluster collection (ID = " << stripSourceProdID << ")\n.";
00225 collectedRegStrip[hit->geographicalId()].insert(clusterReg);
00226 }
00227
00228
00229
00230 void HLTTrackClusterRemover::process(const TrackingRecHit *hit, float chi2) {
00231 DetId detid = hit->geographicalId();
00232 uint32_t subdet = detid.subdetId();
00233
00234 assert ((subdet > 0) && (subdet <= NumberOfParamBlocks));
00235
00236
00237 if (chi2 > pblocks_[subdet-1].maxChi2_) return;
00238
00239 if ((subdet == PixelSubdetector::PixelBarrel) || (subdet == PixelSubdetector::PixelEndcap)) {
00240
00241 if (!doPixel_) return;
00242
00243 const SiPixelRecHit *pixelHit = static_cast<const SiPixelRecHit *>(hit);
00244
00245 SiPixelRecHit::ClusterRef cluster = pixelHit->cluster();
00246
00247 if (cluster.id() != pixelSourceProdID) throw cms::Exception("Inconsistent Data") <<
00248 "HLTTrackClusterRemover: pixel cluster ref from Product ID = " << cluster.id() <<
00249 " does not match with source cluster collection (ID = " << pixelSourceProdID << ")\n.";
00250
00251 assert(cluster.id() == pixelSourceProdID);
00252
00253
00254
00255 if (pblocks_[subdet-1].usesSize_ && (cluster->pixels().size() > pblocks_[subdet-1].maxSize_)) return;
00256
00257
00258
00259
00260 collectedPixel[detid.rawId()].insert(cluster);
00261 } else {
00262 if (!doStrip_) return;
00263 const type_info &hitType = typeid(*hit);
00264 if (hitType == typeid(SiStripRecHit2D)) {
00265 const SiStripRecHit2D *stripHit = static_cast<const SiStripRecHit2D *>(hit);
00266
00267 process(stripHit,subdet);}
00268 else if (hitType == typeid(SiStripRecHit1D)) {
00269 const SiStripRecHit1D *hit1D = static_cast<const SiStripRecHit1D *>(hit);
00270 process(hit1D,subdet);
00271 } else if (hitType == typeid(SiStripMatchedRecHit2D)) {
00272 const SiStripMatchedRecHit2D *matchHit = static_cast<const SiStripMatchedRecHit2D *>(hit);
00273
00274 process(matchHit->monoHit(),subdet);
00275 process(matchHit->stereoHit(),subdet);
00276 } else if (hitType == typeid(ProjectedSiStripRecHit2D)) {
00277 const ProjectedSiStripRecHit2D *projHit = static_cast<const ProjectedSiStripRecHit2D *>(hit);
00278
00279 process(&projHit->originalHit(),subdet);
00280 } else throw cms::Exception("NOT IMPLEMENTED") << "Don't know how to handle " << hitType.name() << " on detid " << detid.rawId() << "\n";
00281 }
00282 }
00283
00284
00285
00286
00287
00288
00289
00290
00291
00292
00293
00294
00295 void
00296 HLTTrackClusterRemover::produce(Event& iEvent, const EventSetup& iSetup)
00297 {
00298 ProductID pixelOldProdID, stripOldProdID;
00299
00300 Handle<edmNew::DetSetVector<SiPixelCluster> > pixelClusters;
00301 if (doPixel_) {
00302 iEvent.getByLabel(pixelClusters_, pixelClusters);
00303 pixelSourceProdID = pixelClusters.id();
00304 }
00305
00306 Handle<edm::LazyGetter<SiStripCluster> > stripClusters;
00307 if (doStrip_) {
00308 iEvent.getByLabel(stripClusters_, stripClusters);
00309 stripSourceProdID = stripClusters.id();
00310 }
00311
00312
00313 Handle<vector<Trajectory> > trajectories;
00314 iEvent.getByLabel(trajectories_, trajectories);
00315
00316
00317
00318
00319 for (vector<Trajectory>::const_iterator it = trajectories->begin(), ed = trajectories->end(); it != ed; ++it) {
00320 const Trajectory &tj = * it;
00321 const vector<TrajectoryMeasurement> &tms = tj.measurements();
00322 vector<TrajectoryMeasurement>::const_iterator itm, endtm;
00323 for (itm = tms.begin(), endtm = tms.end(); itm != endtm; ++itm) {
00324 const TrackingRecHit *hit = itm->recHit()->hit();
00325 if (!hit->isValid()) continue;
00326
00327 process( hit, itm->estimate() );
00328 }
00329 }
00330
00331
00332
00333 auto_ptr<edmNew::DetSetVector<SiPixelClusterRefNew> > removedPixelClsuterRefs(new edmNew::DetSetVector<SiPixelClusterRefNew>());
00334 auto_ptr<edmNew::DetSetVector<SiStripRecHit1D::ClusterRegionalRef> > removedStripClsuterRegRefs(new edmNew::DetSetVector<SiStripRecHit1D::ClusterRegionalRef>());
00335
00336 edm::Handle<edmNew::DetSetVector<SiPixelClusterRefNew> > oldPxlRef;
00337 edm::Handle<edmNew::DetSetVector<SiStripRecHit1D::ClusterRegionalRef> > oldStrRegRef;
00338
00339 if (mergeOld_){
00340 iEvent.getByLabel(oldRemovalInfo_,oldPxlRef);
00341 iEvent.getByLabel(oldRemovalInfo_,oldStrRegRef);
00342 LogDebug("TrackClusterRemover")<<"to merge in, "<<oldStrRegRef->size()<<" strp and "<<oldPxlRef->size()<<" pxl";
00343 }
00344
00345 if (mergeOld_){
00346 for (edmNew::DetSetVector<SiStripRecHit1D::ClusterRegionalRef>::const_iterator itOld=oldStrRegRef->begin();itOld!=oldStrRegRef->end();++itOld){
00347 uint32_t id = itOld->detId();
00348 collectedRegStrip[id].insert(itOld->begin(),itOld->end());
00349 }
00350 }
00351 for (std::map<uint32_t, std::set< SiStripRecHit1D::ClusterRegionalRef > >::iterator itskiped= collectedRegStrip.begin();
00352 itskiped!=collectedRegStrip.end();++itskiped){
00353 edmNew::DetSetVector<SiStripRecHit1D::ClusterRegionalRef>::FastFiller fill(*removedStripClsuterRegRefs, itskiped->first);
00354 for (std::set< SiStripRecHit1D::ClusterRegionalRef >::iterator topush = itskiped->second.begin();
00355 topush!=itskiped->second.end();++topush){
00356 fill.push_back(*topush);
00357 }
00358 if (fill.empty()) fill.abort();
00359 itskiped->second.clear();
00360 }
00361 iEvent.put( removedStripClsuterRegRefs );
00362
00363 if (mergeOld_){
00364 for (edmNew::DetSetVector<SiPixelClusterRefNew>::const_iterator itOld=oldPxlRef->begin();itOld!=oldPxlRef->end();++itOld){
00365 uint32_t id = itOld->detId();
00366 collectedPixel[id].insert(itOld->begin(),itOld->end());
00367 }
00368 }
00369
00370 for (std::map<uint32_t, std::set< SiPixelRecHit::ClusterRef > >::iterator itskiped= collectedPixel.begin();
00371 itskiped!=collectedPixel.end();++itskiped){
00372 edmNew::DetSetVector<SiPixelClusterRefNew>::FastFiller fill(*removedPixelClsuterRefs, itskiped->first);
00373 for (std::set< SiPixelRecHit::ClusterRef >::iterator topush = itskiped->second.begin();
00374 topush!=itskiped->second.end();++topush){
00375 fill.push_back(*topush);
00376 }
00377 if (fill.empty()) fill.abort();
00378 itskiped->second.clear();
00379 }
00380 iEvent.put( removedPixelClsuterRefs );
00381
00382 }
00383
00384 #include "FWCore/PluginManager/interface/ModuleDef.h"
00385 #include "FWCore/Framework/interface/MakerMacros.h"
00386 DEFINE_FWK_MODULE(HLTTrackClusterRemover);