00001
00002
00003
00004
00005
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022 #include <memory>
00023 #include <utility>
00024 #include <vector>
00025 #include <algorithm>
00026
00027
00028 #include "FWCore/Framework/interface/Frameworkfwd.h"
00029 #include "FWCore/Framework/interface/EDProducer.h"
00030 #include "FWCore/Framework/interface/Event.h"
00031 #include "FWCore/Framework/interface/MakerMacros.h"
00032 #include "FWCore/Framework/interface/ESHandle.h"
00033 #include "FWCore/Utilities/interface/InputTag.h"
00034 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00035
00036 #include "DataFormats/Common/interface/Ref.h"
00037 #include "DataFormats/Common/interface/EDProduct.h"
00038 #include "DataFormats/DetId/interface/DetId.h"
00039 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
00040 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00041 #include "DataFormats/TrackReco/interface/Track.h"
00042 #include "DataFormats/TrackReco/interface/TrackExtra.h"
00043 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00044 #include "DataFormats/SiStripDetId/interface/TECDetId.h"
00045 #include "DataFormats/SiStripDetId/interface/TIBDetId.h"
00046 #include "DataFormats/SiStripDetId/interface/TIDDetId.h"
00047 #include "DataFormats/SiStripDetId/interface/TOBDetId.h"
00048 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
00049 #include "DataFormats/SiStripCluster/interface/SiStripClusterCollection.h"
00050 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
00051 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2DCollection.h"
00052 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2DCollection.h"
00053 #include "DataFormats/Candidate/interface/Candidate.h"
00054 #include <DataFormats/SiStripCommon/interface/SiStripEventSummary.h>
00055 #include <DataFormats/SiStripCommon/interface/ConstantsForRunType.h>
00056 #include "DataFormats/SiStripDigi/interface/SiStripRawDigi.h"
00057 #include <DataFormats/SiStripCommon/interface/SiStripFedKey.h>
00058 #include <CondFormats/SiStripObjects/interface/FedChannelConnection.h>
00059 #include <CondFormats/SiStripObjects/interface/SiStripFedCabling.h>
00060 #include "CondFormats/DataRecord/interface/SiStripNoisesRcd.h"
00061 #include <CondFormats/SiStripObjects/interface/SiStripNoises.h>
00062 #include <CondFormats/DataRecord/interface/SiStripFedCablingRcd.h>
00063 #include "CondFormats/DataRecord/interface/SiStripNoisesRcd.h"
00064 #include <CondFormats/SiStripObjects/interface/SiStripNoises.h>
00065
00066
00067 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00068 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
00069 #include "DataFormats/GeometryVector/interface/LocalVector.h"
00070 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00071 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00072 #include "Geometry/CommonDetUnit/interface/GeomDetType.h"
00073 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00074 #include <Geometry/CommonTopologies/interface/Topology.h>
00075 #include <Geometry/CommonTopologies/interface/StripTopology.h>
00076
00077 #include <TrackingTools/PatternTools/interface/Trajectory.h>
00078
00079 #include "DQM/SiStripCommissioningSources/plugins/tracking/SiStripFineDelayHit.h"
00080 #include "DQM/SiStripCommissioningSources/plugins/tracking/SiStripFineDelayTLA.h"
00081 #include "DQM/SiStripCommissioningSources/plugins/tracking/SiStripFineDelayTOF.h"
00082
00083 #include "TMath.h"
00084
00085
00086
00087
00088 SiStripFineDelayHit::SiStripFineDelayHit(const edm::ParameterSet& iConfig):event_(0)
00089 {
00090
00091 produces<edm::DetSetVector<SiStripRawDigi> >("FineDelaySelection");
00092
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 }
00110
00111 SiStripFineDelayHit::~SiStripFineDelayHit()
00112 {
00113
00114
00115 delete anglefinder_;
00116 }
00117
00118
00119
00120
00121 std::pair<uint32_t, uint32_t> SiStripFineDelayHit::deviceMask(const StripSubdetector::SubDetector subdet,const int substructure)
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 }
00153
00154 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,const int substructure)
00155 {
00156 if(substructure==0xff) return detId(tracker,tk,trajVec,0,0);
00157
00158 std::pair<uint32_t, uint32_t> mask = deviceMask(subdet,substructure);
00159
00160 return detId(tracker,tk,trajVec,mask.first,mask.second);
00161 }
00162
00163 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)
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
00169 std::vector<std::pair< std::pair<DetId, LocalPoint> ,float> > hitangle;
00170 if(!cosmic_) {
00171 if(!trajInEvent_) {
00172
00173 hitangle = anglefinder_->findtrackangle(*tk);
00174 }
00175 else {
00176
00177
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
00194 hitangle = anglefinder_->findtrackangle((*(*seedcoll).begin()),*tk);
00195 }
00196 else {
00197
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
00204 for(iter=hitangle.begin();iter!=hitangle.end();iter++){
00205
00206
00207
00208
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
00215 LogDebug("DetId") << "check the angle: " << fabs((iter->second));
00216 if(1-fabs(fabs(iter->second)-1)<cos(maxAngle_/180.*TMath::Pi())) continue;
00217
00218 std::pair<uint32_t,std::pair<double, double> > el;
00219 std::pair<double, double> subel;
00220 el.first = iter->first.first.rawId();
00221
00222
00223
00224
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
00242 result.push_back(el);
00243 usedDetids.push_back(el.first);
00244 }
00245 return result;
00246 }
00247
00248 bool SiStripFineDelayHit::rechit(reco::Track* tk,uint32_t det_id)
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 }
00257
00258 std::pair<const SiStripCluster*,double> SiStripFineDelayHit::closestCluster(const TrackerGeometry& tracker,const reco::Track* tk,const uint32_t& det_id ,const edmNew::DetSetVector<SiStripCluster>& clusters, const edm::DetSetVector<SiStripDigi>& hits)
00259 {
00260 std::pair<const SiStripCluster*,double> result(NULL,0.);
00261 double hitStrip = -1;
00262 int nstrips = -1;
00263
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
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
00278
00279 else if((det_id - (*it)->geographicalId().rawId())==1) {
00280 const SiStripMatchedRecHit2D* hit2D = dynamic_cast<const SiStripMatchedRecHit2D*>(&(**it));
00281 if(!hit2D) continue;
00282 const SiStripRecHit2D* stereo = hit2D->stereoHit();
00283 if(!stereo) continue;
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
00294
00295 else if((det_id - (*it)->geographicalId().rawId())==2) {
00296 const SiStripMatchedRecHit2D* hit2D = dynamic_cast<const SiStripMatchedRecHit2D*>(&(**it));
00297 if(!hit2D) continue;
00298 const SiStripRecHit2D* mono = hit2D->monoHit();
00299 if(!mono) continue;
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
00314 for (edm::DetSetVector<SiStripDigi>::const_iterator DSViter=hits.begin(); DSViter!=hits.end();DSViter++){
00315 if(DSViter->id==det_id) {
00316
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
00330 LogDebug("closestCluster") << "build a fake cluster ";
00331 SiStripCluster* newCluster = new SiStripCluster(det_id,SiStripCluster::SiStripDigiRange(rangeStart,rangeStop));
00332 result.first = newCluster;
00333 result.second = fabs(newCluster->barycenter()-hitStrip);
00334 }
00335 break;
00336 }
00337 }
00338 } else {
00339
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
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 }
00356
00357
00358 void
00359 SiStripFineDelayHit::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00360 {
00361 using namespace edm;
00362
00363 edm::Handle<SiStripEventSummary> runsummary;
00364 iEvent.getByLabel( inputModuleLabel_, runsummary );
00365 if(runsummary->runType()==sistrip::APV_LATENCY) mode_ = 2;
00366 else if(runsummary->runType()==sistrip::FINE_DELAY) mode_ = 1;
00367 else {
00368 mode_ = 0;
00369 return;
00370 }
00371
00372 if(noTracking_) {
00373 produceNoTracking(iEvent,iSetup);
00374 return;
00375 }
00376 event_ = &iEvent;
00377
00378 std::vector< edm::DetSet<SiStripRawDigi> > output;
00379 output.reserve(100);
00380
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
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
00397 edm::Handle<edmNew::DetSetVector<SiStripCluster> > clusters;
00398 iEvent.getByLabel(clusterLabel_, clusters);
00399 const edmNew::DetSetVector<SiStripCluster>* clusterSet = clusters.product();
00400
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
00408 for(reco::TrackCollection::const_iterator itrack = tracks->begin(); itrack<tracks->end(); itrack++) {
00409
00410 if((itrack->px()*itrack->px()+itrack->py()*itrack->py()+itrack->pz()*itrack->pz())<minTrackP2_) continue;
00411
00412 std::vector< std::pair<uint32_t,std::pair<double,double> > > intersections;
00413 if(mode_==1) {
00414
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
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
00435 if(candidateCluster.second>maxClusterDistance_) continue;
00436 LogDebug("produce") << " The cluster is close enough."<< std::endl;
00437
00438
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
00450
00451 std::vector< edm::DetSet<SiStripRawDigi> >::iterator newdsit = output.begin();
00452 for(;newdsit!=output.end()&&newdsit->detId()!=connectionMap_[it->first];++newdsit) {}
00453
00454 if(newdsit==output.end()) {
00455 edm::DetSet<SiStripRawDigi> newds(connectionMap_[it->first]);
00456 output.push_back(newds);
00457 newdsit = output.end()-1;
00458 }
00459
00460 LogDebug("produce") << " New Hit... TOF:" << it->second.first << ", charge: " << int(leadingCharge)
00461 << " at " << int(leadingPosition) << "." << std::endl
00462 << "Angular correction: " << it->second.second
00463 << " giving a final value of " << int(leadingCharge*fabs(it->second.second))
00464 << " for fed key = " << connectionMap_[it->first] << " (detid=" << it->first << ")" ;
00465
00466 if(leadingCharge<255) {
00467
00468 leadingCharge = uint8_t(leadingCharge*fabs(it->second.second));
00469
00470 if((((it->first>>25)&0x7f)==0xd) ||
00471 ((((it->first>>25)&0x7f)==0xe) && (((it->first>>5)&0x7)>4)))
00472 leadingCharge = uint8_t((leadingCharge*0.64));
00473 }
00474
00475 unsigned int tof = abs(int(round(it->second.first*10)));
00476 tof = tof>255 ? 255 : tof;
00477 SiStripRawDigi newSiStrip(leadingCharge + (tof<<8));
00478 newdsit->push_back(newSiStrip);
00479 LogDebug("produce") << "New edm::DetSet<SiStripRawDigi> added.";
00480 }
00481 if(homeMadeClusters_) delete candidateCluster.first;
00482 }
00483 }
00484 }
00485
00486 LogDebug("produce") << "Putting " << output.size() << " new hits in the event.";
00487 std::auto_ptr< edm::DetSetVector<SiStripRawDigi> > formatedOutput(new edm::DetSetVector<SiStripRawDigi>(output) );
00488 iEvent.put(formatedOutput,"FineDelaySelection");
00489 }
00490
00491
00492 void
00493 SiStripFineDelayHit::produceNoTracking(edm::Event& iEvent, const edm::EventSetup& iSetup)
00494 {
00495 event_ = &iEvent;
00496
00497 std::vector< edm::DetSet<SiStripRawDigi> > output;
00498 output.reserve(100);
00499
00500 edm::Handle<SiStripEventSummary> summary;
00501 iEvent.getByLabel( inputModuleLabel_, summary );
00502 uint32_t layerCode = (const_cast<SiStripEventSummary*>(summary.product())->layerScanned())>>16;
00503 StripSubdetector::SubDetector subdet = StripSubdetector::TIB;
00504 if(((layerCode>>6)&0x3)==0) subdet = StripSubdetector::TIB;
00505 else if(((layerCode>>6)&0x3)==1) subdet = StripSubdetector::TOB;
00506 else if(((layerCode>>6)&0x3)==2) subdet = StripSubdetector::TID;
00507 else if(((layerCode>>6)&0x3)==3) subdet = StripSubdetector::TEC;
00508 int32_t layerIdx = (layerCode&0xF)*(((layerCode>>4)&0x3) ? -1 : 1);
00509 std::pair<uint32_t, uint32_t> mask = deviceMask(subdet,layerIdx);
00510
00511 edm::Handle<edmNew::DetSetVector<SiStripCluster> > clusters;
00512 iEvent.getByLabel(clusterLabel_,clusters);
00513 for (edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter=clusters->begin(); DSViter!=clusters->end();DSViter++ ) {
00514
00515 if(mode_==1 && ((DSViter->id() & mask.first) != mask.second) ) continue;
00516
00517 edmNew::DetSet<SiStripCluster>::const_iterator begin=DSViter->begin();
00518 edmNew::DetSet<SiStripCluster>::const_iterator end =DSViter->end();
00519 edm::DetSet<SiStripRawDigi> newds(connectionMap_[DSViter->id()]);
00520 for(edmNew::DetSet<SiStripCluster>::const_iterator iter=begin;iter!=end;++iter) {
00521
00522
00523 const std::vector< uint8_t >& amplitudes = iter->amplitudes();
00524 uint8_t leadingCharge = 0;
00525 uint8_t leadingStrip = iter->firstStrip();
00526 uint8_t leadingPosition = 0;
00527 for(std::vector< uint8_t >::const_iterator amplit = amplitudes.begin();amplit<amplitudes.end();amplit++,leadingStrip++) {
00528 if(leadingCharge<*amplit) {
00529 leadingCharge = *amplit;
00530 leadingPosition = leadingStrip;
00531 }
00532 }
00533
00534
00535
00536
00537 edm::ESHandle<SiStripNoises> noiseHandle_;
00538 iSetup.get<SiStripNoisesRcd>().get(noiseHandle_);
00539 SiStripNoises::Range detNoiseRange = noiseHandle_->getRange(DSViter->id());
00540 float noise=noiseHandle_->getNoise(leadingPosition, detNoiseRange);
00541 if( noise<1.5 ) continue;
00542 if( leadingCharge>=250 || noise>=8 || leadingCharge/noise>50 || leadingCharge/noise<10 ) continue;
00543
00544 if(leadingCharge<255) {
00545
00546 if((((((DSViter->id())>>25)&0x7f)==0xd)||((((DSViter->id())>>25)&0x7f)==0xe))&&((((DSViter->id())>>5)&0x7)>4))
00547 leadingCharge = uint8_t((leadingCharge*0.64));
00548 }
00549
00550 SiStripRawDigi newSiStrip(leadingCharge);
00551 newds.push_back(newSiStrip);
00552 }
00553
00554 output.push_back(newds);
00555 LogDebug("produce") << "New edm::DetSet<SiStripRawDigi> added with fedkey = "
00556 << std::hex << std::setfill('0') << std::setw(8)
00557 << connectionMap_[DSViter->id()] << std::dec;
00558 }
00559
00560 LogDebug("produce") << "Putting " << output.size() << " new hits in the event.";
00561 std::auto_ptr< edm::DetSetVector<SiStripRawDigi> > formatedOutput(new edm::DetSetVector<SiStripRawDigi>(output) );
00562 iEvent.put(formatedOutput,"FineDelaySelection");
00563 }
00564
00565
00566 void
00567 SiStripFineDelayHit::beginRun(edm::Run & run, const edm::EventSetup & iSetup)
00568 {
00569
00570 edm::ESHandle<SiStripFedCabling> cabling;
00571 iSetup.get<SiStripFedCablingRcd>().get( cabling );
00572 const std::vector< uint16_t > & feds = cabling->feds() ;
00573 for(std::vector< uint16_t >::const_iterator fedid = feds.begin();fedid<feds.end();++fedid) {
00574 const std::vector< FedChannelConnection > & connections = cabling->connections(*fedid);
00575 for(std::vector< FedChannelConnection >::const_iterator conn=connections.begin();conn<connections.end();++conn) {
00576
00577
00578
00579
00580
00581
00582
00583 connectionMap_[conn->detId()] = ( ( conn->fedId() & sistrip::invalid_ ) << 16 ) | ( conn->fedCh() & sistrip::invalid_ );
00584 }
00585 }
00586 }