CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SiPixelLorentzAngleReader.cc
Go to the documentation of this file.
1 // system include files
2 #include <cstdio>
3 #include <iostream>
4 #include <sys/time.h>
5 
6 // user include files
25 
26 // ROOT includes
27 #include "TROOT.h"
28 #include "TFile.h"
29 #include "TH2F.h"
30 
31 //
32 //
33 // class decleration
34 //
35 class SiPixelLorentzAngleReader : public edm::one::EDAnalyzer<edm::one::SharedResources> {
36 public:
38  ~SiPixelLorentzAngleReader() override;
39 
40  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
41  void analyze(const edm::Event&, const edm::EventSetup&) override;
42 
43 private:
48  const uint32_t printdebug_;
49  const bool useSimRcd_;
52 };
53 
54 using namespace cms;
55 
57  : siPixelLAToken_(esConsumes(edm::ESInputTag("", iConfig.getParameter<std::string>("recoLabel")))),
58  siPixelSimLAToken_(esConsumes((edm::ESInputTag("", iConfig.getParameter<std::string>("simLabel"))))),
59  printdebug_(iConfig.getUntrackedParameter<uint32_t>("printDebug", 10)),
60  useSimRcd_(iConfig.getParameter<bool>("useSimRcd")) {
61  usesResource(TFileService::kSharedResource);
62 }
63 
65 
67  const SiPixelLorentzAngle* SiPixelLorentzAngle_;
68  if (useSimRcd_) {
69  SiPixelLorentzAngle_ = &iSetup.getData(siPixelSimLAToken_);
70  } else {
71  SiPixelLorentzAngle_ = &iSetup.getData(siPixelLAToken_);
72  }
73 
74  edm::LogInfo("SiPixelLorentzAngleReader")
75  << "[SiPixelLorentzAngleReader::analyze] End Reading SiPixelLorentzAngle" << std::endl;
77  LorentzAngleBarrel_ = fs->make<TH1F>("LorentzAngleBarrelPixel", "LorentzAngleBarrelPixel", 150, 0, 0.15);
78  LorentzAngleForward_ = fs->make<TH1F>("LorentzAngleForwardPixel", "LorentzAngleForwardPixel", 150, 0, 0.15);
79  std::map<unsigned int, float> detid_la = SiPixelLorentzAngle_->getLorentzAngles();
80  std::map<unsigned int, float>::const_iterator it;
81  unsigned int count = 0;
82  for (it = detid_la.begin(); it != detid_la.end(); it++) {
83  count++;
84  if (count <= printdebug_) {
85  edm::LogPrint("SiPixelLorentzAngleReader") << "detid " << it->first << " \t"
86  << " Lorentz angle " << it->second << std::endl;
87  }
88  unsigned int subdet = DetId(it->first).subdetId();
89  if (subdet == static_cast<int>(PixelSubdetector::PixelBarrel)) {
90  LorentzAngleBarrel_->Fill(it->second);
91  } else if (subdet == static_cast<int>(PixelSubdetector::PixelEndcap)) {
92  LorentzAngleForward_->Fill(it->second);
93  }
94  }
95  edm::LogPrint("SiPixelLorentzAngleReader")
96  << "SiPixelLorentzAngleReader::" << __FUNCTION__ << "(...) :examined " << count << " DetIds" << std::endl;
97 }
98 
101  desc.setComment("EDAnalyzer to read per-module SiPixelLorentzAngle payloads in the EventSetup");
102  desc.add<std::string>("recoLabel", "")->setComment("label for the reconstruction tags");
103  desc.add<std::string>("simLabel", "")->setComment("label for the simulation tags");
104  desc.addUntracked<unsigned int>("printDebug", 10);
105  desc.add<bool>("useSimRcd", false);
106  descriptions.addWithDefaultLabel(desc);
107 }
108 
static const std::string kSharedResource
Definition: TFileService.h:76
const std::map< unsigned int, float > & getLorentzAngles() const
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
~SiPixelLorentzAngleReader() override
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
SiPixelLorentzAngleReader(const edm::ParameterSet &)
bool getData(T &iHolder) const
Definition: EventSetup.h:122
void setComment(std::string const &value)
const edm::ESGetToken< SiPixelLorentzAngle, SiPixelLorentzAngleSimRcd > siPixelSimLAToken_
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Log< level::Warning, true > LogPrint
Log< level::Info, false > LogInfo
Definition: DetId.h:17
void analyze(const edm::Event &, const edm::EventSetup &) override
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const edm::ESGetToken< SiPixelLorentzAngle, SiPixelLorentzAngleRcd > siPixelLAToken_
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283