CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
EcalIntercalibHandler.cc
Go to the documentation of this file.
5 
6 #include <iostream>
7 
8 const Int_t kEBChannels = 61200, kEEChannels = 14648;
9 
11  : m_name(ps.getUntrackedParameter<std::string>("name", "EcalIntercalibHandler")) {
12  edm::LogInfo("EcalIntercalib Source handler constructor\n");
13  m_firstRun = static_cast<unsigned int>(atoi(ps.getParameter<std::string>("firstRun").c_str()));
14  m_file_type = ps.getParameter<std::string>("type"); // xml/txt
15  m_file_name = ps.getParameter<std::string>("fileName");
16 }
17 
19 
21  // std::cout << "------- Ecal - > getNewObjects\n";
22  std::ostringstream ss;
23  ss << "ECAL ";
24 
25  unsigned long long irun;
26  std::string file_ = m_file_name;
27  edm::LogInfo("going to open file ") << file_;
28 
29  // EcalCondHeader header;
31  if (m_file_type == "xml")
32  readXML(file_, *payload);
33  else
34  readTXT(file_, *payload);
35  irun = m_firstRun;
36  Time_t snc = (Time_t)irun;
37 
38  popcon::PopConSourceHandler<EcalIntercalibConstants>::m_to_transfer.push_back(std::make_pair(payload, snc));
39 }
40 
42  std::string dummyLine, bid;
43  std::ifstream fxml;
44  fxml.open(file_);
45  if (!fxml.is_open()) {
46  edm::LogInfo("ERROR : cannot open file ") << file_;
47  exit(1);
48  }
49  // header
50  for (int i = 0; i < 6; i++) {
51  getline(fxml, dummyLine); // skip first lines
52  // std::cout << dummyLine << std::endl;
53  }
54  fxml >> bid;
55  std::string stt = bid.substr(7, 5);
56  std::istringstream iEB(stt);
57  int nEB;
58  iEB >> nEB;
59  if (nEB != kEBChannels) {
60  edm::LogInfo("strange number of EB channels ") << nEB;
61  exit(-1);
62  }
63  fxml >> bid; // <item_version>0</item_version>
64  for (int iChannel = 0; iChannel < kEBChannels; iChannel++) {
65  EBDetId myEBDetId = EBDetId::unhashIndex(iChannel);
66  fxml >> bid;
67  std::size_t found = bid.find("</");
68  stt = bid.substr(6, found - 6);
69  float val = std::stof(stt);
70  record[myEBDetId] = val;
71  }
72  for (int i = 0; i < 5; i++) {
73  getline(fxml, dummyLine); // skip first lines
74  // std::cout << dummyLine << std::endl;
75  }
76  fxml >> bid;
77  stt = bid.substr(7, 5);
78  std::istringstream iEE(stt);
79  int nEE;
80  iEE >> nEE;
81  if (nEE != kEEChannels) {
82  edm::LogInfo("strange number of EE channels ") << nEE;
83  exit(-1);
84  }
85  fxml >> bid; // <item_version>0</item_version>
86  // now endcaps
87  for (int iChannel = 0; iChannel < kEEChannels; iChannel++) {
88  EEDetId myEEDetId = EEDetId::unhashIndex(iChannel);
89  fxml >> bid;
90  std::size_t found = bid.find("</");
91  stt = bid.substr(6, found - 6);
92  float val = std::stof(stt);
93  record[myEEDetId] = val;
94  }
95 }
96 
98  std::ifstream ftxt;
99  ftxt.open(file_);
100  if (!ftxt.is_open()) {
101  edm::LogInfo("ERROR : cannot open file ") << file_;
102  exit(1);
103  }
104  int number_of_lines = 0, eta, phi, x, y, z;
105  float val;
107  while (std::getline(ftxt, line)) {
108  if (number_of_lines < kEBChannels) { // barrel
109  sscanf(line.c_str(), "%i %i %i %f", &eta, &phi, &z, &val);
110  EBDetId ebdetid(eta, phi, EBDetId::ETAPHIMODE);
111  record[ebdetid] = val;
112  } else { // endcaps
113  sscanf(line.c_str(), "%i %i %i %f", &x, &y, &z, &val);
114  EEDetId eedetid(x, y, z, EEDetId::XYMODE);
115  record[eedetid] = val;
116  }
117  number_of_lines++;
118  }
119  edm::LogInfo("Number of lines in text file: ") << number_of_lines;
121  if (number_of_lines != kChannels)
122  edm::LogInfo("wrong number of channels! Please check ");
123 }
static const int XYMODE
Definition: EEDetId.h:335
static EEDetId unhashIndex(int hi)
Definition: EEDetId.cc:65
void readTXT(const std::string &filename, EcalFloatCondObjectContainer &record)
void readXML(const std::string &filename, EcalFloatCondObjectContainer &record)
EcalIntercalibHandler(edm::ParameterSet const &)
static const int ETAPHIMODE
Definition: EBDetId.h:158
Log< level::Info, false > LogInfo
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
static EBDetId unhashIndex(int hi)
get a DetId from a compact index for arrays
Definition: EBDetId.h:110
cond::Time_t Time_t
Definition: Time.h:18
EcalIntercalibConstantMap EcalIntercalibConstants
float x
const Int_t kChannels