ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/CMSSW/Alignment/CommonAlignmentAlgorithm/plugins/SiPixelLorentzAngleCalibration.cc
Revision: 1.3
Committed: Wed Sep 19 14:11:49 2012 UTC (12 years, 7 months ago) by flucke
Content type: text/plain
Branch: MAIN
CVS Tags: V04-00-07
Changes since 1.2: +3 -3 lines
Log Message:
fix possible crash introduced in previous version

File Contents

# User Rev Content
1 flucke 1.1 /// \class SiPixelLorentzAngleCalibration
2     ///
3     /// Calibration of Lorentz angle for the pixel,
4     /// integrated in the alignment algorithms.
5     ///
6     /// Note that not all algorithms support this...
7     ///
8     /// \author : Gero Flucke
9     /// date : September 2012
10 flucke 1.3 /// $Revision: 1.2 $
11     /// $Date: 2012/09/14 11:40:43 $
12 flucke 1.1 /// (last update by $Author: flucke $)
13    
14     #include "Alignment/CommonAlignmentAlgorithm/interface/IntegratedCalibrationBase.h"
15    
16     #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
17     #include "CondFormats/DataRecord/interface/SiPixelLorentzAngleRcd.h"
18     #include "CondFormats/SiPixelObjects/interface/SiPixelLorentzAngle.h"
19    
20     #include "DataFormats/DetId/interface/DetId.h"
21     #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
22    
23     #include "FWCore/Framework/interface/ESHandle.h"
24     #include "FWCore/Framework/interface/ESWatcher.h"
25     #include "FWCore/MessageLogger/interface/MessageLogger.h"
26     #include "FWCore/ServiceRegistry/interface/Service.h"
27    
28     #include "Geometry/CommonDetUnit/interface/GeomDet.h"
29     #include "MagneticField/Engine/interface/MagneticField.h"
30     #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
31     #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
32    
33     #include "TTree.h"
34     #include "TFile.h"
35     #include "TString.h"
36    
37     // #include <iostream>
38     #include <vector>
39     #include <sstream>
40    
41     class SiPixelLorentzAngleCalibration : public IntegratedCalibrationBase
42     {
43     public:
44     /// Constructor
45     explicit SiPixelLorentzAngleCalibration(const edm::ParameterSet &cfg);
46    
47     /// Destructor
48     virtual ~SiPixelLorentzAngleCalibration();
49    
50     /// How many parameters does this calibration define?
51     virtual unsigned int numParameters() const;
52    
53     // /// Return all derivatives,
54     // /// default implementation uses other derivatives(..) method,
55     // /// but can be overwritten in derived class for efficiency.
56     // virtual std::vector<double> derivatives(const TransientTrackingRecHit &hit,
57     // const TrajectoryStateOnSurface &tsos,
58     // const edm::EventSetup &setup,
59     // const EventInfo &eventInfo) const;
60    
61     /// Return non-zero derivatives for x- and y-measurements with their indices by reference.
62     /// Return value is their number.
63     virtual unsigned int derivatives(std::vector<ValuesIndexPair> &outDerivInds,
64     const TransientTrackingRecHit &hit,
65     const TrajectoryStateOnSurface &tsos,
66     const edm::EventSetup &setup,
67     const EventInfo &eventInfo) const;
68    
69     /// Setting the determined parameter identified by index,
70     /// should return false if out-of-bounds, true otherwise.
71     virtual bool setParameter(unsigned int index, double value);
72    
73     /// Setting the determined parameter uncertainty identified by index,
74     /// returns false if out-of-bounds, true otherwise.
75     virtual bool setParameterError(unsigned int index, double error);
76    
77     /// Return current value of parameter identified by index.
78     /// Should return 0. if index out-of-bounds.
79     virtual double getParameter(unsigned int index) const;
80    
81     /// Return current value of parameter identified by index.
82     /// Returns 0. if index out-of-bounds or if errors undetermined.
83     virtual double getParameterError(unsigned int index) const;
84    
85     // /// Call at beginning of job:
86     // virtual void beginOfJob(const AlignableTracker *tracker,
87     // const AlignableMuon *muon,
88     // const AlignableExtras *extras);
89    
90     /// Called at end of a the job of the AlignmentProducer.
91     /// Write out determined parameters.
92     virtual void endOfJob();
93    
94     private:
95     /// If called the first time, fill 'siPixelLorentzAngleInput_',
96     /// later check that LorentzAngle has not changed.
97     bool checkLorentzAngleInput(const edm::EventSetup &setup, const EventInfo &eventInfo);
98     /// Input LorentzAngle values:
99     /// - either from EventSetup of first call to derivatives(..)
100     /// - or created from files of passed by configuration (i.e. from parallel processing)
101     const SiPixelLorentzAngle* getLorentzAnglesInput();
102    
103     /// Determined parameter value for this detId (detId not treated => 0.)
104     /// and the given run.
105     double getParameterForDetId(unsigned int detId, edm::RunNumber_t run) const;
106     /// Index of parameter for given detId (detId not treated => < 0)
107     /// and the given run.
108     int getParameterIndexFromDetId(unsigned int detId, edm::RunNumber_t run) const;
109     /// Total number of IOVs.
110     unsigned int numIovs() const;
111     /// First run of iov (0 if iovNum not treated).
112     edm::RunNumber_t firstRunOfIOV(unsigned int iovNum) const;
113    
114     void writeTree(const SiPixelLorentzAngle *lorentzAngle, const char *treeName) const;
115     SiPixelLorentzAngle* createFromTree(const char *fileName, const char *treeName) const;
116    
117     const bool saveToDB_;
118     const std::string recordNameDBwrite_;
119     const std::string outFileName_;
120     const std::vector<std::string> mergeFileNames_;
121    
122     edm::ESWatcher<SiPixelLorentzAngleRcd> watchLorentzAngleRcd_;
123    
124     // const AlignableTracker *alignableTracker_;
125     SiPixelLorentzAngle *siPixelLorentzAngleInput_;
126     std::vector<double> parameters_;
127     std::vector<double> paramUncertainties_;
128     };
129    
130     //======================================================================
131     //======================================================================
132     //======================================================================
133    
134     SiPixelLorentzAngleCalibration::SiPixelLorentzAngleCalibration(const edm::ParameterSet &cfg)
135     : IntegratedCalibrationBase(cfg),
136     saveToDB_(cfg.getParameter<bool>("saveToDB")),
137     recordNameDBwrite_(cfg.getParameter<std::string>("recordNameDBwrite")),
138     outFileName_(cfg.getParameter<std::string>("treeFile")),
139     mergeFileNames_(cfg.getParameter<std::vector<std::string> >("mergeTreeFiles")),
140     // alignableTracker_(0),
141     siPixelLorentzAngleInput_(0)
142     {
143     // FIXME: Which granularity, leading to how many parameters?
144     parameters_.resize(2, 0.); // currently two parameters (ring1-4, 5-8), start value 0.
145     paramUncertainties_.resize(2, 0.); // dito for errors
146    
147     edm::LogInfo("Alignment") << "@SUB=SiPixelLorentzAngleCalibration" << "Created with name "
148     << this->name() << "',\n" << this->numParameters() << " parameters to be determined,"
149     << "\nsaveToDB = " << saveToDB_
150     << "\n outFileName = " << outFileName_
151     << "\n N(merge files) = " << mergeFileNames_.size();
152     if (mergeFileNames_.size()) {
153     edm::LogInfo("Alignment") << "@SUB=SiPixelLorentzAngleCalibration"
154     << "First file to merge: " << mergeFileNames_[0];
155     }
156     }
157    
158     //======================================================================
159     SiPixelLorentzAngleCalibration::~SiPixelLorentzAngleCalibration()
160     {
161     // std::cout << "Destroy SiPixelLorentzAngleCalibration named " << this->name() << std::endl;
162     delete siPixelLorentzAngleInput_;
163     }
164    
165     //======================================================================
166     unsigned int SiPixelLorentzAngleCalibration::numParameters() const
167     {
168     return parameters_.size();
169     }
170    
171     //======================================================================
172     unsigned int
173     SiPixelLorentzAngleCalibration::derivatives(std::vector<ValuesIndexPair> &outDerivInds,
174     const TransientTrackingRecHit &hit,
175     const TrajectoryStateOnSurface &tsos,
176     const edm::EventSetup &setup,
177     const EventInfo &eventInfo) const
178     {
179     // ugly const-cast:
180     // But it is either only first initialisation or throwing an exception...
181     const_cast<SiPixelLorentzAngleCalibration*>(this)->checkLorentzAngleInput(setup, eventInfo);
182    
183     outDerivInds.clear();
184    
185     if (hit.det()) { // otherwise 'constraint hit' or whatever
186    
187     const int index = this->getParameterIndexFromDetId(hit.det()->geographicalId(),
188     eventInfo.eventId_.run());
189     if (index >= 0) { // otherwise not treated
190     edm::ESHandle<MagneticField> magneticField;
191     setup.get<IdealMagneticFieldRecord>().get(magneticField);
192     const GlobalVector bField(magneticField->inTesla(hit.det()->surface().position()));
193     const LocalVector bFieldLocal(hit.det()->surface().toLocal(bField));
194     const double dZ = hit.det()->surface().bounds().thickness(); // it is a float only...
195     // shift due to LA: dx = tan(LA) * dz = mobility * B_y * dz
196     const double xDerivative = bFieldLocal.y() * dZ; // parameter is mobility!
197     if (xDerivative) { // If field is zero, this is zero: do not return it
198     const Values derivs(xDerivative, 0.); // yDerivative = 0.
199     outDerivInds.push_back(ValuesIndexPair(derivs, index));
200     }
201     }
202     } else {
203     edm::LogWarning("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::derivatives2"
204     << "Hit without GeomDet, skip!";
205     }
206    
207     return outDerivInds.size();
208     }
209    
210     //======================================================================
211     bool SiPixelLorentzAngleCalibration::setParameter(unsigned int index, double value)
212     {
213     if (index >= parameters_.size()) {
214     return false;
215     } else {
216     parameters_[index] = value;
217     return true;
218     }
219     }
220    
221     //======================================================================
222     bool SiPixelLorentzAngleCalibration::setParameterError(unsigned int index, double error)
223     {
224     if (index >= paramUncertainties_.size()) {
225     return false;
226     } else {
227     paramUncertainties_[index] = error;
228     return true;
229     }
230     }
231    
232     //======================================================================
233     double SiPixelLorentzAngleCalibration::getParameter(unsigned int index) const
234     {
235     // if (index >= parameters_.size()) {
236     // return 0.;
237     // } else {
238     // return parameters_[index];
239     // }
240     return (index >= parameters_.size() ? 0. : parameters_[index]);
241     }
242    
243     //======================================================================
244     double SiPixelLorentzAngleCalibration::getParameterError(unsigned int index) const
245     {
246     // if (index >= paramUncertainties_.size()) {
247     // return 0.;
248     // } else {
249     // return paramUncertainties_[index];
250     // }
251     return (index >= paramUncertainties_.size() ? 0. : paramUncertainties_[index]);
252     }
253    
254     // //======================================================================
255     // void SiPixelLorentzAngleCalibration::beginOfJob(const AlignableTracker *tracker,
256     // const AlignableMuon */*muon*/,
257     // const AlignableExtras */*extras*/)
258     // {
259     // alignableTracker_ = tracker;
260     // }
261    
262    
263     //======================================================================
264     void SiPixelLorentzAngleCalibration::endOfJob()
265     {
266     // loginfo output
267     std::ostringstream out;
268     out << "Parameter results\n";
269     for (unsigned int iPar = 0; iPar < parameters_.size(); ++iPar) {
270     out << iPar << ": " << parameters_[iPar] << " +- " << paramUncertainties_[iPar] << "\n";
271     }
272     edm::LogInfo("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::endOfJob" << out.str();
273    
274     // now write 'input' tree
275     const SiPixelLorentzAngle *input = this->getLorentzAnglesInput(); // never NULL
276     const std::string treeName(this->name() + '_');
277     this->writeTree(input, (treeName + "input").c_str());
278 flucke 1.2 if (input->getLorentzAngles().empty()) {
279     edm::LogError("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::endOfJob"
280     << "Input Lorentz angle map is empty, skip writing output!";
281     return;
282     }
283 flucke 1.1
284     for (unsigned int iIOV = 0; iIOV < this->numIovs(); ++iIOV) {
285     cond::Time_t firstRunOfIOV = this->firstRunOfIOV(iIOV);
286     SiPixelLorentzAngle *output = new SiPixelLorentzAngle;
287     // Loop on map of values from input and add (possible) parameter results
288     for (auto iterIdValue = input->getLorentzAngles().begin();
289     iterIdValue != input->getLorentzAngles().end(); ++iterIdValue) {
290     // type of (*iterIdValue) is pair<unsigned int, float>
291     const unsigned int detId = iterIdValue->first; // key of map is DetId
292     // Nasty: putLorentzAngle(..) takes float by reference - not even const reference!
293     float value = iterIdValue->second + this->getParameterForDetId(detId, firstRunOfIOV);
294     output->putLorentzAngle(detId, value); // put result in output
295     }
296    
297     // Write this even for mille jobs?
298     this->writeTree(output, (treeName + Form("result_%lld", firstRunOfIOV)).c_str());
299    
300     if (saveToDB_) { // If requested, write out to DB
301     edm::Service<cond::service::PoolDBOutputService> dbService;
302     if (dbService.isAvailable()) {
303     dbService->writeOne(output, firstRunOfIOV, recordNameDBwrite_.c_str());
304     // no 'delete output;': writeOne(..) took over ownership
305     } else {
306     delete output;
307     edm::LogError("BadConfig") << "@SUB=SiPixelLorentzAngleCalibration::endOfJob"
308     << "No PoolDBOutputService available, but saveToDB true!";
309     }
310     } else {
311     delete output;
312     }
313     } // end loop on IOVs
314    
315     }
316    
317     //======================================================================
318     bool SiPixelLorentzAngleCalibration::checkLorentzAngleInput(const edm::EventSetup &setup,
319     const EventInfo &eventInfo)
320     {
321     edm::ESHandle<SiPixelLorentzAngle> lorentzAngleHandle;
322     if (!siPixelLorentzAngleInput_) {
323     setup.get<SiPixelLorentzAngleRcd>().get(lorentzAngleHandle);
324     siPixelLorentzAngleInput_ = new SiPixelLorentzAngle(*lorentzAngleHandle);
325     } else {
326     if (watchLorentzAngleRcd_.check(setup)) { // new IOV of input
327     setup.get<SiPixelLorentzAngleRcd>().get(lorentzAngleHandle);
328     if (lorentzAngleHandle->getLorentzAngles() // but only bad if non-identical values
329     != siPixelLorentzAngleInput_->getLorentzAngles()) { // (comparing maps)
330     // FIXME: Could different maps have same content, but different order?
331     // Or 'floating point comparison' problem?
332     throw cms::Exception("BadInput")
333     << "SiPixelLorentzAngleCalibration::checkLorentzAngleInput:\n"
334     << "Content of SiPixelLorentzAngle changed at run " << eventInfo.eventId_.run()
335     << ", but algorithm expects constant input!\n";
336     return false; // not reached...
337     }
338     }
339     }
340    
341     return true;
342     }
343    
344     //======================================================================
345     const SiPixelLorentzAngle* SiPixelLorentzAngleCalibration::getLorentzAnglesInput()
346     {
347     // For parallel processing in Millepede II, create SiPixelLorentzAngle
348     // from info stored in files of parallel jobs and check that they are identical.
349     // If this job has run on data, still check that LA is identical to the ones
350     // from mergeFileNames_.
351     const std::string treeName(this->name() + "_input");
352     for (auto iFile = mergeFileNames_.begin(); iFile != mergeFileNames_.end(); ++iFile) {
353     SiPixelLorentzAngle* la = this->createFromTree(iFile->c_str(), treeName.c_str());
354     // siPixelLorentzAngleInput_ could be non-null from previous file of this loop
355     // or from checkLorentzAngleInput(..) when running on data in this job as well
356     if (!siPixelLorentzAngleInput_ || siPixelLorentzAngleInput_->getLorentzAngles().empty()) {
357     delete siPixelLorentzAngleInput_; // NULL or empty
358     siPixelLorentzAngleInput_ = la;
359     } else {
360     // FIXME: about comparison of maps see comments in checkLorentzAngleInput
361 flucke 1.3 if (la && !la->getLorentzAngles().empty() && // single job might not have got events
362 flucke 1.1 la->getLorentzAngles() != siPixelLorentzAngleInput_->getLorentzAngles()) {
363     // Throw exception instead of error?
364     edm::LogError("NoInput") << "@SUB=SiPixelLorentzAngleCalibration::getLorentzAnglesInput"
365     << "Different input values from tree " << treeName
366     << " in file " << *iFile << ".";
367    
368     }
369     delete la;
370     }
371     }
372    
373     if (!siPixelLorentzAngleInput_) { // no files nor ran on events
374     siPixelLorentzAngleInput_ = new SiPixelLorentzAngle;
375     edm::LogError("NoInput") << "@SUB=SiPixelLorentzAngleCalibration::getLorentzAnglesInput"
376 flucke 1.2 << "No input, create an empty one!";
377 flucke 1.1 } else if (siPixelLorentzAngleInput_->getLorentzAngles().empty()) {
378     edm::LogError("NoInput") << "@SUB=SiPixelLorentzAngleCalibration::getLorentzAnglesInput"
379 flucke 1.2 << "Empty result!";
380 flucke 1.1 }
381    
382     return siPixelLorentzAngleInput_;
383     }
384    
385     //======================================================================
386     double SiPixelLorentzAngleCalibration::getParameterForDetId(unsigned int detId,
387     edm::RunNumber_t run) const
388     {
389     const int index = this->getParameterIndexFromDetId(detId, run);
390    
391     return (index < 0 ? 0. : parameters_[index]);
392     }
393    
394     //======================================================================
395     int SiPixelLorentzAngleCalibration::getParameterIndexFromDetId(unsigned int detId,
396     edm::RunNumber_t run) const
397     {
398     // Return the index of the parameter that is used for this DetId.
399     // If this DetId is not treated, return values < 0.
400    
401     // FIXME: Extend to configurable granularity?
402     // Including treatment of run dependence?
403     const PXBDetId id(detId);
404     if (id.det() == DetId::Tracker && id.subdetId() == PixelSubdetector::PixelBarrel) {
405     if (id.module() >= 1 && id.module() <= 4) return 0;
406     if (id.module() >= 5 && id.module() <= 8) return 1;
407     edm::LogWarning("Alignment")
408     << "@SUB=SiPixelLorentzAngleCalibration::getParameterIndexFromDetId"
409     << "Module should be 1-8, but is " << id.module() << " => skip!";
410     }
411    
412     return -1;
413     }
414    
415     //======================================================================
416     unsigned int SiPixelLorentzAngleCalibration::numIovs() const
417     {
418     // FIXME: Needed to include treatment of run dependence!
419     return 1;
420     }
421    
422     //======================================================================
423     edm::RunNumber_t SiPixelLorentzAngleCalibration::firstRunOfIOV(unsigned int iovNum) const
424     {
425     // FIXME: Needed to include treatment of run dependence!
426     if (iovNum < this->numIovs()) return 1;
427     else return 0;
428     }
429    
430    
431     //======================================================================
432     void SiPixelLorentzAngleCalibration::writeTree(const SiPixelLorentzAngle *lorentzAngle,
433     const char *treeName) const
434     {
435     if (!lorentzAngle) return;
436    
437     TFile* file = TFile::Open(outFileName_.c_str(), "UPDATE");
438     if (!file) {
439     edm::LogError("BadConfig") << "@SUB=SiPixelLorentzAngleCalibration::writeTree"
440     << "Could not open file '" << outFileName_ << "'.";
441     return;
442     }
443    
444     TTree *tree = new TTree(treeName, treeName);
445     unsigned int id = 0;
446     float value = 0.;
447     tree->Branch("detId", &id, "detId/i");
448     tree->Branch("value", &value, "value/F");
449    
450     for (auto iterIdValue = lorentzAngle->getLorentzAngles().begin();
451     iterIdValue != lorentzAngle->getLorentzAngles().end(); ++iterIdValue) {
452     // type of (*iterIdValue) is pair<unsigned int, float>
453     id = iterIdValue->first; // key of map is DetId
454     value = iterIdValue->second;
455     tree->Fill();
456     }
457     tree->Write();
458     delete file; // tree vanishes with the file... (?)
459    
460     }
461    
462     //======================================================================
463     SiPixelLorentzAngle*
464     SiPixelLorentzAngleCalibration::createFromTree(const char *fileName, const char *treeName) const
465     {
466 flucke 1.2 // Check for file existence on your own to work around
467     // https://hypernews.cern.ch/HyperNews/CMS/get/swDevelopment/2715.html:
468     TFile* file = 0;
469     FILE* testFile = fopen(fileName,"r");
470     if (testFile) {
471     fclose(testFile);
472     file = TFile::Open(fileName, "READ");
473     } // else not existing, see error below
474    
475 flucke 1.1 TTree *tree = 0;
476     if (file) file->GetObject(treeName, tree);
477    
478     SiPixelLorentzAngle *result = 0;
479     if (tree) {
480     unsigned int id = 0;
481     float value = 0.;
482     tree->SetBranchAddress("detId", &id);
483     tree->SetBranchAddress("value", &value);
484    
485     result = new SiPixelLorentzAngle;
486     const Long64_t nEntries = tree->GetEntries();
487     for (Long64_t iEntry = 0; iEntry < nEntries; ++iEntry) {
488     tree->GetEntry(iEntry);
489     result->putLorentzAngle(id, value);
490     }
491 flucke 1.2 } else { // Warning only since could be parallel job on no events.
492     edm::LogWarning("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::createFromTree"
493     << "Could not get TTree '" << treeName << "' from file '"
494     << fileName << (file ? "'." : "' (file does not exist).");
495 flucke 1.1 }
496    
497     delete file; // tree will vanish with file
498     return result;
499     }
500    
501    
502     //======================================================================
503     //======================================================================
504     // Plugin definition
505    
506     #include "Alignment/CommonAlignmentAlgorithm/interface/IntegratedCalibrationPluginFactory.h"
507    
508     DEFINE_EDM_PLUGIN(IntegratedCalibrationPluginFactory,
509     SiPixelLorentzAngleCalibration, "SiPixelLorentzAngleCalibration");