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 TrackClusterRemover : public edm::EDProducer {
00034 public:
00035 TrackClusterRemover(const edm::ParameterSet& iConfig) ;
00036 ~TrackClusterRemover() ;
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 clusterWasteSolution_;
00080 bool filterTracks_;
00081 reco::TrackBase::TrackQuality trackQuality_;
00082 std::map<uint32_t, std::set< SiStripRecHit1D::ClusterRef > > collectedStrip;
00083 std::map<uint32_t, std::set< SiPixelRecHit::ClusterRef > > collectedPixel;
00084 };
00085
00086
00087 using namespace std;
00088 using namespace edm;
00089 using namespace reco;
00090
00091 void
00092 TrackClusterRemover::readPSet(const edm::ParameterSet& iConfig, const std::string &name,
00093 int id1, int id2, int id3, int id4, int id5, int id6)
00094 {
00095 if (iConfig.exists(name)) {
00096 ParamBlock pblock(iConfig.getParameter<ParameterSet>(name));
00097 if (id1 == -1) {
00098 fill(pblocks_, pblocks_+NumberOfParamBlocks, pblock);
00099 } else {
00100 pblocks_[id1] = pblock;
00101 if (id2 != -1) pblocks_[id2] = pblock;
00102 if (id3 != -1) pblocks_[id3] = pblock;
00103 if (id4 != -1) pblocks_[id4] = pblock;
00104 if (id5 != -1) pblocks_[id5] = pblock;
00105 if (id6 != -1) pblocks_[id6] = pblock;
00106 }
00107 }
00108 }
00109
00110 TrackClusterRemover::TrackClusterRemover(const ParameterSet& iConfig):
00111 trajectories_(iConfig.getParameter<InputTag>("trajectories")),
00112 doStrip_(iConfig.existsAs<bool>("doStrip") ? iConfig.getParameter<bool>("doStrip") : true),
00113 doPixel_(iConfig.existsAs<bool>("doPixel") ? iConfig.getParameter<bool>("doPixel") : true),
00114 stripClusters_(doStrip_ ? iConfig.getParameter<InputTag>("stripClusters") : InputTag("NONE")),
00115 pixelClusters_(doPixel_ ? iConfig.getParameter<InputTag>("pixelClusters") : InputTag("NONE")),
00116 mergeOld_(iConfig.exists("oldClusterRemovalInfo")),
00117 oldRemovalInfo_(mergeOld_ ? iConfig.getParameter<InputTag>("oldClusterRemovalInfo") : InputTag("NONE")),
00118 clusterWasteSolution_(true)
00119 {
00120 if (iConfig.exists("clusterLessSolution"))
00121 clusterWasteSolution_=!iConfig.getParameter<bool>("clusterLessSolution");
00122 if (doPixel_ && clusterWasteSolution_) produces< edmNew::DetSetVector<SiPixelCluster> >();
00123 if (doStrip_ && clusterWasteSolution_) produces< edmNew::DetSetVector<SiStripCluster> >();
00124 if (clusterWasteSolution_) produces< ClusterRemovalInfo >();
00125
00126 fill(pblocks_, pblocks_+NumberOfParamBlocks, ParamBlock());
00127 readPSet(iConfig, "Common",-1);
00128 if (doPixel_) {
00129 readPSet(iConfig, "Pixel" ,0,1);
00130 readPSet(iConfig, "PXB" ,0);
00131 readPSet(iConfig, "PXE" ,1);
00132 }
00133 if (doStrip_) {
00134 readPSet(iConfig, "Strip" ,2,3,4,5);
00135 readPSet(iConfig, "StripInner" ,2,3);
00136 readPSet(iConfig, "StripOuter" ,4,5);
00137 readPSet(iConfig, "TIB" ,2);
00138 readPSet(iConfig, "TID" ,3);
00139 readPSet(iConfig, "TOB" ,4);
00140 readPSet(iConfig, "TEC" ,5);
00141 }
00142
00143 bool usingCharge = false;
00144 for (size_t i = 0; i < NumberOfParamBlocks; ++i) {
00145 if (!pblocks_[i].isSet_) throw cms::Exception("Configuration Error") << "TrackClusterRemover: Missing configuration for detector with subDetID = " << (i+1);
00146 if (pblocks_[i].usesCharge_ && !usingCharge) {
00147 throw cms::Exception("Configuration Error") << "TrackClusterRemover: Configuration for subDetID = " << (i+1) << " uses cluster charge, which is not enabled.";
00148 }
00149 }
00150
00151 if (!clusterWasteSolution_){
00152 produces<edmNew::DetSetVector<SiPixelClusterRefNew> >();
00153 produces<edmNew::DetSetVector<SiStripRecHit1D::ClusterRef> >();
00154 }
00155 trackQuality_=reco::TrackBase::undefQuality;
00156 filterTracks_=false;
00157 if (iConfig.exists("TrackQuality")){
00158 filterTracks_=true;
00159 trackQuality_=reco::TrackBase::qualityByName(iConfig.getParameter<std::string>("TrackQuality"));
00160 }
00161
00162 }
00163
00164
00165 TrackClusterRemover::~TrackClusterRemover()
00166 {
00167 }
00168
00169 void TrackClusterRemover::mergeOld(ClusterRemovalInfo::Indices &refs,
00170 const ClusterRemovalInfo::Indices &oldRefs)
00171 {
00172 for (size_t i = 0, n = refs.size(); i < n; ++i) {
00173 refs[i] = oldRefs[refs[i]];
00174 }
00175 }
00176
00177
00178 template<typename T>
00179 auto_ptr<edmNew::DetSetVector<T> >
00180 TrackClusterRemover::cleanup(const edmNew::DetSetVector<T> &oldClusters, const std::vector<uint8_t> &isGood,
00181 reco::ClusterRemovalInfo::Indices &refs, const reco::ClusterRemovalInfo::Indices *oldRefs){
00182 typedef typename edmNew::DetSetVector<T> DSV;
00183 typedef typename edmNew::DetSetVector<T>::FastFiller DSF;
00184 typedef typename edmNew::DetSet<T> DS;
00185 auto_ptr<DSV> output(new DSV());
00186 output->reserve(oldClusters.size(), oldClusters.dataSize());
00187
00188 unsigned int countOld=0;
00189 unsigned int countNew=0;
00190
00191
00192 const T * firstOffset = & oldClusters.data().front();
00193 for (typename DSV::const_iterator itdet = oldClusters.begin(), enddet = oldClusters.end(); itdet != enddet; ++itdet) {
00194 DS oldDS = *itdet;
00195
00196 if (oldDS.empty()) continue;
00197
00198 uint32_t id = oldDS.detId();
00199 DSF outds(*output, id);
00200
00201 for (typename DS::const_iterator it = oldDS.begin(), ed = oldDS.end(); it != ed; ++it) {
00202 uint32_t index = ((&*it) - firstOffset);
00203 countOld++;
00204 if (isGood[index]) {
00205 outds.push_back(*it);
00206 countNew++;
00207 refs.push_back(index);
00208
00209 }
00210
00211 }
00212 if (outds.empty()) outds.abort();
00213 }
00214
00215 if (oldRefs != 0) mergeOld(refs, *oldRefs);
00216 return output;
00217 }
00218
00219
00220 void TrackClusterRemover::process(const SiStripRecHit2D *hit, uint32_t subdet) {
00221 SiStripRecHit2D::ClusterRef cluster = hit->cluster();
00222 if (cluster.id() != stripSourceProdID) throw cms::Exception("Inconsistent Data") <<
00223 "TrackClusterRemover: strip cluster ref from Product ID = " << cluster.id() <<
00224 " does not match with source cluster collection (ID = " << stripSourceProdID << ")\n.";
00225
00226 assert(cluster.id() == stripSourceProdID);
00227 if (pblocks_[subdet-1].usesSize_ && (cluster->amplitudes().size() > pblocks_[subdet-1].maxSize_)) return;
00228
00229 strips[cluster.key()] = false;
00230 if (!clusterWasteSolution_) collectedStrip[hit->geographicalId()].insert(cluster);
00231 }
00232
00233 void TrackClusterRemover::process(const SiStripRecHit1D *hit, uint32_t subdet) {
00234 SiStripRecHit1D::ClusterRef cluster = hit->cluster();
00235 if (cluster.id() != stripSourceProdID) throw cms::Exception("Inconsistent Data") <<
00236 "TrackClusterRemover: strip cluster ref from Product ID = " << cluster.id() <<
00237 " does not match with source cluster collection (ID = " << stripSourceProdID << ")\n.";
00238
00239 assert(cluster.id() == stripSourceProdID);
00240 if (pblocks_[subdet-1].usesSize_ && (cluster->amplitudes().size() > pblocks_[subdet-1].maxSize_)) return;
00241
00242
00243 strips[cluster.key()] = false;
00244 if (!clusterWasteSolution_) collectedStrip[hit->geographicalId()].insert(cluster);
00245 }
00246
00247
00248
00249
00250 void TrackClusterRemover::process(const TrackingRecHit *hit, float chi2) {
00251 DetId detid = hit->geographicalId();
00252 uint32_t subdet = detid.subdetId();
00253
00254 assert ((subdet > 0) && (subdet <= NumberOfParamBlocks));
00255
00256
00257 if (chi2 > pblocks_[subdet-1].maxChi2_) return;
00258
00259 if ((subdet == PixelSubdetector::PixelBarrel) || (subdet == PixelSubdetector::PixelEndcap)) {
00260 if (!doPixel_) return;
00261
00262 const SiPixelRecHit *pixelHit = static_cast<const SiPixelRecHit *>(hit);
00263
00264 SiPixelRecHit::ClusterRef cluster = pixelHit->cluster();
00265
00266 if (cluster.id() != pixelSourceProdID) throw cms::Exception("Inconsistent Data") <<
00267 "TrackClusterRemover: pixel cluster ref from Product ID = " << cluster.id() <<
00268 " does not match with source cluster collection (ID = " << pixelSourceProdID << ")\n.";
00269
00270 assert(cluster.id() == pixelSourceProdID);
00271
00272
00273
00274 if (pblocks_[subdet-1].usesSize_ && (cluster->pixels().size() > pblocks_[subdet-1].maxSize_)) return;
00275
00276
00277 pixels[cluster.key()] = false;
00278
00279 if(!clusterWasteSolution_) collectedPixel[detid.rawId()].insert(cluster);
00280 } else {
00281 if (!doStrip_) return;
00282 const type_info &hitType = typeid(*hit);
00283 if (hitType == typeid(SiStripRecHit2D)) {
00284 const SiStripRecHit2D *stripHit = static_cast<const SiStripRecHit2D *>(hit);
00285
00286 process(stripHit,subdet);}
00287 else if (hitType == typeid(SiStripRecHit1D)) {
00288 const SiStripRecHit1D *hit1D = static_cast<const SiStripRecHit1D *>(hit);
00289 process(hit1D,subdet);
00290 } else if (hitType == typeid(SiStripMatchedRecHit2D)) {
00291 const SiStripMatchedRecHit2D *matchHit = static_cast<const SiStripMatchedRecHit2D *>(hit);
00292
00293 process(matchHit->monoHit(),subdet);
00294 process(matchHit->stereoHit(),subdet);
00295 } else if (hitType == typeid(ProjectedSiStripRecHit2D)) {
00296 const ProjectedSiStripRecHit2D *projHit = static_cast<const ProjectedSiStripRecHit2D *>(hit);
00297
00298 process(&projHit->originalHit(),subdet);
00299 } else throw cms::Exception("NOT IMPLEMENTED") << "Don't know how to handle " << hitType.name() << " on detid " << detid.rawId() << "\n";
00300 }
00301 }
00302
00303
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313
00314 void
00315 TrackClusterRemover::produce(Event& iEvent, const EventSetup& iSetup)
00316 {
00317 ProductID pixelOldProdID, stripOldProdID;
00318
00319 Handle<edmNew::DetSetVector<SiPixelCluster> > pixelClusters;
00320 if (doPixel_) {
00321 iEvent.getByLabel(pixelClusters_, pixelClusters);
00322 pixelSourceProdID = pixelClusters.id();
00323 }
00324
00325
00326 Handle<edmNew::DetSetVector<SiStripCluster> > stripClusters;
00327 if (doStrip_) {
00328 iEvent.getByLabel(stripClusters_, stripClusters);
00329 stripSourceProdID = stripClusters.id();
00330 }
00331
00332
00333 auto_ptr<ClusterRemovalInfo> cri;
00334 if (clusterWasteSolution_){
00335 if (doStrip_ && doPixel_) cri.reset(new ClusterRemovalInfo(pixelClusters, stripClusters));
00336 else if (doStrip_) cri.reset(new ClusterRemovalInfo(stripClusters));
00337 else if (doPixel_) cri.reset(new ClusterRemovalInfo(pixelClusters));
00338 }
00339
00340 Handle<ClusterRemovalInfo> oldRemovalInfo;
00341 if (mergeOld_ && clusterWasteSolution_) {
00342 iEvent.getByLabel(oldRemovalInfo_, oldRemovalInfo);
00343
00344 if ( (oldRemovalInfo->stripNewRefProd().id() == stripClusters.id()) &&
00345 (oldRemovalInfo->pixelNewRefProd().id() == pixelClusters.id()) ) {
00346
00347 cri->getOldClustersFrom(*oldRemovalInfo);
00348
00349 pixelOldProdID = oldRemovalInfo->pixelRefProd().id();
00350 stripOldProdID = oldRemovalInfo->stripRefProd().id();
00351
00352 } else {
00353 throw cms::Exception("Inconsistent Data") << "TrackClusterRemover: " <<
00354 "Input collection product IDs are [pixel: " << pixelClusters.id() << ", strip: " << stripClusters.id() << "] \n" <<
00355 "\t but the *old* ClusterRemovalInfo " << oldRemovalInfo_.encode() << " refers as 'new product ids' to " <<
00356 "[pixel: " << oldRemovalInfo->pixelNewRefProd().id() << ", strip: " << oldRemovalInfo->stripNewRefProd().id() << "]\n" <<
00357 "NOTA BENE: when running TrackClusterRemover with an old ClusterRemovalInfo the hits in the trajectory MUST be already re-keyed.\n";
00358 }
00359 } else {
00360 pixelOldProdID = pixelSourceProdID;
00361 stripOldProdID = stripSourceProdID;
00362 }
00363
00364 Handle<TrajTrackAssociationCollection> trajectories_totrack;
00365 iEvent.getByLabel(trajectories_,trajectories_totrack);
00366
00367 if (doStrip_) {
00368 strips.resize(stripClusters->dataSize()); fill(strips.begin(), strips.end(), true);
00369 }
00370 if (doPixel_) {
00371 pixels.resize(pixelClusters->dataSize()); fill(pixels.begin(), pixels.end(), true);
00372 }
00373
00374 TrajTrackAssociationCollection::const_iterator asst=trajectories_totrack->begin();
00375 for (;asst!=trajectories_totrack->end();++asst){
00376 const Track & track = *(asst->val);
00377 bool goodTk = (track.quality(trackQuality_));
00378 if (filterTracks_ && !goodTk) continue;
00379 const Trajectory &tj = *(asst->key);
00380 const vector<TrajectoryMeasurement> &tms = tj.measurements();
00381 vector<TrajectoryMeasurement>::const_iterator itm, endtm;
00382 for (itm = tms.begin(), endtm = tms.end(); itm != endtm; ++itm) {
00383 const TrackingRecHit *hit = itm->recHit()->hit();
00384 if (!hit->isValid()) continue;
00385 process( hit, itm->estimate() );
00386 }
00387 }
00388
00389
00390 if (doPixel_ && clusterWasteSolution_) {
00391 auto_ptr<edmNew::DetSetVector<SiPixelCluster> > newPixelClusters = cleanup(*pixelClusters, pixels,
00392 cri->pixelIndices(), mergeOld_ ? &oldRemovalInfo->pixelIndices() : 0);
00393 OrphanHandle<edmNew::DetSetVector<SiPixelCluster> > newPixels = iEvent.put(newPixelClusters);
00394
00395 cri->setNewPixelClusters(newPixels);
00396 }
00397 if (doStrip_ && clusterWasteSolution_) {
00398 auto_ptr<edmNew::DetSetVector<SiStripCluster> > newStripClusters = cleanup(*stripClusters, strips,
00399 cri->stripIndices(), mergeOld_ ? &oldRemovalInfo->stripIndices() : 0);
00400 OrphanHandle<edmNew::DetSetVector<SiStripCluster> > newStrips = iEvent.put(newStripClusters);
00401
00402 cri->setNewStripClusters(newStrips);
00403 }
00404
00405
00406 if (clusterWasteSolution_) {
00407
00408
00409
00410 iEvent.put(cri);
00411 }
00412
00413 pixels.clear(); strips.clear();
00414
00415 if (!clusterWasteSolution_){
00416 auto_ptr<edmNew::DetSetVector<SiPixelClusterRefNew> > removedPixelClsuterRefs(new edmNew::DetSetVector<SiPixelClusterRefNew>());
00417 auto_ptr<edmNew::DetSetVector<SiStripRecHit1D::ClusterRef> > removedStripClsuterRefs(new edmNew::DetSetVector<SiStripRecHit1D::ClusterRef>());
00418
00419 edm::Handle<edmNew::DetSetVector<SiPixelClusterRefNew> > oldPxlRef;
00420 edm::Handle<edmNew::DetSetVector<SiStripRecHit1D::ClusterRef> > oldStrRef;
00421 if (mergeOld_){
00422 iEvent.getByLabel(oldRemovalInfo_,oldPxlRef);
00423 iEvent.getByLabel(oldRemovalInfo_,oldStrRef);
00424 LogDebug("TrackClusterRemover")<<"to merge in, "<<oldStrRef->size()<<" strp and "<<oldPxlRef->size()<<" pxl";
00425 }
00426
00427 if (mergeOld_){
00428 for (edmNew::DetSetVector<SiStripRecHit1D::ClusterRef>::const_iterator itOld=oldStrRef->begin();itOld!=oldStrRef->end();++itOld){
00429 uint32_t id = itOld->detId();
00430 collectedStrip[id].insert(itOld->begin(),itOld->end());
00431 }
00432 }
00433 for (std::map<uint32_t, std::set< SiStripRecHit1D::ClusterRef > >::iterator itskiped= collectedStrip.begin();
00434 itskiped!=collectedStrip.end();++itskiped){
00435 edmNew::DetSetVector<SiStripRecHit1D::ClusterRef>::FastFiller fill(*removedStripClsuterRefs, itskiped->first);
00436 for (std::set< SiStripRecHit1D::ClusterRef >::iterator topush = itskiped->second.begin();
00437 topush!=itskiped->second.end();++topush){
00438 fill.push_back(*topush);
00439 LogDebug("TrackClusterRemover")<<"registering strp ref to be skip on: "<<itskiped->first<<" key: "<<topush->key();
00440 }
00441 if (fill.empty()) fill.abort();
00442 itskiped->second.clear();
00443 }
00444 LogDebug("TrackClusterRemover")<<"total strip to skip: "<<removedStripClsuterRefs->size();
00445 iEvent.put( removedStripClsuterRefs );
00446
00447 if (mergeOld_){
00448 for (edmNew::DetSetVector<SiPixelClusterRefNew>::const_iterator itOld=oldPxlRef->begin();itOld!=oldPxlRef->end();++itOld){
00449 uint32_t id = itOld->detId();
00450 collectedPixel[id].insert(itOld->begin(),itOld->end());
00451 }
00452 }
00453
00454 for (std::map<uint32_t, std::set< SiPixelRecHit::ClusterRef > >::iterator itskiped= collectedPixel.begin();
00455 itskiped!=collectedPixel.end();++itskiped){
00456 edmNew::DetSetVector<SiPixelClusterRefNew>::FastFiller fill(*removedPixelClsuterRefs, itskiped->first);
00457 for (std::set< SiPixelRecHit::ClusterRef >::iterator topush = itskiped->second.begin();
00458 topush!=itskiped->second.end();++topush){
00459 fill.push_back(*topush);
00460 LogDebug("TrackClusterRemover")<<"registering pxk ref to be skip on: "<<itskiped->first<<" key: "<<topush->key();
00461 }
00462 if (fill.empty()) fill.abort();
00463 itskiped->second.clear();
00464 }
00465 LogDebug("TrackClusterRemover")<<"total pxl to skip: "<<removedPixelClsuterRefs->size();
00466 iEvent.put( removedPixelClsuterRefs );
00467
00468 }
00469
00470
00471 }
00472
00473 #include "FWCore/PluginManager/interface/ModuleDef.h"
00474 #include "FWCore/Framework/interface/MakerMacros.h"
00475 DEFINE_FWK_MODULE(TrackClusterRemover);