Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:56:08

0001 /// \class SiPixelLorentzAngleCalibration
0002 ///
0003 /// Calibration of Lorentz angle for the pixel,
0004 /// integrated in the alignment algorithms.
0005 ///
0006 /// Note that not all algorithms support this...
0007 ///
0008 ///  \author    : Gero Flucke
0009 ///  date       : September 2012
0010 ///  $Revision: 1.4.2.21 $
0011 ///  $Date: 2013/05/31 08:37:12 $
0012 ///  (last update by $Author: flucke $)
0013 
0014 #include "Alignment/CommonAlignmentAlgorithm/interface/IntegratedCalibrationBase.h"
0015 #include "Alignment/CommonAlignmentAlgorithm/interface/TkModuleGroupSelector.h"
0016 // include 'locally':
0017 #include "TreeStruct.h"
0018 
0019 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
0020 #include "CondFormats/DataRecord/interface/SiPixelLorentzAngleRcd.h"
0021 #include "CondFormats/SiPixelObjects/interface/SiPixelLorentzAngle.h"
0022 
0023 #include "DataFormats/DetId/interface/DetId.h"
0024 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0025 
0026 #include "FWCore/Framework/interface/ESWatcher.h"
0027 #include "FWCore/Framework/interface/Run.h"
0028 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0029 #include "Alignment/CommonAlignment/interface/Alignable.h"
0030 
0031 #include "FWCore/ServiceRegistry/interface/Service.h"
0032 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0033 #include "Alignment/CommonAlignmentAlgorithm/interface/AlignmentParameterSelector.h"
0034 
0035 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
0036 #include "MagneticField/Engine/interface/MagneticField.h"
0037 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0038 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
0039 #include "FWCore/Framework/interface/ConsumesCollector.h"
0040 
0041 #include "TTree.h"
0042 #include "TFile.h"
0043 #include "TString.h"
0044 
0045 #include <vector>
0046 #include <map>
0047 #include <sstream>
0048 #include <list>
0049 #include <utility>
0050 #include <set>
0051 #include <memory>
0052 #include <functional>
0053 
0054 class SiPixelLorentzAngleCalibration : public IntegratedCalibrationBase {
0055 public:
0056   /// Constructor
0057   explicit SiPixelLorentzAngleCalibration(const edm::ParameterSet &cfg, edm::ConsumesCollector &iC);
0058 
0059   /// Destructor
0060   ~SiPixelLorentzAngleCalibration() override = default;
0061 
0062   /// How many parameters does this calibration define?
0063   unsigned int numParameters() const override;
0064 
0065   /// Return non-zero derivatives for x- and y-measurements with their indices by reference.
0066   /// Return value is their number.
0067   unsigned int derivatives(std::vector<ValuesIndexPair> &outDerivInds,
0068                            const TransientTrackingRecHit &hit,
0069                            const TrajectoryStateOnSurface &tsos,
0070                            const edm::EventSetup &setup,
0071                            const EventInfo &eventInfo) const override;
0072 
0073   /// Setting the determined parameter identified by index,
0074   /// returns false if out-of-bounds, true otherwise.
0075   bool setParameter(unsigned int index, double value) override;
0076 
0077   /// Setting the determined parameter uncertainty identified by index,
0078   /// returns false if out-of-bounds, true otherwise.
0079   bool setParameterError(unsigned int index, double error) override;
0080 
0081   /// Return current value of parameter identified by index.
0082   /// Return 0. if index out-of-bounds.
0083   double getParameter(unsigned int index) const override;
0084 
0085   /// Return current value of parameter identified by index.
0086   /// Returns 0. if index out-of-bounds or if errors undetermined.
0087   double getParameterError(unsigned int index) const override;
0088 
0089   /// Call at beginning of job:
0090   void beginOfJob(AlignableTracker *tracker, AlignableMuon *muon, AlignableExtras *extras) override;
0091 
0092   /// Call at beginning of run:
0093   void beginRun(const edm::Run &, const edm::EventSetup &) override;
0094 
0095   /// Called at end of a the job of the AlignmentProducer.
0096   /// Write out determined parameters.
0097   void endOfJob() override;
0098 
0099 private:
0100   /// Input LorentzAngle values:
0101   /// - either from EventSetup of first call to derivatives(..)
0102   /// - or created from files of passed by configuration (i.e. from parallel processing)
0103   const SiPixelLorentzAngle *getLorentzAnglesInput(const align::RunNumber & = 0);
0104 
0105   /// Determined parameter value for this detId (detId not treated => 0.)
0106   /// and the given run.
0107   double getParameterForDetId(unsigned int detId, edm::RunNumber_t run) const;
0108 
0109   void writeTree(const SiPixelLorentzAngle *lorentzAngle,
0110                  const std::map<unsigned int, TreeStruct> &treeInfo,
0111                  const char *treeName) const;
0112   SiPixelLorentzAngle createFromTree(const char *fileName, const char *treeName) const;
0113 
0114   const bool saveToDB_;
0115   const std::string recordNameDBwrite_;
0116   const std::string outFileName_;
0117   const std::vector<std::string> mergeFileNames_;
0118   const std::string lorentzAngleLabel_;
0119   const edm::ESGetToken<SiPixelLorentzAngle, SiPixelLorentzAngleRcd> lorentzAngleToken_;
0120   const edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> magFieldToken_;
0121   edm::ESWatcher<SiPixelLorentzAngleRcd> watchLorentzAngleRcd_;
0122 
0123   // const AlignableTracker *alignableTracker_;
0124   std::map<align::RunNumber, SiPixelLorentzAngle> cachedLorentzAngleInputs_;
0125   SiPixelLorentzAngle *siPixelLorentzAngleInput_{nullptr};
0126   align::RunNumber currentIOV_{0};
0127   std::vector<double> parameters_;
0128   std::vector<double> paramUncertainties_;
0129 
0130   std::unique_ptr<TkModuleGroupSelector> moduleGroupSelector_;
0131   const edm::ParameterSet moduleGroupSelCfg_;
0132 };
0133 
0134 //======================================================================
0135 //======================================================================
0136 //======================================================================
0137 
0138 SiPixelLorentzAngleCalibration::SiPixelLorentzAngleCalibration(const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
0139     : IntegratedCalibrationBase(cfg),
0140       saveToDB_(cfg.getParameter<bool>("saveToDB")),
0141       recordNameDBwrite_(cfg.getParameter<std::string>("recordNameDBwrite")),
0142       outFileName_(cfg.getParameter<std::string>("treeFile")),
0143       mergeFileNames_(cfg.getParameter<std::vector<std::string> >("mergeTreeFiles")),
0144       lorentzAngleLabel_(cfg.getParameter<std::string>("lorentzAngleLabel")),
0145       lorentzAngleToken_(iC.esConsumes<edm::Transition::BeginRun>(edm::ESInputTag("", lorentzAngleLabel_))),
0146       magFieldToken_(iC.esConsumes()),
0147       moduleGroupSelCfg_(cfg.getParameter<edm::ParameterSet>("LorentzAngleModuleGroups")) {}
0148 
0149 //======================================================================
0150 unsigned int SiPixelLorentzAngleCalibration::numParameters() const { return parameters_.size(); }
0151 
0152 //======================================================================
0153 void SiPixelLorentzAngleCalibration::beginRun(const edm::Run &run, const edm::EventSetup &setup) {
0154   // no action needed if the LA record didn't change
0155   if (!(watchLorentzAngleRcd_.check(setup)))
0156     return;
0157 
0158   const auto runNumber = run.run();
0159   auto firstRun = cond::timeTypeSpecs[cond::runnumber].beginValue;
0160 
0161   // avoid infinite loop due to wrap-around of unsigned variable 'i' including
0162   // arrow from i to zero and a nice smiley ;)
0163   for (unsigned int i = moduleGroupSelector_->numIovs(); i-- > 0;) {
0164     const auto firstRunOfIOV = moduleGroupSelector_->firstRunOfIOV(i);
0165     if (runNumber >= firstRunOfIOV) {
0166       firstRun = firstRunOfIOV;
0167       break;
0168     }
0169   }
0170 
0171   const SiPixelLorentzAngle *lorentzAngleHandle = &setup.getData(lorentzAngleToken_);
0172   const auto &lorentzAngleRcd = setup.get<SiPixelLorentzAngleRcd>();
0173 
0174   if (cachedLorentzAngleInputs_.find(firstRun) == cachedLorentzAngleInputs_.end()) {
0175     cachedLorentzAngleInputs_.emplace(firstRun, SiPixelLorentzAngle(*lorentzAngleHandle));
0176   } else {
0177     if (lorentzAngleRcd.validityInterval().first().eventID().run() > firstRun &&
0178         lorentzAngleHandle->getLorentzAngles()                            // only bad if non-identical values
0179             != cachedLorentzAngleInputs_[firstRun].getLorentzAngles()) {  // (comparing maps)
0180       // Maps are containers sorted by key, but comparison problems may arise from
0181       // 'floating point comparison' problems (FIXME?)
0182       throw cms::Exception("BadInput") << "Trying to cache SiPixelLorentzAngle payload for a run (" << runNumber
0183                                        << ") in an IOV (" << firstRun << ") that was already cached.\n"
0184                                        << "The following record in your input database tag has an IOV "
0185                                        << "boundary that does not match your IOV definition:\n"
0186                                        << " - SiPixelLorentzAngleRcd '" << lorentzAngleRcd.key().name() << "' (since "
0187                                        << lorentzAngleRcd.validityInterval().first().eventID().run() << ")\n";
0188     }
0189   }
0190 
0191   siPixelLorentzAngleInput_ = &(cachedLorentzAngleInputs_[firstRun]);
0192   currentIOV_ = firstRun;
0193 }
0194 
0195 //======================================================================
0196 unsigned int SiPixelLorentzAngleCalibration::derivatives(std::vector<ValuesIndexPair> &outDerivInds,
0197                                                          const TransientTrackingRecHit &hit,
0198                                                          const TrajectoryStateOnSurface &tsos,
0199                                                          const edm::EventSetup &setup,
0200                                                          const EventInfo &eventInfo) const {
0201   outDerivInds.clear();
0202 
0203   if (hit.det()) {  // otherwise 'constraint hit' or whatever
0204 
0205     const int index =
0206         moduleGroupSelector_->getParameterIndexFromDetId(hit.det()->geographicalId(), eventInfo.eventId().run());
0207     if (index >= 0) {  // otherwise not treated
0208 
0209       const MagneticField *magneticField = &setup.getData(magFieldToken_);
0210       const GlobalVector bField(magneticField->inTesla(hit.det()->surface().position()));
0211       const LocalVector bFieldLocal(hit.det()->surface().toLocal(bField));
0212       const double dZ = hit.det()->surface().bounds().thickness();  // it is a float only...
0213       // shift due to LA: dx = tan(LA) * dz/2 = mobility * B_y * dz/2,
0214       // shift due to LA: dy = - mobility * B_x * dz/2,
0215       // '-' since we have derivative of the residual r = trk -hit
0216       const double xDerivative = bFieldLocal.y() * dZ * -0.5;  // parameter is mobility!
0217       const double yDerivative = bFieldLocal.x() * dZ * 0.5;   // parameter is mobility!
0218       if (xDerivative || yDerivative) {                        // If field is zero, this is zero: do not return it
0219         const Values derivs{xDerivative, yDerivative};
0220         outDerivInds.push_back(ValuesIndexPair(derivs, index));
0221       }
0222     }
0223   } else {
0224     edm::LogWarning("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::derivatives2"
0225                                  << "Hit without GeomDet, skip!";
0226   }
0227 
0228   return outDerivInds.size();
0229 }
0230 
0231 //======================================================================
0232 bool SiPixelLorentzAngleCalibration::setParameter(unsigned int index, double value) {
0233   if (index >= parameters_.size()) {
0234     return false;
0235   } else {
0236     parameters_[index] = value;
0237     return true;
0238   }
0239 }
0240 
0241 //======================================================================
0242 bool SiPixelLorentzAngleCalibration::setParameterError(unsigned int index, double error) {
0243   if (index >= paramUncertainties_.size()) {
0244     return false;
0245   } else {
0246     paramUncertainties_[index] = error;
0247     return true;
0248   }
0249 }
0250 
0251 //======================================================================
0252 double SiPixelLorentzAngleCalibration::getParameter(unsigned int index) const {
0253   return (index >= parameters_.size() ? 0. : parameters_[index]);
0254 }
0255 
0256 //======================================================================
0257 double SiPixelLorentzAngleCalibration::getParameterError(unsigned int index) const {
0258   return (index >= paramUncertainties_.size() ? 0. : paramUncertainties_[index]);
0259 }
0260 
0261 //======================================================================
0262 void SiPixelLorentzAngleCalibration::beginOfJob(AlignableTracker *aliTracker,
0263                                                 AlignableMuon * /*aliMuon*/,
0264                                                 AlignableExtras * /*aliExtras*/) {
0265   //specify the sub-detectors for which the LA is determined
0266   const std::vector<int> sdets = {PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap};
0267 
0268   moduleGroupSelector_ = std::make_unique<TkModuleGroupSelector>(aliTracker, moduleGroupSelCfg_, sdets);
0269 
0270   parameters_.resize(moduleGroupSelector_->getNumberOfParameters(), 0.);
0271   paramUncertainties_.resize(moduleGroupSelector_->getNumberOfParameters(), 0.);
0272 
0273   edm::LogInfo("Alignment") << "@SUB=SiPixelLorentzAngleCalibration"
0274                             << "Created with name " << this->name() << "',\n"
0275                             << this->numParameters() << " parameters to be determined,"
0276                             << "\n saveToDB = " << saveToDB_ << "\n outFileName = " << outFileName_
0277                             << "\n N(merge files) = " << mergeFileNames_.size()
0278                             << "\n number of IOVs = " << moduleGroupSelector_->numIovs();
0279 
0280   if (!mergeFileNames_.empty()) {
0281     edm::LogInfo("Alignment") << "@SUB=SiPixelLorentzAngleCalibration"
0282                               << "First file to merge: " << mergeFileNames_[0];
0283   }
0284 }
0285 
0286 //======================================================================
0287 void SiPixelLorentzAngleCalibration::endOfJob() {
0288   // loginfo output
0289   std::ostringstream out;
0290   out << "Parameter results\n";
0291   for (unsigned int iPar = 0; iPar < parameters_.size(); ++iPar) {
0292     out << iPar << ": " << parameters_[iPar] << " +- " << paramUncertainties_[iPar] << "\n";
0293   }
0294   edm::LogInfo("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::endOfJob" << out.str();
0295 
0296   std::map<unsigned int, TreeStruct> treeInfo;  // map of TreeStruct for each detId
0297 
0298   // now write 'input' tree(s)
0299   const std::string treeName{this->name() + '_'};
0300   std::vector<const SiPixelLorentzAngle *> inputs{};
0301   inputs.reserve(moduleGroupSelector_->numIovs());
0302   for (unsigned int iIOV = 0; iIOV < moduleGroupSelector_->numIovs(); ++iIOV) {
0303     const auto firstRunOfIOV = moduleGroupSelector_->firstRunOfIOV(iIOV);
0304     inputs.push_back(this->getLorentzAnglesInput(firstRunOfIOV));  // never NULL
0305     this->writeTree(inputs.back(),
0306                     treeInfo,
0307                     (treeName + "input_" + std::to_string(firstRunOfIOV)).c_str());  // empty treeInfo for input...
0308 
0309     if (inputs.back()->getLorentzAngles().empty()) {
0310       edm::LogError("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::endOfJob"
0311                                  << "Input Lorentz angle map is empty, skip writing output!";
0312       return;
0313     }
0314   }
0315 
0316   const unsigned int nonZeroParamsOrErrors =  // Any determined value?
0317       count_if(parameters_.begin(), parameters_.end(), [](auto c) { return c != 0.; }) +
0318       count_if(paramUncertainties_.begin(), paramUncertainties_.end(), [](auto c) { return c != 0.; });
0319 
0320   for (unsigned int iIOV = 0; iIOV < moduleGroupSelector_->numIovs(); ++iIOV) {
0321     auto firstRunOfIOV = static_cast<cond::Time_t>(moduleGroupSelector_->firstRunOfIOV(iIOV));
0322     SiPixelLorentzAngle output{};
0323     // Loop on map of values from input and add (possible) parameter results
0324     for (const auto &iterIdValue : inputs[iIOV]->getLorentzAngles()) {
0325       // type of 'iterIdValue' is pair<unsigned int, float>
0326       const auto detId = iterIdValue.first;  // key of map is DetId
0327       const auto param = this->getParameterForDetId(detId, firstRunOfIOV);
0328       // put result in output, i.e. sum of input and determined parameter, but the nasty
0329       // putLorentzAngle(..) takes float by reference - not even const reference:
0330       auto value = iterIdValue.second + static_cast<float>(param);
0331       output.putLorentzAngle(detId, value);
0332       const int paramIndex = moduleGroupSelector_->getParameterIndexFromDetId(detId, firstRunOfIOV);
0333       treeInfo[detId] = TreeStruct(param, this->getParameterError(paramIndex), paramIndex);
0334     }
0335 
0336     if (saveToDB_ || nonZeroParamsOrErrors != 0) {  // Skip writing mille jobs...
0337       this->writeTree(&output, treeInfo, (treeName + Form("result_%lld", firstRunOfIOV)).c_str());
0338     }
0339 
0340     if (saveToDB_) {  // If requested, write out to DB
0341       edm::Service<cond::service::PoolDBOutputService> dbService;
0342       if (dbService.isAvailable()) {
0343         dbService->writeOneIOV(output, firstRunOfIOV, recordNameDBwrite_);
0344       } else {
0345         edm::LogError("BadConfig") << "@SUB=SiPixelLorentzAngleCalibration::endOfJob"
0346                                    << "No PoolDBOutputService available, but saveToDB true!";
0347       }
0348     }
0349   }  // end loop on IOVs
0350 }
0351 
0352 //======================================================================
0353 const SiPixelLorentzAngle *SiPixelLorentzAngleCalibration::getLorentzAnglesInput(const align::RunNumber &run) {
0354   const auto &resolvedRun = run > 0 ? run : currentIOV_;
0355   // For parallel processing in Millepede II, create SiPixelLorentzAngle
0356   // from info stored in files of parallel jobs and check that they are identical.
0357   // If this job has run on data, still check that LA is identical to the ones
0358   // from mergeFileNames_.
0359   const std::string treeName{this->name() + "_input_" + std::to_string(resolvedRun)};
0360   for (const auto &iFile : mergeFileNames_) {
0361     auto la = this->createFromTree(iFile.c_str(), treeName.c_str());
0362     // siPixelLorentzAngleInput_ could be non-null from previous file of this loop
0363     // or from checkLorentzAngleInput(..) when running on data in this job as well
0364     if (!siPixelLorentzAngleInput_ || siPixelLorentzAngleInput_->getLorentzAngles().empty()) {
0365       cachedLorentzAngleInputs_[resolvedRun] = la;
0366       siPixelLorentzAngleInput_ = &(cachedLorentzAngleInputs_[resolvedRun]);
0367       currentIOV_ = resolvedRun;
0368     } else {
0369       // FIXME: about comparison of maps see comments in checkLorentzAngleInput
0370       if (!la.getLorentzAngles().empty() &&  // single job might not have got events
0371           la.getLorentzAngles() != siPixelLorentzAngleInput_->getLorentzAngles()) {
0372         // Throw exception instead of error?
0373         edm::LogError("NoInput") << "@SUB=SiPixelLorentzAngleCalibration::getLorentzAnglesInput"
0374                                  << "Different input values from tree " << treeName << " in file " << iFile << ".";
0375       }
0376     }
0377   }
0378 
0379   if (!siPixelLorentzAngleInput_) {  // no files nor ran on events
0380     // [] operator default-constructs an empty SiPixelLorentzAngle object in place:
0381     siPixelLorentzAngleInput_ = &(cachedLorentzAngleInputs_[resolvedRun]);
0382     currentIOV_ = resolvedRun;
0383     edm::LogError("NoInput") << "@SUB=SiPixelLorentzAngleCalibration::getLorentzAnglesInput"
0384                              << "No input, create an empty one!";
0385   } else if (siPixelLorentzAngleInput_->getLorentzAngles().empty()) {
0386     edm::LogError("NoInput") << "@SUB=SiPixelLorentzAngleCalibration::getLorentzAnglesInput"
0387                              << "Empty result!";
0388   }
0389 
0390   return siPixelLorentzAngleInput_;
0391 }
0392 
0393 //======================================================================
0394 double SiPixelLorentzAngleCalibration::getParameterForDetId(unsigned int detId, edm::RunNumber_t run) const {
0395   const int index = moduleGroupSelector_->getParameterIndexFromDetId(detId, run);
0396   return (index < 0 ? 0. : parameters_[index]);
0397 }
0398 
0399 //======================================================================
0400 void SiPixelLorentzAngleCalibration::writeTree(const SiPixelLorentzAngle *lorentzAngle,
0401                                                const std::map<unsigned int, TreeStruct> &treeInfo,
0402                                                const char *treeName) const {
0403   if (!lorentzAngle)
0404     return;
0405 
0406   TFile *file = TFile::Open(outFileName_.c_str(), "UPDATE");
0407   if (!file) {
0408     edm::LogError("BadConfig") << "@SUB=SiPixelLorentzAngleCalibration::writeTree"
0409                                << "Could not open file '" << outFileName_ << "'.";
0410     return;
0411   }
0412 
0413   TTree *tree = new TTree(treeName, treeName);
0414   unsigned int id = 0;
0415   float value = 0.;
0416   TreeStruct treeStruct;
0417   tree->Branch("detId", &id, "detId/i");
0418   tree->Branch("value", &value, "value/F");
0419   tree->Branch("treeStruct", &treeStruct, TreeStruct::LeafList());
0420 
0421   for (auto iterIdValue = lorentzAngle->getLorentzAngles().begin();
0422        iterIdValue != lorentzAngle->getLorentzAngles().end();
0423        ++iterIdValue) {
0424     // type of (*iterIdValue) is pair<unsigned int, float>
0425     id = iterIdValue->first;  // key of map is DetId
0426     value = iterIdValue->second;
0427     // type of (*treeStructIter) is pair<unsigned int, TreeStruct>
0428     auto treeStructIter = treeInfo.find(id);  // find info for this id
0429     if (treeStructIter != treeInfo.end()) {
0430       treeStruct = treeStructIter->second;  // info from input map
0431     } else {                                // if none found, fill at least parameter index (using 1st IOV...)
0432       const cond::Time_t run1of1stIov = moduleGroupSelector_->firstRunOfIOV(0);
0433       const int ind = moduleGroupSelector_->getParameterIndexFromDetId(id, run1of1stIov);
0434       treeStruct = TreeStruct(ind);
0435     }
0436     tree->Fill();
0437   }
0438   tree->Write();
0439   delete file;  // tree vanishes with the file... (?)
0440 }
0441 
0442 //======================================================================
0443 SiPixelLorentzAngle SiPixelLorentzAngleCalibration::createFromTree(const char *fileName, const char *treeName) const {
0444   // Check for file existence on your own to work around
0445   // https://hypernews.cern.ch/HyperNews/CMS/get/swDevelopment/2715.html:
0446   TFile *file = nullptr;
0447   FILE *testFile = fopen(fileName, "r");
0448   if (testFile) {
0449     fclose(testFile);
0450     file = TFile::Open(fileName, "READ");
0451   }  // else not existing, see error below
0452 
0453   TTree *tree = nullptr;
0454   if (file)
0455     file->GetObject(treeName, tree);
0456 
0457   SiPixelLorentzAngle result{};
0458   if (tree) {
0459     unsigned int id = 0;
0460     float value = 0.;
0461     tree->SetBranchAddress("detId", &id);
0462     tree->SetBranchAddress("value", &value);
0463 
0464     const Long64_t nEntries = tree->GetEntries();
0465     for (Long64_t iEntry = 0; iEntry < nEntries; ++iEntry) {
0466       tree->GetEntry(iEntry);
0467       result.putLorentzAngle(id, value);
0468     }
0469   } else {  // Warning only since could be parallel job on no events.
0470     edm::LogWarning("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::createFromTree"
0471                                  << "Could not get TTree '" << treeName << "' from file '" << fileName
0472                                  << (file ? "'." : "' (file does not exist).");
0473   }
0474 
0475   delete file;  // tree will vanish with file
0476   return result;
0477 }
0478 
0479 //======================================================================
0480 //======================================================================
0481 // Plugin definition
0482 
0483 #include "Alignment/CommonAlignmentAlgorithm/interface/IntegratedCalibrationPluginFactory.h"
0484 
0485 DEFINE_EDM_PLUGIN(IntegratedCalibrationPluginFactory, SiPixelLorentzAngleCalibration, "SiPixelLorentzAngleCalibration");