Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-02-14 23:30:43

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 
0058 #include "DataFormats/Common/interface/Handle.h"
0059 #include "FWCore/Framework/interface/ESHandle.h"
0060 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0061 #include "FWCore/Utilities/interface/InputTag.h"
0062 #include "FWCore/Utilities/interface/EDPutToken.h"
0063 #include "FWCore/Utilities/interface/ESGetToken.h"
0064 
0065 // Geometry
0066 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0067 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0068 
0069 // Data Formats
0070 #include "DataFormats/DetId/interface/DetId.h"
0071 #include "DataFormats/Common/interface/Ref.h"
0072 #include "DataFormats/Common/interface/DetSet2RangeMap.h"
0073 
0074 // STL
0075 #include <vector>
0076 #include <memory>
0077 #include <string>
0078 #include <iostream>
0079 
0080 // MessageLogger
0081 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0082 
0083 #include "RecoLocalTracker/Records/interface/TkPixelCPERecord.h"
0084 
0085 // Make heterogeneous framework happy
0086 #include "CUDADataFormats/SiPixelCluster/interface/gpuClusteringConstants.h"
0087 #include "CUDADataFormats/Common/interface/HostProduct.h"
0088 
0089 using namespace std;
0090 
0091 namespace cms {
0092 
0093   class SiPixelRecHitConverter : public edm::stream::EDProducer<> {
0094   public:
0095     //--- Constructor, virtual destructor (just in case)
0096     explicit SiPixelRecHitConverter(const edm::ParameterSet& conf);
0097     ~SiPixelRecHitConverter() override;
0098 
0099     //--- Factory method to make CPE's depending on the ParameterSet
0100     //--- Not sure if we need to make more than one CPE to run concurrently
0101     //--- on different parts of the detector (e.g., one for the barrel and the
0102     //--- one for the forward).  The way the CPE's are written now, it's
0103     //--- likely we can use one (and they will switch internally), or
0104     //--- make two of the same but configure them differently.  We need a more
0105     //--- realistic use case...
0106 
0107     //--- The top-level event method.
0108     void produce(edm::Event& e, const edm::EventSetup& c) override;
0109 
0110     //--- Execute the position estimator algorithm(s).
0111     void run(edm::Event& e,
0112              edm::Handle<edmNew::DetSetVector<SiPixelCluster>> inputhandle,
0113              SiPixelRecHitCollectionNew& output,
0114              TrackerGeometry const& geom);
0115 
0116   private:
0117     using HMSstorage = HostProduct<uint32_t[]>;
0118 
0119     // TO DO: maybe allow a map of pointers?
0120     /// const PixelClusterParameterEstimator * cpe_;  // what we got (for now, one ptr to base class)
0121     PixelCPEBase const* cpe_ = nullptr;  // What we got (for now, one ptr to base class)
0122     edm::InputTag const src_;
0123     edm::EDGetTokenT<edmNew::DetSetVector<SiPixelCluster>> const tPixelCluster_;
0124     edm::EDPutTokenT<SiPixelRecHitCollection> const tPut_;
0125     edm::EDPutTokenT<HMSstorage> const tHost_;
0126     edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> const tTrackerGeom_;
0127     edm::ESGetToken<PixelClusterParameterEstimator, TkPixelCPERecord> const tCPE_;
0128     bool m_newCont;  // save also in emdNew::DetSetVector
0129   };
0130 
0131   //---------------------------------------------------------------------------
0132   //!  Constructor: set the ParameterSet and defer all thinking to setupCPE().
0133   //---------------------------------------------------------------------------
0134   SiPixelRecHitConverter::SiPixelRecHitConverter(edm::ParameterSet const& conf)
0135       : src_(conf.getParameter<edm::InputTag>("src")),
0136         tPixelCluster_(consumes<edmNew::DetSetVector<SiPixelCluster>>(src_)),
0137         tPut_(produces<SiPixelRecHitCollection>()),
0138         tHost_(produces<HMSstorage>()),
0139         tTrackerGeom_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>()),
0140         tCPE_(esConsumes<PixelClusterParameterEstimator, TkPixelCPERecord>(
0141             edm::ESInputTag("", conf.getParameter<std::string>("CPE")))) {}
0142 
0143   // Destructor
0144   SiPixelRecHitConverter::~SiPixelRecHitConverter() {}
0145 
0146   //---------------------------------------------------------------------------
0147   //! The "Event" entrypoint: gets called by framework for every event
0148   //---------------------------------------------------------------------------
0149   void SiPixelRecHitConverter::produce(edm::Event& e, const edm::EventSetup& es) {
0150     // Step A.1: get input data
0151     edm::Handle<edmNew::DetSetVector<SiPixelCluster>> input;
0152     e.getByToken(tPixelCluster_, input);
0153 
0154     // Step A.2: get event setup
0155     auto const& geom = es.getData(tTrackerGeom_);
0156 
0157     // Step B: create empty output collection
0158     SiPixelRecHitCollectionNew output;
0159 
0160     // Step B*: create CPE
0161     cpe_ = dynamic_cast<const PixelCPEBase*>(&es.getData(tCPE_));
0162 
0163     // Step C: Iterate over DetIds and invoke the strip CPE algorithm
0164     // on each DetUnit
0165 
0166     run(e, input, output, geom);
0167 
0168     output.shrink_to_fit();
0169     e.emplace(tPut_, std::move(output));
0170   }
0171 
0172   //---------------------------------------------------------------------------
0173   //!  Iterate over DetUnits, then over Clusters and invoke the CPE on each,
0174   //!  and make a RecHit to store the result.
0175   //!  New interface reading DetSetVector by V.Chiochia (May 30th, 2006)
0176   //---------------------------------------------------------------------------
0177   void SiPixelRecHitConverter::run(edm::Event& iEvent,
0178                                    edm::Handle<edmNew::DetSetVector<SiPixelCluster>> inputhandle,
0179                                    SiPixelRecHitCollectionNew& output,
0180                                    TrackerGeometry const& geom) {
0181     if (!cpe_) {
0182       edm::LogError("SiPixelRecHitConverter") << " at least one CPE is not ready -- can't run!";
0183       // TO DO: throw an exception here?  The user may want to know...
0184       assert(0);
0185       return;  // clusterizer is invalid, bail out
0186     }
0187 
0188     int numberOfDetUnits = 0;
0189     int numberOfClusters = 0;
0190 
0191     const edmNew::DetSetVector<SiPixelCluster>& input = *inputhandle;
0192 
0193     // allocate a buffer for the indices of the clusters
0194     auto hmsp = std::make_unique<uint32_t[]>(gpuClustering::maxNumModules + 1);
0195     // hitsModuleStart is a non-owning pointer to the buffer
0196     auto hitsModuleStart = hmsp.get();
0197     // fill cluster arrays
0198     std::array<uint32_t, gpuClustering::maxNumModules + 1> clusInModule{};
0199     for (auto const& dsv : input) {
0200       unsigned int detid = dsv.detId();
0201       DetId detIdObject(detid);
0202       const GeomDetUnit* genericDet = geom.idToDetUnit(detIdObject);
0203       auto gind = genericDet->index();
0204       // FIXME to be changed to support Phase2
0205       if (gind >= int(gpuClustering::maxNumModules))
0206         continue;
0207       auto const nclus = dsv.size();
0208       assert(nclus > 0);
0209       clusInModule[gind] = nclus;
0210       numberOfClusters += nclus;
0211     }
0212     hitsModuleStart[0] = 0;
0213     assert(clusInModule.size() > gpuClustering::maxNumModules);
0214     for (int i = 1, n = clusInModule.size(); i < n; ++i)
0215       hitsModuleStart[i] = hitsModuleStart[i - 1] + clusInModule[i - 1];
0216     assert(numberOfClusters == int(hitsModuleStart[gpuClustering::maxNumModules]));
0217 
0218     // wrap the buffer in a HostProduct, and move it to the Event, without reallocating the buffer or affecting hitsModuleStart
0219     iEvent.emplace(tHost_, std::move(hmsp));
0220 
0221     numberOfClusters = 0;
0222     for (auto const& dsv : input) {
0223       numberOfDetUnits++;
0224       unsigned int detid = dsv.detId();
0225       DetId detIdObject(detid);
0226       const GeomDetUnit* genericDet = geom.idToDetUnit(detIdObject);
0227       const PixelGeomDetUnit* pixDet = dynamic_cast<const PixelGeomDetUnit*>(genericDet);
0228       assert(pixDet);
0229       SiPixelRecHitCollectionNew::FastFiller recHitsOnDetUnit(output, detid);
0230 
0231       edmNew::DetSet<SiPixelCluster>::const_iterator clustIt = dsv.begin(), clustEnd = dsv.end();
0232 
0233       for (; clustIt != clustEnd; clustIt++) {
0234         numberOfClusters++;
0235         std::tuple<LocalPoint, LocalError, SiPixelRecHitQuality::QualWordType> tuple =
0236             cpe_->getParameters(*clustIt, *genericDet);
0237         LocalPoint lp(std::get<0>(tuple));
0238         LocalError le(std::get<1>(tuple));
0239         SiPixelRecHitQuality::QualWordType rqw(std::get<2>(tuple));
0240         // Create a persistent edm::Ref to the cluster
0241         edm::Ref<edmNew::DetSetVector<SiPixelCluster>, SiPixelCluster> cluster =
0242             edmNew::makeRefTo(inputhandle, clustIt);
0243         // Make a RecHit and add it to the DetSet
0244         // old : recHitsOnDetUnit.push_back( new SiPixelRecHit( lp, le, detIdObject, &*clustIt) );
0245         SiPixelRecHit hit(lp, le, rqw, *genericDet, cluster);
0246         //
0247         // Now save it =================
0248         recHitsOnDetUnit.push_back(hit);
0249         // =============================
0250 
0251         // std::cout << "SiPixelRecHitConverterVI " << numberOfClusters << ' '<< lp << " " << le << std::endl;
0252       }  //  <-- End loop on Clusters
0253 
0254       //  LogDebug("SiPixelRecHitConverter")
0255       //std::cout << "SiPixelRecHitConverterVI "
0256       //    << " Found " << recHitsOnDetUnit.size() << " RecHits on " << detid //;
0257       //    << std::endl;
0258 
0259     }  //    <-- End loop on DetUnits
0260 
0261     //    LogDebug ("SiPixelRecHitConverter")
0262     //  std::cout << "SiPixelRecHitConverterVI "
0263     //  << cpeName_ << " converted " << numberOfClusters
0264     //  << " SiPixelClusters into SiPixelRecHits, in "
0265     //  << numberOfDetUnits << " DetUnits." //;
0266     //  << std::endl;
0267   }
0268 }  // end of namespace cms
0269 
0270 using cms::SiPixelRecHitConverter;
0271 
0272 DEFINE_FWK_MODULE(SiPixelRecHitConverter);