CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SiStripNoisesFromDBMiscalibrator.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: CondTools/SiStrip
4 // Class: SiStripNoisesFromDBMiscalibrator
5 //
14 //
15 // Original Author: Marco Musich
16 // Created: Tue, 03 Oct 2017 12:57:34 GMT
17 //
18 //
19 
20 // system include files
21 #include <memory>
22 #include <iostream>
23 #include <fstream>
24 
25 // user include files
26 #include "CLHEP/Random/RandGauss.h"
46 
47 //
48 // class declaration
49 //
50 
52 public:
54  ~SiStripNoisesFromDBMiscalibrator() override = default;
55 
56  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
57 
58 private:
59  void analyze(const edm::Event&, const edm::EventSetup&) override;
60  SiStripNoises getNewObject(const std::map<std::pair<uint32_t, int>, float>& theMap);
61  SiStripNoises getNewObject_withDefaults(const std::map<std::pair<uint32_t, int>, float>& theMap,
62  const float theDefault);
63  void endJob() override;
64 
65  // ----------member data ---------------------------
66  const uint32_t m_printdebug;
67  const bool m_fillDefaults;
68  const bool m_saveMaps;
69  const std::vector<edm::ParameterSet> m_parameters;
73 
74  std::unique_ptr<TrackerMap> scale_map;
75  std::unique_ptr<TrackerMap> smear_map;
76  std::unique_ptr<TrackerMap> ratio_map;
77  std::unique_ptr<TrackerMap> old_payload_map;
78  std::unique_ptr<TrackerMap> new_payload_map;
79  std::unique_ptr<TrackerMap> missing_map;
80 };
81 
82 //
83 // constructors and destructor
84 //
86  : m_printdebug{iConfig.getUntrackedParameter<uint32_t>("printDebug", 1)},
87  m_fillDefaults{iConfig.getUntrackedParameter<bool>("fillDefaults", false)},
88  m_saveMaps{iConfig.getUntrackedParameter<bool>("saveMaps", true)},
89  m_parameters{iConfig.getParameter<std::vector<edm::ParameterSet> >("params")},
94  //now do what ever initialization is needed
95 
96  scale_map = std::make_unique<TrackerMap>("scale");
97  scale_map->setTitle("Tracker Map of Scale factor averaged by module");
98  scale_map->setPalette(1);
99 
100  smear_map = std::make_unique<TrackerMap>("smear");
101  smear_map->setTitle("Tracker Map of Smear factor averaged by module");
102  smear_map->setPalette(1);
103 
104  old_payload_map = std::make_unique<TrackerMap>("old_payload");
105  old_payload_map->setTitle("Tracker Map of Starting Noise Payload averaged by module");
106  old_payload_map->setPalette(1);
107 
108  new_payload_map = std::make_unique<TrackerMap>("new_payload");
109  new_payload_map->setTitle("Tracker Map of Modified Noise Payload averaged by module");
110  new_payload_map->setPalette(1);
111 
112  ratio_map = std::make_unique<TrackerMap>("ratio");
113  ratio_map->setTitle("Tracker Map of Average by module of the payload ratio (new/old)");
114  ratio_map->setPalette(1);
115 
116  if (m_fillDefaults) {
117  missing_map = std::make_unique<TrackerMap>("uncabled");
118  missing_map->setTitle("Tracker Map of uncabled modules");
119  missing_map->setPalette(1);
120  }
121 }
122 
123 //
124 // member functions
125 //
126 
127 // ------------ method called for each event ------------
129  using namespace edm;
130 
131  const auto tTopo = &iSetup.getData(m_tTopoToken);
132 
133  std::vector<std::string> partitions;
134 
135  // fill the list of partitions
136  for (auto& thePSet : m_parameters) {
137  const std::string partition(thePSet.getParameter<std::string>("partition"));
138  // only if it is not yet in the list
139  if (std::find(partitions.begin(), partitions.end(), partition) == partitions.end()) {
140  partitions.push_back(partition);
141  }
142  }
143 
144  std::map<sistripsummary::TrackerRegion, SiStripMiscalibrate::Smearings> mapOfSmearings;
145 
146  for (auto& thePSet : m_parameters) {
147  const std::string partition(thePSet.getParameter<std::string>("partition"));
149 
150  bool m_doScale(thePSet.getParameter<bool>("doScale"));
151  bool m_doSmear(thePSet.getParameter<bool>("doSmear"));
152  double m_scaleFactor(thePSet.getParameter<double>("scaleFactor"));
153  double m_smearFactor(thePSet.getParameter<double>("smearFactor"));
154 
156  params.setSmearing(m_doScale, m_doSmear, m_scaleFactor, m_smearFactor);
157  mapOfSmearings[region] = params;
158  }
159 
160  const auto& stripNoises = iSetup.getData(m_noiseToken);
161 
162  std::map<std::pair<uint32_t, int>, float> theMap, oldPayloadMap;
163 
164  std::vector<uint32_t> detid;
165  stripNoises.getDetIds(detid);
166  for (const auto& d : detid) {
167  SiStripNoises::Range range = stripNoises.getRange(d);
168 
170 
171  // sort by largest to smallest
172  std::sort(regions.rbegin(), regions.rend());
173 
175 
176  for (unsigned int j = 0; j < regions.size(); j++) {
177  bool checkRegion = (mapOfSmearings.count(regions[j]) != 0);
178 
179  if (!checkRegion) {
180  // if the subdetector is not in the list and there's no indication for the whole tracker, just use the default
181  // i.e. no change
182  continue;
183  } else {
184  params = mapOfSmearings[regions[j]];
185  break;
186  }
187  }
188 
189  scale_map->fill(d, params.m_scaleFactor);
190  smear_map->fill(d, params.m_smearFactor);
191 
192  int nStrips = 0;
193  for (int it = 0; it < (range.second - range.first) * 8 / 9; ++it) {
194  auto noise = stripNoises.getNoise(it, range);
195  std::pair<uint32_t, int> index = std::make_pair(d, nStrips);
196 
197  oldPayloadMap[index] = noise;
198 
199  if (params.m_doScale) {
200  noise *= params.m_scaleFactor;
201  }
202 
203  if (params.m_doSmear) {
204  float smearedNoise = CLHEP::RandGauss::shoot(noise, params.m_smearFactor);
205  noise = smearedNoise;
206  }
207 
208  theMap[index] = noise;
209 
210  nStrips += 1;
211 
212  } // loop over APVs
213  } // loop over DetIds
214 
215  SiStripNoises theSiStripNoises{};
216  if (!m_fillDefaults) {
217  theSiStripNoises = this->getNewObject(theMap);
218  } else {
219  theSiStripNoises = this->getNewObject_withDefaults(theMap, -1.);
220  }
221 
222  // make the payload ratio map
223  uint32_t cachedId(0);
224  SiStripMiscalibrate::Entry noise_ratio;
227  unsigned int countDetIds(0); // count DetIds to print
228  for (const auto& element : theMap) {
229  uint32_t DetId = element.first.first;
230  int nstrip = element.first.second;
231  float new_noise = element.second;
232  float old_noise = oldPayloadMap[std::make_pair(DetId, nstrip)];
233 
234  // flush the counters
235  if (cachedId != 0 && DetId != cachedId) {
236  ratio_map->fill(cachedId, noise_ratio.mean());
237  old_payload_map->fill(cachedId, o_noise.mean());
238  new_payload_map->fill(cachedId, n_noise.mean());
239 
240  //auto test = new_payload_map.get()->smoduleMap;
241 
242  noise_ratio.reset();
243  o_noise.reset();
244  n_noise.reset();
245  countDetIds++;
246  }
247 
248  // printout for debug
249  if (countDetIds < m_printdebug) {
250  edm::LogPrint("SiStripNoisesFromDBMiscalibrator")
251  << "SiStripNoisesFromDBMiscalibrator"
252  << "::" << __FUNCTION__ << " detid " << DetId << " \t"
253  << " strip " << nstrip << " \t new noise: " << std::setw(5) << std::setprecision(2) << new_noise
254  << " \t old noise: " << old_noise << " \t" << std::endl;
255  }
256 
257  cachedId = DetId;
258  noise_ratio.add(new_noise / old_noise);
259  o_noise.add(old_noise);
260  n_noise.add(new_noise);
261  }
262 
263  // write out the SiStripNoises record
265 
266  if (poolDbService.isAvailable()) {
267  if (poolDbService->isNewTagRequest("SiStripNoisesRcd")) {
268  poolDbService->createOneIOV(theSiStripNoises, poolDbService->currentTime(), "SiStripNoisesRcd");
269  } else {
270  poolDbService->appendOneIOV(theSiStripNoises, poolDbService->currentTime(), "SiStripNoisesRcd");
271  }
272  } else {
273  throw std::runtime_error("PoolDBService required.");
274  }
275 }
276 
277 // ------------ method called once each job just after ending the event loop ------------
279  if (m_saveMaps) {
280  scale_map->save(true, 0, 0, "noise_scale_map.pdf");
281  scale_map->save(true, 0, 0, "noise_scale_map.png");
282 
283  smear_map->save(true, 0, 0, "noise_smear_map.pdf");
284  smear_map->save(true, 0, 0, "noise_smear_map.png");
285 
286  ratio_map->save(true, 0, 0, "noise_ratio_map.pdf");
287  ratio_map->save(true, 0, 0, "noise_ratio_map.png");
288 
290 
291  old_payload_map->save(true, range.first, range.second, "starting_noise_payload_map.pdf");
292  old_payload_map->save(true, range.first, range.second, "starting_noise_payload_map.png");
293 
295 
296  new_payload_map->save(true, range.first, range.second, "new_noise_payload_map.pdf");
297  new_payload_map->save(true, range.first, range.second, "new_noise_payload_map.png");
298 
299  if (m_fillDefaults) {
300  missing_map->save(true, 0, 0, "missing_map.pdf");
301  missing_map->save(true, 0, 0, "missing_map.png");
302  }
303  }
304 }
305 
306 //********************************************************************************//
308  const std::map<std::pair<uint32_t, int>, float>& theMap, const float theDefault) {
309  SiStripNoises obj{};
310 
311  std::vector<uint32_t> missingDetIds;
312 
314  const auto& DetInfos = reader.getAllData();
315 
316  for (const auto& it : DetInfos) {
317  const auto& nAPVs = it.second.nApvs;
318  //Generate Noise for det detid
319  bool isMissing(false);
320  SiStripNoises::InputVector theSiStripVector;
321  for (int t_strip = 0; t_strip < 128 * nAPVs; ++t_strip) {
322  std::pair<uint32_t, int> index = std::make_pair(it.first, t_strip);
323 
324  if (theMap.find(index) == theMap.end()) {
325  LogDebug("SiStripNoisesFromDBMiscalibrator") << "detid " << it.first << " \t"
326  << " strip " << t_strip << " \t"
327  << " not found" << std::endl;
328 
329  isMissing = true;
330  obj.setData(theDefault, theSiStripVector);
331 
332  } else {
333  float noise = theMap.at(index);
334  obj.setData(noise, theSiStripVector);
335  }
336  }
337 
338  if (isMissing)
339  missingDetIds.push_back(it.first);
340 
341  if (!obj.put(it.first, theSiStripVector)) {
342  edm::LogError("SiStripNoisesFromDBMiscalibrator")
343  << "[SiStripNoisesFromDBMiscalibrator::analyze] detid already exists" << std::endl;
344  }
345  }
346 
347  if (!missingDetIds.empty()) {
348  // open output file
349  std::stringstream name;
350  name << "missing_modules.txt";
351  std::ofstream* ofile = new std::ofstream(name.str(), std::ofstream::trunc);
352  if (!ofile->is_open())
353  throw "cannot open output file!";
354  for (const auto& missing : missingDetIds) {
355  edm::LogVerbatim("SiStripNoisesFromDBMiscalibrator") << missing << " " << 1 << std::endl;
356  (*ofile) << missing << " " << 1 << std::endl;
357  missing_map->fill(missing, 1);
358  }
359 
360  ofile->close();
361  delete ofile;
362  }
363 
364  return obj;
365 }
366 
367 //********************************************************************************//
368 SiStripNoises SiStripNoisesFromDBMiscalibrator::getNewObject(const std::map<std::pair<uint32_t, int>, float>& theMap) {
369  SiStripNoises obj{};
370 
371  uint32_t PreviousDetId = 0;
372  SiStripNoises::InputVector theSiStripVector;
373  for (const auto& element : theMap) {
374  uint32_t DetId = element.first.first;
375  float noise = element.second;
376 
377  if (DetId != PreviousDetId) {
378  if (!theSiStripVector.empty()) {
379  if (!obj.put(PreviousDetId, theSiStripVector)) {
380  edm::LogError("SiStripNoisesFromDBMiscalibrator")
381  << "[SiStripNoisesFromDBMiscalibrator::analyze] detid already exists" << std::endl;
382  }
383  }
384 
385  theSiStripVector.clear();
386  PreviousDetId = DetId;
387  }
388  obj.setData(noise, theSiStripVector);
389  }
390  return obj;
391 }
392 
393 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
396 
397  desc.setComment(
398  "Creates rescaled / smeared SiStrip Noise payload."
399  "PoolDBOutputService must be set up for 'SiSiStripNoisesRcd'.");
400 
401  edm::ParameterSetDescription descScaler;
402  descScaler.setComment(
403  "ParameterSet specifying the Strip tracker partition to be scaled / smeared "
404  "by a given factor.");
405 
406  descScaler.add<std::string>("partition", "Tracker");
407  descScaler.add<bool>("doScale", true);
408  descScaler.add<bool>("doSmear", true);
409  descScaler.add<double>("scaleFactor", 1.0);
410  descScaler.add<double>("smearFactor", 1.0);
411  desc.addVPSet("params", descScaler, std::vector<edm::ParameterSet>(1));
412 
413  desc.addUntracked<unsigned int>("printDebug", 1);
414  desc.addUntracked<bool>("fillDefaults", false);
415  desc.addUntracked<bool>("saveMaps", true);
416 
417  descriptions.add("scaleAndSmearSiStripNoises", desc);
418 }
419 
420 //define this as a plug-in
Log< level::Info, true > LogVerbatim
dictionary missing
Definition: combine.py:5
T getUntrackedParameter(std::string const &, T const &) const
ParameterDescriptionBase * addVPSet(U const &iLabel, ParameterSetDescription const &validator, std::vector< ParameterSet > const &defaults)
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
std::vector< sistripsummary::TrackerRegion > getRegionsFromDetId(const TrackerTopology *m_trackerTopo, DetId detid)
sistripsummary::TrackerRegion getRegionFromString(std::string region)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< uint16_t > InputVector
Definition: SiStripNoises.h:50
Log< level::Error, false > LogError
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
tuple nStrips
1.2 is to make the matching window safely the two nearest strips 0.35 is the size of an ME0 chamber i...
void createOneIOV(const T &payload, cond::Time_t firstSinceTime, const std::string &recordName)
const uint16_t range(const Frame &aFrame)
SiStripNoises getNewObject_withDefaults(const std::map< std::pair< uint32_t, int >, float > &theMap, const float theDefault)
bool getData(T &iHolder) const
Definition: EventSetup.h:128
tuple d
Definition: ztail.py:151
void setComment(std::string const &value)
int iEvent
Definition: GenABIO.cc:224
~SiStripNoisesFromDBMiscalibrator() override=default
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void analyze(const edm::Event &, const edm::EventSetup &) override
tuple reader
Definition: DQM.py:105
SiStripDetInfo read(std::string filePath)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Log< level::Warning, true > LogPrint
__shared__ int noise
Definition: DetId.h:17
void setSmearing(bool doScale, bool doSmear, double the_scaleFactor, double the_smearFactor)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::pair< float, float > getTruncatedRange(const TrackerMap *theMap)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const std::vector< edm::ParameterSet > m_parameters
SiStripNoises getNewObject(const std::map< std::pair< uint32_t, int >, float > &theMap)
SiStripNoisesFromDBMiscalibrator(const edm::ParameterSet &)
std::string fullPath() const
Definition: FileInPath.cc:161
static constexpr char const *const kDefaultFile
const edm::ESGetToken< SiStripNoises, SiStripNoisesRcd > m_noiseToken
std::pair< ContainerIterator, ContainerIterator > Range
Definition: SiStripNoises.h:47
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > m_tTopoToken
#define LogDebug(id)