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
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
|
#include <string>
#include <iostream>
#include <map>
#include "FWCore/Framework/interface/one/EDAnalyzer.h"
#include "FWCore/Framework/interface/ESHandle.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "CondFormats/SiPixelObjects/interface/SiPixelQuality.h"
#include "CondFormats/DataRecord/interface/SiPixelQualityFromDbRcd.h"
#include "CondFormats/SiPixelObjects/interface/SiPixelFEDChannelContainer.h"
#include "CondFormats/DataRecord/interface/SiPixelStatusScenariosRcd.h"
#include "CondFormats/SiPixelObjects/interface/SiPixelQualityProbabilities.h"
#include "CondFormats/DataRecord/interface/SiPixelStatusScenarioProbabilityRcd.h"
class SiPixelBadFEDChannelSimulationSanityChecker : public edm::one::EDAnalyzer<> {
public:
explicit SiPixelBadFEDChannelSimulationSanityChecker(edm::ParameterSet const& p);
~SiPixelBadFEDChannelSimulationSanityChecker() override;
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
private:
void analyze(const edm::Event& e, const edm::EventSetup& c) override;
// ----------member data ---------------------------
const edm::ESGetToken<SiPixelFEDChannelContainer, SiPixelStatusScenariosRcd> siPixelBadFEDChToken_;
const edm::ESGetToken<SiPixelQualityProbabilities, SiPixelStatusScenarioProbabilityRcd> siPixelQPToken_;
const bool printdebug_;
};
SiPixelBadFEDChannelSimulationSanityChecker::SiPixelBadFEDChannelSimulationSanityChecker(edm::ParameterSet const& p)
: siPixelBadFEDChToken_(esConsumes()),
siPixelQPToken_(esConsumes()),
printdebug_(p.getUntrackedParameter<bool>("printDebug", true)) {
edm::LogInfo("SiPixelBadFEDChannelSimulationSanityChecker")
<< "SiPixelBadFEDChannelSimulationSanityChecker" << std::endl;
}
SiPixelBadFEDChannelSimulationSanityChecker::~SiPixelBadFEDChannelSimulationSanityChecker() {
edm::LogInfo("SiPixelBadFEDChannelSimulationSanityChecker")
<< "~SiPixelBadFEDChannelSimulationSanityChecker " << std::endl;
}
void SiPixelBadFEDChannelSimulationSanityChecker::analyze(const edm::Event& e, const edm::EventSetup& context) {
edm::LogInfo("SiPixelBadFEDChannelSimulationSanityChecker")
<< "### SiPixelBadFEDChannelSimulationSanityChecker::analyze ###" << std::endl;
edm::LogInfo("SiPixelBadFEDChannelSimulationSanityChecker") << " I AM IN RUN NUMBER " << e.id().run() << std::endl;
edm::LogInfo("SiPixelBadFEDChannelSimulationSanityChecker") << " ---EVENT NUMBER " << e.id().event() << std::endl;
edm::eventsetup::EventSetupRecordKey recordKey(
edm::eventsetup::EventSetupRecordKey::TypeTag::findType("SiPixelStatusScenariosRcd"));
if (recordKey.type() == edm::eventsetup::EventSetupRecordKey::TypeTag()) {
//record not found
edm::LogInfo("SiPixelBadFEDChannelSimulationSanityChecker") << "Record \"SiPixelStatusScenariosRcd"
<< "\" does not exist " << std::endl;
}
//this part gets the handle of the event source and the record (i.e. the Database)
const SiPixelFEDChannelContainer* quality_map = &context.getData(siPixelBadFEDChToken_);
edm::eventsetup::EventSetupRecordKey recordKey2(
edm::eventsetup::EventSetupRecordKey::TypeTag::findType("SiPixelStatusScenarioProbabilityRcd>"));
if (recordKey2.type() == edm::eventsetup::EventSetupRecordKey::TypeTag()) {
//record not found
edm::LogWarning("SiPixelQualityProbabilitiesTestReader") << "Record \"SiPixelStatusScenarioProbabilityRcd>"
<< "\" does not exist " << std::endl;
}
//this part gets the handle of the event source and the record (i.e. the Database)
const SiPixelQualityProbabilities* myProbabilities = &context.getData(siPixelQPToken_);
SiPixelFEDChannelContainer::SiPixelBadFEDChannelsScenarioMap m_qualities = quality_map->getScenarioMap();
SiPixelQualityProbabilities::probabilityMap m_probabilities = myProbabilities->getProbability_Map();
std::vector<std::string> allScenarios = quality_map->getScenarioList();
std::vector<std::string> allScenariosInProb;
for (auto it = m_probabilities.begin(); it != m_probabilities.end(); ++it) {
//int PUbin = it->first;
for (const auto& entry : it->second) {
auto scenario = entry.first;
auto probability = entry.second;
if (probability != 0) {
if (std::find(allScenariosInProb.begin(), allScenariosInProb.end(), scenario) == allScenariosInProb.end()) {
allScenariosInProb.push_back(scenario);
}
// if(m_qualities.find(scenario) == m_qualities.end()){
// edm::LogWarning("SiPixelBadFEDChannelSimulationSanityChecker") <<"Pretty worrying! the scenario: " << scenario << " (prob:" << probability << ") is not found in the map!!"<<std::endl;
// } else {
// edm::LogInfo("SiPixelBadFEDChannelSimulationSanityChecker") << "scenario: "<< scenario << " is in the map! (all is good)"<< std::endl;
// }
} // if prob!=0
} // loop on the scenarios for that PU bin
} // loop on PU bins
std::vector<std::string> notFound;
std::copy_if(allScenariosInProb.begin(),
allScenariosInProb.end(),
std::back_inserter(notFound),
[&allScenarios](const std::string& arg) {
return (std::find(allScenarios.begin(), allScenarios.end(), arg) == allScenarios.end());
});
if (!notFound.empty()) {
for (const auto& entry : notFound) {
edm::LogWarning("SiPixelBadFEDChannelSimulationSanityChecker")
<< "Pretty worrying! the scenario: " << entry << " is not found in the map!!" << std::endl;
if (printdebug_) {
edm::LogVerbatim("SiPixelBadFEDChannelSimulationSanityChecker") << " This scenario is found in: " << std::endl;
for (auto it = m_probabilities.begin(); it != m_probabilities.end(); ++it) {
int PUbin = it->first;
for (const auto& pair : it->second) {
if (pair.first == entry) {
edm::LogVerbatim("SiPixelBadFEDChannelSimulationSanityChecker")
<< " - PU bin " << PUbin << " with probability: " << pair.second << std::endl;
} // if the scenario matches
} // loop on scenarios
} // loop on PU
} // if printdebug
edm::LogVerbatim("SiPixelBadFEDChannelSimulationSanityChecker")
<< "==============================================" << std::endl;
} // loop on scenarios not found
edm::LogWarning("SiPixelBadFEDChannelSimulationSanityChecker")
<< " ====> A total of " << notFound.size() << " scenarios are not found in the map!" << std::endl;
} else {
edm::LogVerbatim("SiPixelBadFEDChannelSimulationSanityChecker")
<< "=================================================================================" << std::endl;
edm::LogInfo("SiPixelBadFEDChannelSimulationSanityChecker")
<< " All scenarios in probability record are found in the scenario map, (all is good)!" << std::endl;
edm::LogVerbatim("SiPixelBadFEDChannelSimulationSanityChecker")
<< "=================================================================================" << std::endl;
}
}
void SiPixelBadFEDChannelSimulationSanityChecker::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.setComment("Tries sanity of Pixel Stuck TBM simulation");
desc.addUntracked<bool>("printDebug", true);
descriptions.add("SiPixelBadFEDChannelSimulationSanityChecker", desc);
}
DEFINE_FWK_MODULE(SiPixelBadFEDChannelSimulationSanityChecker);
|