ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/CMSSW/Alignment/CommonAlignmentAlgorithm/plugins/SiPixelLorentzAngleCalibration.cc
Revision: 1.4
Committed: Thu Sep 20 13:17:23 2012 UTC (12 years, 7 months ago) by flucke
Content type: text/plain
Branch: MAIN
CVS Tags: CMSSW_6_1_2_SLHC6_patch1, CMSSW_6_1_2_SLHC6, CMSSW_6_1_2_SLHC5, CMSSW_6_1_2_SLHC4_patch1, CMSSW_6_1_2_SLHC4, CMSSW_6_1_2_SLHC2_patch3, CMSSW_6_1_2_SLHC2_patch2, CMSSW_6_1_2_SLHC3, CMSSW_6_1_2_SLHC2_patch1, CMSSW_6_1_2_SLHC2, CMSSW_6_1_2_SLHC1, CMSSW_6_1_X_2012-12-19-0200, CMSSW_6_1_2, CMSSW_6_1_1_SLHCphase2tk1, CMSSW_6_1_1_SLHCphase1tk1, CMSSW_6_1_1, CMSSW_6_1_0, CMSSW_6_1_0_pre8, V04-00-10, CMSSW_6_1_0_pre7_TS127013, CMSSW_6_1_0_pre7, V04-00-09, V04-00-08
Branch point for: branch53X_calibration
Changes since 1.3: +5 -4 lines
Log Message:
Correct derivative by factor -0.5:
* 0.5 since LA shift is only half of maximum drift distance
* -1. since we need the derivative of the residual r = x_trk - x_hit,
  so everything affecting x_hit gets the sign flip

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.4 /// $Revision: 1.3 $
11     /// $Date: 2012/09/19 14:11:49 $
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 flucke 1.4 // shift due to LA: dx = tan(LA) * dz/2 = mobility * B_y * dz/2,
196     // '-' since we have derivative of the residual r = trk -hit
197     const double xDerivative = bFieldLocal.y() * dZ * -0.5; // parameter is mobility!
198 flucke 1.1 if (xDerivative) { // If field is zero, this is zero: do not return it
199     const Values derivs(xDerivative, 0.); // yDerivative = 0.
200     outDerivInds.push_back(ValuesIndexPair(derivs, index));
201     }
202     }
203     } else {
204     edm::LogWarning("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::derivatives2"
205     << "Hit without GeomDet, skip!";
206     }
207    
208     return outDerivInds.size();
209     }
210    
211     //======================================================================
212     bool SiPixelLorentzAngleCalibration::setParameter(unsigned int index, double value)
213     {
214     if (index >= parameters_.size()) {
215     return false;
216     } else {
217     parameters_[index] = value;
218     return true;
219     }
220     }
221    
222     //======================================================================
223     bool SiPixelLorentzAngleCalibration::setParameterError(unsigned int index, double error)
224     {
225     if (index >= paramUncertainties_.size()) {
226     return false;
227     } else {
228     paramUncertainties_[index] = error;
229     return true;
230     }
231     }
232    
233     //======================================================================
234     double SiPixelLorentzAngleCalibration::getParameter(unsigned int index) const
235     {
236     // if (index >= parameters_.size()) {
237     // return 0.;
238     // } else {
239     // return parameters_[index];
240     // }
241     return (index >= parameters_.size() ? 0. : parameters_[index]);
242     }
243    
244     //======================================================================
245     double SiPixelLorentzAngleCalibration::getParameterError(unsigned int index) const
246     {
247     // if (index >= paramUncertainties_.size()) {
248     // return 0.;
249     // } else {
250     // return paramUncertainties_[index];
251     // }
252     return (index >= paramUncertainties_.size() ? 0. : paramUncertainties_[index]);
253     }
254    
255     // //======================================================================
256     // void SiPixelLorentzAngleCalibration::beginOfJob(const AlignableTracker *tracker,
257     // const AlignableMuon */*muon*/,
258     // const AlignableExtras */*extras*/)
259     // {
260     // alignableTracker_ = tracker;
261     // }
262    
263    
264     //======================================================================
265     void SiPixelLorentzAngleCalibration::endOfJob()
266     {
267     // loginfo output
268     std::ostringstream out;
269     out << "Parameter results\n";
270     for (unsigned int iPar = 0; iPar < parameters_.size(); ++iPar) {
271     out << iPar << ": " << parameters_[iPar] << " +- " << paramUncertainties_[iPar] << "\n";
272     }
273     edm::LogInfo("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::endOfJob" << out.str();
274    
275     // now write 'input' tree
276     const SiPixelLorentzAngle *input = this->getLorentzAnglesInput(); // never NULL
277     const std::string treeName(this->name() + '_');
278     this->writeTree(input, (treeName + "input").c_str());
279 flucke 1.2 if (input->getLorentzAngles().empty()) {
280     edm::LogError("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::endOfJob"
281     << "Input Lorentz angle map is empty, skip writing output!";
282     return;
283     }
284 flucke 1.1
285     for (unsigned int iIOV = 0; iIOV < this->numIovs(); ++iIOV) {
286     cond::Time_t firstRunOfIOV = this->firstRunOfIOV(iIOV);
287     SiPixelLorentzAngle *output = new SiPixelLorentzAngle;
288     // Loop on map of values from input and add (possible) parameter results
289     for (auto iterIdValue = input->getLorentzAngles().begin();
290     iterIdValue != input->getLorentzAngles().end(); ++iterIdValue) {
291     // type of (*iterIdValue) is pair<unsigned int, float>
292     const unsigned int detId = iterIdValue->first; // key of map is DetId
293     // Nasty: putLorentzAngle(..) takes float by reference - not even const reference!
294     float value = iterIdValue->second + this->getParameterForDetId(detId, firstRunOfIOV);
295     output->putLorentzAngle(detId, value); // put result in output
296     }
297    
298     // Write this even for mille jobs?
299     this->writeTree(output, (treeName + Form("result_%lld", firstRunOfIOV)).c_str());
300    
301     if (saveToDB_) { // If requested, write out to DB
302     edm::Service<cond::service::PoolDBOutputService> dbService;
303     if (dbService.isAvailable()) {
304     dbService->writeOne(output, firstRunOfIOV, recordNameDBwrite_.c_str());
305     // no 'delete output;': writeOne(..) took over ownership
306     } else {
307     delete output;
308     edm::LogError("BadConfig") << "@SUB=SiPixelLorentzAngleCalibration::endOfJob"
309     << "No PoolDBOutputService available, but saveToDB true!";
310     }
311     } else {
312     delete output;
313     }
314     } // end loop on IOVs
315    
316     }
317    
318     //======================================================================
319     bool SiPixelLorentzAngleCalibration::checkLorentzAngleInput(const edm::EventSetup &setup,
320     const EventInfo &eventInfo)
321     {
322     edm::ESHandle<SiPixelLorentzAngle> lorentzAngleHandle;
323     if (!siPixelLorentzAngleInput_) {
324     setup.get<SiPixelLorentzAngleRcd>().get(lorentzAngleHandle);
325     siPixelLorentzAngleInput_ = new SiPixelLorentzAngle(*lorentzAngleHandle);
326     } else {
327     if (watchLorentzAngleRcd_.check(setup)) { // new IOV of input
328     setup.get<SiPixelLorentzAngleRcd>().get(lorentzAngleHandle);
329     if (lorentzAngleHandle->getLorentzAngles() // but only bad if non-identical values
330     != siPixelLorentzAngleInput_->getLorentzAngles()) { // (comparing maps)
331     // FIXME: Could different maps have same content, but different order?
332     // Or 'floating point comparison' problem?
333     throw cms::Exception("BadInput")
334     << "SiPixelLorentzAngleCalibration::checkLorentzAngleInput:\n"
335     << "Content of SiPixelLorentzAngle changed at run " << eventInfo.eventId_.run()
336     << ", but algorithm expects constant input!\n";
337     return false; // not reached...
338     }
339     }
340     }
341    
342     return true;
343     }
344    
345     //======================================================================
346     const SiPixelLorentzAngle* SiPixelLorentzAngleCalibration::getLorentzAnglesInput()
347     {
348     // For parallel processing in Millepede II, create SiPixelLorentzAngle
349     // from info stored in files of parallel jobs and check that they are identical.
350     // If this job has run on data, still check that LA is identical to the ones
351     // from mergeFileNames_.
352     const std::string treeName(this->name() + "_input");
353     for (auto iFile = mergeFileNames_.begin(); iFile != mergeFileNames_.end(); ++iFile) {
354     SiPixelLorentzAngle* la = this->createFromTree(iFile->c_str(), treeName.c_str());
355     // siPixelLorentzAngleInput_ could be non-null from previous file of this loop
356     // or from checkLorentzAngleInput(..) when running on data in this job as well
357     if (!siPixelLorentzAngleInput_ || siPixelLorentzAngleInput_->getLorentzAngles().empty()) {
358     delete siPixelLorentzAngleInput_; // NULL or empty
359     siPixelLorentzAngleInput_ = la;
360     } else {
361     // FIXME: about comparison of maps see comments in checkLorentzAngleInput
362 flucke 1.3 if (la && !la->getLorentzAngles().empty() && // single job might not have got events
363 flucke 1.1 la->getLorentzAngles() != siPixelLorentzAngleInput_->getLorentzAngles()) {
364     // Throw exception instead of error?
365     edm::LogError("NoInput") << "@SUB=SiPixelLorentzAngleCalibration::getLorentzAnglesInput"
366     << "Different input values from tree " << treeName
367     << " in file " << *iFile << ".";
368    
369     }
370     delete la;
371     }
372     }
373    
374     if (!siPixelLorentzAngleInput_) { // no files nor ran on events
375     siPixelLorentzAngleInput_ = new SiPixelLorentzAngle;
376     edm::LogError("NoInput") << "@SUB=SiPixelLorentzAngleCalibration::getLorentzAnglesInput"
377 flucke 1.2 << "No input, create an empty one!";
378 flucke 1.1 } else if (siPixelLorentzAngleInput_->getLorentzAngles().empty()) {
379     edm::LogError("NoInput") << "@SUB=SiPixelLorentzAngleCalibration::getLorentzAnglesInput"
380 flucke 1.2 << "Empty result!";
381 flucke 1.1 }
382    
383     return siPixelLorentzAngleInput_;
384     }
385    
386     //======================================================================
387     double SiPixelLorentzAngleCalibration::getParameterForDetId(unsigned int detId,
388     edm::RunNumber_t run) const
389     {
390     const int index = this->getParameterIndexFromDetId(detId, run);
391    
392     return (index < 0 ? 0. : parameters_[index]);
393     }
394    
395     //======================================================================
396     int SiPixelLorentzAngleCalibration::getParameterIndexFromDetId(unsigned int detId,
397     edm::RunNumber_t run) const
398     {
399     // Return the index of the parameter that is used for this DetId.
400     // If this DetId is not treated, return values < 0.
401    
402     // FIXME: Extend to configurable granularity?
403     // Including treatment of run dependence?
404     const PXBDetId id(detId);
405     if (id.det() == DetId::Tracker && id.subdetId() == PixelSubdetector::PixelBarrel) {
406     if (id.module() >= 1 && id.module() <= 4) return 0;
407     if (id.module() >= 5 && id.module() <= 8) return 1;
408     edm::LogWarning("Alignment")
409     << "@SUB=SiPixelLorentzAngleCalibration::getParameterIndexFromDetId"
410     << "Module should be 1-8, but is " << id.module() << " => skip!";
411     }
412    
413     return -1;
414     }
415    
416     //======================================================================
417     unsigned int SiPixelLorentzAngleCalibration::numIovs() const
418     {
419     // FIXME: Needed to include treatment of run dependence!
420     return 1;
421     }
422    
423     //======================================================================
424     edm::RunNumber_t SiPixelLorentzAngleCalibration::firstRunOfIOV(unsigned int iovNum) const
425     {
426     // FIXME: Needed to include treatment of run dependence!
427     if (iovNum < this->numIovs()) return 1;
428     else return 0;
429     }
430    
431    
432     //======================================================================
433     void SiPixelLorentzAngleCalibration::writeTree(const SiPixelLorentzAngle *lorentzAngle,
434     const char *treeName) const
435     {
436     if (!lorentzAngle) return;
437    
438     TFile* file = TFile::Open(outFileName_.c_str(), "UPDATE");
439     if (!file) {
440     edm::LogError("BadConfig") << "@SUB=SiPixelLorentzAngleCalibration::writeTree"
441     << "Could not open file '" << outFileName_ << "'.";
442     return;
443     }
444    
445     TTree *tree = new TTree(treeName, treeName);
446     unsigned int id = 0;
447     float value = 0.;
448     tree->Branch("detId", &id, "detId/i");
449     tree->Branch("value", &value, "value/F");
450    
451     for (auto iterIdValue = lorentzAngle->getLorentzAngles().begin();
452     iterIdValue != lorentzAngle->getLorentzAngles().end(); ++iterIdValue) {
453     // type of (*iterIdValue) is pair<unsigned int, float>
454     id = iterIdValue->first; // key of map is DetId
455     value = iterIdValue->second;
456     tree->Fill();
457     }
458     tree->Write();
459     delete file; // tree vanishes with the file... (?)
460    
461     }
462    
463     //======================================================================
464     SiPixelLorentzAngle*
465     SiPixelLorentzAngleCalibration::createFromTree(const char *fileName, const char *treeName) const
466     {
467 flucke 1.2 // Check for file existence on your own to work around
468     // https://hypernews.cern.ch/HyperNews/CMS/get/swDevelopment/2715.html:
469     TFile* file = 0;
470     FILE* testFile = fopen(fileName,"r");
471     if (testFile) {
472     fclose(testFile);
473     file = TFile::Open(fileName, "READ");
474     } // else not existing, see error below
475    
476 flucke 1.1 TTree *tree = 0;
477     if (file) file->GetObject(treeName, tree);
478    
479     SiPixelLorentzAngle *result = 0;
480     if (tree) {
481     unsigned int id = 0;
482     float value = 0.;
483     tree->SetBranchAddress("detId", &id);
484     tree->SetBranchAddress("value", &value);
485    
486     result = new SiPixelLorentzAngle;
487     const Long64_t nEntries = tree->GetEntries();
488     for (Long64_t iEntry = 0; iEntry < nEntries; ++iEntry) {
489     tree->GetEntry(iEntry);
490     result->putLorentzAngle(id, value);
491     }
492 flucke 1.2 } else { // Warning only since could be parallel job on no events.
493     edm::LogWarning("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::createFromTree"
494     << "Could not get TTree '" << treeName << "' from file '"
495     << fileName << (file ? "'." : "' (file does not exist).");
496 flucke 1.1 }
497    
498     delete file; // tree will vanish with file
499     return result;
500     }
501    
502    
503     //======================================================================
504     //======================================================================
505     // Plugin definition
506    
507     #include "Alignment/CommonAlignmentAlgorithm/interface/IntegratedCalibrationPluginFactory.h"
508    
509     DEFINE_EDM_PLUGIN(IntegratedCalibrationPluginFactory,
510     SiPixelLorentzAngleCalibration, "SiPixelLorentzAngleCalibration");