CMS 3D CMS Logo

SiPixelTemplateDBObjectReader.cc
Go to the documentation of this file.
1 // system includes
2 #include <cmath>
3 #include <fstream>
4 #include <iomanip>
5 #include <iostream>
6 #include <memory>
7 
8 // user includes
25 
27 public:
29  ~SiPixelTemplateDBObjectReader() override = default;
30 
31 private:
32  void analyze(const edm::Event&, const edm::EventSetup&) override;
33 
36 
38 
42  const bool testGlobalTag;
46 };
47 
49  : hasTriggeredWatcher(false),
50  theTemplateCalibrationLocation(iConfig.getParameter<std::string>("siPixelTemplateCalibrationLocation")),
51  theDetailedTemplateDBErrorOutput(iConfig.getParameter<bool>("wantDetailedTemplateDBErrorOutput")),
52  theFullTemplateDBOutput(iConfig.getParameter<bool>("wantFullTemplateDBOutput")),
53  testGlobalTag(iConfig.getParameter<bool>("TestGlobalTag")),
54  magneticFieldToken_(esConsumes()),
55  the1DTemplateESProdToken_(esConsumes()),
56  the1DTemplateToken_(esConsumes()) {}
57 
59  //To test with the ESProducer
60  SiPixelTemplateDBObject dbobject;
61  if (testGlobalTag) {
62  // Get magnetic field
63  GlobalPoint center(0.0, 0.0, 0.0);
65  float theMagField = magfield.product()->inTesla(center).mag();
66 
67  edm::LogPrint("SiPixelTemplateDBObjectReader") << "\nTesting global tag at magfield = " << theMagField;
68  if (SiPixTemplDBObjWatcher_.check(iSetup)) {
69  edm::LogPrint("SiPixelTemplateDBObjectReader") << "With record SiPixelTemplateDBObjectESProducerRcd";
70  dbobject = *&iSetup.getData(the1DTemplateESProdToken_);
71  hasTriggeredWatcher = true;
72  }
73  } else {
74  edm::LogPrint("SiPixelTemplateDBObjectReader") << "\nLoading from file " << std::endl;
75  if (SiPixTemplDBObjWatcher_.check(iSetup)) {
76  edm::LogPrint("SiPixelTemplateDBObjectReader") << "With record SiPixelTemplateDBObjectRcd";
77  dbobject = *&iSetup.getData(the1DTemplateToken_);
78  hasTriggeredWatcher = true;
79  }
80  }
81 
82  if (hasTriggeredWatcher) {
83  std::vector<short> tempMapId;
84 
86  edm::LogPrint("SiPixelTemplateDBObjectReader") << "Map info" << std::endl;
87  std::map<unsigned int, short> templMap = dbobject.getTemplateIDs();
88  for (std::map<unsigned int, short>::const_iterator it = templMap.begin(); it != templMap.end(); ++it) {
89  if (tempMapId.empty())
90  tempMapId.push_back(it->second);
91  for (unsigned int i = 0; i < tempMapId.size(); ++i) {
92  if (tempMapId[i] == it->second)
93  continue;
94  else if (i == tempMapId.size() - 1) {
95  tempMapId.push_back(it->second);
96  break;
97  }
98  }
100  edm::LogPrint("SiPixelTemplateDBObjectReader")
101  << "DetId: " << it->first << " TemplateID: " << it->second << "\n";
102  }
103 
104  edm::LogPrint("SiPixelTemplateDBObjectReader") << "\nMap stores template Id(s): ";
105  for (unsigned int vindex = 0; vindex < tempMapId.size(); ++vindex)
106  edm::LogPrint("SiPixelTemplateDBObjectReader") << tempMapId[vindex] << " ";
107  edm::LogPrint("SiPixelTemplateDBObjectReader") << std::endl;
108 
109  //local variables
110  int numOfTempl = dbobject.numOfTempl();
111  int index = 0;
112  float tempnum = 0, diff = 0;
113  float tol = 1.0E-23;
114  bool error = false, givenErrorMsg = false;
115 
116  edm::LogPrint("SiPixelTemplateDBObjectReader")
117  << "\nChecking Template DB object version " << dbobject.version() << " containing " << numOfTempl
118  << " calibration(s) at " << dbobject.sVector()[index + 22] << "T\n";
119 
120  /*
121  for(unsigned int kk=0;kk < dbobject.sVector().size(); kk++){
122  edm::LogPrint("SiPixelTemplateDBObjectReader") << "dbobject.sVector()[" << kk <<"] = " << dbobject.sVector()[kk] << "\n";
123  }
124  */
125 
126  for (int i = 0; i < numOfTempl; ++i) {
127  //Removes header in db object from diff
128  index += 20;
129 
130  //Tell the person viewing the output what the template ID and version are -- note that version is only valid for >=13
131  edm::LogPrint("SiPixelTemplateDBObjectReader")
132  << "Calibration " << i + 1 << " of " << numOfTempl << ", with Template ID " << dbobject.sVector()[index]
133  << "\tand Version " << dbobject.sVector()[index + 1] << "\t-------- ";
134 
135  //Opening the text-based template calibration
136  std::ostringstream tout;
137  tout << theTemplateCalibrationLocation.c_str() << "/data/template_summary_zp" << std::setw(4) << std::setfill('0')
138  << std::right << dbobject.sVector()[index] << ".out" << std::ends;
139 
140  if (testGlobalTag)
141  continue;
142 
143  edm::FileInPath file(tout.str());
144  std::ifstream in_file(file.fullPath(), std::ios::in);
145 
146  if (in_file.is_open()) {
147  //Removes header in textfile from diff
148  //First read in from the text file -- this will be compared with index = 20
149  in_file >> tempnum;
150 
151  //Read until the end of the current text file
152  while (!in_file.eof()) {
153  //Calculate the difference between the text file and the db object
154  diff = std::abs(tempnum - dbobject.sVector()[index]);
155 
156  //Is there a difference?
157  if (diff > tol) {
158  //We have the if statement to output the message only once
159  if (!givenErrorMsg)
160  edm::LogPrint("SiPixelTemplateDBObjectReader") << "does NOT match\n";
161  //If there is an error we want to display a message upon completion
162  error = true;
163  givenErrorMsg = true;
164  //Do we want more detailed output?
166  edm::LogPrint("SiPixelTemplateDBObjectReader")
167  << "from file = " << tempnum << "\t from dbobject = " << dbobject.sVector()[index]
168  << "\tdiff = " << diff << "\t db index = " << index << std::endl;
169  }
170  }
171  //Go to the next entries
172  in_file >> tempnum;
173  ++index;
174  }
175  //There were no errors, the two files match.
176  if (!givenErrorMsg)
177  edm::LogPrint("SiPixelTemplateDBObjectReader") << "MATCHES\n";
178  } //end current file
179  in_file.close();
180  givenErrorMsg = false;
181  } //end loop over all files
182 
184  edm::LogPrint("SiPixelTemplateDBObjectReader")
185  << "\nThe were differences found between the files and the database.\n"
186  << "If you would like more detailed information please set\n"
187  << "wantDetailedOutput = True in the cfg file. If you would like a\n"
188  << "full output of the contents of the database file please set\n"
189  << "wantFullOutput = True. Make sure that you pipe the output to a\n"
190  << "log file. This could take a few minutes.\n\n";
191 
193  edm::LogPrint("SiPixelTemplateDBObjectReader") << dbobject << std::endl;
194  }
195 }
196 
197 std::ostream& operator<<(std::ostream& s, const SiPixelTemplateDBObject& dbobject) {
199  int index = 0;
201  int txsize[4] = {7, 13, 0, 0};
202  int tysize[4] = {21, 21, 0, 0};
204  int entries[4] = {0};
206  int i, j, k, l, m, n, entry_it;
208  int sizeSetter = 0, templateVersion = 0;
209 
210  edm::LogPrint("SiPixelTemplateDBObjectReader") << "\n\nDBobject version: " << dbobject.version() << std::endl;
211 
212  for (m = 0; m < dbobject.numOfTempl(); ++m) {
213  //To change the size of the output based on which template version we are using"
214  templateVersion = (int)dbobject.sVector_[index + 21];
215  if (templateVersion <= 10) {
216  edm::LogPrint("SiPixelTemplateDBObjectReader")
217  << "*****WARNING***** This code will not format this template version properly *****WARNING*****\n";
218  sizeSetter = 0;
219  } else if (templateVersion <= 16)
220  sizeSetter = 1;
221  else
222  edm::LogPrint("SiPixelTemplateDBObjectReader")
223  << "*****WARNING***** This code has not been tested at formatting this version *****WARNING*****\n";
224 
225  edm::LogPrint("SiPixelTemplateDBObjectReader")
226  << "\n\n*********************************************************************************************"
227  << std::endl;
228  edm::LogPrint("SiPixelTemplateDBObjectReader")
229  << "*************** Reading Template ID " << dbobject.sVector_[index + 20] << "\t(" << m + 1
230  << "/" << dbobject.numOfTempl_ << ") ***************" << std::endl;
231  edm::LogPrint("SiPixelTemplateDBObjectReader")
232  << "*********************************************************************************************\n\n"
233  << std::endl;
234 
235  //Header Title
237  for (n = 0; n < 20; ++n) {
238  temp.f = dbobject.sVector_[index];
239  s << temp.c[0] << temp.c[1] << temp.c[2] << temp.c[3];
240  ++index;
241  }
242 
243  entries[0] = (int)dbobject.sVector_[index + 3]; // Y
244  entries[1] = (int)(dbobject.sVector_[index + 4] * dbobject.sVector_[index + 5]); // X
245 
246  //Header
247  s << dbobject.sVector_[index] << "\t" << dbobject.sVector_[index + 1] << "\t" << dbobject.sVector_[index + 2]
248  << "\t" << dbobject.sVector_[index + 3] << "\t" << dbobject.sVector_[index + 4] << "\t"
249  << dbobject.sVector_[index + 5] << "\t" << dbobject.sVector_[index + 6] << "\t" << dbobject.sVector_[index + 7]
250  << "\t" << dbobject.sVector_[index + 8] << "\t" << dbobject.sVector_[index + 9] << "\t"
251  << dbobject.sVector_[index + 10] << "\t" << dbobject.sVector_[index + 11] << "\t" << dbobject.sVector_[index + 12]
252  << "\t" << dbobject.sVector_[index + 13] << "\t" << dbobject.sVector_[index + 14] << "\t"
253  << dbobject.sVector_[index + 15] << "\t" << dbobject.sVector_[index + 16] << std::endl;
254  index += 17;
255 
256  //Loop over By,Bx,Fy,Fx
257  for (entry_it = 0; entry_it < 4; ++entry_it) {
258  //Run,costrk,qavg,...,clslenx
259  for (i = 0; i < entries[entry_it]; ++i) {
260  s << dbobject.sVector_[index] << "\t" << dbobject.sVector_[index + 1] << "\t" << dbobject.sVector_[index + 2]
261  << "\t" << dbobject.sVector_[index + 3] << "\n"
262  << dbobject.sVector_[index + 4] << "\t" << dbobject.sVector_[index + 5] << "\t"
263  << dbobject.sVector_[index + 6] << "\t" << dbobject.sVector_[index + 7] << "\t"
264  << dbobject.sVector_[index + 8] << "\t" << dbobject.sVector_[index + 9] << "\t"
265  << dbobject.sVector_[index + 10] << "\t" << dbobject.sVector_[index + 11] << "\n"
266  << dbobject.sVector_[index + 12] << "\t" << dbobject.sVector_[index + 13] << "\t"
267  << dbobject.sVector_[index + 14] << "\t" << dbobject.sVector_[index + 15] << "\t"
268  << dbobject.sVector_[index + 16] << "\t" << dbobject.sVector_[index + 17] << "\t"
269  << dbobject.sVector_[index + 18] << std::endl;
270  index += 19;
271  //YPar
272  for (j = 0; j < 2; ++j) {
273  for (k = 0; k < 5; ++k) {
274  s << dbobject.sVector_[index] << "\t";
275  ++index;
276  }
277  s << std::endl;
278  }
279  //YTemp
280  for (j = 0; j < 9; ++j) {
281  for (k = 0; k < tysize[sizeSetter]; ++k) {
282  s << dbobject.sVector_[index] << "\t";
283  ++index;
284  }
285  s << std::endl;
286  }
287  //XPar
288  for (j = 0; j < 2; ++j) {
289  for (k = 0; k < 5; ++k) {
290  s << dbobject.sVector_[index] << "\t";
291  ++index;
292  }
293  s << std::endl;
294  }
295  //XTemp
296  for (j = 0; j < 9; ++j) {
297  for (k = 0; k < txsize[sizeSetter]; ++k) {
298  s << dbobject.sVector_[index] << "\t";
299  ++index;
300  }
301  s << std::endl;
302  }
303  //Y average reco params
304  for (j = 0; j < 4; ++j) {
305  for (k = 0; k < 4; ++k) {
306  s << dbobject.sVector_[index] << "\t";
307  ++index;
308  }
309  s << std::endl;
310  }
311  //Yflpar
312  for (j = 0; j < 4; ++j) {
313  for (k = 0; k < 6; ++k) {
314  s << dbobject.sVector_[index] << "\t";
315  ++index;
316  }
317  s << std::endl;
318  }
319  //X average reco params
320  for (j = 0; j < 4; ++j) {
321  for (k = 0; k < 4; ++k) {
322  s << dbobject.sVector_[index] << "\t";
323  ++index;
324  }
325  s << std::endl;
326  }
327  //Xflpar
328  for (j = 0; j < 4; ++j) {
329  for (k = 0; k < 6; ++k) {
330  s << dbobject.sVector_[index] << "\t";
331  ++index;
332  }
333  s << std::endl;
334  }
335  //Chi2X,Y
336  for (j = 0; j < 4; ++j) {
337  for (k = 0; k < 2; ++k) {
338  for (l = 0; l < 2; ++l) {
339  s << dbobject.sVector_[index] << "\t";
340  ++index;
341  }
342  }
343  s << std::endl;
344  }
345  //Y average Chi2 params
346  for (j = 0; j < 4; ++j) {
347  for (k = 0; k < 4; ++k) {
348  s << dbobject.sVector_[index] << "\t";
349  ++index;
350  }
351  s << std::endl;
352  }
353  //X average Chi2 params
354  for (j = 0; j < 4; ++j) {
355  for (k = 0; k < 4; ++k) {
356  s << dbobject.sVector_[index] << "\t";
357  ++index;
358  }
359  s << std::endl;
360  }
361  //Y average reco params for CPE Generic
362  for (j = 0; j < 4; ++j) {
363  for (k = 0; k < 4; ++k) {
364  s << dbobject.sVector_[index] << "\t";
365  ++index;
366  }
367  s << std::endl;
368  }
369  //X average reco params for CPE Generic
370  for (j = 0; j < 4; ++j) {
371  for (k = 0; k < 4; ++k) {
372  s << dbobject.sVector_[index] << "\t";
373  ++index;
374  }
375  s << std::endl;
376  }
377  //SpareX,Y
378  for (j = 0; j < 20; ++j) {
379  s << dbobject.sVector_[index] << "\t";
380  ++index;
381  if (j == 9 || j == 19)
382  s << std::endl;
383  }
384  }
385  }
386  }
387  return s;
388 }
389 
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
const edm::ESGetToken< SiPixelTemplateDBObject, SiPixelTemplateDBObjectRcd > the1DTemplateToken_
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
SiPixelTemplateDBObjectReader(const edm::ParameterSet &)
edm::ESWatcher< SiPixelTemplateDBObjectESProducerRcd > SiPixTemplDBObjectWatcher_
int iEvent
Definition: GenABIO.cc:224
const edm::ESGetToken< SiPixelTemplateDBObject, SiPixelTemplateDBObjectESProducerRcd > the1DTemplateESProdToken_
std::ostream & operator<<(std::ostream &s, const SiPixelTemplateDBObject &dbobject)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::ESWatcher< SiPixelTemplateDBObjectRcd > SiPixTemplDBObjWatcher_
std::vector< float > const & sVector() const
Log< level::Warning, true > LogPrint
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
const std::map< unsigned int, short > & getTemplateIDs() const
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magneticFieldToken_
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
~SiPixelTemplateDBObjectReader() override=default
void analyze(const edm::Event &, const edm::EventSetup &) override
if(threadIdxLocalY==0 &&threadIdxLocalX==0)