ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/CMSSW/Alignment/CommonAlignmentAlgorithm/plugins/SiStripLorentzAngleCalibration.cc
Revision: 1.4
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.3: +3 -4 lines
Log Message:
fix possible crash introduced in previous version

File Contents

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