CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiPixelLorentzAngleCalibration.cc
Go to the documentation of this file.
1 
16 // include 'locally':
17 #include "TreeStruct.h"
18 
22 
25 
30 
34 
39 
40 #include "TTree.h"
41 #include "TFile.h"
42 #include "TString.h"
43 
44 #include <boost/assign/list_of.hpp>
45 #include <vector>
46 #include <map>
47 #include <sstream>
48 #include <list>
49 #include <utility>
50 #include <set>
51 
52 #include <functional>
53 
55 {
56 public:
59 
62 
64  virtual unsigned int numParameters() const override;
65 
66  // /// Return all derivatives,
67  // /// default implementation uses other derivatives(..) method,
68  // /// but can be overwritten in derived class for efficiency.
69  // virtual std::vector<double> derivatives(const TransientTrackingRecHit &hit,
70  // const TrajectoryStateOnSurface &tsos,
71  // const edm::EventSetup &setup,
72  // const EventInfo &eventInfo) const;
73 
76  virtual unsigned int derivatives(std::vector<ValuesIndexPair> &outDerivInds,
78  const TrajectoryStateOnSurface &tsos,
79  const edm::EventSetup &setup,
80  const EventInfo &eventInfo) const override;
81 
84  virtual bool setParameter(unsigned int index, double value) override;
85 
88  virtual bool setParameterError(unsigned int index, double error) override;
89 
92  virtual double getParameter(unsigned int index) const override;
93 
96  virtual double getParameterError(unsigned int index) const override;
97 
98  // /// Call at beginning of job:
99  virtual void beginOfJob(AlignableTracker *tracker,
101  AlignableExtras *extras) override;
102 
103 
104 
107  virtual void endOfJob() override;
108 
109 private:
112  bool checkLorentzAngleInput(const edm::EventSetup &setup, const EventInfo &eventInfo);
117 
120  double getParameterForDetId(unsigned int detId, edm::RunNumber_t run) const;
121 
122  void writeTree(const SiPixelLorentzAngle *lorentzAngle,
123  const std::map<unsigned int,TreeStruct>& treeInfo, const char *treeName) const;
124  SiPixelLorentzAngle* createFromTree(const char *fileName, const char *treeName) const;
125 
126  const bool saveToDB_;
129  const std::vector<std::string> mergeFileNames_;
130 
132 
133  // const AlignableTracker *alignableTracker_;
135  std::vector<double> parameters_;
136  std::vector<double> paramUncertainties_;
137 
140 };
141 
142 //======================================================================
143 //======================================================================
144 //======================================================================
145 
148  saveToDB_(cfg.getParameter<bool>("saveToDB")),
149  recordNameDBwrite_(cfg.getParameter<std::string>("recordNameDBwrite")),
150  outFileName_(cfg.getParameter<std::string>("treeFile")),
151  mergeFileNames_(cfg.getParameter<std::vector<std::string> >("mergeTreeFiles")),
152  siPixelLorentzAngleInput_(0),
153  moduleGroupSelector_(0),
154  moduleGroupSelCfg_(cfg.getParameter<edm::ParameterSet>("LorentzAngleModuleGroups"))
155 {
156 
157 }
158 
159 //======================================================================
161 {
162  delete moduleGroupSelector_;
163  // std::cout << "Destroy SiPixelLorentzAngleCalibration named " << this->name() << std::endl;
165 }
166 
167 //======================================================================
169 {
170  return parameters_.size();
171 }
172 
173 //======================================================================
174 unsigned int
175 SiPixelLorentzAngleCalibration::derivatives(std::vector<ValuesIndexPair> &outDerivInds,
177  const TrajectoryStateOnSurface &tsos,
178  const edm::EventSetup &setup,
179  const EventInfo &eventInfo) const
180 {
181  // ugly const-cast:
182  // But it is either only first initialisation or throwing an exception...
183  const_cast<SiPixelLorentzAngleCalibration*>(this)->checkLorentzAngleInput(setup, eventInfo);
184 
185  outDerivInds.clear();
186 
187  if (hit.det()) { // otherwise 'constraint hit' or whatever
188 
190  eventInfo.eventId().run());
191  if (index >= 0) { // otherwise not treated
193  setup.get<IdealMagneticFieldRecord>().get(magneticField);
194  const GlobalVector bField(magneticField->inTesla(hit.det()->surface().position()));
195  const LocalVector bFieldLocal(hit.det()->surface().toLocal(bField));
196  const double dZ = hit.det()->surface().bounds().thickness(); // it is a float only...
197  // shift due to LA: dx = tan(LA) * dz/2 = mobility * B_y * dz/2,
198  // '-' since we have derivative of the residual r = trk -hit
199  const double xDerivative = bFieldLocal.y() * dZ * -0.5; // parameter is mobility!
200  // FIXME: Have to treat that for FPIX yDerivative != 0., due to higher order effects!
201  if (xDerivative) { // If field is zero, this is zero: do not return it
202  const Values derivs(xDerivative, 0.); // yDerivative = 0.
203  outDerivInds.push_back(ValuesIndexPair(derivs, index));
204  }
205  }
206  } else {
207  edm::LogWarning("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::derivatives2"
208  << "Hit without GeomDet, skip!";
209  }
210 
211  return outDerivInds.size();
212 }
213 
214 //======================================================================
216 {
217  if (index >= parameters_.size()) {
218  return false;
219  } else {
221  return true;
222  }
223 }
224 
225 //======================================================================
227 {
228  if (index >= paramUncertainties_.size()) {
229  return false;
230  } else {
232  return true;
233  }
234 }
235 
236 //======================================================================
238 {
239  return (index >= parameters_.size() ? 0. : parameters_[index]);
240 }
241 
242 //======================================================================
244 {
245  return (index >= paramUncertainties_.size() ? 0. : paramUncertainties_[index]);
246 }
247 
248 
249 
250 
251 //======================================================================
253  AlignableMuon * /*aliMuon*/,
254  AlignableExtras * /*aliExtras*/)
255 {
256  //specify the sub-detectors for which the LA is determined
257  const std::vector<int> sdets = boost::assign::list_of(PixelSubdetector::PixelBarrel)(PixelSubdetector::PixelEndcap);
258 
260 
263 
264  edm::LogInfo("Alignment") << "@SUB=SiPixelLorentzAngleCalibration" << "Created with name "
265  << this->name() << "',\n" << this->numParameters() << " parameters to be determined,"
266  << "\n saveToDB = " << saveToDB_
267  << "\n outFileName = " << outFileName_
268  << "\n N(merge files) = " << mergeFileNames_.size()
269  << "\n number of IOVs = " << moduleGroupSelector_->numIovs();
270 
271  if (mergeFileNames_.size()) {
272  edm::LogInfo("Alignment") << "@SUB=SiPixelLorentzAngleCalibration"
273  << "First file to merge: " << mergeFileNames_[0];
274  }
275 }
276 
277 
278 //======================================================================
280 {
281  // loginfo output
282  std::ostringstream out;
283  out << "Parameter results\n";
284  for (unsigned int iPar = 0; iPar < parameters_.size(); ++iPar) {
285  out << iPar << ": " << parameters_[iPar] << " +- " << paramUncertainties_[iPar] << "\n";
286  }
287  edm::LogInfo("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::endOfJob" << out.str();
288 
289  std::map<unsigned int, TreeStruct> treeInfo; // map of TreeStruct for each detId
290 
291  // now write 'input' tree
292  const SiPixelLorentzAngle *input = this->getLorentzAnglesInput(); // never NULL
293  const std::string treeName(this->name() + '_');
294  this->writeTree(input, treeInfo, (treeName + "input").c_str()); // empty treeInfo for input...
295 
296  if (input->getLorentzAngles().empty()) {
297  edm::LogError("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::endOfJob"
298  << "Input Lorentz angle map is empty, skip writing output!";
299  return;
300  }
301 
302  const unsigned int nonZeroParamsOrErrors = // Any determined value?
303  count_if (parameters_.begin(), parameters_.end(), std::bind2nd(std::not_equal_to<double>(),0.))
304  + count_if(paramUncertainties_.begin(), paramUncertainties_.end(),
305  std::bind2nd(std::not_equal_to<double>(), 0.));
306 
307  for (unsigned int iIOV = 0; iIOV < moduleGroupSelector_->numIovs(); ++iIOV) {
308  // for (unsigned int iIOV = 0; iIOV < 1; ++iIOV) { // For writing out the modified values
309  cond::Time_t firstRunOfIOV = moduleGroupSelector_->firstRunOfIOV(iIOV);
311  // Loop on map of values from input and add (possible) parameter results
312  for (auto iterIdValue = input->getLorentzAngles().begin();
313  iterIdValue != input->getLorentzAngles().end(); ++iterIdValue) {
314  // type of (*iterIdValue) is pair<unsigned int, float>
315  const unsigned int detId = iterIdValue->first; // key of map is DetId
316  const double param = this->getParameterForDetId(detId, firstRunOfIOV);
317  // put result in output, i.e. sum of input and determined parameter, but the nasty
318  // putLorentzAngle(..) takes float by reference - not even const reference:
319  float value = iterIdValue->second + param;
320  output->putLorentzAngle(detId, value);
321  const int paramIndex = moduleGroupSelector_->getParameterIndexFromDetId(detId,firstRunOfIOV);
322  treeInfo[detId] = TreeStruct(param, this->getParameterError(paramIndex), paramIndex);
323  }
324 
325  if (saveToDB_ || nonZeroParamsOrErrors != 0) { // Skip writing mille jobs...
326  this->writeTree(output, treeInfo, (treeName + Form("result_%lld", firstRunOfIOV)).c_str());
327  }
328 
329  if (saveToDB_) { // If requested, write out to DB
331  if (dbService.isAvailable()) {
332  dbService->writeOne(output, firstRunOfIOV, recordNameDBwrite_.c_str());
333  // no 'delete output;': writeOne(..) took over ownership
334  } else {
335  delete output;
336  edm::LogError("BadConfig") << "@SUB=SiPixelLorentzAngleCalibration::endOfJob"
337  << "No PoolDBOutputService available, but saveToDB true!";
338  }
339  } else {
340  delete output;
341  }
342  } // end loop on IOVs
343 
344 }
345 
346 //======================================================================
348  const EventInfo &eventInfo)
349 {
350  edm::ESHandle<SiPixelLorentzAngle> lorentzAngleHandle;
352  setup.get<SiPixelLorentzAngleRcd>().get(lorentzAngleHandle);
353  siPixelLorentzAngleInput_ = new SiPixelLorentzAngle(*lorentzAngleHandle);
354  } else {
355  if (watchLorentzAngleRcd_.check(setup)) { // new IOV of input
356  setup.get<SiPixelLorentzAngleRcd>().get(lorentzAngleHandle);
357  if (lorentzAngleHandle->getLorentzAngles() // but only bad if non-identical values
358  != siPixelLorentzAngleInput_->getLorentzAngles()) { // (comparing maps)
359  // Maps are containers sorted by key, but comparison problems may arise from
360  // 'floating point comparison' problems (FIXME?)
361  throw cms::Exception("BadInput")
362  << "SiPixelLorentzAngleCalibration::checkLorentzAngleInput:\n"
363  << "Content of SiPixelLorentzAngle changed at run " << eventInfo.eventId().run()
364  << ", but algorithm expects constant input!\n";
365  return false; // not reached...
366  }
367  }
368  }
369 
370  return true;
371 }
372 
373 //======================================================================
375 {
376  // For parallel processing in Millepede II, create SiPixelLorentzAngle
377  // from info stored in files of parallel jobs and check that they are identical.
378  // If this job has run on data, still check that LA is identical to the ones
379  // from mergeFileNames_.
380  const std::string treeName(this->name() + "_input");
381  for (auto iFile = mergeFileNames_.begin(); iFile != mergeFileNames_.end(); ++iFile) {
382  SiPixelLorentzAngle* la = this->createFromTree(iFile->c_str(), treeName.c_str());
383  // siPixelLorentzAngleInput_ could be non-null from previous file of this loop
384  // or from checkLorentzAngleInput(..) when running on data in this job as well
386  delete siPixelLorentzAngleInput_; // NULL or empty
388  } else {
389  // FIXME: about comparison of maps see comments in checkLorentzAngleInput
390  if (la && !la->getLorentzAngles().empty() && // single job might not have got events
392  // Throw exception instead of error?
393  edm::LogError("NoInput") << "@SUB=SiPixelLorentzAngleCalibration::getLorentzAnglesInput"
394  << "Different input values from tree " << treeName
395  << " in file " << *iFile << ".";
396 
397  }
398  delete la;
399  }
400  }
401 
402  if (!siPixelLorentzAngleInput_) { // no files nor ran on events
404  edm::LogError("NoInput") << "@SUB=SiPixelLorentzAngleCalibration::getLorentzAnglesInput"
405  << "No input, create an empty one!";
406  } else if (siPixelLorentzAngleInput_->getLorentzAngles().empty()) {
407  edm::LogError("NoInput") << "@SUB=SiPixelLorentzAngleCalibration::getLorentzAnglesInput"
408  << "Empty result!";
409  }
410 
412 }
413 
414 //======================================================================
416  edm::RunNumber_t run) const
417 {
418  const int index = moduleGroupSelector_->getParameterIndexFromDetId(detId, run);
419  return (index < 0 ? 0. : parameters_[index]);
420 }
421 
422 //======================================================================
424  const std::map<unsigned int,TreeStruct> &treeInfo,
425  const char *treeName) const
426 {
427  if (!lorentzAngle) return;
428 
429  TFile* file = TFile::Open(outFileName_.c_str(), "UPDATE");
430  if (!file) {
431  edm::LogError("BadConfig") << "@SUB=SiPixelLorentzAngleCalibration::writeTree"
432  << "Could not open file '" << outFileName_ << "'.";
433  return;
434  }
435 
436  TTree *tree = new TTree(treeName, treeName);
437  unsigned int id = 0;
438  float value = 0.;
439  TreeStruct treeStruct;
440  tree->Branch("detId", &id, "detId/i");
441  tree->Branch("value", &value, "value/F");
442  tree->Branch("treeStruct", &treeStruct, TreeStruct::LeafList());
443 
444  for (auto iterIdValue = lorentzAngle->getLorentzAngles().begin();
445  iterIdValue != lorentzAngle->getLorentzAngles().end(); ++iterIdValue) {
446  // type of (*iterIdValue) is pair<unsigned int, float>
447  id = iterIdValue->first; // key of map is DetId
448  value = iterIdValue->second;
449  // type of (*treeStructIter) is pair<unsigned int, TreeStruct>
450  auto treeStructIter = treeInfo.find(id); // find info for this id
451  if (treeStructIter != treeInfo.end()) {
452  treeStruct = treeStructIter->second; // info from input map
453  } else { // if none found, fill at least parameter index (using 1st IOV...)
454  const cond::Time_t run1of1stIov = moduleGroupSelector_->firstRunOfIOV(0);
455  const int ind = moduleGroupSelector_->getParameterIndexFromDetId(id, run1of1stIov);
456  treeStruct = TreeStruct(ind);
457  }
458  tree->Fill();
459  }
460  tree->Write();
461  delete file; // tree vanishes with the file... (?)
462 
463 }
464 
465 //======================================================================
467 SiPixelLorentzAngleCalibration::createFromTree(const char *fileName, const char *treeName) const
468 {
469  // Check for file existence on your own to work around
470  // https://hypernews.cern.ch/HyperNews/CMS/get/swDevelopment/2715.html:
471  TFile* file = 0;
472  FILE* testFile = fopen(fileName,"r");
473  if (testFile) {
474  fclose(testFile);
475  file = TFile::Open(fileName, "READ");
476  } // else not existing, see error below
477 
478  TTree *tree = 0;
479  if (file) file->GetObject(treeName, tree);
480 
482  if (tree) {
483  unsigned int id = 0;
484  float value = 0.;
485  tree->SetBranchAddress("detId", &id);
486  tree->SetBranchAddress("value", &value);
487 
488  result = new SiPixelLorentzAngle;
489  const Long64_t nEntries = tree->GetEntries();
490  for (Long64_t iEntry = 0; iEntry < nEntries; ++iEntry) {
491  tree->GetEntry(iEntry);
492  result->putLorentzAngle(id, value);
493  }
494  } else { // Warning only since could be parallel job on no events.
495  edm::LogWarning("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::createFromTree"
496  << "Could not get TTree '" << treeName << "' from file '"
497  << fileName << (file ? "'." : "' (file does not exist).");
498  }
499 
500  delete file; // tree will vanish with file
501  return result;
502 }
503 
504 
505 //======================================================================
506 //======================================================================
507 // Plugin definition
508 
510 
512  SiPixelLorentzAngleCalibration, "SiPixelLorentzAngleCalibration");
RunNumber_t run() const
Definition: EventID.h:39
const std::map< unsigned int, float > & getLorentzAngles() const
virtual bool setParameterError(unsigned int index, double error) override
tuple cfg
Definition: looper.py:293
edm::RunNumber_t firstRunOfIOV(unsigned int iovNum) const
First run of iov (0 if iovNum not treated).
virtual bool setParameter(unsigned int index, double value) override
virtual unsigned int numParameters() const override
How many parameters does this calibration define?
const std::string & name() const
name of this calibration
tuple magneticField
const Bounds & bounds() const
Definition: Surface.h:120
bool checkLorentzAngleInput(const edm::EventSetup &setup, const EventInfo &eventInfo)
int getParameterIndexFromDetId(unsigned int detId, edm::RunNumber_t run) const
virtual double getParameter(unsigned int index) const override
const edm::EventID eventId() const
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:40
virtual unsigned int derivatives(std::vector< ValuesIndexPair > &outDerivInds, const TransientTrackingRecHit &hit, const TrajectoryStateOnSurface &tsos, const edm::EventSetup &setup, const EventInfo &eventInfo) const override
define event information passed to algorithms
static std::string const input
Definition: EdmProvDump.cc:44
void writeTree(const SiPixelLorentzAngle *lorentzAngle, const std::map< unsigned int, TreeStruct > &treeInfo, const char *treeName) const
virtual float thickness() const =0
structure to store algorithm results in a TTree
Definition: TreeStruct.h:6
edm::ESWatcher< SiPixelLorentzAngleRcd > watchLorentzAngleRcd_
unsigned long long Time_t
Definition: Time.h:16
SiPixelLorentzAngleCalibration(const edm::ParameterSet &cfg)
Constructor.
const GeomDet * det() const
LocalPoint toLocal(const GlobalPoint &gp) const
tuple result
Definition: query.py:137
std::pair< Values, unsigned int > ValuesIndexPair
x- and y-values
virtual double getParameterError(unsigned int index) const override
virtual void beginOfJob(AlignableTracker *tracker, AlignableMuon *muon, AlignableExtras *extras) override
bool isAvailable() const
Definition: Service.h:46
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:77
void writeOne(T *payload, Time_t time, const std::string &recordName, bool withlogging=false)
bool putLorentzAngle(const uint32_t &, float &)
const SiPixelLorentzAngle * getLorentzAnglesInput()
double getParameterForDetId(unsigned int detId, edm::RunNumber_t run) const
std::pair< double, double > Values
SiPixelLorentzAngle * createFromTree(const char *fileName, const char *treeName) const
const T & get() const
Definition: EventSetup.h:56
static const char * LeafList()
Definition: TreeStruct.h:17
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
unsigned int getNumberOfParameters() const
unsigned int numIovs() const
Total number of IOVs.
const std::vector< std::string > mergeFileNames_
#define DEFINE_EDM_PLUGIN(factory, type, name)
unsigned int RunNumber_t
Constructor of the full muon geometry.
Definition: AlignableMuon.h:36
const PositionType & position() const
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")