CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/DQM/SiStripCommissioningSources/plugins/tracking/SiStripFineDelayHit.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    SiStripFineDelayHit
00004 // Class:      SiStripFineDelayHit
00005 // 
00013 //
00014 // Original Author:  Christophe DELAERE
00015 //         Created:  Fri Nov 17 10:52:42 CET 2006
00016 // $Id: SiStripFineDelayHit.cc,v 1.17 2012/01/17 10:04:37 innocent Exp $
00017 //
00018 //
00019 
00020 
00021 // system include files
00022 #include <memory>
00023 #include <utility>
00024 #include <vector>
00025 #include <algorithm>
00026 
00027 // user include files
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/DetId/interface/DetId.h"
00038 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
00039 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00040 #include "DataFormats/TrackReco/interface/Track.h"
00041 #include "DataFormats/TrackReco/interface/TrackExtra.h"
00042 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00043 #include "DataFormats/SiStripDetId/interface/TECDetId.h"
00044 #include "DataFormats/SiStripDetId/interface/TIBDetId.h"
00045 #include "DataFormats/SiStripDetId/interface/TIDDetId.h"
00046 #include "DataFormats/SiStripDetId/interface/TOBDetId.h"
00047 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
00048 #include "DataFormats/SiStripCluster/interface/SiStripClusterCollection.h"
00049 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
00050 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2DCollection.h"
00051 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2DCollection.h"
00052 #include "DataFormats/Candidate/interface/Candidate.h"
00053 #include <DataFormats/SiStripCommon/interface/SiStripEventSummary.h>
00054 #include <DataFormats/SiStripCommon/interface/ConstantsForRunType.h>
00055 #include "DataFormats/SiStripDigi/interface/SiStripRawDigi.h"
00056 #include <DataFormats/SiStripCommon/interface/SiStripFedKey.h>
00057 #include <CondFormats/SiStripObjects/interface/FedChannelConnection.h>
00058 #include <CondFormats/SiStripObjects/interface/SiStripFedCabling.h>
00059 #include "CondFormats/DataRecord/interface/SiStripNoisesRcd.h"
00060 #include <CondFormats/SiStripObjects/interface/SiStripNoises.h>
00061 #include <CondFormats/DataRecord/interface/SiStripFedCablingRcd.h>
00062 #include "CondFormats/DataRecord/interface/SiStripNoisesRcd.h"
00063 #include <CondFormats/SiStripObjects/interface/SiStripNoises.h>
00064 
00065 
00066 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00067 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
00068 #include "DataFormats/GeometryVector/interface/LocalVector.h"
00069 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00070 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00071 #include "Geometry/CommonDetUnit/interface/GeomDetType.h"
00072 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00073 #include <Geometry/CommonTopologies/interface/Topology.h>
00074 #include <Geometry/CommonTopologies/interface/StripTopology.h>
00075 
00076 #include <TrackingTools/PatternTools/interface/Trajectory.h>
00077 
00078 #include "DQM/SiStripCommissioningSources/plugins/tracking/SiStripFineDelayHit.h"
00079 #include "DQM/SiStripCommissioningSources/plugins/tracking/SiStripFineDelayTLA.h"
00080 #include "DQM/SiStripCommissioningSources/plugins/tracking/SiStripFineDelayTOF.h"
00081 
00082 #include "TMath.h"
00083 
00084 //
00085 // constructors and destructor
00086 //
00087 SiStripFineDelayHit::SiStripFineDelayHit(const edm::ParameterSet& iConfig):event_(0)
00088 {
00089    //register your products
00090    produces<edm::DetSetVector<SiStripRawDigi> >("FineDelaySelection");
00091    //now do what ever other initialization is needed
00092    anglefinder_=new SiStripFineDelayTLA(iConfig);
00093    cosmic_ = iConfig.getParameter<bool>("cosmic");
00094    field_ = iConfig.getParameter<bool>("MagneticField");
00095    maxAngle_ = iConfig.getParameter<double>("MaxTrackAngle");
00096    minTrackP2_ = iConfig.getParameter<double>("MinTrackMomentum")*iConfig.getParameter<double>("MinTrackMomentum");
00097    maxClusterDistance_ = iConfig.getParameter<double>("MaxClusterDistance");
00098    clusterLabel_ = iConfig.getParameter<edm::InputTag>("ClustersLabel");
00099    trackLabel_ = iConfig.getParameter<edm::InputTag>("TracksLabel");
00100    seedLabel_  = iConfig.getParameter<edm::InputTag>("SeedsLabel");
00101    inputModuleLabel_ = iConfig.getParameter<edm::InputTag>( "InputModuleLabel" ) ;
00102    digiLabel_ = iConfig.getParameter<edm::InputTag>("DigiLabel");
00103    homeMadeClusters_ = iConfig.getParameter<bool>("NoClustering");
00104    explorationWindow_ = iConfig.getParameter<uint32_t>("ExplorationWindow");
00105    noTracking_ = iConfig.getParameter<bool>("NoTracking");
00106    mode_=0;
00107 }
00108 
00109 SiStripFineDelayHit::~SiStripFineDelayHit()
00110 {
00111    // do anything here that needs to be done at desctruction time
00112    // (e.g. close files, deallocate resources etc.)
00113    delete anglefinder_;
00114 }
00115 
00116 //
00117 // member functions
00118 //
00119 std::pair<uint32_t, uint32_t> SiStripFineDelayHit::deviceMask(const StripSubdetector::SubDetector subdet,const int substructure)
00120 {
00121   uint32_t rootDetId = 0;
00122   uint32_t maskDetId = 0;
00123   switch(subdet){
00124     case (int)StripSubdetector::TIB :
00125     {
00126       rootDetId = TIBDetId(substructure,0,0,0,0,0).rawId();
00127       maskDetId = TIBDetId(15,0,0,0,0,0).rawId();
00128       break;
00129     }
00130     case (int)StripSubdetector::TID :
00131     {
00132       rootDetId = TIDDetId(substructure>0 ? 2 : 1,abs(substructure),0,0,0,0).rawId();
00133       maskDetId = TIDDetId(3,15,0,0,0,0).rawId();
00134       break;
00135     }
00136     case (int)StripSubdetector::TOB :
00137     {
00138       rootDetId = TOBDetId(substructure,0,0,0,0).rawId();
00139       maskDetId = TOBDetId(15,0,0,0,0).rawId();
00140       break;
00141     }
00142     case (int)StripSubdetector::TEC :
00143     {
00144       rootDetId = TECDetId(substructure>0 ? 2 : 1,abs(substructure),0,0,0,0,0).rawId();
00145       maskDetId = TECDetId(3,15,0,0,0,0,0).rawId();
00146       break;
00147     }
00148   }
00149   return std::make_pair(maskDetId,rootDetId);
00150 }
00151 
00152 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)
00153 {
00154   if(substructure==0xff) return detId(tracker,tk,trajVec,0,0);
00155   // first determine the root detId we are looking for
00156   std::pair<uint32_t, uint32_t> mask = deviceMask(subdet,substructure);
00157   // then call the method that loops on recHits
00158   return detId(tracker,tk,trajVec,mask.first,mask.second);
00159 }
00160 
00161 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)
00162 {
00163   bool onDisk = ((maskDetId==TIDDetId(3,15,0,0,0,0).rawId())||(maskDetId==TECDetId(3,15,0,0,0,0,0).rawId())) ;
00164   std::vector< std::pair<uint32_t,std::pair<double, double> > > result;
00165   std::vector<uint32_t> usedDetids;
00166   // now loop on recHits to find the right detId plus the track local angle
00167   std::vector<std::pair< std::pair<DetId, LocalPoint> ,float> > hitangle;
00168   if(!cosmic_) {
00169     // use trajectories in event.
00170     // we have first to find the right trajectory for the considered track.
00171     for(std::vector<Trajectory>::const_iterator traj = trajVec.begin(); traj< trajVec.end(); ++traj) {
00172       if(
00173         ((traj->lastMeasurement().recHit()->geographicalId().rawId() == (*(tk->recHitsEnd()-1))->geographicalId().rawId()) &&
00174         ( traj->lastMeasurement().recHit()->localPosition().x() == (*(tk->recHitsEnd()-1))->localPosition().x())               ) ||
00175         ((traj->firstMeasurement().recHit()->geographicalId().rawId() == (*(tk->recHitsEnd()-1))->geographicalId().rawId()) &&
00176         ( traj->firstMeasurement().recHit()->localPosition().x() == (*(tk->recHitsEnd()-1))->localPosition().x())              )   ) {
00177           hitangle = anglefinder_->findtrackangle(*traj);
00178           break;
00179       }
00180     }
00181   } else {
00182     edm::Handle<TrajectorySeedCollection> seedcoll;
00183     event_->getByLabel(seedLabel_,seedcoll);
00184     // use trajectories in event.
00185     hitangle = anglefinder_->findtrackangle(trajVec);
00186   }
00187   LogDebug("DetId") << "number of hits for the track: " << hitangle.size();
00188   std::vector<std::pair< std::pair<DetId, LocalPoint> ,float> >::iterator iter;
00189   // select the interesting DetIds, based on the ID and TLA
00190   for(iter=hitangle.begin();iter!=hitangle.end();iter++){
00191     // check the detId.
00192     // if substructure was 0xff, then maskDetId and rootDetId == 0 
00193     // this implies all detids are accepted. (also if maskDetId=rootDetId=0 explicitely).
00194     // That "unusual" mode of operation allows to analyze also Latency scans
00195     LogDebug("DetId") << "check the detid: " << std::hex << (iter->first.first.rawId()) << " vs " << rootDetId
00196                       << " with a mask of "  << maskDetId << std::dec << std::endl;
00197 
00198     if(((iter->first.first.rawId() & maskDetId) != rootDetId)) continue;
00199     if(std::find(usedDetids.begin(),usedDetids.end(),iter->first.first.rawId())!=usedDetids.end()) continue;
00200     // check the local angle (extended to the equivalent angle correction)
00201     LogDebug("DetId") << "check the angle: " << fabs((iter->second));
00202     if(1-fabs(fabs(iter->second)-1)<cos(maxAngle_/180.*TMath::Pi())) continue;
00203     // returns the detid + the time of flight to there
00204     std::pair<uint32_t,std::pair<double, double> > el;
00205     std::pair<double, double> subel;
00206     el.first = iter->first.first.rawId();
00207     // here, we compute the TOF.
00208     // For cosmics, some track parameters are missing. Parameters are recomputed.
00209     // for our calculation, the track momemtum at any point is enough:
00210     // only used without B field or for the sign of Pz.
00211     double trackParameters[5];
00212     for(int i=0;i<5;i++) trackParameters[i] = tk->parameters()[i];
00213     if(cosmic_) SiStripFineDelayTOF::trackParameters(*tk,trackParameters);
00214     double hit[3];
00215     const GeomDetUnit* det(tracker.idToDetUnit(iter->first.first));
00216     Surface::GlobalPoint gp = det->surface().toGlobal(iter->first.second);
00217     hit[0]=gp.x();
00218     hit[1]=gp.y();
00219     hit[2]=gp.z();
00220     double phit[3];
00221     phit[0] = tk->momentum().x();
00222     phit[1] = tk->momentum().y();
00223     phit[2] = tk->momentum().z();
00224     subel.first = SiStripFineDelayTOF::timeOfFlight(cosmic_,field_,trackParameters,hit,phit,onDisk);
00225     subel.second = iter->second;
00226     el.second = subel;
00227     // returns the detid + TOF
00228     result.push_back(el);
00229     usedDetids.push_back(el.first);
00230   }
00231   return result;
00232 }
00233 
00234 bool SiStripFineDelayHit::rechit(reco::Track* tk,uint32_t det_id)
00235 {
00236   for(trackingRecHit_iterator it = tk->recHitsBegin(); it != tk->recHitsEnd(); it++) 
00237     if((*it)->geographicalId().rawId() == det_id) {
00238       return (*it)->isValid();
00239       break;
00240     }
00241   return false;
00242 }
00243 
00244 // VI January 2012: FIXME
00245 // do not understand what is going on here: each hit has a cluster: by definition will be the closest!
00246 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)
00247 {
00248   std::pair<const SiStripCluster*,double> result(NULL,0.);
00249   double hitStrip = -1;
00250   int nstrips = -1;
00251   // localize the crossing point of the track on the module
00252   for(trackingRecHit_iterator it = tk->recHitsBegin(); it != tk->recHitsEnd(); it++) {
00253     LogDebug("closestCluster") << "(*it)->geographicalId().rawId() vs det_id" << (*it)->geographicalId().rawId() << " " <<  det_id;
00254     //handle the mono rechits
00255     if((*it)->geographicalId().rawId() == det_id) {
00256       if(!(*it)->isValid()) continue;
00257       LogDebug("closestCluster") << " using the single mono hit";
00258       LocalPoint lp = (*it)->localPosition();
00259       const GeomDetUnit* gdu = static_cast<const GeomDetUnit*>(tracker.idToDet((*it)->geographicalId()));
00260       MeasurementPoint p = gdu->topology().measurementPosition(lp);
00261       hitStrip = p.x();
00262       nstrips = (dynamic_cast<const StripTopology*>(&(gdu->topology())))->nstrips();
00263       break;
00264     }
00265     /* FIXME: local position is not there anymore...
00266     //handle stereo part of matched hits
00267     //one could try to cast to SiStripMatchedRecHit2D but it is faster to look at the detid
00268     else if((det_id - (*it)->geographicalId().rawId())==1) {
00269       const SiStripMatchedRecHit2D* hit2D = dynamic_cast<const SiStripMatchedRecHit2D*>(&(**it));
00270       if(!hit2D) continue; // this is a security that should never trigger
00271       const SiStripRecHit2D* stereo = hit2D->stereoHit();
00272       if(!stereo) continue; // this is a security that should never trigger
00273       if(!stereo->isValid()) continue;
00274       LogDebug("closestCluster") << " using the stereo hit";
00275       LocalPoint lp = stereo->localPosition();
00276       const GeomDetUnit* gdu = static_cast<const GeomDetUnit*>(tracker.idToDet(stereo->geographicalId()));
00277       MeasurementPoint p = gdu->topology().measurementPosition(lp);
00278       hitStrip = p.x();
00279       nstrips = (dynamic_cast<const StripTopology*>(&(gdu->topology())))->nstrips();
00280       break;
00281     }
00282     //handle mono part of matched hits
00283     //one could try to cast to SiStripMatchedRecHit2D but it is faster to look at the detid
00284     else if((det_id - (*it)->geographicalId().rawId())==2) {
00285       const SiStripMatchedRecHit2D* hit2D = dynamic_cast<const SiStripMatchedRecHit2D*>(&(**it));
00286       if(!hit2D) continue; // this is a security that should never trigger
00287       const SiStripRecHit2D* mono = hit2D->monoHit();
00288       if(!mono) continue; // this is a security that should never trigger
00289       if(!mono->isValid()) continue;
00290       LogDebug("closestCluster") << " using the mono hit";
00291       LocalPoint lp = mono->localPosition();
00292       const GeomDetUnit* gdu = static_cast<const GeomDetUnit*>(tracker.idToDet(mono->geographicalId()));
00293       MeasurementPoint p = gdu->topology().measurementPosition(lp);
00294       hitStrip = p.x();
00295       nstrips = (dynamic_cast<const StripTopology*>(&(gdu->topology())))->nstrips();
00296       break;
00297     }
00298     */
00299   }
00300   LogDebug("closestCluster") << " hit strip = " << hitStrip;
00301   if(hitStrip<0) return result;
00302   if(homeMadeClusters_) {
00303     // take the list of digis on the module
00304     for (edm::DetSetVector<SiStripDigi>::const_iterator DSViter=hits.begin(); DSViter!=hits.end();DSViter++){
00305       if(DSViter->id==det_id)  {
00306         // loop from hitstrip-n to hitstrip+n (explorationWindow_) and select the highest strip
00307         int minStrip = int(round(hitStrip))- explorationWindow_;
00308         minStrip = minStrip<0 ? 0 : minStrip;
00309         int maxStrip = int(round(hitStrip)) + explorationWindow_ + 1;
00310         maxStrip = maxStrip>=nstrips ? nstrips-1 : maxStrip;
00311         edm::DetSet<SiStripDigi>::const_iterator rangeStart = DSViter->end();
00312         edm::DetSet<SiStripDigi>::const_iterator rangeStop  = DSViter->end();
00313         for(edm::DetSet<SiStripDigi>::const_iterator digiIt = DSViter->begin(); digiIt!=DSViter->end(); ++digiIt) {
00314           if(digiIt->strip()>=minStrip && rangeStart == DSViter->end()) rangeStart = digiIt;
00315           if(digiIt->strip()<=maxStrip) rangeStop = digiIt;
00316         }
00317         if(rangeStart != DSViter->end()) {
00318           if(rangeStop !=DSViter->end()) ++rangeStop;
00319           // build a fake cluster 
00320           LogDebug("closestCluster") << "build a fake cluster ";
00321           SiStripCluster* newCluster = new SiStripCluster(det_id,SiStripCluster::SiStripDigiRange(rangeStart,rangeStop)); // /!\ ownership transfered
00322           result.first = newCluster;
00323           result.second = fabs(newCluster->barycenter()-hitStrip);
00324         }
00325         break;
00326       }
00327     }
00328   } else {
00329   // loop on the detsetvector<cluster> to find the right one
00330    for (edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter=clusters.begin(); DSViter!=clusters.end();DSViter++ ) 
00331      if(DSViter->id()==det_id)  {
00332         LogDebug("closestCluster") << " detset with the right detid. ";
00333         edmNew::DetSet<SiStripCluster>::const_iterator begin=DSViter->begin();
00334         edmNew::DetSet<SiStripCluster>::const_iterator end  =DSViter->end();
00335         //find the cluster close to the hitStrip
00336         result.second = 1000.;
00337         for(edmNew::DetSet<SiStripCluster>::const_iterator iter=begin;iter!=end;++iter) {
00338           double dist = fabs(iter->barycenter()-hitStrip);
00339           if(dist<result.second) { result.second = dist; result.first = &(*iter); }
00340         }
00341         break;
00342      }
00343   }
00344   return result;
00345 }
00346 
00347 // ------------ method called to produce the data  ------------
00348 void
00349 SiStripFineDelayHit::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00350 {
00351    using namespace edm;
00352    // Retrieve commissioning information from "event summary"
00353    edm::Handle<SiStripEventSummary> runsummary;
00354    iEvent.getByLabel( inputModuleLabel_, runsummary );
00355    if(runsummary->runType()==sistrip::APV_LATENCY) mode_ = 2; // LatencyScan
00356    else if(runsummary->runType()==sistrip::FINE_DELAY) mode_ = 1; // DelayScan
00357    else { 
00358     mode_ = 0; //unknown
00359     return;
00360    }
00361 
00362    if(noTracking_) {
00363       produceNoTracking(iEvent,iSetup);
00364       return;
00365    }
00366    event_ = &iEvent;
00367    // container for the selected hits
00368    std::vector< edm::DetSet<SiStripRawDigi> > output;
00369    output.reserve(100);
00370    // access the tracks
00371    edm::Handle<reco::TrackCollection> trackCollection;
00372    iEvent.getByLabel(trackLabel_,trackCollection);  
00373    const reco::TrackCollection *tracks=trackCollection.product();
00374    edm::ESHandle<TrackerGeometry> tracker;
00375    iSetup.get<TrackerDigiGeometryRecord>().get(tracker);
00376    if (tracks->size()) {
00377      anglefinder_->init(iEvent,iSetup);
00378      LogDebug("produce") << "Found " << tracks->size() << " tracks.";
00379      // look at the hits if one needs them
00380      edm::Handle< edm::DetSetVector<SiStripDigi> > hits;
00381      const edm::DetSetVector<SiStripDigi>* hitSet = NULL;
00382      if(homeMadeClusters_) {
00383        iEvent.getByLabel(digiLabel_,hits);
00384        hitSet = hits.product();
00385      }
00386      // look at the clusters 
00387      edm::Handle<edmNew::DetSetVector<SiStripCluster> > clusters;
00388      iEvent.getByLabel(clusterLabel_, clusters);
00389      const edmNew::DetSetVector<SiStripCluster>* clusterSet = clusters.product();
00390      // look at the trajectories if they are in the event
00391      std::vector<Trajectory> trajVec;
00392      edm::Handle<std::vector<Trajectory> > TrajectoryCollection;
00393      iEvent.getByLabel(trackLabel_,TrajectoryCollection);
00394      trajVec = *(TrajectoryCollection.product());
00395      // loop on tracks
00396      for(reco::TrackCollection::const_iterator itrack = tracks->begin(); itrack<tracks->end(); itrack++) {
00397        // first check the track Pt
00398        if((itrack->px()*itrack->px()+itrack->py()*itrack->py()+itrack->pz()*itrack->pz())<minTrackP2_) continue;
00399        // check that we have something in the layer we are interested in
00400        std::vector< std::pair<uint32_t,std::pair<double,double> > > intersections;
00401        if(mode_==1) {
00402          // Retrieve and decode commissioning information from "event summary"
00403          edm::Handle<SiStripEventSummary> summary;
00404          iEvent.getByLabel( inputModuleLabel_, summary );
00405          uint32_t layerCode = (const_cast<SiStripEventSummary*>(summary.product())->layerScanned())>>16;
00406          StripSubdetector::SubDetector subdet = StripSubdetector::TIB;
00407          if(((layerCode>>6)&0x3)==0) subdet = StripSubdetector::TIB;
00408          else if(((layerCode>>6)&0x3)==1) subdet = StripSubdetector::TOB;
00409          else if(((layerCode>>6)&0x3)==2) subdet = StripSubdetector::TID;
00410          else if(((layerCode>>6)&0x3)==3) subdet = StripSubdetector::TEC;
00411          int32_t layerIdx = (layerCode&0xF)*(((layerCode>>4)&0x3) ? -1 : 1);
00412          intersections = detId(*tracker,&(*itrack),trajVec,subdet,layerIdx);
00413        } else {
00414          // for latency scans, no layer is specified -> no cut on detid
00415          intersections = detId(*tracker,&(*itrack),trajVec);
00416        }
00417        LogDebug("produce") << "  Found " << intersections.size() << " interesting intersections." << std::endl;
00418        for(std::vector< std::pair<uint32_t,std::pair<double,double> > >::iterator it = intersections.begin();it<intersections.end();it++) {
00419          std::pair<const SiStripCluster*,double> candidateCluster = closestCluster(*tracker,&(*itrack),it->first,*clusterSet,*hitSet);
00420          if(candidateCluster.first) {
00421            LogDebug("produce") << "    Found a cluster."<< std::endl;
00422            // cut on the distance 
00423          if(candidateCluster.second>maxClusterDistance_) continue; 
00424            LogDebug("produce") << "    The cluster is close enough."<< std::endl;
00425          // build the rawdigi corresponding to the leading strip and save it
00426          // here, only the leading strip is retained. All other rawdigis in the module are set to 0.
00427          const std::vector< uint8_t >& amplitudes = candidateCluster.first->amplitudes();
00428          uint8_t leadingCharge = 0;
00429          uint8_t leadingStrip = candidateCluster.first->firstStrip();
00430          uint8_t leadingPosition = 0;
00431          for(std::vector< uint8_t >::const_iterator amplit = amplitudes.begin();amplit<amplitudes.end();amplit++,leadingStrip++) {
00432            if(leadingCharge<*amplit) {
00433              leadingCharge = *amplit;
00434              leadingPosition = leadingStrip;
00435            }
00436          }
00437 
00438          // look for an existing detset
00439          std::vector< edm::DetSet<SiStripRawDigi> >::iterator newdsit = output.begin();
00440          for(;newdsit!=output.end()&&newdsit->detId()!=connectionMap_[it->first];++newdsit) {}
00441          // if there is no detset yet, create it.
00442          if(newdsit==output.end()) {
00443            edm::DetSet<SiStripRawDigi> newds(connectionMap_[it->first]);
00444            output.push_back(newds);
00445            newdsit = output.end()-1;
00446          }
00447 
00448          LogDebug("produce") << " New Hit...   TOF:" << it->second.first << ", charge: " << int(leadingCharge) 
00449                              << " at " << int(leadingPosition) << "." << std::endl
00450                              << "Angular correction: " << it->second.second 
00451                              << " giving a final value of " << int(leadingCharge*fabs(it->second.second)) 
00452                              << " for fed key = " << connectionMap_[it->first] << " (detid=" << it->first << ")" ;
00453          // apply corrections to the leading charge, but only if it has not saturated.
00454          if(leadingCharge<255) {
00455            // correct the leading charge for the crossing angle
00456            leadingCharge = uint8_t(leadingCharge*fabs(it->second.second));
00457            // correct for module thickness for TEC and TOB
00458            if((((it->first>>25)&0x7f)==0xd) ||
00459               ((((it->first>>25)&0x7f)==0xe) && (((it->first>>5)&0x7)>4)))
00460              leadingCharge = uint8_t((leadingCharge*0.64));
00461          }
00462          //code the time of flight in the digi
00463          unsigned int tof = abs(int(round(it->second.first*10)));
00464          tof = tof>255 ? 255 : tof;
00465          SiStripRawDigi newSiStrip(leadingCharge + (tof<<8));
00466          newdsit->push_back(newSiStrip);
00467          LogDebug("produce") << "New edm::DetSet<SiStripRawDigi> added.";
00468          }
00469          if(homeMadeClusters_) delete candidateCluster.first; // we are owner of home-made clusters
00470        }
00471      }
00472    }
00473    // add the selected hits to the event.
00474    LogDebug("produce") << "Putting " << output.size() << " new hits in the event.";
00475    std::auto_ptr< edm::DetSetVector<SiStripRawDigi> > formatedOutput(new edm::DetSetVector<SiStripRawDigi>(output) );
00476    iEvent.put(formatedOutput,"FineDelaySelection");
00477 }
00478 
00479 // Simple solution when tracking is not available/ not working
00480 void
00481 SiStripFineDelayHit::produceNoTracking(edm::Event& iEvent, const edm::EventSetup& iSetup)
00482 {
00483    event_ = &iEvent;
00484    // container for the selected hits
00485    std::vector< edm::DetSet<SiStripRawDigi> > output;
00486    output.reserve(100);
00487    // Retrieve and decode commissioning information from "event summary"
00488    edm::Handle<SiStripEventSummary> summary;
00489    iEvent.getByLabel( inputModuleLabel_, summary );
00490    uint32_t layerCode = (const_cast<SiStripEventSummary*>(summary.product())->layerScanned())>>16;
00491    StripSubdetector::SubDetector subdet = StripSubdetector::TIB;
00492    if(((layerCode>>6)&0x3)==0) subdet = StripSubdetector::TIB;
00493    else if(((layerCode>>6)&0x3)==1) subdet = StripSubdetector::TOB;
00494    else if(((layerCode>>6)&0x3)==2) subdet = StripSubdetector::TID;
00495    else if(((layerCode>>6)&0x3)==3) subdet = StripSubdetector::TEC;
00496    int32_t layerIdx = (layerCode&0xF)*(((layerCode>>4)&0x3) ? -1 : 1);
00497    std::pair<uint32_t, uint32_t> mask = deviceMask(subdet,layerIdx);
00498    // look at the clusters 
00499    edm::Handle<edmNew::DetSetVector<SiStripCluster> > clusters;
00500    iEvent.getByLabel(clusterLabel_,clusters);
00501    for (edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter=clusters->begin(); DSViter!=clusters->end();DSViter++ ) {
00502      // check that we are in the layer of interest
00503      if(mode_==1 && ((DSViter->id() & mask.first) != mask.second) ) continue;
00504      // iterate over clusters
00505      edmNew::DetSet<SiStripCluster>::const_iterator begin=DSViter->begin();
00506      edmNew::DetSet<SiStripCluster>::const_iterator end  =DSViter->end();
00507      edm::DetSet<SiStripRawDigi> newds(connectionMap_[DSViter->id()]);
00508      for(edmNew::DetSet<SiStripCluster>::const_iterator iter=begin;iter!=end;++iter) {
00509          // build the rawdigi corresponding to the leading strip and save it
00510          // here, only the leading strip is retained. All other rawdigis in the module are set to 0.
00511          const std::vector< uint8_t >& amplitudes = iter->amplitudes();
00512          uint8_t leadingCharge = 0;
00513          uint8_t leadingStrip = iter->firstStrip();
00514          uint8_t leadingPosition = 0;
00515          for(std::vector< uint8_t >::const_iterator amplit = amplitudes.begin();amplit<amplitudes.end();amplit++,leadingStrip++) {
00516            if(leadingCharge<*amplit) {
00517              leadingCharge = *amplit;
00518              leadingPosition = leadingStrip;
00519            }
00520          }
00521          // apply some sanity cuts. This is needed since we don't use tracking to clean clusters
00522          // 1.5< noise <8
00523          // charge<250
00524          // 50 > s/n > 10
00525          edm::ESHandle<SiStripNoises> noiseHandle_;
00526          iSetup.get<SiStripNoisesRcd>().get(noiseHandle_);
00527          SiStripNoises::Range detNoiseRange = noiseHandle_->getRange(DSViter->id());  
00528          float noise=noiseHandle_->getNoise(leadingPosition, detNoiseRange);   
00529          if( noise<1.5 ) continue;
00530          if( leadingCharge>=250 || noise>=8 || leadingCharge/noise>50 || leadingCharge/noise<10 ) continue;
00531          // apply some correction to the leading charge, but only if it has not saturated.
00532          if(leadingCharge<255) {
00533            // correct for modulethickness for TEC and TOB
00534            if((((((DSViter->id())>>25)&0x7f)==0xd)||((((DSViter->id())>>25)&0x7f)==0xe))&&((((DSViter->id())>>5)&0x7)>4)) 
00535               leadingCharge = uint8_t((leadingCharge*0.64));
00536          }
00537          //code the time of flight == 0 in the digi
00538          SiStripRawDigi newSiStrip(leadingCharge);
00539          newds.push_back(newSiStrip);
00540      }
00541      //store into the detsetvector
00542      output.push_back(newds);
00543      LogDebug("produce") << "New edm::DetSet<SiStripRawDigi> added with fedkey = " 
00544                          << std::hex << std::setfill('0') << std::setw(8) 
00545                          << connectionMap_[DSViter->id()] << std::dec;
00546    }
00547    // add the selected hits to the event.
00548    LogDebug("produce") << "Putting " << output.size() << " new hits in the event.";
00549    std::auto_ptr< edm::DetSetVector<SiStripRawDigi> > formatedOutput(new edm::DetSetVector<SiStripRawDigi>(output) );
00550    iEvent.put(formatedOutput,"FineDelaySelection");
00551 }
00552 
00553 // ------------ method called once each job just before starting event loop  ------------
00554 void 
00555 SiStripFineDelayHit::beginRun(edm::Run & run, const edm::EventSetup & iSetup)
00556 {
00557    // Retrieve FED cabling object
00558    edm::ESHandle<SiStripFedCabling> cabling;
00559    iSetup.get<SiStripFedCablingRcd>().get( cabling );
00560    const std::vector< uint16_t > & feds = cabling->feds() ;
00561    for(std::vector< uint16_t >::const_iterator fedid = feds.begin();fedid<feds.end();++fedid) {
00562      const std::vector< FedChannelConnection > & connections = cabling->connections(*fedid);
00563      for(std::vector< FedChannelConnection >::const_iterator conn=connections.begin();conn<connections.end();++conn) {
00564      /*
00565        SiStripFedKey key(conn->fedId(),
00566                          SiStripFedKey::feUnit(conn->fedCh()),
00567                          SiStripFedKey::feChan(conn->fedCh()));
00568        connectionMap_[conn->detId()] = key.key();
00569      */
00570      // the key is computed using an alternate formula for performance reasons.
00571      connectionMap_[conn->detId()] = ( ( conn->fedId() & sistrip::invalid_ ) << 16 ) | ( conn->fedCh() & sistrip::invalid_ );
00572      }
00573    }
00574 }