Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-01-08 03:36:28

0001 /** SiPixelRecHitConverter.cc
0002  * ------------------------------------------------------
0003  * Description:  see SiPixelRecHitConverter.h
0004  * Authors:  P. Maksimovic (JHU), V.Chiochia (Uni Zurich)
0005  * History: Feb 27, 2006 -  initial version
0006  *          May 30, 2006 -  edm::DetSetVector and edm::Ref
0007  *          Aug 30, 2007 -  edmNew::DetSetVector
0008  *          Jan 31, 2008 -  change to use Lorentz angle from DB (Lotte Wilke)
0009  * ------------------------------------------------------
0010  */
0011 
0012 //---------------------------------------------------------------------------
0013 //! \class SiPixelRecHitConverter
0014 //!
0015 //! \brief EDProducer to covert SiPixelClusters into SiPixelRecHits
0016 //!
0017 //! SiPixelRecHitConverter is an EDProducer subclass (i.e., a module)
0018 //! which orchestrates the conversion of SiPixelClusters into SiPixelRecHits.
0019 //! Consequently, the input is a edm::DetSetVector<SiPixelCluster> and the output is
0020 //! SiPixelRecHitCollection.
0021 //!
0022 //! SiPixelRecHitConverter invokes one of descendents from
0023 //! ClusterParameterEstimator (templated on SiPixelCluster), e.g.
0024 //! CPEFromDetPosition (which is the only available option
0025 //! right now).  SiPixelRecHitConverter loads the SiPixelClusterCollection,
0026 //! and then iterates over DetIds, invoking the chosen CPE's methods
0027 //! localPosition() and localError() to perform the correction (some of which
0028 //! may be rather involved).  A RecHit is made on the spot, and appended
0029 //! to the output collection.
0030 //!
0031 //! The calibrations are not loaded at the moment,
0032 //! although that is being planned for the near future.
0033 //!
0034 //! \author Porting from ORCA by Petar Maksimovic (JHU). Implementation of the
0035 //!         DetSetVector by V.Chiochia (Zurich University).
0036 //!
0037 //! \version v2, May 30, 2006
0038 //! change to use Lorentz angle from DB Lotte Wilke, Jan. 31st, 2008
0039 //!
0040 //---------------------------------------------------------------------------
0041 
0042 //--- Base class for CPEs:
0043 
0044 #include "RecoLocalTracker/SiPixelRecHits/interface/PixelCPEBase.h"
0045 
0046 //--- Geometry + DataFormats
0047 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0048 #include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
0049 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHitCollection.h"
0050 #include "DataFormats/Common/interface/DetSetVector.h"
0051 
0052 //--- Framework
0053 #include "FWCore/Framework/interface/stream/EDProducer.h"
0054 #include "FWCore/Framework/interface/Event.h"
0055 #include "FWCore/Framework/interface/EventSetup.h"
0056 #include "FWCore/Framework/interface/MakerMacros.h"
0057 #include "DataFormats/Common/interface/Handle.h"
0058 #include "FWCore/Framework/interface/ESHandle.h"
0059 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0060 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0061 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0062 #include "FWCore/Utilities/interface/InputTag.h"
0063 #include "FWCore/Utilities/interface/EDPutToken.h"
0064 #include "FWCore/Utilities/interface/ESGetToken.h"
0065 
0066 // Geometry
0067 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0068 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0069 
0070 // Data Formats
0071 #include "DataFormats/DetId/interface/DetId.h"
0072 #include "DataFormats/Common/interface/Ref.h"
0073 #include "DataFormats/Common/interface/DetSet2RangeMap.h"
0074 
0075 // STL
0076 #include <vector>
0077 #include <memory>
0078 #include <string>
0079 #include <iostream>
0080 
0081 // MessageLogger
0082 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0083 
0084 #include "RecoLocalTracker/Records/interface/TkPixelCPERecord.h"
0085 
0086 // Make heterogeneous framework happy
0087 #include "CUDADataFormats/SiPixelCluster/interface/gpuClusteringConstants.h"
0088 #include "CUDADataFormats/Common/interface/HostProduct.h"
0089 
0090 using namespace std;
0091 
0092 namespace cms {
0093 
0094   class SiPixelRecHitConverter : public edm::stream::EDProducer<> {
0095   public:
0096     //--- Constructor, virtual destructor (just in case)
0097     explicit SiPixelRecHitConverter(const edm::ParameterSet& conf);
0098     ~SiPixelRecHitConverter() override;
0099 
0100     static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0101 
0102     //--- Factory method to make CPE's depending on the ParameterSet
0103     //--- Not sure if we need to make more than one CPE to run concurrently
0104     //--- on different parts of the detector (e.g., one for the barrel and the
0105     //--- one for the forward).  The way the CPE's are written now, it's
0106     //--- likely we can use one (and they will switch internally), or
0107     //--- make two of the same but configure them differently.  We need a more
0108     //--- realistic use case...
0109 
0110     //--- The top-level event method.
0111     void produce(edm::Event& e, const edm::EventSetup& c) override;
0112 
0113     //--- Execute the position estimator algorithm(s).
0114     void run(edm::Event& e,
0115              edm::Handle<SiPixelClusterCollectionNew> inputhandle,
0116              SiPixelRecHitCollectionNew& output,
0117              TrackerGeometry const& geom);
0118 
0119   private:
0120     using HMSstorage = HostProduct<uint32_t[]>;
0121 
0122     // TO DO: maybe allow a map of pointers?
0123     PixelCPEBase const* cpe_ = nullptr;  // What we got (for now, one ptr to base class)
0124     edm::InputTag const src_;
0125     std::string const cpeName_;
0126     edm::EDGetTokenT<SiPixelClusterCollectionNew> const tPixelCluster_;
0127     edm::EDPutTokenT<SiPixelRecHitCollection> const tPut_;
0128     edm::EDPutTokenT<HMSstorage> const tHost_;
0129     edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> const tTrackerGeom_;
0130     edm::ESGetToken<PixelClusterParameterEstimator, TkPixelCPERecord> const tCPE_;
0131     bool m_newCont;  // save also in emdNew::DetSetVector
0132   };
0133 
0134   //---------------------------------------------------------------------------
0135   //!  Constructor: set the ParameterSet and defer all thinking to setupCPE().
0136   //---------------------------------------------------------------------------
0137   SiPixelRecHitConverter::SiPixelRecHitConverter(edm::ParameterSet const& conf)
0138       : src_(conf.getParameter<edm::InputTag>("src")),
0139         cpeName_(conf.getParameter<std::string>("CPE")),
0140         tPixelCluster_(consumes<SiPixelClusterCollectionNew>(src_)),
0141         tPut_(produces<SiPixelRecHitCollection>()),
0142         tHost_(produces<HMSstorage>()),
0143         tTrackerGeom_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>()),
0144         tCPE_(esConsumes<PixelClusterParameterEstimator, TkPixelCPERecord>(edm::ESInputTag("", cpeName_))) {}
0145 
0146   // Destructor
0147   SiPixelRecHitConverter::~SiPixelRecHitConverter() {}
0148 
0149   //---------------------------------------------------------------------------
0150   //! The "Event" entrypoint: gets called by framework for every event
0151   //---------------------------------------------------------------------------
0152   void SiPixelRecHitConverter::produce(edm::Event& e, const edm::EventSetup& es) {
0153     // Step A.1: get input data
0154     edm::Handle<SiPixelClusterCollectionNew> input;
0155     e.getByToken(tPixelCluster_, input);
0156 
0157     // Step A.2: get event setup
0158     auto const& geom = es.getData(tTrackerGeom_);
0159 
0160     // Step B: create empty output collection
0161     SiPixelRecHitCollectionNew output;
0162 
0163     // Step B*: create CPE
0164     cpe_ = dynamic_cast<const PixelCPEBase*>(&es.getData(tCPE_));
0165 
0166     // Step C: Iterate over DetIds and invoke the strip CPE algorithm
0167     // on each DetUnit
0168 
0169     run(e, input, output, geom);
0170 
0171     output.shrink_to_fit();
0172     e.emplace(tPut_, std::move(output));
0173   }
0174 
0175   //---------------------------------------------------------------------------
0176   //!  Iterate over DetUnits, then over Clusters and invoke the CPE on each,
0177   //!  and make a RecHit to store the result.
0178   //!  New interface reading DetSetVector by V.Chiochia (May 30th, 2006)
0179   //---------------------------------------------------------------------------
0180   void SiPixelRecHitConverter::run(edm::Event& iEvent,
0181                                    edm::Handle<SiPixelClusterCollectionNew> inputhandle,
0182                                    SiPixelRecHitCollectionNew& output,
0183                                    TrackerGeometry const& geom) {
0184     if (!cpe_) {
0185       edm::LogError("SiPixelRecHitConverter") << " at least one CPE is not ready -- can't run!";
0186       // TO DO: throw an exception here?  The user may want to know...
0187       assert(0);
0188       return;  // clusterizer is invalid, bail out
0189     }
0190 
0191     int numberOfDetUnits = 0;
0192     int numberOfClusters = 0;
0193 
0194     const SiPixelClusterCollectionNew& input = *inputhandle;
0195 
0196     // allocate a buffer for the indices of the clusters
0197     auto hmsp = std::make_unique<uint32_t[]>(gpuClustering::maxNumModules + 1);
0198     // hitsModuleStart is a non-owning pointer to the buffer
0199     auto hitsModuleStart = hmsp.get();
0200     // fill cluster arrays
0201     std::array<uint32_t, gpuClustering::maxNumModules + 1> clusInModule{};
0202     for (auto const& dsv : input) {
0203       unsigned int detid = dsv.detId();
0204       DetId detIdObject(detid);
0205       const GeomDetUnit* genericDet = geom.idToDetUnit(detIdObject);
0206       auto gind = genericDet->index();
0207       // FIXME to be changed to support Phase2
0208       if (gind >= int(gpuClustering::maxNumModules))
0209         continue;
0210       auto const nclus = dsv.size();
0211       assert(nclus > 0);
0212       clusInModule[gind] = nclus;
0213       numberOfClusters += nclus;
0214     }
0215     hitsModuleStart[0] = 0;
0216     assert(clusInModule.size() > gpuClustering::maxNumModules);
0217     for (int i = 1, n = clusInModule.size(); i < n; ++i)
0218       hitsModuleStart[i] = hitsModuleStart[i - 1] + clusInModule[i - 1];
0219     assert(numberOfClusters == int(hitsModuleStart[gpuClustering::maxNumModules]));
0220 
0221     // wrap the buffer in a HostProduct, and move it to the Event, without reallocating the buffer or affecting hitsModuleStart
0222     iEvent.emplace(tHost_, std::move(hmsp));
0223 
0224     numberOfClusters = 0;
0225     for (auto const& dsv : input) {
0226       numberOfDetUnits++;
0227       unsigned int detid = dsv.detId();
0228       DetId detIdObject(detid);
0229       const GeomDetUnit* genericDet = geom.idToDetUnit(detIdObject);
0230       const PixelGeomDetUnit* pixDet = dynamic_cast<const PixelGeomDetUnit*>(genericDet);
0231       assert(pixDet);
0232       SiPixelRecHitCollectionNew::FastFiller recHitsOnDetUnit(output, detid);
0233 
0234       edmNew::DetSet<SiPixelCluster>::const_iterator clustIt = dsv.begin(), clustEnd = dsv.end();
0235 
0236       for (; clustIt != clustEnd; clustIt++) {
0237         numberOfClusters++;
0238         std::tuple<LocalPoint, LocalError, SiPixelRecHitQuality::QualWordType> tuple =
0239             cpe_->getParameters(*clustIt, *genericDet);
0240         LocalPoint lp(std::get<0>(tuple));
0241         LocalError le(std::get<1>(tuple));
0242         SiPixelRecHitQuality::QualWordType rqw(std::get<2>(tuple));
0243         // Create a persistent edm::Ref to the cluster
0244         SiPixelClusterRefNew cluster = edmNew::makeRefTo(inputhandle, clustIt);
0245         // Make a RecHit and add it to the DetSet
0246         SiPixelRecHit hit(lp, le, rqw, *genericDet, cluster);
0247         recHitsOnDetUnit.push_back(hit);
0248 
0249         LogDebug("SiPixelRecHitConverter") << "RecHit " << (numberOfClusters - 1)  //
0250                                            << " with local position " << lp << " and local error " << le;
0251       }  //  <-- End loop on Clusters
0252 
0253       LogDebug("SiPixelRecHitConverter") << "Found " << recHitsOnDetUnit.size() << " RecHits on " << detid;
0254 
0255     }  //    <-- End loop on DetUnits
0256 
0257     LogDebug("SiPixelRecHitConverter") << cpeName_ << " converted " << numberOfClusters
0258                                        << " SiPixelClusters into SiPixelRecHits, in " << numberOfDetUnits
0259                                        << " DetUnits.";
0260   }
0261 }  // end of namespace cms
0262 
0263 using cms::SiPixelRecHitConverter;
0264 
0265 void SiPixelRecHitConverter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0266   edm::ParameterSetDescription desc;
0267   desc.add<edm::InputTag>("src", edm::InputTag("siPixelClusters"));
0268   desc.add<std::string>("CPE", "PixelCPEGeneric");
0269   descriptions.addWithDefaultLabel(desc);
0270 }
0271 
0272 DEFINE_FWK_MODULE(SiPixelRecHitConverter);