#include <DQM/SiStripCommissioningSources/plugins/tracking/SiStripFineDelayHit.cc>
Public Member Functions | |
SiStripFineDelayHit (const edm::ParameterSet &) | |
virtual | ~SiStripFineDelayHit () |
Private Member Functions | |
virtual void | beginJob (const edm::EventSetup &) |
std::pair< const SiStripCluster *, double > | closestCluster (const TrackerGeometry &tracker, const reco::Track *tk, const uint32_t &detId, const edmNew::DetSetVector< SiStripCluster > &clusters, const edm::DetSetVector< SiStripDigi > &hits) |
std::vector< std::pair < uint32_t, std::pair< double, double > > > | detId (const TrackerGeometry &tracker, const reco::Track *tk, const std::vector< Trajectory > &trajVec, const uint32_t &maskDetId, const uint32_t &rootDetId) |
std::vector< std::pair < uint32_t, std::pair< double, double > > > | detId (const TrackerGeometry &tracker, const reco::Track *tk, const std::vector< Trajectory > &trajVec, const StripSubdetector::SubDetector subdet=StripSubdetector::TIB, const int substructure=0xff) |
std::pair< uint32_t, uint32_t > | deviceMask (const StripSubdetector::SubDetector subdet, const int substructure) |
virtual void | endJob () |
virtual void | produce (edm::Event &, const edm::EventSetup &) |
virtual void | produceNoTracking (edm::Event &, const edm::EventSetup &) |
bool | rechit (reco::Track *tk, uint32_t detId) |
Private Attributes | |
SiStripFineDelayTLA * | anglefinder_ |
edm::InputTag | clusterLabel_ |
std::map< uint32_t, uint32_t > | connectionMap_ |
bool | cosmic_ |
edm::InputTag | digiLabel_ |
const edm::Event * | event_ |
int | explorationWindow_ |
bool | field_ |
bool | homeMadeClusters_ |
edm::InputTag | inputModuleLabel_ |
double | maxAngle_ |
double | maxClusterDistance_ |
double | minTrackP2_ |
int | mode_ |
bool | noTracking_ |
edm::InputTag | seedLabel_ |
edm::InputTag | trackLabel_ |
bool | trajInEvent_ |
Implementation: <Notes on="" implementation>="">
Definition at line 32 of file SiStripFineDelayHit.h.
SiStripFineDelayHit::SiStripFineDelayHit | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 88 of file SiStripFineDelayHit.cc.
References anglefinder_, clusterLabel_, cosmic_, digiLabel_, explorationWindow_, field_, edm::ParameterSet::getParameter(), homeMadeClusters_, inputModuleLabel_, maxAngle_, maxClusterDistance_, minTrackP2_, mode_, noTracking_, seedLabel_, trackLabel_, and trajInEvent_.
00088 :event_(0) 00089 { 00090 //register your products 00091 produces<edm::DetSetVector<SiStripRawDigi> >("FineDelaySelection"); 00092 //now do what ever other initialization is needed 00093 anglefinder_=new SiStripFineDelayTLA(iConfig); 00094 cosmic_ = iConfig.getParameter<bool>("cosmic"); 00095 field_ = iConfig.getParameter<bool>("MagneticField"); 00096 trajInEvent_ = iConfig.getParameter<bool>("TrajInEvent"); 00097 maxAngle_ = iConfig.getParameter<double>("MaxTrackAngle"); 00098 minTrackP2_ = iConfig.getParameter<double>("MinTrackMomentum")*iConfig.getParameter<double>("MinTrackMomentum"); 00099 maxClusterDistance_ = iConfig.getParameter<double>("MaxClusterDistance"); 00100 clusterLabel_ = iConfig.getParameter<edm::InputTag>("ClustersLabel"); 00101 trackLabel_ = iConfig.getParameter<edm::InputTag>("TracksLabel"); 00102 seedLabel_ = iConfig.getParameter<edm::InputTag>("SeedsLabel"); 00103 inputModuleLabel_ = iConfig.getParameter<edm::InputTag>( "InputModuleLabel" ) ; 00104 digiLabel_ = iConfig.getParameter<edm::InputTag>("DigiLabel"); 00105 homeMadeClusters_ = iConfig.getParameter<bool>("NoClustering"); 00106 explorationWindow_ = iConfig.getParameter<uint32_t>("ExplorationWindow"); 00107 noTracking_ = iConfig.getParameter<bool>("NoTracking"); 00108 mode_=0; 00109 }
SiStripFineDelayHit::~SiStripFineDelayHit | ( | ) | [virtual] |
Definition at line 111 of file SiStripFineDelayHit.cc.
References anglefinder_.
00112 { 00113 // do anything here that needs to be done at desctruction time 00114 // (e.g. close files, deallocate resources etc.) 00115 delete anglefinder_; 00116 }
void SiStripFineDelayHit::beginJob | ( | const edm::EventSetup & | iSetup | ) | [private, virtual] |
Reimplemented from edm::EDProducer.
Definition at line 558 of file SiStripFineDelayHit.cc.
References connectionMap_, edm::EventSetup::get(), and sistrip::invalid_.
00559 { 00560 // Retrieve FED cabling object 00561 edm::ESHandle<SiStripFedCabling> cabling; 00562 iSetup.get<SiStripFedCablingRcd>().get( cabling ); 00563 const std::vector< uint16_t > & feds = cabling->feds() ; 00564 for(std::vector< uint16_t >::const_iterator fedid = feds.begin();fedid<feds.end();++fedid) { 00565 const std::vector< FedChannelConnection > & connections = cabling->connections(*fedid); 00566 for(std::vector< FedChannelConnection >::const_iterator conn=connections.begin();conn<connections.end();++conn) { 00567 /* 00568 SiStripFedKey key(conn->fedId(), 00569 SiStripFedKey::feUnit(conn->fedCh()), 00570 SiStripFedKey::feChan(conn->fedCh())); 00571 connectionMap_[conn->detId()] = key.key(); 00572 */ 00573 // the key is computed using an alternate formula for performance reasons. 00574 connectionMap_[conn->detId()] = ( ( conn->fedId() & sistrip::invalid_ ) << 16 ) | ( conn->fedCh() & sistrip::invalid_ ); 00575 } 00576 } 00577 }
std::pair< const SiStripCluster *, double > SiStripFineDelayHit::closestCluster | ( | const TrackerGeometry & | tracker, | |
const reco::Track * | tk, | |||
const uint32_t & | detId, | |||
const edmNew::DetSetVector< SiStripCluster > & | clusters, | |||
const edm::DetSetVector< SiStripDigi > & | hits | |||
) | [private] |
Definition at line 258 of file SiStripFineDelayHit.cc.
References SiStripCluster::barycenter(), edmNew::DetSetVector< T >::begin(), edm::DetSetVector< T >::begin(), begin, dist(), edmNew::DetSetVector< T >::end(), edm::DetSet< T >::end(), end, edm::DetSetVector< T >::end(), explorationWindow_, TrackingRecHit::geographicalId(), homeMadeClusters_, TrackerGeometry::idToDet(), int, TrackingRecHit::isValid(), it, iter, BaseSiTrackerRecHit2DLocalPos::localPosition(), LogDebug, lp, Topology::measurementPosition(), SiStripMatchedRecHit2D::monoHit(), NULL, p, reco::Track::recHitsBegin(), reco::Track::recHitsEnd(), HLT_VtxMuL3::result, SiStripMatchedRecHit2D::stereoHit(), GeomDetUnit::topology(), and PV2DBase< T, PVType, FrameType >::x().
Referenced by produce().
00259 { 00260 std::pair<const SiStripCluster*,double> result(NULL,0.); 00261 double hitStrip = -1; 00262 int nstrips = -1; 00263 // localize the crossing point of the track on the module 00264 for(trackingRecHit_iterator it = tk->recHitsBegin(); it != tk->recHitsEnd(); it++) { 00265 LogDebug("closestCluster") << "(*it)->geographicalId().rawId() vs det_id" << (*it)->geographicalId().rawId() << " " << det_id; 00266 //handle the mono rechits 00267 if((*it)->geographicalId().rawId() == det_id) { 00268 if(!(*it)->isValid()) continue; 00269 LogDebug("closestCluster") << " using the single mono hit"; 00270 LocalPoint lp = (*it)->localPosition(); 00271 const GeomDetUnit* gdu = static_cast<const GeomDetUnit*>(tracker.idToDet((*it)->geographicalId())); 00272 MeasurementPoint p = gdu->topology().measurementPosition(lp); 00273 hitStrip = p.x(); 00274 nstrips = (dynamic_cast<const StripTopology*>(&(gdu->topology())))->nstrips(); 00275 break; 00276 } 00277 //handle stereo part of matched hits 00278 //one could try to cast to SiStripMatchedRecHit2D but it is faster to look at the detid 00279 else if((det_id - (*it)->geographicalId().rawId())==1) { 00280 const SiStripMatchedRecHit2D* hit2D = dynamic_cast<const SiStripMatchedRecHit2D*>(&(**it)); 00281 if(!hit2D) continue; // this is a security that should never trigger 00282 const SiStripRecHit2D* stereo = hit2D->stereoHit(); 00283 if(!stereo) continue; // this is a security that should never trigger 00284 if(!stereo->isValid()) continue; 00285 LogDebug("closestCluster") << " using the stereo hit"; 00286 LocalPoint lp = stereo->localPosition(); 00287 const GeomDetUnit* gdu = static_cast<const GeomDetUnit*>(tracker.idToDet(stereo->geographicalId())); 00288 MeasurementPoint p = gdu->topology().measurementPosition(lp); 00289 hitStrip = p.x(); 00290 nstrips = (dynamic_cast<const StripTopology*>(&(gdu->topology())))->nstrips(); 00291 break; 00292 } 00293 //handle mono part of matched hits 00294 //one could try to cast to SiStripMatchedRecHit2D but it is faster to look at the detid 00295 else if((det_id - (*it)->geographicalId().rawId())==2) { 00296 const SiStripMatchedRecHit2D* hit2D = dynamic_cast<const SiStripMatchedRecHit2D*>(&(**it)); 00297 if(!hit2D) continue; // this is a security that should never trigger 00298 const SiStripRecHit2D* mono = hit2D->monoHit(); 00299 if(!mono) continue; // this is a security that should never trigger 00300 if(!mono->isValid()) continue; 00301 LogDebug("closestCluster") << " using the mono hit"; 00302 LocalPoint lp = mono->localPosition(); 00303 const GeomDetUnit* gdu = static_cast<const GeomDetUnit*>(tracker.idToDet(mono->geographicalId())); 00304 MeasurementPoint p = gdu->topology().measurementPosition(lp); 00305 hitStrip = p.x(); 00306 nstrips = (dynamic_cast<const StripTopology*>(&(gdu->topology())))->nstrips(); 00307 break; 00308 } 00309 } 00310 LogDebug("closestCluster") << " hit strip = " << hitStrip; 00311 if(hitStrip<0) return result; 00312 if(homeMadeClusters_) { 00313 // take the list of digis on the module 00314 for (edm::DetSetVector<SiStripDigi>::const_iterator DSViter=hits.begin(); DSViter!=hits.end();DSViter++){ 00315 if(DSViter->id==det_id) { 00316 // loop from hitstrip-n to hitstrip+n (explorationWindow_) and select the highest strip 00317 int minStrip = int(round(hitStrip))- explorationWindow_; 00318 minStrip = minStrip<0 ? 0 : minStrip; 00319 int maxStrip = int(round(hitStrip)) + explorationWindow_ + 1; 00320 maxStrip = maxStrip>=nstrips ? nstrips-1 : maxStrip; 00321 edm::DetSet<SiStripDigi>::const_iterator rangeStart = DSViter->end(); 00322 edm::DetSet<SiStripDigi>::const_iterator rangeStop = DSViter->end(); 00323 for(edm::DetSet<SiStripDigi>::const_iterator digiIt = DSViter->begin(); digiIt!=DSViter->end(); ++digiIt) { 00324 if(digiIt->strip()>=minStrip && rangeStart == DSViter->end()) rangeStart = digiIt; 00325 if(digiIt->strip()<=maxStrip) rangeStop = digiIt; 00326 } 00327 if(rangeStart != DSViter->end()) { 00328 if(rangeStop !=DSViter->end()) ++rangeStop; 00329 // build a fake cluster 00330 LogDebug("closestCluster") << "build a fake cluster "; 00331 SiStripCluster* newCluster = new SiStripCluster(det_id,SiStripCluster::SiStripDigiRange(rangeStart,rangeStop)); // /!\ ownership transfered 00332 result.first = newCluster; 00333 result.second = fabs(newCluster->barycenter()-hitStrip); 00334 } 00335 break; 00336 } 00337 } 00338 } else { 00339 // loop on the detsetvector<cluster> to find the right one 00340 for (edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter=clusters.begin(); DSViter!=clusters.end();DSViter++ ) 00341 if(DSViter->id()==det_id) { 00342 LogDebug("closestCluster") << " detset with the right detid. "; 00343 edmNew::DetSet<SiStripCluster>::const_iterator begin=DSViter->begin(); 00344 edmNew::DetSet<SiStripCluster>::const_iterator end =DSViter->end(); 00345 //find the cluster close to the hitStrip 00346 result.second = 1000.; 00347 for(edmNew::DetSet<SiStripCluster>::const_iterator iter=begin;iter!=end;++iter) { 00348 double dist = fabs(iter->barycenter()-hitStrip); 00349 if(dist<result.second) { result.second = dist; result.first = &(*iter); } 00350 } 00351 break; 00352 } 00353 } 00354 return result; 00355 }
std::vector< std::pair< uint32_t, std::pair< double, double > > > SiStripFineDelayHit::detId | ( | const TrackerGeometry & | tracker, | |
const reco::Track * | tk, | |||
const std::vector< Trajectory > & | trajVec, | |||
const uint32_t & | maskDetId, | |||
const uint32_t & | rootDetId | |||
) | [private] |
Definition at line 163 of file SiStripFineDelayHit.cc.
References anglefinder_, funct::cos(), cosmic_, lat::endl(), event_, field_, find(), SiStripFineDelayTLA::findtrackangle(), edm::Event::getByLabel(), i, TrackerGeometry::idToDetUnit(), iter, LogDebug, maxAngle_, reco::TrackBase::momentum(), reco::TrackBase::parameters(), Pi, reco::Track::recHitsEnd(), HLT_VtxMuL3::result, seedLabel_, SiStripFineDelayTOF::timeOfFlight(), SiStripFineDelayTOF::trackParameters(), and trajInEvent_.
00164 { 00165 bool onDisk = ((maskDetId==TIDDetId(3,15,0,0,0,0).rawId())||(maskDetId==TECDetId(3,15,0,0,0,0,0).rawId())) ; 00166 std::vector< std::pair<uint32_t,std::pair<double, double> > > result; 00167 std::vector<uint32_t> usedDetids; 00168 // now loop on recHits to find the right detId plus the track local angle 00169 std::vector<std::pair< std::pair<DetId, LocalPoint> ,float> > hitangle; 00170 if(!cosmic_) { 00171 if(!trajInEvent_) { 00172 //use the track. It will be refitted by the angleFinder 00173 hitangle = anglefinder_->findtrackangle(*tk); 00174 } 00175 else { 00176 // use trajectories in event. 00177 // we have first to find the right trajectory for the considered track. 00178 for(std::vector<Trajectory>::const_iterator traj = trajVec.begin(); traj< trajVec.end(); ++traj) { 00179 if( 00180 ((traj->lastMeasurement().recHit()->geographicalId().rawId() == (*(tk->recHitsEnd()-1))->geographicalId().rawId()) && 00181 ( traj->lastMeasurement().recHit()->localPosition().x() == (*(tk->recHitsEnd()-1))->localPosition().x()) ) || 00182 ((traj->firstMeasurement().recHit()->geographicalId().rawId() == (*(tk->recHitsEnd()-1))->geographicalId().rawId()) && 00183 ( traj->firstMeasurement().recHit()->localPosition().x() == (*(tk->recHitsEnd()-1))->localPosition().x()) ) ) { 00184 hitangle = anglefinder_->findtrackangle(*traj); 00185 break; 00186 } 00187 } 00188 } 00189 } else { 00190 edm::Handle<TrajectorySeedCollection> seedcoll; 00191 event_->getByLabel(seedLabel_,seedcoll); 00192 if(!trajInEvent_) { 00193 //use the track. It will be refitted by the angleFinder 00194 hitangle = anglefinder_->findtrackangle((*(*seedcoll).begin()),*tk); 00195 } 00196 else { 00197 // use trajectories in event. 00198 hitangle = anglefinder_->findtrackangle(trajVec); 00199 } 00200 } 00201 LogDebug("DetId") << "number of hits for the track: " << hitangle.size(); 00202 std::vector<std::pair< std::pair<DetId, LocalPoint> ,float> >::iterator iter; 00203 // select the interesting DetIds, based on the ID and TLA 00204 for(iter=hitangle.begin();iter!=hitangle.end();iter++){ 00205 // check the detId. 00206 // if substructure was 0xff, then maskDetId and rootDetId == 0 00207 // this implies all detids are accepted. (also if maskDetId=rootDetId=0 explicitely). 00208 // That "unusual" mode of operation allows to analyze also Latency scans 00209 LogDebug("DetId") << "check the detid: " << std::hex << (iter->first.first.rawId()) << " vs " << rootDetId 00210 << " with a mask of " << maskDetId << std::dec << std::endl; 00211 00212 if(((iter->first.first.rawId() & maskDetId) != rootDetId)) continue; 00213 if(std::find(usedDetids.begin(),usedDetids.end(),iter->first.first.rawId())!=usedDetids.end()) continue; 00214 // check the local angle (extended to the equivalent angle correction) 00215 LogDebug("DetId") << "check the angle: " << fabs((iter->second)); 00216 if(1-fabs(fabs(iter->second)-1)<cos(maxAngle_/180.*TMath::Pi())) continue; 00217 // returns the detid + the time of flight to there 00218 std::pair<uint32_t,std::pair<double, double> > el; 00219 std::pair<double, double> subel; 00220 el.first = iter->first.first.rawId(); 00221 // here, we compute the TOF. 00222 // For cosmics, some track parameters are missing. Parameters are recomputed. 00223 // for our calculation, the track momemtum at any point is enough: 00224 // only used without B field or for the sign of Pz. 00225 double trackParameters[5]; 00226 for(int i=0;i<5;i++) trackParameters[i] = tk->parameters()[i]; 00227 if(cosmic_) SiStripFineDelayTOF::trackParameters(*tk,trackParameters); 00228 double hit[3]; 00229 const GeomDetUnit* det(tracker.idToDetUnit(iter->first.first)); 00230 Surface::GlobalPoint gp = det->surface().toGlobal(iter->first.second); 00231 hit[0]=gp.x(); 00232 hit[1]=gp.y(); 00233 hit[2]=gp.z(); 00234 double phit[3]; 00235 phit[0] = tk->momentum().x(); 00236 phit[1] = tk->momentum().y(); 00237 phit[2] = tk->momentum().z(); 00238 subel.first = SiStripFineDelayTOF::timeOfFlight(cosmic_,field_,trackParameters,hit,phit,onDisk); 00239 subel.second = iter->second; 00240 el.second = subel; 00241 // returns the detid + TOF 00242 result.push_back(el); 00243 usedDetids.push_back(el.first); 00244 } 00245 return result; 00246 }
std::vector< std::pair< uint32_t, std::pair< double, double > > > SiStripFineDelayHit::detId | ( | const TrackerGeometry & | tracker, | |
const reco::Track * | tk, | |||
const std::vector< Trajectory > & | trajVec, | |||
const StripSubdetector::SubDetector | subdet = StripSubdetector::TIB , |
|||
const int | substructure = 0xff | |||
) | [private] |
Definition at line 154 of file SiStripFineDelayHit.cc.
References deviceMask().
Referenced by produce().
00155 { 00156 if(substructure==0xff) return detId(tracker,tk,trajVec,0,0); 00157 // first determine the root detId we are looking for 00158 std::pair<uint32_t, uint32_t> mask = deviceMask(subdet,substructure); 00159 // then call the method that loops on recHits 00160 return detId(tracker,tk,trajVec,mask.first,mask.second); 00161 }
std::pair< uint32_t, uint32_t > SiStripFineDelayHit::deviceMask | ( | const StripSubdetector::SubDetector | subdet, | |
const int | substructure | |||
) | [private] |
Definition at line 121 of file SiStripFineDelayHit.cc.
References funct::abs(), StripSubdetector::TEC, StripSubdetector::TIB, StripSubdetector::TID, and StripSubdetector::TOB.
Referenced by detId(), and produceNoTracking().
00122 { 00123 uint32_t rootDetId = 0; 00124 uint32_t maskDetId = 0; 00125 switch(subdet){ 00126 case (int)StripSubdetector::TIB : 00127 { 00128 rootDetId = TIBDetId(substructure,0,0,0,0,0).rawId(); 00129 maskDetId = TIBDetId(15,0,0,0,0,0).rawId(); 00130 break; 00131 } 00132 case (int)StripSubdetector::TID : 00133 { 00134 rootDetId = TIDDetId(substructure>0 ? 2 : 1,abs(substructure),0,0,0,0).rawId(); 00135 maskDetId = TIDDetId(3,15,0,0,0,0).rawId(); 00136 break; 00137 } 00138 case (int)StripSubdetector::TOB : 00139 { 00140 rootDetId = TOBDetId(substructure,0,0,0,0).rawId(); 00141 maskDetId = TOBDetId(15,0,0,0,0).rawId(); 00142 break; 00143 } 00144 case (int)StripSubdetector::TEC : 00145 { 00146 rootDetId = TECDetId(substructure>0 ? 2 : 1,abs(substructure),0,0,0,0,0).rawId(); 00147 maskDetId = TECDetId(3,15,0,0,0,0,0).rawId(); 00148 break; 00149 } 00150 } 00151 return std::make_pair(maskDetId,rootDetId); 00152 }
void SiStripFineDelayHit::produce | ( | edm::Event & | iEvent, | |
const edm::EventSetup & | iSetup | |||
) | [private, virtual] |
Implements edm::EDProducer.
Definition at line 359 of file SiStripFineDelayHit.cc.
References funct::abs(), anglefinder_, sistrip::APV_LATENCY, closestCluster(), clusterLabel_, connectionMap_, detId(), digiLabel_, lat::endl(), event_, sistrip::FINE_DELAY, edm::EventSetup::get(), edm::Event::getByLabel(), homeMadeClusters_, SiStripFineDelayTLA::init(), inputModuleLabel_, int, it, LogDebug, maxClusterDistance_, minTrackP2_, mode_, noTracking_, NULL, output(), produceNoTracking(), edm::Handle< T >::product(), edm::DetSet< T >::push_back(), edm::Event::put(), summary, StripSubdetector::TEC, StripSubdetector::TIB, StripSubdetector::TID, StripSubdetector::TOB, trackLabel_, tracks, and trajInEvent_.
00360 { 00361 using namespace edm; 00362 // Retrieve commissioning information from "event summary" 00363 edm::Handle<SiStripEventSummary> runsummary; 00364 iEvent.getByLabel( inputModuleLabel_, runsummary ); 00365 if(runsummary->runType()==sistrip::APV_LATENCY) mode_ = 2; // LatencyScan 00366 else if(runsummary->runType()==sistrip::FINE_DELAY) mode_ = 1; // DelayScan 00367 else { 00368 mode_ = 0; //unknown 00369 return; 00370 } 00371 00372 if(noTracking_) { 00373 produceNoTracking(iEvent,iSetup); 00374 return; 00375 } 00376 event_ = &iEvent; 00377 // container for the selected hits 00378 std::vector< edm::DetSet<SiStripRawDigi> > output; 00379 output.reserve(100); 00380 // access the tracks 00381 edm::Handle<reco::TrackCollection> trackCollection; 00382 iEvent.getByLabel(trackLabel_,trackCollection); 00383 const reco::TrackCollection *tracks=trackCollection.product(); 00384 edm::ESHandle<TrackerGeometry> tracker; 00385 iSetup.get<TrackerDigiGeometryRecord>().get(tracker); 00386 if (tracks->size()) { 00387 anglefinder_->init(iEvent,iSetup); 00388 LogDebug("produce") << "Found " << tracks->size() << " tracks."; 00389 // look at the hits if one needs them 00390 edm::Handle< edm::DetSetVector<SiStripDigi> > hits; 00391 const edm::DetSetVector<SiStripDigi>* hitSet = NULL; 00392 if(homeMadeClusters_) { 00393 iEvent.getByLabel(digiLabel_,hits); 00394 hitSet = hits.product(); 00395 } 00396 // look at the clusters 00397 edm::Handle<edmNew::DetSetVector<SiStripCluster> > clusters; 00398 iEvent.getByLabel(clusterLabel_, clusters); 00399 const edmNew::DetSetVector<SiStripCluster>* clusterSet = clusters.product(); 00400 // look at the trajectories if they are in the event 00401 std::vector<Trajectory> trajVec; 00402 if(trajInEvent_) { 00403 edm::Handle<std::vector<Trajectory> > TrajectoryCollection; 00404 iEvent.getByLabel(trackLabel_,TrajectoryCollection); 00405 trajVec = *(TrajectoryCollection.product()); 00406 } 00407 // loop on tracks 00408 for(reco::TrackCollection::const_iterator itrack = tracks->begin(); itrack<tracks->end(); itrack++) { 00409 // first check the track Pt 00410 if((itrack->px()*itrack->px()+itrack->py()*itrack->py()+itrack->pz()*itrack->pz())<minTrackP2_) continue; 00411 // check that we have something in the layer we are interested in 00412 std::vector< std::pair<uint32_t,std::pair<double,double> > > intersections; 00413 if(mode_==1) { 00414 // Retrieve and decode commissioning information from "event summary" 00415 edm::Handle<SiStripEventSummary> summary; 00416 iEvent.getByLabel( inputModuleLabel_, summary ); 00417 uint32_t layerCode = (const_cast<SiStripEventSummary*>(summary.product())->layerScanned())>>16; 00418 StripSubdetector::SubDetector subdet = StripSubdetector::TIB; 00419 if(((layerCode>>6)&0x3)==0) subdet = StripSubdetector::TIB; 00420 else if(((layerCode>>6)&0x3)==1) subdet = StripSubdetector::TOB; 00421 else if(((layerCode>>6)&0x3)==2) subdet = StripSubdetector::TID; 00422 else if(((layerCode>>6)&0x3)==3) subdet = StripSubdetector::TEC; 00423 int32_t layerIdx = (layerCode&0xF)*(((layerCode>>4)&0x3) ? -1 : 1); 00424 intersections = detId(*tracker,&(*itrack),trajVec,subdet,layerIdx); 00425 } else { 00426 // for latency scans, no layer is specified -> no cut on detid 00427 intersections = detId(*tracker,&(*itrack),trajVec); 00428 } 00429 LogDebug("produce") << " Found " << intersections.size() << " interesting intersections." << std::endl; 00430 for(std::vector< std::pair<uint32_t,std::pair<double,double> > >::iterator it = intersections.begin();it<intersections.end();it++) { 00431 std::pair<const SiStripCluster*,double> candidateCluster = closestCluster(*tracker,&(*itrack),it->first,*clusterSet,*hitSet); 00432 if(candidateCluster.first) { 00433 LogDebug("produce") << " Found a cluster."<< std::endl; 00434 // cut on the distance 00435 if(candidateCluster.second>maxClusterDistance_) continue; 00436 LogDebug("produce") << " The cluster is close enough."<< std::endl; 00437 // build the rawdigi corresponding to the leading strip and save it 00438 // here, only the leading strip is retained. All other rawdigis in the module are set to 0. 00439 const std::vector< uint8_t >& amplitudes = candidateCluster.first->amplitudes(); 00440 uint8_t leadingCharge = 0; 00441 uint8_t leadingStrip = candidateCluster.first->firstStrip(); 00442 uint8_t leadingPosition = 0; 00443 for(std::vector< uint8_t >::const_iterator amplit = amplitudes.begin();amplit<amplitudes.end();amplit++,leadingStrip++) { 00444 if(leadingCharge<*amplit) { 00445 leadingCharge = *amplit; 00446 leadingPosition = leadingStrip; 00447 } 00448 } 00449 edm::DetSet<SiStripRawDigi> newds(connectionMap_[it->first]); 00450 LogDebug("produce") << " New Hit... TOF:" << it->second.first << ", charge: " << int(leadingCharge) 00451 << " at " << int(leadingPosition) << "." << std::endl 00452 << "Angular correction: " << it->second.second 00453 << " giving a final value of " << int(leadingCharge*fabs(it->second.second)); 00454 // apply corrections to the leading charge, but only if it has not saturated. 00455 if(leadingCharge<255) { 00456 // correct the leading charge for the crossing angle 00457 leadingCharge = uint8_t(leadingCharge*fabs(it->second.second)); 00458 // correct for module thickness for TEC and TOB 00459 if((((it->first>>25)&0x7f)==0xd) || 00460 ((((it->first>>25)&0x7f)==0xe) && (((it->first>>5)&0x7)>4))) 00461 leadingCharge = uint8_t((leadingCharge*0.64)); 00462 } 00463 //code the time of flight in the digi 00464 unsigned int tof = abs(int(round(it->second.first*10))); 00465 tof = tof>255 ? 255 : tof; 00466 SiStripRawDigi newSiStrip(leadingCharge + (tof<<8)); 00467 newds.push_back(newSiStrip); 00468 //store into the detsetvector 00469 output.push_back(newds); 00470 LogDebug("produce") << "New edm::DetSet<SiStripRawDigi> added."; 00471 } 00472 if(homeMadeClusters_) delete candidateCluster.first; // we are owner of home-made clusters 00473 } 00474 } 00475 } 00476 // add the selected hits to the event. 00477 LogDebug("produce") << "Putting " << output.size() << " new hits in the event."; 00478 std::auto_ptr< edm::DetSetVector<SiStripRawDigi> > formatedOutput(new edm::DetSetVector<SiStripRawDigi>(output) ); 00479 iEvent.put(formatedOutput,"FineDelaySelection"); 00480 }
void SiStripFineDelayHit::produceNoTracking | ( | edm::Event & | iEvent, | |
const edm::EventSetup & | iSetup | |||
) | [private, virtual] |
Definition at line 484 of file SiStripFineDelayHit.cc.
References begin, clusterLabel_, connectionMap_, deviceMask(), end, event_, edm::EventSetup::get(), edm::Event::getByLabel(), inputModuleLabel_, iter, LogDebug, mode_, output(), edm::Handle< T >::product(), edm::Event::put(), summary, StripSubdetector::TEC, StripSubdetector::TIB, StripSubdetector::TID, and StripSubdetector::TOB.
Referenced by produce().
00485 { 00486 event_ = &iEvent; 00487 // container for the selected hits 00488 std::vector< edm::DetSet<SiStripRawDigi> > output; 00489 output.reserve(100); 00490 // Retrieve and decode commissioning information from "event summary" 00491 edm::Handle<SiStripEventSummary> summary; 00492 iEvent.getByLabel( inputModuleLabel_, summary ); 00493 uint32_t layerCode = (const_cast<SiStripEventSummary*>(summary.product())->layerScanned())>>16; 00494 StripSubdetector::SubDetector subdet = StripSubdetector::TIB; 00495 if(((layerCode>>6)&0x3)==0) subdet = StripSubdetector::TIB; 00496 else if(((layerCode>>6)&0x3)==1) subdet = StripSubdetector::TOB; 00497 else if(((layerCode>>6)&0x3)==2) subdet = StripSubdetector::TID; 00498 else if(((layerCode>>6)&0x3)==3) subdet = StripSubdetector::TEC; 00499 int32_t layerIdx = (layerCode&0xF)*(((layerCode>>4)&0x3) ? -1 : 1); 00500 std::pair<uint32_t, uint32_t> mask = deviceMask(subdet,layerIdx); 00501 // look at the clusters 00502 edm::Handle<edmNew::DetSetVector<SiStripCluster> > clusters; 00503 iEvent.getByLabel(clusterLabel_,clusters); 00504 for (edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter=clusters->begin(); DSViter!=clusters->end();DSViter++ ) { 00505 // check that we are in the layer of interest 00506 if(mode_==1 && ((DSViter->id() & mask.first) != mask.second) ) continue; 00507 // iterate over clusters 00508 edmNew::DetSet<SiStripCluster>::const_iterator begin=DSViter->begin(); 00509 edmNew::DetSet<SiStripCluster>::const_iterator end =DSViter->end(); 00510 edm::DetSet<SiStripRawDigi> newds(connectionMap_[DSViter->id()]); 00511 for(edmNew::DetSet<SiStripCluster>::const_iterator iter=begin;iter!=end;++iter) { 00512 // build the rawdigi corresponding to the leading strip and save it 00513 // here, only the leading strip is retained. All other rawdigis in the module are set to 0. 00514 const std::vector< uint8_t >& amplitudes = iter->amplitudes(); 00515 uint8_t leadingCharge = 0; 00516 uint8_t leadingStrip = iter->firstStrip(); 00517 uint8_t leadingPosition = 0; 00518 for(std::vector< uint8_t >::const_iterator amplit = amplitudes.begin();amplit<amplitudes.end();amplit++,leadingStrip++) { 00519 if(leadingCharge<*amplit) { 00520 leadingCharge = *amplit; 00521 leadingPosition = leadingStrip; 00522 } 00523 } 00524 // apply some sanity cuts. This is needed since we don't use tracking to clean clusters 00525 // 1.5< noise <8 00526 // charge<250 00527 // 50 > s/n > 10 00528 edm::ESHandle<SiStripNoises> noiseHandle_; 00529 iSetup.get<SiStripNoisesRcd>().get(noiseHandle_); 00530 SiStripNoises::Range detNoiseRange = noiseHandle_->getRange(DSViter->id()); 00531 float noise=noiseHandle_->getNoise(leadingPosition, detNoiseRange); 00532 if( noise<1.5 ) continue; 00533 if( leadingCharge>=250 || noise>=8 || leadingCharge/noise>50 || leadingCharge/noise<10 ) continue; 00534 // apply some correction to the leading charge, but only if it has not saturated. 00535 if(leadingCharge<255) { 00536 // correct for modulethickness for TEC and TOB 00537 if((((((DSViter->id())>>25)&0x7f)==0xd)||((((DSViter->id())>>25)&0x7f)==0xe))&&((((DSViter->id())>>5)&0x7)>4)) 00538 leadingCharge = uint8_t((leadingCharge*0.64)); 00539 } 00540 //code the time of flight == 0 in the digi 00541 SiStripRawDigi newSiStrip(leadingCharge); 00542 newds.push_back(newSiStrip); 00543 } 00544 //store into the detsetvector 00545 output.push_back(newds); 00546 LogDebug("produce") << "New edm::DetSet<SiStripRawDigi> added with fedkey = " 00547 << std::hex << std::setfill('0') << std::setw(8) 00548 << connectionMap_[DSViter->id()] << std::dec; 00549 } 00550 // add the selected hits to the event. 00551 LogDebug("produce") << "Putting " << output.size() << " new hits in the event."; 00552 std::auto_ptr< edm::DetSetVector<SiStripRawDigi> > formatedOutput(new edm::DetSetVector<SiStripRawDigi>(output) ); 00553 iEvent.put(formatedOutput,"FineDelaySelection"); 00554 }
bool SiStripFineDelayHit::rechit | ( | reco::Track * | tk, | |
uint32_t | detId | |||
) | [private] |
Definition at line 248 of file SiStripFineDelayHit.cc.
References it, reco::Track::recHitsBegin(), and reco::Track::recHitsEnd().
00249 { 00250 for(trackingRecHit_iterator it = tk->recHitsBegin(); it != tk->recHitsEnd(); it++) 00251 if((*it)->geographicalId().rawId() == det_id) { 00252 return (*it)->isValid(); 00253 break; 00254 } 00255 return false; 00256 }
Definition at line 53 of file SiStripFineDelayHit.h.
Referenced by detId(), produce(), SiStripFineDelayHit(), and ~SiStripFineDelayHit().
Definition at line 59 of file SiStripFineDelayHit.h.
Referenced by produce(), produceNoTracking(), and SiStripFineDelayHit().
std::map<uint32_t,uint32_t> SiStripFineDelayHit::connectionMap_ [private] |
Definition at line 60 of file SiStripFineDelayHit.h.
Referenced by beginJob(), produce(), and produceNoTracking().
bool SiStripFineDelayHit::cosmic_ [private] |
Definition at line 55 of file SiStripFineDelayHit.h.
Referenced by detId(), and SiStripFineDelayHit().
edm::InputTag SiStripFineDelayHit::digiLabel_ [private] |
Definition at line 59 of file SiStripFineDelayHit.h.
Referenced by produce(), and SiStripFineDelayHit().
const edm::Event* SiStripFineDelayHit::event_ [private] |
Definition at line 54 of file SiStripFineDelayHit.h.
Referenced by detId(), produce(), and produceNoTracking().
int SiStripFineDelayHit::explorationWindow_ [private] |
Definition at line 58 of file SiStripFineDelayHit.h.
Referenced by closestCluster(), and SiStripFineDelayHit().
bool SiStripFineDelayHit::field_ [private] |
Definition at line 55 of file SiStripFineDelayHit.h.
Referenced by detId(), and SiStripFineDelayHit().
bool SiStripFineDelayHit::homeMadeClusters_ [private] |
Definition at line 55 of file SiStripFineDelayHit.h.
Referenced by closestCluster(), produce(), and SiStripFineDelayHit().
Definition at line 59 of file SiStripFineDelayHit.h.
Referenced by produce(), produceNoTracking(), and SiStripFineDelayHit().
double SiStripFineDelayHit::maxAngle_ [private] |
Definition at line 56 of file SiStripFineDelayHit.h.
Referenced by detId(), and SiStripFineDelayHit().
double SiStripFineDelayHit::maxClusterDistance_ [private] |
Definition at line 56 of file SiStripFineDelayHit.h.
Referenced by produce(), and SiStripFineDelayHit().
double SiStripFineDelayHit::minTrackP2_ [private] |
Definition at line 56 of file SiStripFineDelayHit.h.
Referenced by produce(), and SiStripFineDelayHit().
int SiStripFineDelayHit::mode_ [private] |
Definition at line 57 of file SiStripFineDelayHit.h.
Referenced by produce(), produceNoTracking(), and SiStripFineDelayHit().
bool SiStripFineDelayHit::noTracking_ [private] |
Definition at line 55 of file SiStripFineDelayHit.h.
Referenced by produce(), and SiStripFineDelayHit().
edm::InputTag SiStripFineDelayHit::seedLabel_ [private] |
Definition at line 59 of file SiStripFineDelayHit.h.
Referenced by detId(), and SiStripFineDelayHit().
Definition at line 59 of file SiStripFineDelayHit.h.
Referenced by produce(), and SiStripFineDelayHit().
bool SiStripFineDelayHit::trajInEvent_ [private] |
Definition at line 55 of file SiStripFineDelayHit.h.
Referenced by detId(), produce(), and SiStripFineDelayHit().