CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch2/src/FastSimulation/TrackingRecHitProducer/src/SiClusterTranslator.cc

Go to the documentation of this file.
00001 
00010 // SiTracker Gaussian Smearing
00011 #include "FastSimulation/TrackingRecHitProducer/interface/SiClusterTranslator.h"
00012 
00013 // Geometry
00014 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00015 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00016 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00017 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetUnit.h"
00018 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetType.h"
00019 #include "Geometry/CommonTopologies/interface/RectangularStripTopology.h"
00020 #include "Geometry/CommonTopologies/interface/RadialStripTopology.h"
00021 #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetUnit.h"
00022 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
00023 
00024 //CPEs
00025 #include "RecoLocalTracker/Records/interface/TkPixelCPERecord.h"
00026 #include "RecoLocalTracker/Records/interface/TkStripCPERecord.h"
00027 #include "FastSimulation/TrackingRecHitProducer/interface/FastPixelCPE.h"
00028 #include "FastSimulation/TrackingRecHitProducer/interface/FastStripCPE.h"
00029 
00030 // Framework
00031 #include "FWCore/Framework/interface/Event.h"
00032 #include "FWCore/Framework/interface/EventSetup.h"
00033 #include "FWCore/Framework/interface/ESHandle.h"
00034 #include "FWCore/Framework/interface/Frameworkfwd.h"
00035 #include "FWCore/Framework/interface/Event.h"
00036 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00037 
00038 // Data Formats
00039 #include "DataFormats/Common/interface/Handle.h"
00040 #include "DataFormats/DetId/interface/DetId.h"
00041 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
00042 #include "DataFormats/GeometrySurface/interface/LocalError.h"
00043 #include "DataFormats/Common/interface/DetSetVector.h"
00044 #include "DataFormats/SiPixelDetId/interface/PixelChannelIdentifier.h"
00045 
00046 //for the SimHit
00047 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
00048 #include "SimDataFormats/CrossingFrame/interface/CrossingFrame.h"
00049 #include "SimDataFormats/CrossingFrame/interface/MixCollection.h"
00050 #include "SimDataFormats/EncodedEventId/interface/EncodedEventId.h"
00051 
00052 // STL
00053 #include <memory>
00054 #include <string>
00055 
00056 SiClusterTranslator::SiClusterTranslator(edm::ParameterSet const& conf) 
00057 {
00058   produces<edmNew::DetSetVector<SiStripCluster> >();
00059   produces<edmNew::DetSetVector<SiPixelCluster> >();
00060   produces<edm::DetSetVector<StripDigiSimLink> >();
00061   produces<edm::DetSetVector<PixelDigiSimLink> >();
00062 }
00063 
00064 // Destructor
00065 SiClusterTranslator::~SiClusterTranslator() {}  
00066 
00067 void 
00068 SiClusterTranslator::beginRun(edm::Run & run, const edm::EventSetup & es) {
00069 
00070   // Initialize the Tracker Geometry
00071   edm::ESHandle<TrackerGeometry> theGeometry;
00072   es.get<TrackerDigiGeometryRecord> ().get (theGeometry);
00073   geometry = &(*theGeometry);
00074 }
00075 
00076 void 
00077 SiClusterTranslator::produce(edm::Event& e, const edm::EventSetup& es) 
00078 {
00079   // Step A: Get Inputs (FastGSRecHit's)
00080   edm::Handle<FastTrackerClusterCollection> theFastClusters; 
00081   e.getByType(theFastClusters);
00082   
00083   edm::ESHandle<TrackerGeometry> tkgeom;
00084   es.get<TrackerDigiGeometryRecord>().get( tkgeom ); 
00085   const TrackerGeometry &tracker(*tkgeom);
00086   
00087   edm::ESHandle<PixelClusterParameterEstimator> pixelCPE;
00088   es.get<TkPixelCPERecord>().get("FastPixelCPE",pixelCPE);
00089   const PixelClusterParameterEstimator &pixelcpe(*pixelCPE);
00090   
00091   edm::ESHandle<StripClusterParameterEstimator> stripCPE;
00092   es.get<TkStripCPERecord>().get("FastStripCPE", stripCPE); 
00093   const StripClusterParameterEstimator &stripcpe(*stripCPE);
00094     
00095   edm::Handle<CrossingFrame<PSimHit> > cf_simhit;
00096   std::vector<const CrossingFrame<PSimHit> *> cf_simhitvec;
00097   e.getByLabel("mix","famosSimHitsTrackerHits", cf_simhit);
00098   cf_simhitvec.push_back(cf_simhit.product());
00099   
00100   std::auto_ptr<MixCollection<PSimHit> > allTrackerHits(new MixCollection<PSimHit>(cf_simhitvec));
00101   int counter =0;
00102   
00103   //Clearing vector to hopefully make it run faster.
00104   theNewSimHitList.clear();
00105   thePixelDigiLinkVector.clear();
00106   theStripDigiLinkVector.clear();
00107 
00108   for(MixCollection<PSimHit>::iterator it = allTrackerHits->begin(); it!= allTrackerHits->end();it++){
00109     counter++;
00110     theNewSimHitList.push_back(std::make_pair((*it), counter));
00111   }
00112 
00113   // Step B: fill a temporary full Cluster collection from the fast Cluster collection
00114   FastTrackerClusterCollection::const_iterator aCluster = theFastClusters->begin();
00115   FastTrackerClusterCollection::const_iterator theLastHit = theFastClusters->end();
00116   std::map< DetId, std::vector<SiPixelCluster> > temporaryPixelClusters;
00117   std::map< DetId, std::vector<SiStripCluster> > temporaryStripClusters;
00118   
00119   //Clearing CPE maps from previous event.
00120   stripcpe.clearParameters();
00121   pixelcpe.clearParameters();
00122   
00123   int ClusterNum = 0;
00124   
00125   // loop on Fast GS Hits
00126   for ( ; aCluster != theLastHit; ++aCluster ) {
00127     ClusterNum++;
00128     
00129     //Finding SubDet Id of cluster: 1 & 2 = Pixel, 3 & 4 & 5 & 6 = Strip, >6 = ?
00130     DetId det = aCluster->id();
00131     unsigned int subdet   = det.subdetId();
00132     int sim_counter = 0;
00133     for (std::vector<std::pair<PSimHit,int> >::const_iterator 
00134            simcount = theNewSimHitList.begin() ; simcount != theNewSimHitList.end(); simcount ++){
00135       if((aCluster->simtrackId() == (int)(*simcount).first.trackId())&&(det.rawId() == (*simcount).first.detUnitId())&&(aCluster->eeId() == (*simcount).first.eventId().rawId()))
00136         sim_counter = (*simcount).second;
00137     }
00138     if (sim_counter == 0)  throw cms::Exception("SiClusterTranslator") << "No Matching SimHit found.";
00139     
00140     //Separating into Pixels and Strips
00141     
00142     //Pixel
00143     if (subdet < 3) {
00144       //Here is the hard part. From the position of the FastSim Cluster I need to figure out the Pixel location for the Cluster.
00145       LocalPoint position = aCluster->localPosition();
00146       LocalError error = aCluster->localPositionError();
00147       //std::cout << "The pixel charge is " << aCluster->charge() << std::endl;
00148       int charge = (int)(aCluster->charge() + 0.5);
00149       //std::cout << "The pixel charge after integer conversion is " << charge << std::endl;
00150 
00151       //std::vector<int> digi_vec;
00152       //while (charge > 255) { 
00153       //digi_vec.push_back(charge);
00154       //charge -= 256;
00155       //}
00156       //digi_vec.push_back(charge);
00157       
00158       const GeomDetUnit *  geoDet = tracker.idToDetUnit(det);
00159       const PixelGeomDetUnit * pixelDet=(const PixelGeomDetUnit*)(geoDet);
00160       const PixelTopology& topol=(PixelTopology&)pixelDet->topology();
00161       
00162       //Out of pixel is float, but input of pixel is int. Hopeful it works...
00163       std::pair<float,float> pixelPos_out = topol.pixel(position);
00164       SiPixelCluster::PixelPos pixelPos((int)pixelPos_out.first, (int)pixelPos_out.second);
00165       
00166       //Filling Pixel CPE with information.
00167       std::pair<int,int> row_col((int)pixelPos_out.first,(int)pixelPos_out.second);
00168       pixelcpe.enterLocalParameters((unsigned int) det.rawId() , row_col, std::make_pair(position,error));
00169       
00170       unsigned int ch = PixelChannelIdentifier::pixelToChannel((int)pixelPos_out.first, (int)pixelPos_out.second);
00171       
00172       //Creating a new pixel cluster.
00173       SiPixelCluster temporaryPixelCluster(pixelPos, charge);
00174       temporaryPixelClusters[det].push_back(temporaryPixelCluster);
00175       
00176       //Making a PixelDigiSimLink.
00177       edm::DetSet<PixelDigiSimLink> pixelDetSet;
00178       pixelDetSet.id = det.rawId();
00179       pixelDetSet.data.push_back(PixelDigiSimLink(ch,
00180                                                   aCluster->simtrackId(),
00181                                                   EncodedEventId(aCluster->eeId()),
00182                                                   1.0));
00183       thePixelDigiLinkVector.push_back(pixelDetSet);
00184     } 
00185     
00186     //Strips
00187     else if ((subdet > 2) && (subdet < 7)) {
00188       //Getting pos/err info from the cluster
00189       LocalPoint position = aCluster->localPosition();
00190       LocalError error = aCluster->localPositionError();
00191       
00192       //Will have to make charge into ADC eventually...
00193       uint16_t charge = (uint16_t)(aCluster->charge() + 0.5);
00194       
00195       //std::cout << "The charge is " << charge << std::endl;
00196       
00197       uint16_t strip_num = 0;
00198       std::vector<uint16_t> digi_vec;
00199       while (charge > 255) { 
00200         digi_vec.push_back(255);
00201         charge -= 255;
00202       }
00203       if (charge > 0) digi_vec.push_back(charge);
00204       //std::cout << "The digi_vec size is " << digi_vec.size() << std::endl;
00205       //int totcharge = 0;
00206       //for(int i = 0; i < digi_vec.size(); ++i) {
00207       //totcharge += digi_vec[i];
00208       //} 
00209       const GeomDetUnit *  geoDet = tracker.idToDetUnit(det);
00210       const StripGeomDetUnit * stripDet = (const StripGeomDetUnit*)(geoDet);
00211       
00212       //3 = TIB, 4 = TID, 5 = TOB, 6 = TEC
00213       if((subdet == 3) || (subdet == 5)) {
00214         const RectangularStripTopology& topol=(RectangularStripTopology&)stripDet->type().topology();
00215         strip_num = (uint16_t)topol.strip(position);
00216       } else if ((subdet == 4) || (subdet == 6)) {
00217         const RadialStripTopology& topol=(RadialStripTopology&)stripDet->type().topology();
00218         strip_num = (uint16_t)topol.strip(position);
00219       }
00220       
00221       //Filling Strip CPE with info.
00222       stripcpe.enterLocalParameters(det.rawId(), strip_num, std::make_pair(position,error));
00223       
00224       //Creating a new strip cluster
00225       SiStripCluster temporaryStripCluster(det.rawId(), strip_num, digi_vec.begin(), digi_vec.end());
00226       temporaryStripClusters[det].push_back(temporaryStripCluster);
00227       
00228       //Making a StripDigiSimLink
00229       edm::DetSet<StripDigiSimLink> stripDetSet;
00230       stripDetSet.id = det.rawId();
00231       stripDetSet.data.push_back(StripDigiSimLink(strip_num,
00232                                                   aCluster->simtrackId(),
00233                                                   sim_counter,
00234                                                   EncodedEventId(aCluster->eeId()),
00235                                                   1.0));
00236       theStripDigiLinkVector.push_back(stripDetSet);
00237     } 
00238     
00239     //?????
00240     else {
00241       throw cms::Exception("SiClusterTranslator") <<
00242         "Trying to build a cluster that is not in the SiStripTracker or Pixels.\n";
00243     }
00244     
00245   }//Cluster loop
00246 
00247   // Step C: from the temporary Cluster collections, create the real ones.
00248   
00249   //Pixels
00250   std::auto_ptr<edmNew::DetSetVector<SiPixelCluster> >
00251     siPixelClusterCollection(new edmNew::DetSetVector<SiPixelCluster>);
00252   loadPixelClusters(temporaryPixelClusters, *siPixelClusterCollection);
00253   std::auto_ptr<edm::DetSetVector<StripDigiSimLink> > stripoutputlink(new edm::DetSetVector<StripDigiSimLink>(theStripDigiLinkVector) );
00254   
00255   
00256   //Strips
00257   std::auto_ptr<edmNew::DetSetVector<SiStripCluster> > 
00258     siStripClusterCollection(new edmNew::DetSetVector<SiStripCluster>);
00259   loadStripClusters(temporaryStripClusters, *siStripClusterCollection);
00260   std::auto_ptr<edm::DetSetVector<PixelDigiSimLink> > pixeloutputlink(new edm::DetSetVector<PixelDigiSimLink>(thePixelDigiLinkVector) );
00261   
00262   // Step D: write output to file
00263   e.put(siPixelClusterCollection);
00264   e.put(siStripClusterCollection);
00265   e.put(stripoutputlink);
00266   e.put(pixeloutputlink);
00267 }
00268 
00269 void 
00270 SiClusterTranslator::loadStripClusters(
00271                                        std::map<DetId,std::vector<SiStripCluster> >& theClusters,
00272                                        edmNew::DetSetVector<SiStripCluster>& theClusterCollection) const
00273 {
00274   std::map<DetId,std::vector<SiStripCluster> >::const_iterator 
00275     it = theClusters.begin();
00276   std::map<DetId,std::vector<SiStripCluster> >::const_iterator 
00277     lastDet = theClusters.end();
00278   for( ; it != lastDet ; ++it ) { 
00279     edmNew::DetSetVector<SiStripCluster>::FastFiller cluster_col(theClusterCollection, it->first);
00280     
00281     std::vector<SiStripCluster>::const_iterator clust_it = it->second.begin();
00282     std::vector<SiStripCluster>::const_iterator clust_end = it->second.end();
00283     for( ; clust_it != clust_end ; ++clust_it)  {
00284       cluster_col.push_back(*clust_it);
00285     }
00286   }
00287 }
00288 
00289 void 
00290 SiClusterTranslator::loadPixelClusters(
00291                                        std::map<DetId,std::vector<SiPixelCluster> >& theClusters,
00292                                        edmNew::DetSetVector<SiPixelCluster>& theClusterCollection) const
00293 {
00294   
00295   std::map<DetId,std::vector<SiPixelCluster> >::const_iterator 
00296     it = theClusters.begin();
00297   std::map<DetId,std::vector<SiPixelCluster> >::const_iterator 
00298     lastCluster = theClusters.end();
00299   for( ; it != lastCluster ; ++it ) { 
00300     edmNew::DetSetVector<SiPixelCluster>::FastFiller spc(theClusterCollection, it->first);
00301     
00302     std::vector<SiPixelCluster>::const_iterator clust_it = it->second.begin();
00303     std::vector<SiPixelCluster>::const_iterator clust_end = it->second.end();
00304     for( ; clust_it != clust_end ; ++clust_it)  {
00305       spc.push_back(*clust_it);
00306     }
00307   }
00308 }