Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:19:39

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