00001 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h" 00002 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h" 00003 #include "Geometry/CommonDetUnit/interface/GeomDetType.h" 00004 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h" 00005 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h" 00006 #include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h" 00007 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHitCollection.h" 00008 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2DCollection.h" 00009 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2DCollection.h" 00010 #include "DataFormats/TrackingRecHit/interface/TrackingRecHitFwd.h" 00011 #include "DataFormats/Common/interface/Handle.h" 00012 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h" 00013 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h" 00014 #include "DataFormats/Common/interface/DetSetVector.h" 00015 #include "RecoLocalTracker/SubCollectionProducers/interface/RemainingClusterProducer.h" 00016 // 00017 // class decleration 00018 // 00019 using namespace std; 00020 using namespace edm; 00021 00022 00023 00024 RemainingClusterProducer::RemainingClusterProducer(const ParameterSet& iConfig): 00025 conf_(iConfig) 00026 { 00027 produces< DetSetVector<SiPixelCluster> >(); 00028 produces< DetSetVector<SiStripCluster> > (); 00029 00030 00031 } 00032 00033 00034 RemainingClusterProducer::~RemainingClusterProducer() 00035 { 00036 00037 00038 00039 } 00040 00041 00042 00043 void 00044 RemainingClusterProducer::produce(Event& iEvent, const EventSetup& iSetup) 00045 { 00046 00047 00048 using namespace edm; 00049 InputTag rphirecHitsTag = conf_.getParameter<InputTag>("rphirecHits"); 00050 InputTag stereorecHitsTag = conf_.getParameter<InputTag>("stereorecHits"); 00051 InputTag pixelTag = conf_.getParameter<InputTag>("pixelHits"); 00052 InputTag tkTag = conf_.getParameter<InputTag>("recTracks"); 00053 InputTag matchedrecHitsTag = conf_.getParameter<InputTag>("matchedRecHits"); 00054 00055 //HANDLE 00056 Handle<SiStripMatchedRecHit2DCollection> matchedrecHits; 00057 Handle<SiStripRecHit2DCollection> rphirecHits; 00058 Handle<SiStripRecHit2DCollection> stereorecHits; 00059 Handle<SiPixelRecHitCollection> pixelHits; 00060 Handle<TrackingRecHitCollection> trackingHits; 00061 00062 //GET COLLECTIONS 00063 iEvent.getByLabel( matchedrecHitsTag , matchedrecHits); 00064 iEvent.getByLabel( rphirecHitsTag , rphirecHits); 00065 iEvent.getByLabel( stereorecHitsTag , stereorecHits); 00066 iEvent.getByLabel( pixelTag , pixelHits); 00067 iEvent.getByLabel( tkTag , trackingHits); 00068 00069 // 00070 thePixelClusterVector.clear(); 00071 theStripClusterVector.clear(); 00072 00073 00074 00075 //DETSETVECTOR OF TRACKING HITS 00076 DetSetVector<const SiPixelRecHit *> pixelInTrack; 00077 DetSetVector<const SiStripRecHit2D*> monoInTrack; 00078 DetSetVector<const SiStripRecHit2D*> stereoInTrack; 00079 00080 for(TrackerGeometry::DetContainer::const_iterator it = 00081 pDD->dets().begin(); it != pDD->dets().end(); it++){ 00082 DetId detid = ((*it)->geographicalId()); 00083 unsigned int isub=detid.subdetId(); 00084 //PIXEL 00085 if ((isub== PixelSubdetector::PixelBarrel) || (isub== PixelSubdetector::PixelEndcap)) { 00086 SiPixelRecHitCollection::range pixelrechitRange = (pixelHits.product())->get((detid)); 00087 if (pixelrechitRange.second-pixelrechitRange.first!=0){ 00088 DetSet<const SiPixelRecHit *> Pix_tk(detid.rawId()); 00089 pixelInTrack.insert(Pix_tk); 00090 } 00091 continue; // se e' pixel non e' tracker, no? 00092 } 00093 00094 //std::cout << "detOD: " << detid.rawId() << std::endl; 00095 //STRIP 00096 if ((uint(detid.subdetId())==StripSubdetector::TIB) || 00097 (uint(detid.subdetId())==StripSubdetector::TOB) || 00098 (uint(detid.subdetId())==StripSubdetector::TID) || 00099 (uint(detid.subdetId())==StripSubdetector::TEC)){ 00100 //rphi 00101 SiStripRecHit2DCollection::range monorechitRange = (rphirecHits.product())->get((detid)); 00102 if (monorechitRange.second-monorechitRange.first!=0){ 00103 //std::cout << "++ IN MONO: " << detid.rawId() << std::endl; 00104 DetSet<const SiStripRecHit2D*> Mono_tk(detid.rawId()); 00105 monoInTrack.insert(Mono_tk); 00106 } 00107 //stereo 00108 SiStripRecHit2DCollection::range stereorechitRange = (stereorecHits.product())->get((detid)); 00109 if (stereorechitRange.second-stereorechitRange.first!=0){ 00110 //std::cout << "++ IN STEREO: " << detid.rawId() << std::endl; 00111 DetSet<const SiStripRecHit2D*> Stereo_tk(detid.rawId()); 00112 stereoInTrack.insert(Stereo_tk); 00113 } 00114 } 00115 } 00116 00117 //FILL THE DETSETVECTOR 00118 TrackingRecHitCollection::const_iterator hit; 00119 for(hit=trackingHits.product()->begin();hit!=trackingHits.product()->end();hit++){ 00120 if ((*hit).isValid()){ 00121 const SiPixelRecHit* pixelhit=dynamic_cast<const SiPixelRecHit*>(&(*hit)); 00122 if (pixelhit!=0) 00123 (*pixelInTrack.find(pixelhit->geographicalId().rawId())).push_back(pixelhit); 00124 00125 00126 00127 const SiStripMatchedRecHit2D* matchedhit=dynamic_cast<const SiStripMatchedRecHit2D*>(&(*hit)); 00128 if (matchedhit!=0) { 00129 (*stereoInTrack.find(matchedhit->stereoHit()->geographicalId().rawId())). 00130 push_back(matchedhit->stereoHit()); 00131 (*monoInTrack.find(matchedhit->monoHit()->geographicalId().rawId())). 00132 push_back(matchedhit->monoHit()); 00133 } 00134 00135 const SiStripRecHit2D* monohit=dynamic_cast<const SiStripRecHit2D*>(&(*hit)); 00136 if (monohit!=0) { 00137 //std::cout << "-- MONO HIT IN: " << monohit->geographicalId().rawId() << std::endl; 00138 uint32_t detid = (monohit->geographicalId().rawId()); 00139 DetSetVector<const SiStripRecHit2D*>::iterator point = monoInTrack.find(detid); 00140 if (point != monoInTrack.end()) { 00141 point->push_back(monohit); 00142 } else { 00143 point = stereoInTrack.find(detid); 00144 if (point == stereoInTrack.end()) { 00145 edm::LogError("RemainingClusterProducer") << "Don't know how to handle detid " << detid; 00146 continue; 00147 } else { 00148 point->push_back(monohit); 00149 } 00150 } 00151 } 00152 } 00153 } 00154 vector<const SiPixelRecHit*>::const_iterator ipix; 00155 vector<const SiStripRecHit2D*>::const_iterator istrip; 00156 00157 00158 //LOOP OVER ALL THE HITS 00159 for(TrackerGeometry::DetContainer::const_iterator it = 00160 pDD->dets().begin(); it != pDD->dets().end(); it++){ 00161 00162 DetId detid = ((*it)->geographicalId()); 00163 unsigned int isub=detid.subdetId(); 00164 00165 //PIXEL 00166 if ((isub== PixelSubdetector::PixelBarrel) || (isub== PixelSubdetector::PixelEndcap)) { 00167 00168 SiPixelRecHitCollection::range pixelrechitRange = (pixelHits.product())->get((detid)); 00169 SiPixelRecHitCollection::const_iterator pixeliter =pixelrechitRange.first; 00170 SiPixelRecHitCollection::const_iterator pixelrechitRangeIteratorEnd = pixelrechitRange.second; 00171 if ((pixelrechitRangeIteratorEnd-pixeliter)>0){ 00172 00173 DetSet<const SiPixelRecHit*> sipix= *(pixelInTrack.find(detid.rawId())); 00174 DetSet<SiPixelCluster> Pix_second(detid.rawId()); 00175 00176 for ( ; pixeliter != pixelrechitRangeIteratorEnd; ++pixeliter) { 00177 int itk=0; 00178 00179 for (ipix=sipix.begin();ipix!=sipix.end();ipix++){ 00180 if ((*ipix)->sharesInput(&(*pixeliter) , TrackingRecHit::all )) itk++; 00181 } 00182 00183 if (itk==0) Pix_second.push_back((*(*pixeliter).cluster())); 00184 } 00185 if(!Pix_second.empty()) thePixelClusterVector.push_back(Pix_second); 00186 } 00187 } //end pixel detid 00188 00189 00190 //STRIP 00191 if ((uint(detid.subdetId())==StripSubdetector::TIB) || 00192 (uint(detid.subdetId())==StripSubdetector::TOB) || 00193 (uint(detid.subdetId())==StripSubdetector::TID) || 00194 (uint(detid.subdetId())==StripSubdetector::TEC)){ 00195 00196 //rphi hits 00197 SiStripRecHit2DCollection::range monorechitRange = (rphirecHits.product())->get((detid)); 00198 SiStripRecHit2DCollection::const_iterator monoiter =monorechitRange.first; 00199 SiStripRecHit2DCollection::const_iterator monorechitRangeIteratorEnd = monorechitRange.second; 00200 if ((monorechitRangeIteratorEnd-monoiter)>0){ 00201 DetSet<SiStripCluster> Strip_second(detid.rawId()); 00202 DetSet<const SiStripRecHit2D*> simono= *(monoInTrack.find(detid.rawId())); 00203 for ( ; monoiter != monorechitRangeIteratorEnd; ++monoiter) { 00204 00205 int itk=0; 00206 for (istrip=simono.begin();istrip!=simono.end();istrip++){ 00207 00208 if ((*istrip)->sharesInput(&(*monoiter) , TrackingRecHit::all ))itk++; 00209 } 00210 if (itk==0) Strip_second.push_back((*(*monoiter).cluster())); 00211 } 00212 00213 if(!Strip_second.empty()) theStripClusterVector.push_back(Strip_second); 00214 } 00215 00216 00217 //stereo hits 00218 SiStripRecHit2DCollection::range stereorechitRange = (stereorecHits.product())->get((detid)); 00219 SiStripRecHit2DCollection::const_iterator stereoiter =stereorechitRange.first; 00220 SiStripRecHit2DCollection::const_iterator stereorechitRangeIteratorEnd =stereorechitRange.second; 00221 00222 if ((stereorechitRangeIteratorEnd-stereoiter)>0){ 00223 DetSet<SiStripCluster> Strip_second(detid.rawId()); 00224 DetSet<const SiStripRecHit2D*> sistereo= *(stereoInTrack.find(detid.rawId())); 00225 for ( ; stereoiter != stereorechitRangeIteratorEnd; ++stereoiter) { 00226 int itk=0; 00227 for (istrip=sistereo.begin();istrip!=sistereo.end();istrip++){ 00228 if ((*istrip)->sharesInput(&(*stereoiter) , TrackingRecHit::all ))itk++; 00229 } 00230 if (itk==0) Strip_second.push_back((*(*stereoiter).cluster())); 00231 } 00232 if(!Strip_second.empty()) theStripClusterVector.push_back(Strip_second); 00233 } 00234 } // end strips 00235 00236 00237 } //end loop over tracker detid 00238 00239 auto_ptr<DetSetVector<SiPixelCluster> > 00240 output_pixel(new DetSetVector<SiPixelCluster>(thePixelClusterVector) ); 00241 00242 auto_ptr<DetSetVector<SiStripCluster> > 00243 output_strip(new DetSetVector<SiStripCluster>(theStripClusterVector) ); 00244 00245 iEvent.put(output_pixel); 00246 iEvent.put(output_strip); 00247 } 00248 00249 // ------------ method called once each job just before starting event loop ------------ 00250 void 00251 RemainingClusterProducer::beginJob(const EventSetup& iSetup) 00252 { 00253 iSetup.get<TrackerDigiGeometryRecord> ().get (pDD); 00254 } 00255 00256 // ------------ method called once each job just after ending the event loop ------------ 00257 void 00258 RemainingClusterProducer::endJob() { 00259 } 00260 00261