CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/FastSimulation/TrackingRecHitProducer/src/SiClusterTranslator.cc

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