CMS 3D CMS Logo

PPSDirectSimulationData.cc
Go to the documentation of this file.
2 
9 
10 #include "TFile.h"
11 
12 #include <iostream>
13 #include <string>
14 
16  : empiricalAperture45_(""),
17  empiricalAperture56_(""),
18 
19  timeResolutionDiamonds45_(""),
20  timeResolutionDiamonds56_("") {}
21 
23 
24 // Getters
27 
30 
31 std::map<unsigned int, PPSDirectSimulationData::FileObject> &PPSDirectSimulationData::getEfficienciesPerRP() {
32  return efficienciesPerRP_;
33 }
34 std::map<unsigned int, PPSDirectSimulationData::FileObject> &PPSDirectSimulationData::getEfficienciesPerPlane() {
35  return efficienciesPerPlane_;
36 };
37 
38 // Setters
41 
44 
45 std::map<unsigned int, std::unique_ptr<TH2F>> PPSDirectSimulationData::loadEffeciencyHistogramsPerRP() const {
46  std::map<unsigned int, std::unique_ptr<TH2F>> result;
47 
48  for (const auto &it : efficienciesPerRP_)
49  result[it.first] = loadObject(it.second.first, it.second.second);
50 
51  return result;
52 }
53 
54 std::map<unsigned int, std::unique_ptr<TH2F>> PPSDirectSimulationData::loadEffeciencyHistogramsPerPlane() const {
55  std::map<unsigned int, std::unique_ptr<TH2F>> result;
56 
57  for (const auto &it : efficienciesPerPlane_) {
58  CTPPSDetId rpId(it.first);
59 
60  if (rpId.subdetId() == CTPPSDetId::sdTrackingStrip) {
61  for (unsigned int pl = 0; pl < 10; ++pl) {
62  TotemRPDetId plId(rpId.arm(), rpId.station(), rpId.rp(), pl);
63  result[plId] = loadObject(it.second.first, replace(it.second.second, "<detid>", std::to_string(pl)));
64  }
65  }
66 
67  if (rpId.subdetId() == CTPPSDetId::sdTrackingPixel) {
68  for (unsigned int pl = 0; pl < 6; ++pl) {
69  CTPPSPixelDetId plId(rpId.arm(), rpId.station(), rpId.rp(), pl);
70  result[plId] = loadObject(it.second.first, replace(it.second.second, "<detid>", std::to_string(pl)));
71  }
72  }
73 
74  if (rpId.subdetId() == CTPPSDetId::sdTimingDiamond) {
75  for (unsigned int pl = 0; pl < 4; ++pl) {
76  CTPPSDiamondDetId plId(rpId.arm(), rpId.station(), rpId.rp(), pl);
77  result[plId] = loadObject(it.second.first, replace(it.second.second, "<detid>", std::to_string(pl)));
78  }
79  }
80  }
81 
82  return result;
83 }
84 
85 std::unique_ptr<TH2F> PPSDirectSimulationData::loadObject(const std::string &file, const std::string &object) {
86  edm::FileInPath fip(file.c_str());
87  TFile *f_in = TFile::Open(fip.fullPath().c_str());
88  if (!f_in)
89  throw cms::Exception("PPS") << "Cannot open file '" << fip.fullPath() << "'.";
90 
91  TH2F *o_in = (TH2F *)f_in->Get(object.c_str());
92  if (!o_in)
93  throw cms::Exception("PPS") << "Cannot load object '" << object << "' from file '" << fip.fullPath() << "'.";
94 
95  // disassociate histogram from the file
96  o_in->SetDirectory(nullptr);
97 
98  delete f_in;
99 
100  return std::unique_ptr<TH2F>(o_in);
101 }
102 
104  size_t start_pos = 0;
105  while ((start_pos = input.find(from, start_pos)) != std::string::npos) {
106  input.replace(start_pos, from.length(), to);
107  start_pos += to.length();
108  }
109  return input;
110 }
Detector ID class for TOTEM Si strip detectors.
Definition: TotemRPDetId.h:30
const std::string & getEmpiricalAperture56() const
const std::string & getEmpiricalAperture45() const
static std::string replace(std::string input, const std::string &from, const std::string &to)
const std::string & getTimeResolutionDiamonds45() const
const std::string & getTimeResolutionDiamonds56() const
std::map< unsigned int, FileObject > efficienciesPerRP_
void setTimeResolutionDiamonds56(std::string s)
static std::unique_ptr< TH2F > loadObject(const std::string &file, const std::string &object)
static std::string to_string(const XMLCh *ch)
static std::string const input
Definition: EdmProvDump.cc:50
std::map< unsigned int, std::unique_ptr< TH2F > > loadEffeciencyHistogramsPerRP() const
void setEmpiricalAperture45(std::string s)
std::map< unsigned int, FileObject > & getEfficienciesPerRP()
void setEmpiricalAperture56(std::string s)
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
Detector ID class for CTPPS Timing Diamond detectors. Bits [19:31] : Assigend in CTPPSDetId Calss Bit...
std::map< unsigned int, std::unique_ptr< TH2F > > loadEffeciencyHistogramsPerPlane() const
void setTimeResolutionDiamonds45(std::string s)
std::map< unsigned int, FileObject > efficienciesPerPlane_
std::map< unsigned int, FileObject > & getEfficienciesPerPlane()