CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SiStripApvGainReader.cc
Go to the documentation of this file.
1 // system include files
2 #include <iostream>
3 #include <cstdio>
4 #include <sys/time.h>
5 
6 // user include files
19 
20 // root objects
21 #include "TROOT.h"
22 #include "TSystem.h"
23 #include "TFile.h"
24 #include "TDirectory.h"
25 #include "TTree.h"
26 
27 class SiStripApvGainReader : public edm::one::EDAnalyzer<edm::one::SharedResources> {
28 public:
29  explicit SiStripApvGainReader(const edm::ParameterSet&);
30  ~SiStripApvGainReader() override = default;
31 
32  void analyze(const edm::Event&, const edm::EventSetup&) override;
33 
34 private:
35  // initializers list
36  const bool printdebug_;
38  const uint32_t gainType_;
41  TTree* tree_ = nullptr;
42  int id_ = 0, detId_ = 0, apvId_ = 0;
43  double gain_ = 0;
44 };
45 
46 using namespace cms;
47 
49  : printdebug_(iConfig.getUntrackedParameter<bool>("printDebug", true)),
50  formatedOutput_(iConfig.getUntrackedParameter<std::string>("outputFile", "")),
51  gainType_(iConfig.getUntrackedParameter<uint32_t>("gainType", 1)),
53  usesResource(TFileService::kSharedResource);
54 
55  if (fs_.isAvailable()) {
56  tree_ = fs_->make<TTree>("Gains", "Gains");
57 
58  tree_->Branch("Index", &id_, "Index/I");
59  tree_->Branch("DetId", &detId_, "DetId/I");
60  tree_->Branch("APVId", &apvId_, "APVId/I");
61  tree_->Branch("Gain", &gain_, "Gain/D");
62  }
63 }
64 
66  const auto& stripApvGain = iSetup.getData(gainToken_);
67  edm::LogInfo("SiStripApvGainReader") << "[SiStripApvGainReader::analyze] End Reading SiStripApvGain" << std::endl;
68  std::vector<uint32_t> detid;
69  stripApvGain.getDetIds(detid);
70  edm::LogInfo("Number of detids ") << detid.size() << std::endl;
71 
72  FILE* pFile = nullptr;
73  if (!formatedOutput_.empty())
74  pFile = fopen(formatedOutput_.c_str(), "w");
75  for (size_t id = 0; id < detid.size(); id++) {
76  SiStripApvGain::Range range = stripApvGain.getRange(detid[id], gainType_);
77  if (printdebug_) {
78  int apv = 0;
79  for (int it = 0; it < range.second - range.first; it++) {
80  edm::LogInfo("SiStripApvGainReader")
81  << "detid " << detid[id] << " \t " << apv++ << " \t " << stripApvGain.getApvGain(it, range) << std::endl;
82  id_++;
83 
84  if (tree_) {
85  detId_ = detid[id];
86  apvId_ = apv;
87  gain_ = stripApvGain.getApvGain(it, range);
88  tree_->Fill();
89  }
90  }
91  }
92 
93  if (pFile) {
94  fprintf(pFile, "%i ", detid[id]);
95  for (int it = 0; it < range.second - range.first; it++) {
96  fprintf(pFile, "%f ", stripApvGain.getApvGain(it, range));
97  }
98  fprintf(pFile, "\n");
99  }
100  }
101 
102  if (pFile)
103  fclose(pFile);
104 }
105 
108 
static const std::string kSharedResource
Definition: TFileService.h:76
uint16_t *__restrict__ id
gainToken_(esConsumes())
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
~SiStripApvGainReader() override=default
const uint16_t range(const Frame &aFrame)
bool getData(T &iHolder) const
Definition: EventSetup.h:128
const std::string formatedOutput_
bool isAvailable() const
Definition: Service.h:40
std::pair< ContainerIterator, ContainerIterator > Range
Log< level::Info, false > LogInfo
edm::Service< TFileService > fs_
void analyze(const edm::Event &, const edm::EventSetup &) override
SiStripApvGainReader(const edm::ParameterSet &)
const edm::ESGetToken< SiStripGain, SiStripGainRcd > gainToken_
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283