Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-04-29 02:32:48

0001 #include "RecoLocalTracker/SiPixelRecHits/interface/PixelCPETemplateReco.h"
0002 #include "RecoLocalTracker/Records/interface/TkPixelCPERecord.h"
0003 #include "RecoLocalTracker/ClusterParameterEstimator/interface/PixelClusterParameterEstimator.h"
0004 #include "MagneticField/Engine/interface/MagneticField.h"
0005 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0006 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0007 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0008 #include "Geometry/Records/interface/TrackerTopologyRcd.h"
0009 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0010 
0011 #include "FWCore/Framework/interface/EventSetup.h"
0012 #include "FWCore/Framework/interface/ESHandle.h"
0013 #include "FWCore/Framework/interface/ModuleFactory.h"
0014 #include "FWCore/Framework/interface/ESProducer.h"
0015 
0016 #include <string>
0017 #include <memory>
0018 
0019 class PixelCPETemplateRecoESProducer : public edm::ESProducer {
0020 public:
0021   PixelCPETemplateRecoESProducer(const edm::ParameterSet& p);
0022   std::unique_ptr<PixelClusterParameterEstimator> produce(const TkPixelCPERecord&);
0023   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0024 
0025 private:
0026   edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> magfieldToken_;
0027   edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> pDDToken_;
0028   edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> hTTToken_;
0029   edm::ESGetToken<SiPixelLorentzAngle, SiPixelLorentzAngleRcd> lorentzAngleToken_;
0030   edm::ESGetToken<SiPixelTemplateDBObject, SiPixelTemplateDBObjectESProducerRcd> templateDBobjectToken_;
0031 
0032   edm::ParameterSet pset_;
0033   bool doLorentzFromAlignment_;
0034   bool useLAFromDB_;
0035 };
0036 
0037 using namespace edm;
0038 
0039 PixelCPETemplateRecoESProducer::PixelCPETemplateRecoESProducer(const edm::ParameterSet& p) {
0040   std::string myname = p.getParameter<std::string>("ComponentName");
0041 
0042   useLAFromDB_ = p.getParameter<bool>("useLAFromDB");
0043   doLorentzFromAlignment_ = p.getParameter<bool>("doLorentzFromAlignment");
0044 
0045   pset_ = p;
0046   auto c = setWhatProduced(this, myname);
0047   magfieldToken_ = c.consumes();
0048   pDDToken_ = c.consumes();
0049   hTTToken_ = c.consumes();
0050   templateDBobjectToken_ = c.consumes();
0051   if (useLAFromDB_ || doLorentzFromAlignment_) {
0052     char const* laLabel = doLorentzFromAlignment_ ? "fromAlignment" : "";
0053     lorentzAngleToken_ = c.consumes(edm::ESInputTag("", laLabel));
0054   }
0055 }
0056 
0057 std::unique_ptr<PixelClusterParameterEstimator> PixelCPETemplateRecoESProducer::produce(
0058     const TkPixelCPERecord& iRecord) {
0059   // Normal, default LA is used in case of template failure, load it unless
0060   // turned off
0061   // if turned off, null is ok, becomes zero
0062   const SiPixelLorentzAngle* lorentzAngleProduct = nullptr;
0063   if (useLAFromDB_ || doLorentzFromAlignment_) {
0064     lorentzAngleProduct = &iRecord.get(lorentzAngleToken_);
0065   }
0066 
0067   return std::make_unique<PixelCPETemplateReco>(pset_,
0068                                                 &iRecord.get(magfieldToken_),
0069                                                 iRecord.get(pDDToken_),
0070                                                 iRecord.get(hTTToken_),
0071                                                 lorentzAngleProduct,
0072                                                 &iRecord.get(templateDBobjectToken_));
0073 }
0074 
0075 void PixelCPETemplateRecoESProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0076   edm::ParameterSetDescription desc;
0077 
0078   // from PixelCPEBase
0079   PixelCPEBase::fillPSetDescription(desc);
0080 
0081   // from PixelCPETemplateReco
0082   PixelCPETemplateReco::fillPSetDescription(desc);
0083 
0084   // specific to PixelCPETemplateRecoESProducer
0085   desc.add<std::string>("ComponentName", "PixelCPETemplateReco");
0086   descriptions.add("_templates_default", desc);
0087 }
0088 
0089 DEFINE_FWK_EVENTSETUP_MODULE(PixelCPETemplateRecoESProducer);