File indexing completed on 2024-04-06 12:03:24
0001
0002 #include <cstdio>
0003 #include <iostream>
0004 #include <sys/time.h>
0005
0006
0007 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0008 #include "CondFormats/DataRecord/interface/SiPixelLorentzAngleRcd.h"
0009 #include "CondFormats/DataRecord/interface/SiPixelLorentzAngleSimRcd.h"
0010 #include "CondFormats/SiPixelObjects/interface/SiPixelLorentzAngle.h"
0011 #include "DataFormats/DetId/interface/DetId.h"
0012 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0013 #include "DataFormats/TrackerCommon/interface/PixelBarrelName.h"
0014 #include "DataFormats/TrackerCommon/interface/PixelEndcapName.h"
0015 #include "FWCore/Framework/interface/Event.h"
0016 #include "FWCore/Framework/interface/EventSetup.h"
0017 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0019 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0020 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0021 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0022 #include "FWCore/ServiceRegistry/interface/Service.h"
0023 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0024 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0025
0026
0027 #include "TROOT.h"
0028 #include "TFile.h"
0029 #include "TH2F.h"
0030
0031
0032
0033
0034
0035 class SiPixelLorentzAngleReader : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0036 public:
0037 explicit SiPixelLorentzAngleReader(const edm::ParameterSet&);
0038 ~SiPixelLorentzAngleReader() override;
0039
0040 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0041 void analyze(const edm::Event&, const edm::EventSetup&) override;
0042
0043 private:
0044 const edm::ESGetToken<SiPixelLorentzAngle, SiPixelLorentzAngleRcd> siPixelLAToken_;
0045 const edm::ESGetToken<SiPixelLorentzAngle, SiPixelLorentzAngleSimRcd> siPixelSimLAToken_;
0046 const std::string siPixelLALabel_;
0047 const std::string siPixelSimLALabel_;
0048 const uint32_t printdebug_;
0049 const bool useSimRcd_;
0050 TH1F* LorentzAngleBarrel_;
0051 TH1F* LorentzAngleForward_;
0052 };
0053
0054 using namespace cms;
0055
0056 SiPixelLorentzAngleReader::SiPixelLorentzAngleReader(const edm::ParameterSet& iConfig)
0057 : siPixelLAToken_(esConsumes(edm::ESInputTag("", iConfig.getParameter<std::string>("recoLabel")))),
0058 siPixelSimLAToken_(esConsumes((edm::ESInputTag("", iConfig.getParameter<std::string>("simLabel"))))),
0059 printdebug_(iConfig.getUntrackedParameter<uint32_t>("printDebug", 10)),
0060 useSimRcd_(iConfig.getParameter<bool>("useSimRcd")) {
0061 usesResource(TFileService::kSharedResource);
0062 }
0063
0064 SiPixelLorentzAngleReader::~SiPixelLorentzAngleReader() = default;
0065
0066 void SiPixelLorentzAngleReader::analyze(const edm::Event& e, const edm::EventSetup& iSetup) {
0067 const SiPixelLorentzAngle* SiPixelLorentzAngle_;
0068 if (useSimRcd_) {
0069 SiPixelLorentzAngle_ = &iSetup.getData(siPixelSimLAToken_);
0070 } else {
0071 SiPixelLorentzAngle_ = &iSetup.getData(siPixelLAToken_);
0072 }
0073
0074 edm::LogInfo("SiPixelLorentzAngleReader")
0075 << "[SiPixelLorentzAngleReader::analyze] End Reading SiPixelLorentzAngle" << std::endl;
0076 edm::Service<TFileService> fs;
0077 LorentzAngleBarrel_ = fs->make<TH1F>("LorentzAngleBarrelPixel", "LorentzAngleBarrelPixel", 150, 0, 0.15);
0078 LorentzAngleForward_ = fs->make<TH1F>("LorentzAngleForwardPixel", "LorentzAngleForwardPixel", 150, 0, 0.15);
0079 std::map<unsigned int, float> detid_la = SiPixelLorentzAngle_->getLorentzAngles();
0080 std::map<unsigned int, float>::const_iterator it;
0081 unsigned int count = 0;
0082 for (it = detid_la.begin(); it != detid_la.end(); it++) {
0083 count++;
0084 if (count <= printdebug_) {
0085 edm::LogPrint("SiPixelLorentzAngleReader") << "detid " << it->first << " \t"
0086 << " Lorentz angle " << it->second << std::endl;
0087 }
0088 unsigned int subdet = DetId(it->first).subdetId();
0089 if (subdet == static_cast<int>(PixelSubdetector::PixelBarrel)) {
0090 LorentzAngleBarrel_->Fill(it->second);
0091 } else if (subdet == static_cast<int>(PixelSubdetector::PixelEndcap)) {
0092 LorentzAngleForward_->Fill(it->second);
0093 }
0094 }
0095 edm::LogPrint("SiPixelLorentzAngleReader")
0096 << "SiPixelLorentzAngleReader::" << __FUNCTION__ << "(...) :examined " << count << " DetIds" << std::endl;
0097 }
0098
0099 void SiPixelLorentzAngleReader::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0100 edm::ParameterSetDescription desc;
0101 desc.setComment("EDAnalyzer to read per-module SiPixelLorentzAngle payloads in the EventSetup");
0102 desc.add<std::string>("recoLabel", "")->setComment("label for the reconstruction tags");
0103 desc.add<std::string>("simLabel", "")->setComment("label for the simulation tags");
0104 desc.addUntracked<unsigned int>("printDebug", 10);
0105 desc.add<bool>("useSimRcd", false);
0106 descriptions.addWithDefaultLabel(desc);
0107 }
0108
0109 #include "FWCore/Framework/interface/MakerMacros.h"
0110 #include "FWCore/Framework/interface/Frameworkfwd.h"
0111 DEFINE_FWK_MODULE(SiPixelLorentzAngleReader);