ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/CMSSW/Alignment/CommonAlignmentAlgorithm/plugins/SiPixelLorentzAngleCalibration.cc
Revision: 1.6
Committed: Fri May 31 12:13:40 2013 UTC (11 years, 11 months ago) by flucke
Content type: text/plain
Branch: MAIN
CVS Tags: CMSSW_6_2_0, CMSSW_6_2_0_pre8, V04-00-14, V04-00-13, V04-00-09-53X-calib-V01, V04-00-12, HEAD
Changes since 1.5: +101 -120 lines
Log Message:
merging calibration development from branch 'branch53X_calibration' (tag
'br53-00-09') to HEAD:
- add backplane calibration
- configurable granularity in time (i.e. run) and space for lorentz angle
  calibrations
- more info in diagnostics tree

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.6 /// $Revision: 1.4.2.21 $
11     /// $Date: 2013/05/31 08:37:12 $
12 flucke 1.1 /// (last update by $Author: flucke $)
13    
14     #include "Alignment/CommonAlignmentAlgorithm/interface/IntegratedCalibrationBase.h"
15 flucke 1.6 #include "Alignment/CommonAlignmentAlgorithm/interface/TkModuleGroupSelector.h"
16     // include 'locally':
17     #include "TreeStruct.h"
18 flucke 1.1
19     #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
20     #include "CondFormats/DataRecord/interface/SiPixelLorentzAngleRcd.h"
21     #include "CondFormats/SiPixelObjects/interface/SiPixelLorentzAngle.h"
22    
23     #include "DataFormats/DetId/interface/DetId.h"
24 wmtan 1.5 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
25 flucke 1.1
26     #include "FWCore/Framework/interface/ESHandle.h"
27     #include "FWCore/Framework/interface/ESWatcher.h"
28     #include "FWCore/MessageLogger/interface/MessageLogger.h"
29 flucke 1.6 #include "Alignment/CommonAlignment/interface/Alignable.h"
30    
31 flucke 1.1 #include "FWCore/ServiceRegistry/interface/Service.h"
32 flucke 1.6 #include "FWCore/ParameterSet/interface/ParameterSet.h"
33     #include "Alignment/CommonAlignmentAlgorithm/interface/AlignmentParameterSelector.h"
34 flucke 1.1
35     #include "Geometry/CommonDetUnit/interface/GeomDet.h"
36     #include "MagneticField/Engine/interface/MagneticField.h"
37     #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
38     #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
39    
40     #include "TTree.h"
41     #include "TFile.h"
42     #include "TString.h"
43    
44 flucke 1.6 #include <boost/assign/list_of.hpp>
45 flucke 1.1 #include <vector>
46 flucke 1.6 #include <map>
47 flucke 1.1 #include <sstream>
48 flucke 1.6 #include <list>
49     #include <utility>
50     #include <set>
51    
52     #include <functional>
53 flucke 1.1
54     class SiPixelLorentzAngleCalibration : public IntegratedCalibrationBase
55     {
56     public:
57     /// Constructor
58     explicit SiPixelLorentzAngleCalibration(const edm::ParameterSet &cfg);
59    
60     /// Destructor
61     virtual ~SiPixelLorentzAngleCalibration();
62    
63     /// How many parameters does this calibration define?
64     virtual unsigned int numParameters() const;
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    
74     /// Return non-zero derivatives for x- and y-measurements with their indices by reference.
75     /// Return value is their number.
76     virtual unsigned int derivatives(std::vector<ValuesIndexPair> &outDerivInds,
77     const TransientTrackingRecHit &hit,
78     const TrajectoryStateOnSurface &tsos,
79     const edm::EventSetup &setup,
80     const EventInfo &eventInfo) const;
81    
82     /// Setting the determined parameter identified by index,
83 flucke 1.6 /// returns false if out-of-bounds, true otherwise.
84 flucke 1.1 virtual bool setParameter(unsigned int index, double value);
85    
86     /// Setting the determined parameter uncertainty identified by index,
87     /// returns false if out-of-bounds, true otherwise.
88     virtual bool setParameterError(unsigned int index, double error);
89    
90     /// Return current value of parameter identified by index.
91 flucke 1.6 /// Return 0. if index out-of-bounds.
92 flucke 1.1 virtual double getParameter(unsigned int index) const;
93    
94     /// Return current value of parameter identified by index.
95     /// Returns 0. if index out-of-bounds or if errors undetermined.
96     virtual double getParameterError(unsigned int index) const;
97    
98     // /// Call at beginning of job:
99 flucke 1.6 virtual void beginOfJob(AlignableTracker *tracker,
100     AlignableMuon *muon,
101     AlignableExtras *extras);
102    
103    
104 flucke 1.1
105     /// Called at end of a the job of the AlignmentProducer.
106     /// Write out determined parameters.
107     virtual void endOfJob();
108    
109     private:
110     /// If called the first time, fill 'siPixelLorentzAngleInput_',
111     /// later check that LorentzAngle has not changed.
112     bool checkLorentzAngleInput(const edm::EventSetup &setup, const EventInfo &eventInfo);
113     /// Input LorentzAngle values:
114     /// - either from EventSetup of first call to derivatives(..)
115     /// - or created from files of passed by configuration (i.e. from parallel processing)
116     const SiPixelLorentzAngle* getLorentzAnglesInput();
117    
118     /// Determined parameter value for this detId (detId not treated => 0.)
119     /// and the given run.
120 flucke 1.6 double getParameterForDetId(unsigned int detId, edm::RunNumber_t run) const;
121 flucke 1.1
122 flucke 1.6 void writeTree(const SiPixelLorentzAngle *lorentzAngle,
123     const std::map<unsigned int,TreeStruct>& treeInfo, const char *treeName) const;
124 flucke 1.1 SiPixelLorentzAngle* createFromTree(const char *fileName, const char *treeName) const;
125 flucke 1.6
126 flucke 1.1 const bool saveToDB_;
127     const std::string recordNameDBwrite_;
128     const std::string outFileName_;
129     const std::vector<std::string> mergeFileNames_;
130    
131     edm::ESWatcher<SiPixelLorentzAngleRcd> watchLorentzAngleRcd_;
132    
133     // const AlignableTracker *alignableTracker_;
134     SiPixelLorentzAngle *siPixelLorentzAngleInput_;
135     std::vector<double> parameters_;
136     std::vector<double> paramUncertainties_;
137 wmtan 1.5
138 flucke 1.6 TkModuleGroupSelector *moduleGroupSelector_;
139     const edm::ParameterSet moduleGroupSelCfg_;
140 flucke 1.1 };
141    
142     //======================================================================
143     //======================================================================
144     //======================================================================
145    
146     SiPixelLorentzAngleCalibration::SiPixelLorentzAngleCalibration(const edm::ParameterSet &cfg)
147     : IntegratedCalibrationBase(cfg),
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 wmtan 1.5 siPixelLorentzAngleInput_(0),
153 flucke 1.6 moduleGroupSelector_(0),
154     moduleGroupSelCfg_(cfg.getParameter<edm::ParameterSet>("LorentzAngleModuleGroups"))
155 flucke 1.1 {
156    
157     }
158 wmtan 1.5
159 flucke 1.1 //======================================================================
160     SiPixelLorentzAngleCalibration::~SiPixelLorentzAngleCalibration()
161     {
162 flucke 1.6 delete moduleGroupSelector_;
163 flucke 1.1 // std::cout << "Destroy SiPixelLorentzAngleCalibration named " << this->name() << std::endl;
164     delete siPixelLorentzAngleInput_;
165     }
166    
167     //======================================================================
168     unsigned int SiPixelLorentzAngleCalibration::numParameters() const
169     {
170     return parameters_.size();
171     }
172    
173     //======================================================================
174     unsigned int
175     SiPixelLorentzAngleCalibration::derivatives(std::vector<ValuesIndexPair> &outDerivInds,
176     const TransientTrackingRecHit &hit,
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    
189 flucke 1.6 const int index = moduleGroupSelector_->getParameterIndexFromDetId(hit.det()->geographicalId(),
190     eventInfo.eventId_.run());
191 flucke 1.1 if (index >= 0) { // otherwise not treated
192     edm::ESHandle<MagneticField> magneticField;
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 flucke 1.4 // 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 flucke 1.6 // FIXME: Have to treat that for FPIX yDerivative != 0., due to higher order effects!
201 flucke 1.1 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     //======================================================================
215     bool SiPixelLorentzAngleCalibration::setParameter(unsigned int index, double value)
216     {
217     if (index >= parameters_.size()) {
218     return false;
219     } else {
220     parameters_[index] = value;
221     return true;
222     }
223     }
224    
225     //======================================================================
226     bool SiPixelLorentzAngleCalibration::setParameterError(unsigned int index, double error)
227     {
228     if (index >= paramUncertainties_.size()) {
229     return false;
230     } else {
231     paramUncertainties_[index] = error;
232     return true;
233     }
234     }
235    
236     //======================================================================
237     double SiPixelLorentzAngleCalibration::getParameter(unsigned int index) const
238     {
239     return (index >= parameters_.size() ? 0. : parameters_[index]);
240     }
241    
242     //======================================================================
243     double SiPixelLorentzAngleCalibration::getParameterError(unsigned int index) const
244     {
245     return (index >= paramUncertainties_.size() ? 0. : paramUncertainties_[index]);
246     }
247    
248 flucke 1.6
249    
250    
251     //======================================================================
252     void SiPixelLorentzAngleCalibration::beginOfJob(AlignableTracker *aliTracker,
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    
259     moduleGroupSelector_ = new TkModuleGroupSelector(aliTracker, moduleGroupSelCfg_, sdets);
260    
261     parameters_.resize(moduleGroupSelector_->getNumberOfParameters(), 0.);
262     paramUncertainties_.resize(moduleGroupSelector_->getNumberOfParameters(), 0.);
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 flucke 1.1
277    
278     //======================================================================
279     void SiPixelLorentzAngleCalibration::endOfJob()
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 flucke 1.6 std::map<unsigned int, TreeStruct> treeInfo; // map of TreeStruct for each detId
290    
291 flucke 1.1 // now write 'input' tree
292     const SiPixelLorentzAngle *input = this->getLorentzAnglesInput(); // never NULL
293     const std::string treeName(this->name() + '_');
294 flucke 1.6 this->writeTree(input, treeInfo, (treeName + "input").c_str()); // empty treeInfo for input...
295    
296 flucke 1.2 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 flucke 1.1
302 flucke 1.6 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);
310 flucke 1.1 SiPixelLorentzAngle *output = new SiPixelLorentzAngle;
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 flucke 1.6 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 flucke 1.1 }
324    
325 flucke 1.6 if (saveToDB_ || nonZeroParamsOrErrors != 0) { // Skip writing mille jobs...
326     this->writeTree(output, treeInfo, (treeName + Form("result_%lld", firstRunOfIOV)).c_str());
327     }
328 flucke 1.1
329     if (saveToDB_) { // If requested, write out to DB
330     edm::Service<cond::service::PoolDBOutputService> dbService;
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     //======================================================================
347     bool SiPixelLorentzAngleCalibration::checkLorentzAngleInput(const edm::EventSetup &setup,
348     const EventInfo &eventInfo)
349     {
350     edm::ESHandle<SiPixelLorentzAngle> lorentzAngleHandle;
351     if (!siPixelLorentzAngleInput_) {
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 flucke 1.6 // Maps are containers sorted by key, but comparison problems may arise from
360     // 'floating point comparison' problems (FIXME?)
361 flucke 1.1 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     //======================================================================
374     const SiPixelLorentzAngle* SiPixelLorentzAngleCalibration::getLorentzAnglesInput()
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
385     if (!siPixelLorentzAngleInput_ || siPixelLorentzAngleInput_->getLorentzAngles().empty()) {
386     delete siPixelLorentzAngleInput_; // NULL or empty
387     siPixelLorentzAngleInput_ = la;
388     } else {
389     // FIXME: about comparison of maps see comments in checkLorentzAngleInput
390 flucke 1.3 if (la && !la->getLorentzAngles().empty() && // single job might not have got events
391 flucke 1.1 la->getLorentzAngles() != siPixelLorentzAngleInput_->getLorentzAngles()) {
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
403     siPixelLorentzAngleInput_ = new SiPixelLorentzAngle;
404     edm::LogError("NoInput") << "@SUB=SiPixelLorentzAngleCalibration::getLorentzAnglesInput"
405 flucke 1.2 << "No input, create an empty one!";
406 flucke 1.1 } else if (siPixelLorentzAngleInput_->getLorentzAngles().empty()) {
407     edm::LogError("NoInput") << "@SUB=SiPixelLorentzAngleCalibration::getLorentzAnglesInput"
408 flucke 1.2 << "Empty result!";
409 flucke 1.1 }
410    
411     return siPixelLorentzAngleInput_;
412     }
413    
414     //======================================================================
415     double SiPixelLorentzAngleCalibration::getParameterForDetId(unsigned int detId,
416     edm::RunNumber_t run) const
417     {
418 flucke 1.6 const int index = moduleGroupSelector_->getParameterIndexFromDetId(detId, run);
419 flucke 1.1 return (index < 0 ? 0. : parameters_[index]);
420     }
421    
422     //======================================================================
423     void SiPixelLorentzAngleCalibration::writeTree(const SiPixelLorentzAngle *lorentzAngle,
424 flucke 1.6 const std::map<unsigned int,TreeStruct> &treeInfo,
425 flucke 1.1 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 flucke 1.6 TreeStruct treeStruct;
440 flucke 1.1 tree->Branch("detId", &id, "detId/i");
441     tree->Branch("value", &value, "value/F");
442 flucke 1.6 tree->Branch("treeStruct", &treeStruct, TreeStruct::LeafList());
443 flucke 1.1
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 flucke 1.6 // 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 flucke 1.1 tree->Fill();
459     }
460     tree->Write();
461     delete file; // tree vanishes with the file... (?)
462    
463     }
464    
465     //======================================================================
466     SiPixelLorentzAngle*
467     SiPixelLorentzAngleCalibration::createFromTree(const char *fileName, const char *treeName) const
468     {
469 flucke 1.2 // 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 flucke 1.1 TTree *tree = 0;
479     if (file) file->GetObject(treeName, tree);
480    
481     SiPixelLorentzAngle *result = 0;
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 flucke 1.2 } 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 flucke 1.1 }
499    
500     delete file; // tree will vanish with file
501     return result;
502     }
503    
504    
505     //======================================================================
506     //======================================================================
507     // Plugin definition
508    
509     #include "Alignment/CommonAlignmentAlgorithm/interface/IntegratedCalibrationPluginFactory.h"
510    
511     DEFINE_EDM_PLUGIN(IntegratedCalibrationPluginFactory,
512 flucke 1.6 SiPixelLorentzAngleCalibration, "SiPixelLorentzAngleCalibration");