CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/Alignment/CommonAlignmentAlgorithm/plugins/SiPixelLorentzAngleCalibration.cc

Go to the documentation of this file.
00001 
00002 
00003 
00004 
00005 
00006 
00007 
00008 
00009 
00010 
00011 
00012 
00013 
00014 #include "Alignment/CommonAlignmentAlgorithm/interface/IntegratedCalibrationBase.h"
00015 
00016 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
00017 #include "CondFormats/DataRecord/interface/SiPixelLorentzAngleRcd.h"
00018 #include "CondFormats/SiPixelObjects/interface/SiPixelLorentzAngle.h"
00019 
00020 #include "DataFormats/DetId/interface/DetId.h"
00021 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
00022 
00023 #include "FWCore/Framework/interface/ESHandle.h"
00024 #include "FWCore/Framework/interface/ESWatcher.h"
00025 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00026 #include "FWCore/ServiceRegistry/interface/Service.h"
00027 
00028 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
00029 #include "MagneticField/Engine/interface/MagneticField.h"
00030 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00031 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
00032 
00033 #include "TTree.h"
00034 #include "TFile.h"
00035 #include "TString.h"
00036 
00037 // #include <iostream>
00038 #include <vector>
00039 #include <sstream>
00040 
00041 class SiPixelLorentzAngleCalibration : public IntegratedCalibrationBase
00042 {
00043 public:
00045   explicit SiPixelLorentzAngleCalibration(const edm::ParameterSet &cfg);
00046   
00048   virtual ~SiPixelLorentzAngleCalibration();
00049 
00051   virtual unsigned int numParameters() const;
00052 
00053   // /// Return all derivatives,
00054   // /// default implementation uses other derivatives(..) method,
00055   // /// but can be overwritten in derived class for efficiency.
00056   // virtual std::vector<double> derivatives(const TransientTrackingRecHit &hit,
00057   //                                      const TrajectoryStateOnSurface &tsos,
00058   //                                      const edm::EventSetup &setup,
00059   //                                      const EventInfo &eventInfo) const;
00060 
00063   virtual unsigned int derivatives(std::vector<ValuesIndexPair> &outDerivInds,
00064                                    const TransientTrackingRecHit &hit,
00065                                    const TrajectoryStateOnSurface &tsos,
00066                                    const edm::EventSetup &setup,
00067                                    const EventInfo &eventInfo) const;
00068 
00071   virtual bool setParameter(unsigned int index, double value);
00072 
00075   virtual bool setParameterError(unsigned int index, double error);
00076 
00079   virtual double getParameter(unsigned int index) const;
00080 
00083   virtual double getParameterError(unsigned int index) const;
00084 
00085   // /// Call at beginning of job:
00086   // virtual void beginOfJob(const AlignableTracker *tracker,
00087   //                      const AlignableMuon *muon,
00088   //                      const AlignableExtras *extras);
00089 
00092   virtual void endOfJob();
00093 
00094 private:
00097   bool checkLorentzAngleInput(const edm::EventSetup &setup, const EventInfo &eventInfo);
00101   const SiPixelLorentzAngle* getLorentzAnglesInput();
00102 
00105   double getParameterForDetId(unsigned int detId, edm::RunNumber_t run) const;
00108   int getParameterIndexFromDetId(unsigned int detId, edm::RunNumber_t run) const;
00110   unsigned int numIovs() const;
00112   edm::RunNumber_t firstRunOfIOV(unsigned int iovNum) const;
00113 
00114   void writeTree(const SiPixelLorentzAngle *lorentzAngle, const char *treeName) const;
00115   SiPixelLorentzAngle* createFromTree(const char *fileName, const char *treeName) const;
00116 
00117   const bool saveToDB_;
00118   const std::string recordNameDBwrite_;
00119   const std::string outFileName_;
00120   const std::vector<std::string> mergeFileNames_;
00121 
00122   edm::ESWatcher<SiPixelLorentzAngleRcd> watchLorentzAngleRcd_;
00123 
00124   // const AlignableTracker *alignableTracker_;
00125   SiPixelLorentzAngle *siPixelLorentzAngleInput_;
00126   std::vector<double> parameters_;
00127   std::vector<double> paramUncertainties_;
00128 };
00129 
00130 //======================================================================
00131 //======================================================================
00132 //======================================================================
00133 
00134 SiPixelLorentzAngleCalibration::SiPixelLorentzAngleCalibration(const edm::ParameterSet &cfg)
00135   : IntegratedCalibrationBase(cfg),
00136     saveToDB_(cfg.getParameter<bool>("saveToDB")),
00137     recordNameDBwrite_(cfg.getParameter<std::string>("recordNameDBwrite")),
00138     outFileName_(cfg.getParameter<std::string>("treeFile")),
00139     mergeFileNames_(cfg.getParameter<std::vector<std::string> >("mergeTreeFiles")),
00140     //    alignableTracker_(0),
00141     siPixelLorentzAngleInput_(0)
00142 {
00143   // FIXME: Which granularity, leading to how many parameters?
00144   parameters_.resize(2, 0.); // currently two parameters (ring1-4, 5-8), start value 0.
00145   paramUncertainties_.resize(2, 0.); // dito for errors
00146 
00147   edm::LogInfo("Alignment") << "@SUB=SiPixelLorentzAngleCalibration" << "Created with name "
00148                             << this->name() << "',\n" << this->numParameters() << " parameters to be determined,"
00149                             << "\nsaveToDB = " << saveToDB_
00150                             << "\n outFileName = " << outFileName_
00151                             << "\n N(merge files) = " << mergeFileNames_.size();
00152   if (mergeFileNames_.size()) {
00153     edm::LogInfo("Alignment") << "@SUB=SiPixelLorentzAngleCalibration"
00154                               << "First file to merge: " << mergeFileNames_[0];
00155   }
00156 }
00157   
00158 //======================================================================
00159 SiPixelLorentzAngleCalibration::~SiPixelLorentzAngleCalibration()
00160 {
00161   //  std::cout << "Destroy SiPixelLorentzAngleCalibration named " << this->name() << std::endl;
00162   delete siPixelLorentzAngleInput_;
00163 }
00164 
00165 //======================================================================
00166 unsigned int SiPixelLorentzAngleCalibration::numParameters() const
00167 {
00168   return parameters_.size();
00169 }
00170 
00171 //======================================================================
00172 unsigned int
00173 SiPixelLorentzAngleCalibration::derivatives(std::vector<ValuesIndexPair> &outDerivInds,
00174                                             const TransientTrackingRecHit &hit,
00175                                             const TrajectoryStateOnSurface &tsos,
00176                                             const edm::EventSetup &setup,
00177                                             const EventInfo &eventInfo) const
00178 {
00179   // ugly const-cast:
00180   // But it is either only first initialisation or throwing an exception...
00181   const_cast<SiPixelLorentzAngleCalibration*>(this)->checkLorentzAngleInput(setup, eventInfo);
00182 
00183   outDerivInds.clear();
00184 
00185   if (hit.det()) { // otherwise 'constraint hit' or whatever
00186     
00187     const int index = this->getParameterIndexFromDetId(hit.det()->geographicalId(),
00188                                                        eventInfo.eventId_.run());
00189     if (index >= 0) { // otherwise not treated
00190       edm::ESHandle<MagneticField> magneticField;
00191       setup.get<IdealMagneticFieldRecord>().get(magneticField);
00192       const GlobalVector bField(magneticField->inTesla(hit.det()->surface().position()));
00193       const LocalVector bFieldLocal(hit.det()->surface().toLocal(bField));
00194       const double dZ = hit.det()->surface().bounds().thickness(); // it is a float only...
00195       // shift due to LA: dx = tan(LA) * dz/2 = mobility * B_y * dz/2,
00196       // '-' since we have derivative of the residual r = trk -hit
00197       const double xDerivative = bFieldLocal.y() * dZ * -0.5; // parameter is mobility!
00198       if (xDerivative) { // If field is zero, this is zero: do not return it
00199         const Values derivs(xDerivative, 0.); // yDerivative = 0.
00200         outDerivInds.push_back(ValuesIndexPair(derivs, index));
00201       }
00202     }
00203   } else {
00204     edm::LogWarning("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::derivatives2"
00205                                  << "Hit without GeomDet, skip!";
00206   }
00207   
00208   return outDerivInds.size();
00209 }
00210 
00211 //======================================================================
00212 bool SiPixelLorentzAngleCalibration::setParameter(unsigned int index, double value)
00213 {
00214   if (index >= parameters_.size()) {
00215     return false;
00216   } else {
00217     parameters_[index] = value;
00218     return true;
00219   }
00220 }
00221 
00222 //======================================================================
00223 bool SiPixelLorentzAngleCalibration::setParameterError(unsigned int index, double error)
00224 {
00225   if (index >= paramUncertainties_.size()) {
00226     return false;
00227   } else {
00228     paramUncertainties_[index] = error;
00229     return true;
00230   }
00231 }
00232 
00233 //======================================================================
00234 double SiPixelLorentzAngleCalibration::getParameter(unsigned int index) const
00235 {
00236   //   if (index >= parameters_.size()) {
00237   //     return 0.;
00238   //   } else {
00239   //     return parameters_[index];
00240   //   }
00241   return (index >= parameters_.size() ? 0. : parameters_[index]);
00242 }
00243 
00244 //======================================================================
00245 double SiPixelLorentzAngleCalibration::getParameterError(unsigned int index) const
00246 {
00247   //   if (index >= paramUncertainties_.size()) {
00248   //     return 0.;
00249   //   } else {
00250   //     return paramUncertainties_[index];
00251   //   }
00252   return (index >= paramUncertainties_.size() ? 0. : paramUncertainties_[index]);
00253 }
00254 
00255 // //======================================================================
00256 // void SiPixelLorentzAngleCalibration::beginOfJob(const AlignableTracker *tracker,
00257 //                                                 const AlignableMuon */*muon*/,
00258 //                                                 const AlignableExtras */*extras*/)
00259 // {
00260 //   alignableTracker_ = tracker;
00261 // }
00262 
00263 
00264 //======================================================================
00265 void SiPixelLorentzAngleCalibration::endOfJob()
00266 {
00267   // loginfo output
00268   std::ostringstream out;
00269   out << "Parameter results\n";
00270   for (unsigned int iPar = 0; iPar < parameters_.size(); ++iPar) {
00271     out << iPar << ": " << parameters_[iPar] << " +- " << paramUncertainties_[iPar] << "\n";
00272   }
00273   edm::LogInfo("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::endOfJob" << out.str();
00274 
00275   // now write 'input' tree
00276   const SiPixelLorentzAngle *input = this->getLorentzAnglesInput(); // never NULL
00277   const std::string treeName(this->name() + '_');
00278   this->writeTree(input, (treeName + "input").c_str());
00279   if (input->getLorentzAngles().empty()) {
00280     edm::LogError("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::endOfJob"
00281                                << "Input Lorentz angle map is empty, skip writing output!";
00282     return;
00283   }
00284 
00285   for (unsigned int iIOV = 0; iIOV < this->numIovs(); ++iIOV) {
00286     cond::Time_t firstRunOfIOV = this->firstRunOfIOV(iIOV);
00287     SiPixelLorentzAngle *output = new SiPixelLorentzAngle;
00288     // Loop on map of values from input and add (possible) parameter results
00289     for (auto iterIdValue = input->getLorentzAngles().begin();
00290          iterIdValue != input->getLorentzAngles().end(); ++iterIdValue) {
00291       // type of (*iterIdValue) is pair<unsigned int, float>
00292       const unsigned int detId = iterIdValue->first; // key of map is DetId
00293       // Nasty: putLorentzAngle(..) takes float by reference - not even const reference!
00294       float value = iterIdValue->second + this->getParameterForDetId(detId, firstRunOfIOV);
00295       output->putLorentzAngle(detId, value); // put result in output
00296     }
00297 
00298     // Write this even for mille jobs?
00299     this->writeTree(output, (treeName + Form("result_%lld", firstRunOfIOV)).c_str());
00300 
00301     if (saveToDB_) { // If requested, write out to DB 
00302       edm::Service<cond::service::PoolDBOutputService> dbService;
00303       if (dbService.isAvailable()) {
00304         dbService->writeOne(output, firstRunOfIOV, recordNameDBwrite_.c_str());
00305         // no 'delete output;': writeOne(..) took over ownership
00306       } else {
00307         delete output;
00308         edm::LogError("BadConfig") << "@SUB=SiPixelLorentzAngleCalibration::endOfJob"
00309                                    << "No PoolDBOutputService available, but saveToDB true!";
00310       }
00311     } else {
00312       delete output;
00313     }
00314   } // end loop on IOVs
00315 
00316 }
00317 
00318 //======================================================================
00319 bool SiPixelLorentzAngleCalibration::checkLorentzAngleInput(const edm::EventSetup &setup,
00320                                                             const EventInfo &eventInfo)
00321 {
00322   edm::ESHandle<SiPixelLorentzAngle> lorentzAngleHandle;
00323   if (!siPixelLorentzAngleInput_) {
00324     setup.get<SiPixelLorentzAngleRcd>().get(lorentzAngleHandle);
00325     siPixelLorentzAngleInput_ = new SiPixelLorentzAngle(*lorentzAngleHandle);
00326   } else {
00327     if (watchLorentzAngleRcd_.check(setup)) { // new IOV of input
00328       setup.get<SiPixelLorentzAngleRcd>().get(lorentzAngleHandle);
00329       if (lorentzAngleHandle->getLorentzAngles() // but only bad if non-identical values
00330           != siPixelLorentzAngleInput_->getLorentzAngles()) { // (comparing maps)
00331         // FIXME: Could different maps have same content, but different order?
00332         //        Or 'floating point comparison' problem?
00333         throw cms::Exception("BadInput")
00334           << "SiPixelLorentzAngleCalibration::checkLorentzAngleInput:\n"
00335           << "Content of SiPixelLorentzAngle changed at run " << eventInfo.eventId_.run()
00336           << ", but algorithm expects constant input!\n";
00337         return false; // not reached...
00338       }
00339     }
00340   }
00341   
00342   return true;
00343 }
00344 
00345 //======================================================================
00346 const SiPixelLorentzAngle* SiPixelLorentzAngleCalibration::getLorentzAnglesInput()
00347 {
00348   // For parallel processing in Millepede II, create SiPixelLorentzAngle
00349   // from info stored in files of parallel jobs and check that they are identical.
00350   // If this job has run on data, still check that LA is identical to the ones
00351   // from mergeFileNames_.
00352   const std::string treeName(this->name() + "_input");
00353   for (auto iFile = mergeFileNames_.begin(); iFile != mergeFileNames_.end(); ++iFile) {
00354     SiPixelLorentzAngle* la = this->createFromTree(iFile->c_str(), treeName.c_str());
00355     // siPixelLorentzAngleInput_ could be non-null from previous file of this loop
00356     // or from checkLorentzAngleInput(..) when running on data in this job as well
00357     if (!siPixelLorentzAngleInput_ || siPixelLorentzAngleInput_->getLorentzAngles().empty()) {
00358       delete siPixelLorentzAngleInput_; // NULL or empty
00359       siPixelLorentzAngleInput_ = la;
00360     } else {
00361       // FIXME: about comparison of maps see comments in checkLorentzAngleInput
00362       if (la && !la->getLorentzAngles().empty() && // single job might not have got events
00363           la->getLorentzAngles() != siPixelLorentzAngleInput_->getLorentzAngles()) {
00364         // Throw exception instead of error?
00365         edm::LogError("NoInput") << "@SUB=SiPixelLorentzAngleCalibration::getLorentzAnglesInput"
00366                                  << "Different input values from tree " << treeName
00367                                  << " in file " << *iFile << ".";
00368         
00369       }
00370       delete la;
00371     }
00372   }
00373 
00374   if (!siPixelLorentzAngleInput_) { // no files nor ran on events
00375     siPixelLorentzAngleInput_ = new SiPixelLorentzAngle;
00376     edm::LogError("NoInput") << "@SUB=SiPixelLorentzAngleCalibration::getLorentzAnglesInput"
00377                              << "No input, create an empty one!";
00378   } else if (siPixelLorentzAngleInput_->getLorentzAngles().empty()) {
00379     edm::LogError("NoInput") << "@SUB=SiPixelLorentzAngleCalibration::getLorentzAnglesInput"
00380                              << "Empty result!";
00381   }
00382 
00383   return siPixelLorentzAngleInput_;
00384 }
00385 
00386 //======================================================================
00387 double SiPixelLorentzAngleCalibration::getParameterForDetId(unsigned int detId,
00388                                                             edm::RunNumber_t run) const
00389 {
00390   const int index = this->getParameterIndexFromDetId(detId, run);
00391 
00392   return (index < 0 ? 0. : parameters_[index]);
00393 }
00394 
00395 //======================================================================
00396 int SiPixelLorentzAngleCalibration::getParameterIndexFromDetId(unsigned int detId,
00397                                                                edm::RunNumber_t run) const
00398 {
00399   // Return the index of the parameter that is used for this DetId.
00400   // If this DetId is not treated, return values < 0.
00401   
00402   // FIXME: Extend to configurable granularity? 
00403   //        Including treatment of run dependence?
00404   const PXBDetId id(detId);
00405   if (id.det() == DetId::Tracker && id.subdetId() == PixelSubdetector::PixelBarrel) {
00406     if (id.module() >= 1 && id.module() <= 4) return 0;
00407     if (id.module() >= 5 && id.module() <= 8) return 1;
00408     edm::LogWarning("Alignment")
00409         << "@SUB=SiPixelLorentzAngleCalibration::getParameterIndexFromDetId"
00410         << "Module should be 1-8, but is " << id.module() << " => skip!";
00411   }
00412 
00413   return -1;
00414 }
00415 
00416 //======================================================================
00417 unsigned int SiPixelLorentzAngleCalibration::numIovs() const
00418 {
00419   // FIXME: Needed to include treatment of run dependence!
00420   return 1; 
00421 }
00422 
00423 //======================================================================
00424 edm::RunNumber_t SiPixelLorentzAngleCalibration::firstRunOfIOV(unsigned int iovNum) const
00425 {
00426   // FIXME: Needed to include treatment of run dependence!
00427   if (iovNum < this->numIovs()) return 1;
00428   else return 0;
00429 }
00430 
00431 
00432 //======================================================================
00433 void SiPixelLorentzAngleCalibration::writeTree(const SiPixelLorentzAngle *lorentzAngle,
00434                                                const char *treeName) const
00435 {
00436   if (!lorentzAngle) return;
00437 
00438   TFile* file = TFile::Open(outFileName_.c_str(), "UPDATE");
00439   if (!file) {
00440     edm::LogError("BadConfig") << "@SUB=SiPixelLorentzAngleCalibration::writeTree"
00441                                << "Could not open file '" << outFileName_ << "'.";
00442     return;
00443   }
00444 
00445   TTree *tree = new TTree(treeName, treeName);
00446   unsigned int id = 0;
00447   float value = 0.;
00448   tree->Branch("detId", &id, "detId/i");
00449   tree->Branch("value", &value, "value/F");
00450 
00451   for (auto iterIdValue = lorentzAngle->getLorentzAngles().begin();
00452        iterIdValue != lorentzAngle->getLorentzAngles().end(); ++iterIdValue) {
00453     // type of (*iterIdValue) is pair<unsigned int, float>
00454     id = iterIdValue->first; // key of map is DetId
00455     value = iterIdValue->second;
00456     tree->Fill();
00457   }
00458   tree->Write();
00459   delete file; // tree vanishes with the file... (?)
00460 
00461 }
00462 
00463 //======================================================================
00464 SiPixelLorentzAngle* 
00465 SiPixelLorentzAngleCalibration::createFromTree(const char *fileName, const char *treeName) const
00466 {
00467   // Check for file existence on your own to work around
00468   // https://hypernews.cern.ch/HyperNews/CMS/get/swDevelopment/2715.html:
00469   TFile* file = 0;
00470   FILE* testFile = fopen(fileName,"r");
00471   if (testFile) {
00472     fclose(testFile);
00473     file = TFile::Open(fileName, "READ");
00474   } // else not existing, see error below
00475 
00476   TTree *tree = 0;
00477   if (file) file->GetObject(treeName, tree);
00478 
00479   SiPixelLorentzAngle *result = 0;
00480   if (tree) {
00481     unsigned int id = 0;
00482     float value = 0.;
00483     tree->SetBranchAddress("detId", &id);
00484     tree->SetBranchAddress("value", &value);
00485 
00486     result = new SiPixelLorentzAngle;
00487     const Long64_t nEntries = tree->GetEntries();
00488     for (Long64_t iEntry = 0; iEntry < nEntries; ++iEntry) {
00489       tree->GetEntry(iEntry);
00490       result->putLorentzAngle(id, value);
00491     }
00492   } else { // Warning only since could be parallel job on no events.
00493     edm::LogWarning("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::createFromTree"
00494                                  << "Could not get TTree '" << treeName << "' from file '"
00495                                  << fileName << (file ? "'." : "' (file does not exist).");
00496   }
00497 
00498   delete file; // tree will vanish with file
00499   return result;
00500 }
00501 
00502 
00503 //======================================================================
00504 //======================================================================
00505 // Plugin definition
00506 
00507 #include "Alignment/CommonAlignmentAlgorithm/interface/IntegratedCalibrationPluginFactory.h"
00508 
00509 DEFINE_EDM_PLUGIN(IntegratedCalibrationPluginFactory,
00510                    SiPixelLorentzAngleCalibration, "SiPixelLorentzAngleCalibration");