ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/CMSSW/Alignment/CommonAlignmentAlgorithm/plugins/SiPixelLorentzAngleCalibration.cc
Revision: 1.5
Committed: Mon Jan 7 20:56:25 2013 UTC (12 years, 4 months ago) by wmtan
Content type: text/plain
Branch: MAIN
CVS Tags: CMSSW_6_2_0_pre7_TS133806, CMSSW_6_2_0_pre7_TS132947, CMSSW_6_2_0_pre7_g496p02, CMSSW_6_2_0_pre7, CMSSW_6_2_0_pre6_patch1, CMSSW_6_2_0_pre6, CMSSW_6_2_0_pre5slc6, CMSSW_6_2_0_pre5, CMSSW_6_2_0_pre4, CMSSW_6_2_0_pre3, CMSSW_6_2_0_pre2, CMSSW_6_2_0_pre1, V04-00-11
Changes since 1.4: +35 -14 lines
Log Message:
Get geometry from the EventSetup system, rather than having it hard coded in the data formats

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