CMS 3D CMS Logo

AlignPCLThresholdsReader.cc
Go to the documentation of this file.
1 #include <array>
2 #include <string>
3 #include <iostream>
4 #include <map>
14 
15 namespace edmtest {
16  template <typename T, typename R>
18  public:
20  ~AlignPCLThresholdsReader() override;
21 
22  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
23 
24  private:
25  void analyze(const edm::Event& e, const edm::EventSetup& c) override;
26 
27  // ----------member data ---------------------------
29  const bool printdebug_;
31  };
32 
33  template <typename T, typename R>
35  : thresholdToken_(esConsumes()),
36  printdebug_(p.getUntrackedParameter<bool>("printDebug", true)),
37  formatedOutput_(p.getUntrackedParameter<std::string>("outputFile", "")) {
38  edm::LogInfo("AlignPCLThresholdsReader") << "AlignPCLThresholdsReader" << std::endl;
39  }
40 
41  template <typename T, typename R>
43  edm::LogInfo("AlignPCLThresholdsReader") << "~AlignPCLThresholdsReader " << std::endl;
44  }
45 
46  template <typename T, typename R>
48  edm::LogInfo("AlignPCLThresholdsReader") << "### AlignPCLThresholdsReader::analyze ###" << std::endl;
49  edm::LogInfo("AlignPCLThresholdsReader") << " I AM IN RUN NUMBER " << e.id().run() << std::endl;
50  edm::LogInfo("AlignPCLThresholdsReader") << " ---EVENT NUMBER " << e.id().event() << std::endl;
51 
52  edm::eventsetup::EventSetupRecordKey inputKey = edm::eventsetup::EventSetupRecordKey::makeKey<R>();
55 
56  if (recordKey.type() == edm::eventsetup::EventSetupRecordKey::TypeTag()) {
57  //record not found
58  edm::LogInfo("AlignPCLThresholdsReader")
59  << "Record \"" << inputKey.type().name() << "\" does not exist " << std::endl;
60  }
61 
62  //this part gets the handle of the event source and the record (i.e. the Database)
63  edm::ESHandle<T> thresholdHandle = context.getHandle(thresholdToken_);
64  edm::LogInfo("AlignPCLThresholdsReader") << "got eshandle" << std::endl;
65 
66  if (!thresholdHandle.isValid()) {
67  edm::LogError("AlignPCLThresholdsReader") << " Could not get Handle" << std::endl;
68  return;
69  }
70 
71  const T* thresholds = thresholdHandle.product();
72  edm::LogInfo("AlignPCLThresholdsReader") << "got AlignPCLThresholds* " << std::endl;
73  edm::LogInfo("AlignPCLThresholdsReader") << "print pointer address : ";
74  edm::LogInfo("AlignPCLThresholdsReader") << thresholds << std::endl;
75 
76  edm::LogInfo("AlignPCLThresholdsReader") << "Size " << thresholds->size() << std::endl;
77  edm::LogInfo("AlignPCLThresholdsReader") << "Content of myThresholds " << std::endl;
78  // use built-in method in the CondFormat to print the content
79  if (thresholds && printdebug_) {
80  thresholds->printAll();
81  }
82 
83  FILE* pFile = nullptr;
84  if (!formatedOutput_.empty())
85  pFile = fopen(formatedOutput_.c_str(), "w");
86  if (pFile) {
87  fprintf(pFile, "AlignPCLThresholds::printAll() \n");
88  fprintf(pFile,
89  " ======================================================================================================="
90  "============\n");
91  fprintf(pFile, "N records cut: %i \n", thresholds->getNrecords());
92 
93  AlignPCLThresholds::threshold_map m_thresholds = thresholds->getThreshold_Map();
95 
96  if constexpr (std::is_same_v<T, AlignPCLThresholdsHG>) {
97  m_floatMap = thresholds->getFloatMap();
98  }
99 
100  for (auto it = m_thresholds.begin(); it != m_thresholds.end(); ++it) {
101  bool hasFractionCut = (m_floatMap.find(it->first) != m_floatMap.end());
102 
103  fprintf(pFile,
104  " ====================================================================================================="
105  "==============\n");
106  fprintf(pFile, "key : %s \n", (it->first).c_str());
107  fprintf(pFile, "- Xcut : %8.3f um ", (it->second).getXcut());
108  fprintf(pFile, "| sigXcut : %8.3f ", (it->second).getSigXcut());
109  fprintf(pFile, "| maxMoveXcut : %8.3f um ", (it->second).getMaxMoveXcut());
110  fprintf(pFile, "| ErrorXcut : %8.3f um ", (it->second).getErrorXcut());
111  if constexpr (std::is_same_v<T, AlignPCLThresholdsHG>) {
112  if (hasFractionCut) {
113  fprintf(pFile,
114  "| X_fractionCut : %8.3f \n",
115  thresholds->getFractionCut(it->first, AlignPCLThresholds::coordType::X));
116  } else {
117  fprintf(pFile, "\n");
118  }
119  } else {
120  fprintf(pFile, "\n");
121  }
122 
123  fprintf(pFile, "- thetaXcut : %8.3f urad ", (it->second).getThetaXcut());
124  fprintf(pFile, "| sigThetaXcut : %8.3f ", (it->second).getSigThetaXcut());
125  fprintf(pFile, "| maxMoveThetaXcut : %8.3f urad ", (it->second).getMaxMoveThetaXcut());
126  fprintf(pFile, "| ErrorThetaXcut : %8.3f urad ", (it->second).getErrorThetaXcut());
127  if constexpr (std::is_same_v<T, AlignPCLThresholdsHG>) {
128  if (hasFractionCut) {
129  fprintf(pFile,
130  "| thetaX_fractionCut : %8.3f \n",
131  thresholds->getFractionCut(it->first, AlignPCLThresholds::coordType::theta_X));
132  } else {
133  fprintf(pFile, "\n");
134  }
135  } else {
136  fprintf(pFile, "\n");
137  }
138 
139  fprintf(pFile, "- Ycut : %8.3f um ", (it->second).getYcut());
140  fprintf(pFile, "| sigYcut : %8.3f ", (it->second).getSigXcut());
141  fprintf(pFile, "| maxMoveYcut : %8.3f um ", (it->second).getMaxMoveYcut());
142  fprintf(pFile, "| ErrorYcut : %8.3f um ", (it->second).getErrorYcut());
143  if constexpr (std::is_same_v<T, AlignPCLThresholdsHG>) {
144  if (hasFractionCut) {
145  fprintf(pFile,
146  "| Y_fractionCut : %8.3f \n",
147  thresholds->getFractionCut(it->first, AlignPCLThresholds::coordType::Y));
148  } else {
149  fprintf(pFile, "\n");
150  }
151  } else {
152  fprintf(pFile, "\n");
153  }
154 
155  fprintf(pFile, "- thetaYcut : %8.3f urad ", (it->second).getThetaYcut());
156  fprintf(pFile, "| sigThetaYcut : %8.3f ", (it->second).getSigThetaYcut());
157  fprintf(pFile, "| maxMoveThetaYcut : %8.3f urad ", (it->second).getMaxMoveThetaYcut());
158  fprintf(pFile, "| ErrorThetaYcut : %8.3f urad ", (it->second).getErrorThetaYcut());
159  if constexpr (std::is_same_v<T, AlignPCLThresholdsHG>) {
160  if (hasFractionCut) {
161  fprintf(pFile,
162  "| thetaY_fractionCut : %8.3f \n",
163  thresholds->getFractionCut(it->first, AlignPCLThresholds::coordType::theta_Y));
164  } else {
165  fprintf(pFile, "\n");
166  }
167  } else {
168  fprintf(pFile, "\n");
169  }
170 
171  fprintf(pFile, "- Zcut : %8.3f um ", (it->second).getZcut());
172  fprintf(pFile, "| sigZcut : %8.3f ", (it->second).getSigZcut());
173  fprintf(pFile, "| maxMoveZcut : %8.3f um ", (it->second).getMaxMoveZcut());
174  fprintf(pFile, "| ErrorZcut : %8.3f um ", (it->second).getErrorZcut());
175  if constexpr (std::is_same_v<T, AlignPCLThresholdsHG>) {
176  if (hasFractionCut) {
177  fprintf(pFile,
178  "| Z_fractionCut : %8.3f \n",
179  thresholds->getFractionCut(it->first, AlignPCLThresholds::coordType::Z));
180  } else {
181  fprintf(pFile, "\n");
182  }
183  } else {
184  fprintf(pFile, "\n");
185  }
186 
187  fprintf(pFile, "- thetaZcut : %8.3f urad ", (it->second).getThetaZcut());
188  fprintf(pFile, "| sigThetaZcut : %8.3f ", (it->second).getSigThetaZcut());
189  fprintf(pFile, "| maxMoveThetaZcut : %8.3f urad ", (it->second).getMaxMoveThetaZcut());
190  fprintf(pFile, "| ErrorThetaZcut : %8.3f urad ", (it->second).getErrorThetaZcut());
191  if constexpr (std::is_same_v<T, AlignPCLThresholdsHG>) {
192  if (hasFractionCut) {
193  fprintf(pFile,
194  "| thetaZ_fractionCut : %8.3f \n",
195  thresholds->getFractionCut(it->first, AlignPCLThresholds::coordType::theta_Z));
196  } else {
197  fprintf(pFile, "\n");
198  }
199  } else {
200  fprintf(pFile, "\n");
201  }
202 
203  if ((it->second).hasExtraDOF()) {
204  for (unsigned int j = 0; j < (it->second).extraDOFSize(); j++) {
205  std::array<float, 4> extraDOFCuts = thresholds->getExtraDOFCutsForAlignable(it->first, j);
206  fprintf(pFile,
207  "Extra DOF: %i with label %s \n ",
208  j,
209  thresholds->getExtraDOFLabelForAlignable(it->first, j).c_str());
210  fprintf(pFile, "- cut : %8.3f ", extraDOFCuts.at(0));
211  fprintf(pFile, "| sigCut : %8.3f ", extraDOFCuts.at(1));
212  fprintf(pFile, "| maxMoveCut : %8.3f ", extraDOFCuts.at(2));
213  fprintf(pFile, "| maxErrorCut : %8.3f \n ", extraDOFCuts.at(3));
214  }
215  }
216  }
217  }
218  }
219 
220  template <typename T, typename R>
223  desc.setComment("Reads payloads of type AlignPCLThresholds");
224  desc.addUntracked<bool>("printDebug", true);
225  desc.addUntracked<std::string>("outputFile", "");
227  }
228 
231 
234 } // namespace edmtest
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
std::unordered_map< std::string, std::vector< float > > param_map
std::map< std::string, AlignPCLThreshold > threshold_map
const edm::ESGetToken< T, R > thresholdToken_
AlignPCLThresholdsReader(edm::ParameterSet const &p)
#define X(str)
Definition: MuonsGrabber.cc:38
std::string defaultModuleLabel()
Log< level::Error, false > LogError
T const * product() const
Definition: ESHandle.h:86
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
AlignPCLThresholdsReader< AlignPCLThresholdsHG, AlignPCLThresholdsHGRcd > AlignPCLThresholdsHGReader
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool isValid() const
Definition: ESHandle.h:44
Log< level::Info, false > LogInfo
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void analyze(const edm::Event &e, const edm::EventSetup &c) override
AlignPCLThresholdsReader< AlignPCLThresholds, AlignPCLThresholdsRcd > AlignPCLThresholdsLGReader
heterocontainer::HCTypeTag TypeTag
long double T
static HCTypeTag findType(char const *iTypeName)
find a type based on the types name, if not found will return default HCTypeTag
Definition: HCTypeTag.cc:121