1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
|
// system include files
//#include <memory>
#include <iostream>
#include <cstdio>
#include <sys/time.h>
// user include files
#include "CondFormats/DataRecord/interface/SiStripPedestalsRcd.h"
#include "CondFormats/SiStripObjects/interface/SiStripPedestals.h"
#include "FWCore/Framework/interface/ESHandle.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/one/EDAnalyzer.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ServiceRegistry/interface/Service.h"
class SiStripPedestalsReader : public edm::one::EDAnalyzer<> {
public:
explicit SiStripPedestalsReader(const edm::ParameterSet&);
~SiStripPedestalsReader() override = default;
void analyze(const edm::Event&, const edm::EventSetup&) override;
private:
uint32_t printdebug_;
const edm::ESGetToken<SiStripPedestals, SiStripPedestalsRcd> pedestalsToken_;
};
using namespace std;
using namespace cms;
SiStripPedestalsReader::SiStripPedestalsReader(const edm::ParameterSet& iConfig)
: printdebug_(iConfig.getUntrackedParameter<uint32_t>("printDebug", 1)), pedestalsToken_(esConsumes()) {}
void SiStripPedestalsReader::analyze(const edm::Event& e, const edm::EventSetup& iSetup) {
const auto& pedestals = iSetup.getData(pedestalsToken_);
edm::LogInfo("SiStripPedestalsReader") << "[SiStripPedestalsReader::analyze] End Reading SiStripPedestals"
<< std::endl;
std::vector<uint32_t> detid;
pedestals.getDetIds(detid);
edm::LogInfo("Number of detids ") << detid.size() << std::endl;
if (printdebug_)
for (size_t id = 0; id < detid.size() && id < printdebug_; id++) {
SiStripPedestals::Range range = pedestals.getRange(detid[id]);
int strip = 0;
for (int it = 0; it < (range.second - range.first) * 8 / 10; it++) {
edm::LogInfo("SiStripPedestalsReader")
<< "detid " << detid[id] << " \t"
<< " strip " << strip++ << " \t" << pedestals.getPed(it, range) << " \t" << std::endl;
}
}
}
#include "FWCore/PluginManager/interface/ModuleDef.h"
#include "FWCore/Framework/interface/MakerMacros.h"
DEFINE_FWK_MODULE(SiStripPedestalsReader);
|