CMS 3D CMS Logo

SiPixelCondObjBuilder.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiPixelCondObjBuilder
4 // Class: SiPixelCondObjBuilder
5 //
13 //
14 // Original Author: Vincenzo CHIOCHIA
15 // Created: Tue Oct 17 17:40:56 CEST 2006
16 // $Id: SiPixelCondObjBuilder.h,v 1.9 2009/05/28 22:12:54 dlange Exp $
17 //
18 //
19 
20 // system includes
21 #include <memory>
22 #include <iostream>
23 #include <string>
24 
25 // user includes
41 
42 #include "CLHEP/Random/RandGauss.h"
43 
44 namespace cms {
46  public:
47  explicit SiPixelCondObjBuilder(const edm::ParameterSet& iConfig);
48 
49  ~SiPixelCondObjBuilder() override = default;
50  void beginJob() override;
51  void analyze(const edm::Event&, const edm::EventSetup&) override;
52  bool loadFromFile();
53 
54  private:
55  const bool appendMode_;
57  std::unique_ptr<SiPixelGainCalibration> SiPixelGainCalibration_;
60 
61  const double meanPed_;
62  const double rmsPed_;
63  const double meanGain_;
64  const double rmsGain_;
66  const double secondRocRowPedOffset_;
67  const int numberOfModules_;
68  const bool fromFile_;
70 
71  // Internal class
72  class CalParameters {
73  public:
74  float p0;
75  float p1;
76  };
77  // Map for storing calibration constants
78  std::map<int, CalParameters, std::less<int> > calmap_;
79  PixelIndices* pIndexConverter_; // Pointer to the index converter
80  };
81 } // namespace cms
82 
83 namespace cms {
85  : appendMode_(iConfig.getUntrackedParameter<bool>("appendMode", true)),
86  pddToken_(esConsumes()),
87  SiPixelGainCalibration_(nullptr),
88  SiPixelGainCalibrationService_(iConfig, consumesCollector()),
89  recordName_(iConfig.getParameter<std::string>("record")),
90  meanPed_(iConfig.getParameter<double>("meanPed")),
91  rmsPed_(iConfig.getParameter<double>("rmsPed")),
92  meanGain_(iConfig.getParameter<double>("meanGain")),
93  rmsGain_(iConfig.getParameter<double>("rmsGain")),
94  secondRocRowGainOffset_(iConfig.getParameter<double>("secondRocRowGainOffset")),
95  secondRocRowPedOffset_(iConfig.getParameter<double>("secondRocRowPedOffset")),
96  numberOfModules_(iConfig.getParameter<int>("numberOfModules")),
97  fromFile_(iConfig.getParameter<bool>("fromFile")),
98  fileName_(iConfig.getParameter<std::string>("fileName")) {
99  ::putenv((char*)"CORAL_AUTH_USER=me");
100  ::putenv((char*)"CORAL_AUTH_PASSWORD=test");
101  }
102 
104  using namespace edm;
105  unsigned int run = iEvent.id().run();
106  int nmodules = 0;
107  uint32_t nchannels = 0;
108  // int mycol = 415;
109  // int myrow = 159;
110 
111  edm::LogInfo("SiPixelCondObjBuilder")
112  << "... creating dummy SiPixelGainCalibration Data for Run " << run << "\n " << std::endl;
113  //
114  // Instantiate Gain calibration offset and define pedestal/gain range
115  //
116  float mingain = 0;
117  float maxgain = 10;
118  float minped = 0;
119  float maxped = 255;
120  SiPixelGainCalibration_ = std::make_unique<SiPixelGainCalibration>(minped, maxped, mingain, maxgain);
121 
122  const TrackerGeometry* pDD = &iSetup.getData(pddToken_);
123  edm::LogInfo("SiPixelCondObjBuilder") << " There are " << pDD->dets().size() << " detectors" << std::endl;
124 
125  for (TrackerGeometry::DetContainer::const_iterator it = pDD->dets().begin(); it != pDD->dets().end(); it++) {
126  if (dynamic_cast<PixelGeomDetUnit const*>((*it)) != nullptr) {
127  uint32_t detid = ((*it)->geographicalId()).rawId();
128 
129  // Stop if module limit reached
130  nmodules++;
131  if (nmodules > numberOfModules_)
132  break;
133 
134  const PixelGeomDetUnit* pixDet = dynamic_cast<const PixelGeomDetUnit*>((*it));
135  const PixelTopology& topol = pixDet->specificTopology();
136  // Get the module sizes.
137  int nrows = topol.nrows(); // rows in x
138  int ncols = topol.ncolumns(); // cols in y
139  //edm::LogPrint("SiPixelCondObjBuilder") << " ---> PIXEL DETID " << detid << " Cols " << ncols << " Rows " << nrows << std::endl;
140 
141  PixelIndices pIndexConverter(ncols, nrows);
142 
143  std::vector<char> theSiPixelGainCalibration;
144 
145  // Loop over columns and rows of this DetID
146  for (int i = 0; i < ncols; i++) {
147  for (int j = 0; j < nrows; j++) {
148  nchannels++;
149 
150  float ped = 0.0, gain = 0.0;
151 
152  if (fromFile_) {
153  // Use calibration from a file
154  int chipIndex = 0, colROC = 0, rowROC = 0;
155 
156  pIndexConverter.transformToROC(i, j, chipIndex, colROC, rowROC);
157  int chanROC = PixelIndices::pixelToChannelROC(rowROC, colROC); // use ROC coordinates
158  // float pp0=0, pp1=0;
159  std::map<int, CalParameters, std::less<int> >::const_iterator it = calmap_.find(chanROC);
160  CalParameters theCalParameters = (*it).second;
161  ped = theCalParameters.p0;
162  gain = theCalParameters.p1;
163 
164  } else {
165  if (rmsPed_ > 0) {
166  ped = CLHEP::RandGauss::shoot(meanPed_, rmsPed_);
167  while (minped > ped || maxped < ped)
168  ped = CLHEP::RandGauss::shoot(meanPed_, rmsPed_);
169 
170  } else
171  ped = meanPed_;
172  if (rmsGain_ > 0) {
173  gain = CLHEP::RandGauss::shoot(meanGain_, rmsGain_);
174  while (mingain > gain || maxgain < gain)
175  gain = CLHEP::RandGauss::shoot(meanGain_, rmsGain_);
176  } else
177  gain = meanGain_;
178  }
179 
180  //if in the second row of rocs (i.e. a 2xN plaquette) add an offset (if desired) for testing
181  if (j >= 80) {
182  ped += secondRocRowPedOffset_;
184 
185  if (gain > maxgain)
186  gain = maxgain;
187  else if (gain < mingain)
188  gain = mingain;
189 
190  if (ped > maxped)
191  ped = maxped;
192  else if (ped < minped)
193  ped = minped;
194  }
195 
196  // if(i==mycol && j==myrow) {
197  // edm::LogPrint("SiPixelCondObjBuilder") << " Col "<<i<<" Row "<<j<<" Ped "<<ped<<" Gain "<<gain<<std::endl;
198  // }
199 
200  // gain = 2.8;
201  // ped = 28.2;
202 
203  // Insert data in the container
204  SiPixelGainCalibration_->setData(ped, gain, theSiPixelGainCalibration);
205  }
206  }
207 
208  SiPixelGainCalibration::Range range(theSiPixelGainCalibration.begin(), theSiPixelGainCalibration.end());
209  if (!SiPixelGainCalibration_->put(detid, range, ncols))
210  edm::LogError("SiPixelCondObjBuilder")
211  << "[SiPixelCondObjBuilder::analyze] detid already exists" << std::endl;
212  }
213  }
214  edm::LogPrint("SiPixelCondObjBuilder") << " ---> PIXEL Modules " << nmodules << std::endl;
215  edm::LogPrint("SiPixelCondObjBuilder") << " ---> PIXEL Channels " << nchannels << std::endl;
216 
217  // // Try to read object
218  // int mynmodules =0;
219  // for(TrackerGeometry::DetContainer::const_iterator it = pDD->dets().begin(); it != pDD->dets().end(); it++){
220  // if( dynamic_cast<PixelGeomDetUnit const*>((*it))!=0){
221  // uint32_t mydetid=((*it)->geographicalId()).rawId();
222  // mynmodules++;
223  // if( mynmodules > numberOfModules_) break;
224  // SiPixelGainCalibration::Range myrange = SiPixelGainCalibration_->getRange(mydetid);
225  // float mypedestal = SiPixelGainCalibration_->getPed (mycol,myrow,myrange,416);
226  // float mygain = SiPixelGainCalibration_->getGain(mycol,myrow,myrange,416);
227  // //edm::LogPrint("SiPixelCondObjBuilder")<<" PEDESTAL "<< mypedestal<<" GAIN "<<mygain<<std::endl;
228  // }
229  // }
230  // Write into DB
231  edm::LogInfo(" --- writeing to DB!");
233  if (!mydbservice.isAvailable()) {
234  edm::LogPrint("SiPixelCondObjBuilder") << "Didn't get DB" << std::endl;
235  edm::LogError("db service unavailable");
236  return;
237  } else {
238  edm::LogInfo("DB service OK");
239  }
240 
241  try {
242  // size_t callbackToken = mydbservice->callbackToken("SiPixelGainCalibration");
243  // edm::LogInfo("SiPixelCondObjBuilder")<<"CallbackToken SiPixelGainCalibration "
244  // <<callbackToken<<std::endl;
245  // unsigned long long tillTime;
246  // if ( appendMode_)
247  // tillTime = mydbservice->currentTime();
248  // else
249  // tillTime = mydbservice->endOfTime();
250  // edm::LogInfo("SiPixelCondObjBuilder")<<"[SiPixelCondObjBuilder::analyze] tillTime = "
251  // <<tillTime<<std::endl;
252  // mydbservice->newValidityForNewPayload<SiPixelGainCalibration>(
253  // SiPixelGainCalibration_, tillTime , callbackToken);
254 
255  if (mydbservice->isNewTagRequest(recordName_)) {
256  mydbservice->createOneIOV<SiPixelGainCalibration>(
258  } else {
259  mydbservice->appendOneIOV<SiPixelGainCalibration>(
261  }
262  edm::LogInfo(" --- all OK");
263  } catch (const cond::Exception& er) {
264  edm::LogPrint("SiPixelCondObjBuilder") << "Database exception! " << er.what() << std::endl;
265  edm::LogError("SiPixelCondObjBuilder") << er.what() << std::endl;
266  } catch (const std::exception& er) {
267  edm::LogPrint("SiPixelCondObjBuilder") << "Standard exception! " << er.what() << std::endl;
268  edm::LogError("SiPixelCondObjBuilder") << "caught std::exception " << er.what() << std::endl;
269  } catch (...) {
270  edm::LogError("SiPixelCondObjBuilder") << "Funny error" << std::endl;
271  }
272  }
273 
274  // ------------ method called once each job just before starting event loop ------------
276  if (fromFile_) {
277  if (loadFromFile()) {
278  edm::LogInfo("SiPixelCondObjBuilder") << " Calibration loaded: Map size " << calmap_.size() << " max "
279  << calmap_.max_size() << " " << calmap_.empty() << std::endl;
280  }
281  }
282  }
283 
285  float par0, par1; //,par2,par3;
286  int colid, rowid; //rocid
288 
289  std::ifstream in_file; // data file pointer
290  in_file.open(fileName_.c_str(), std::ios::in); // in C++
291  if (in_file.bad()) {
292  edm::LogError("SiPixelCondObjBuilder") << "Input file not found" << std::endl;
293  }
294  if (in_file.eof() != 0) {
295  edm::LogError("SiPixelCondObjBuilder") << in_file.eof() << " " << in_file.gcount() << " " << in_file.fail() << " "
296  << in_file.good() << " end of file " << std::endl;
297  return false;
298  }
299  //Load file header
300  char line[500];
301  for (int i = 0; i < 3; i++) {
302  in_file.getline(line, 500, '\n');
303  edm::LogInfo("SiPixelCondObjBuilder") << line << std::endl;
304  }
305  //Loading calibration constants from file, loop on pixels
306  for (int i = 0; i < (52 * 80); i++) {
307  in_file >> par0 >> par1 >> name >> colid >> rowid;
308 
309  edm::LogPrint("SiPixelCondObjBuilder")
310  << " Col " << colid << " Row " << rowid << " P0 " << par0 << " P1 " << par1 << std::endl;
311 
312  CalParameters onePix;
313  onePix.p0 = par0;
314  onePix.p1 = par1;
315 
316  // Convert ROC pixel index to channel
317  int chan = PixelIndices::pixelToChannelROC(rowid, colid);
318  calmap_.insert(std::pair<int, CalParameters>(chan, onePix));
319  }
320 
321  bool flag;
322  if (calmap_.size() == 4160) {
323  flag = true;
324  } else {
325  flag = false;
326  }
327  return flag;
328  }
329 
330 } // namespace cms
331 
332 using namespace cms;
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > pddToken_
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
SiPixelGainCalibrationService SiPixelGainCalibrationService_
Base exception class for the object to relational access.
Definition: Exception.h:11
std::unique_ptr< SiPixelGainCalibration > SiPixelGainCalibration_
virtual int ncolumns() const =0
virtual int nrows() const =0
Log< level::Error, false > LogError
std::map< int, CalParameters, std::less< int > > calmap_
void createOneIOV(const T &payload, cond::Time_t firstSinceTime, const std::string &recordName)
void appendOneIOV(const T &payload, cond::Time_t sinceTime, const std::string &recordName)
int transformToROC(const int col, const int row, int &rocId, int &colROC, int &rowROC) const
Definition: PixelIndices.h:149
int iEvent
Definition: GenABIO.cc:224
bool isNewTagRequest(const std::string &recordName)
static int pixelToChannelROC(const int rowROC, const int colROC)
Definition: PixelIndices.h:233
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool getData(T &iHolder) const
Definition: EventSetup.h:122
Log< level::Warning, true > LogPrint
Namespace of DDCMS conversion namespace.
Log< level::Info, false > LogInfo
chan
lumi = TPaveText(lowX+0.38, lowY+0.061, lowX+0.45, lowY+0.161, "NDC") lumi.SetBorderSize( 0 ) lumi...
std::pair< ContainerIterator, ContainerIterator > Range
HLT enums.
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
bool isAvailable() const
Definition: Service.h:40
void analyze(const edm::Event &, const edm::EventSetup &) override
SiPixelCondObjBuilder(const edm::ParameterSet &iConfig)
char const * what() const noexcept override
Definition: Exception.cc:103
~SiPixelCondObjBuilder() override=default