00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014 #include "Alignment/CommonAlignmentAlgorithm/interface/IntegratedCalibrationBase.h"
00015 #include "Alignment/CommonAlignmentAlgorithm/interface/TkModuleGroupSelector.h"
00016
00017 #include "TreeStruct.h"
00018
00019 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
00020 #include "CondFormats/DataRecord/interface/SiPixelLorentzAngleRcd.h"
00021 #include "CondFormats/SiPixelObjects/interface/SiPixelLorentzAngle.h"
00022
00023 #include "DataFormats/DetId/interface/DetId.h"
00024 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00025
00026 #include "FWCore/Framework/interface/ESHandle.h"
00027 #include "FWCore/Framework/interface/ESWatcher.h"
00028 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00029 #include "Alignment/CommonAlignment/interface/Alignable.h"
00030
00031 #include "FWCore/ServiceRegistry/interface/Service.h"
00032 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00033 #include "Alignment/CommonAlignmentAlgorithm/interface/AlignmentParameterSelector.h"
00034
00035 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
00036 #include "MagneticField/Engine/interface/MagneticField.h"
00037 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00038 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
00039
00040 #include "TTree.h"
00041 #include "TFile.h"
00042 #include "TString.h"
00043
00044 #include <boost/assign/list_of.hpp>
00045 #include <vector>
00046 #include <map>
00047 #include <sstream>
00048 #include <list>
00049 #include <utility>
00050 #include <set>
00051
00052 #include <functional>
00053
00054 class SiPixelLorentzAngleCalibration : public IntegratedCalibrationBase
00055 {
00056 public:
00058 explicit SiPixelLorentzAngleCalibration(const edm::ParameterSet &cfg);
00059
00061 virtual ~SiPixelLorentzAngleCalibration();
00062
00064 virtual unsigned int numParameters() const;
00065
00066
00067
00068
00069
00070
00071
00072
00073
00076 virtual unsigned int derivatives(std::vector<ValuesIndexPair> &outDerivInds,
00077 const TransientTrackingRecHit &hit,
00078 const TrajectoryStateOnSurface &tsos,
00079 const edm::EventSetup &setup,
00080 const EventInfo &eventInfo) const;
00081
00084 virtual bool setParameter(unsigned int index, double value);
00085
00088 virtual bool setParameterError(unsigned int index, double error);
00089
00092 virtual double getParameter(unsigned int index) const;
00093
00096 virtual double getParameterError(unsigned int index) const;
00097
00098
00099 virtual void beginOfJob(AlignableTracker *tracker,
00100 AlignableMuon *muon,
00101 AlignableExtras *extras);
00102
00103
00104
00107 virtual void endOfJob();
00108
00109 private:
00112 bool checkLorentzAngleInput(const edm::EventSetup &setup, const EventInfo &eventInfo);
00116 const SiPixelLorentzAngle* getLorentzAnglesInput();
00117
00120 double getParameterForDetId(unsigned int detId, edm::RunNumber_t run) const;
00121
00122 void writeTree(const SiPixelLorentzAngle *lorentzAngle,
00123 const std::map<unsigned int,TreeStruct>& treeInfo, const char *treeName) const;
00124 SiPixelLorentzAngle* createFromTree(const char *fileName, const char *treeName) const;
00125
00126 const bool saveToDB_;
00127 const std::string recordNameDBwrite_;
00128 const std::string outFileName_;
00129 const std::vector<std::string> mergeFileNames_;
00130
00131 edm::ESWatcher<SiPixelLorentzAngleRcd> watchLorentzAngleRcd_;
00132
00133
00134 SiPixelLorentzAngle *siPixelLorentzAngleInput_;
00135 std::vector<double> parameters_;
00136 std::vector<double> paramUncertainties_;
00137
00138 TkModuleGroupSelector *moduleGroupSelector_;
00139 const edm::ParameterSet moduleGroupSelCfg_;
00140 };
00141
00142
00143
00144
00145
00146 SiPixelLorentzAngleCalibration::SiPixelLorentzAngleCalibration(const edm::ParameterSet &cfg)
00147 : IntegratedCalibrationBase(cfg),
00148 saveToDB_(cfg.getParameter<bool>("saveToDB")),
00149 recordNameDBwrite_(cfg.getParameter<std::string>("recordNameDBwrite")),
00150 outFileName_(cfg.getParameter<std::string>("treeFile")),
00151 mergeFileNames_(cfg.getParameter<std::vector<std::string> >("mergeTreeFiles")),
00152 siPixelLorentzAngleInput_(0),
00153 moduleGroupSelector_(0),
00154 moduleGroupSelCfg_(cfg.getParameter<edm::ParameterSet>("LorentzAngleModuleGroups"))
00155 {
00156
00157 }
00158
00159
00160 SiPixelLorentzAngleCalibration::~SiPixelLorentzAngleCalibration()
00161 {
00162 delete moduleGroupSelector_;
00163
00164 delete siPixelLorentzAngleInput_;
00165 }
00166
00167
00168 unsigned int SiPixelLorentzAngleCalibration::numParameters() const
00169 {
00170 return parameters_.size();
00171 }
00172
00173
00174 unsigned int
00175 SiPixelLorentzAngleCalibration::derivatives(std::vector<ValuesIndexPair> &outDerivInds,
00176 const TransientTrackingRecHit &hit,
00177 const TrajectoryStateOnSurface &tsos,
00178 const edm::EventSetup &setup,
00179 const EventInfo &eventInfo) const
00180 {
00181
00182
00183 const_cast<SiPixelLorentzAngleCalibration*>(this)->checkLorentzAngleInput(setup, eventInfo);
00184
00185 outDerivInds.clear();
00186
00187 if (hit.det()) {
00188
00189 const int index = moduleGroupSelector_->getParameterIndexFromDetId(hit.det()->geographicalId(),
00190 eventInfo.eventId_.run());
00191 if (index >= 0) {
00192 edm::ESHandle<MagneticField> magneticField;
00193 setup.get<IdealMagneticFieldRecord>().get(magneticField);
00194 const GlobalVector bField(magneticField->inTesla(hit.det()->surface().position()));
00195 const LocalVector bFieldLocal(hit.det()->surface().toLocal(bField));
00196 const double dZ = hit.det()->surface().bounds().thickness();
00197
00198
00199 const double xDerivative = bFieldLocal.y() * dZ * -0.5;
00200
00201 if (xDerivative) {
00202 const Values derivs(xDerivative, 0.);
00203 outDerivInds.push_back(ValuesIndexPair(derivs, index));
00204 }
00205 }
00206 } else {
00207 edm::LogWarning("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::derivatives2"
00208 << "Hit without GeomDet, skip!";
00209 }
00210
00211 return outDerivInds.size();
00212 }
00213
00214
00215 bool SiPixelLorentzAngleCalibration::setParameter(unsigned int index, double value)
00216 {
00217 if (index >= parameters_.size()) {
00218 return false;
00219 } else {
00220 parameters_[index] = value;
00221 return true;
00222 }
00223 }
00224
00225
00226 bool SiPixelLorentzAngleCalibration::setParameterError(unsigned int index, double error)
00227 {
00228 if (index >= paramUncertainties_.size()) {
00229 return false;
00230 } else {
00231 paramUncertainties_[index] = error;
00232 return true;
00233 }
00234 }
00235
00236
00237 double SiPixelLorentzAngleCalibration::getParameter(unsigned int index) const
00238 {
00239 return (index >= parameters_.size() ? 0. : parameters_[index]);
00240 }
00241
00242
00243 double SiPixelLorentzAngleCalibration::getParameterError(unsigned int index) const
00244 {
00245 return (index >= paramUncertainties_.size() ? 0. : paramUncertainties_[index]);
00246 }
00247
00248
00249
00250
00251
00252 void SiPixelLorentzAngleCalibration::beginOfJob(AlignableTracker *aliTracker,
00253 AlignableMuon * ,
00254 AlignableExtras * )
00255 {
00256
00257 const std::vector<int> sdets = boost::assign::list_of(PixelSubdetector::PixelBarrel)(PixelSubdetector::PixelEndcap);
00258
00259 moduleGroupSelector_ = new TkModuleGroupSelector(aliTracker, moduleGroupSelCfg_, sdets);
00260
00261 parameters_.resize(moduleGroupSelector_->getNumberOfParameters(), 0.);
00262 paramUncertainties_.resize(moduleGroupSelector_->getNumberOfParameters(), 0.);
00263
00264 edm::LogInfo("Alignment") << "@SUB=SiPixelLorentzAngleCalibration" << "Created with name "
00265 << this->name() << "',\n" << this->numParameters() << " parameters to be determined,"
00266 << "\n saveToDB = " << saveToDB_
00267 << "\n outFileName = " << outFileName_
00268 << "\n N(merge files) = " << mergeFileNames_.size()
00269 << "\n number of IOVs = " << moduleGroupSelector_->numIovs();
00270
00271 if (mergeFileNames_.size()) {
00272 edm::LogInfo("Alignment") << "@SUB=SiPixelLorentzAngleCalibration"
00273 << "First file to merge: " << mergeFileNames_[0];
00274 }
00275 }
00276
00277
00278
00279 void SiPixelLorentzAngleCalibration::endOfJob()
00280 {
00281
00282 std::ostringstream out;
00283 out << "Parameter results\n";
00284 for (unsigned int iPar = 0; iPar < parameters_.size(); ++iPar) {
00285 out << iPar << ": " << parameters_[iPar] << " +- " << paramUncertainties_[iPar] << "\n";
00286 }
00287 edm::LogInfo("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::endOfJob" << out.str();
00288
00289 std::map<unsigned int, TreeStruct> treeInfo;
00290
00291
00292 const SiPixelLorentzAngle *input = this->getLorentzAnglesInput();
00293 const std::string treeName(this->name() + '_');
00294 this->writeTree(input, treeInfo, (treeName + "input").c_str());
00295
00296 if (input->getLorentzAngles().empty()) {
00297 edm::LogError("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::endOfJob"
00298 << "Input Lorentz angle map is empty, skip writing output!";
00299 return;
00300 }
00301
00302 const unsigned int nonZeroParamsOrErrors =
00303 count_if (parameters_.begin(), parameters_.end(), std::bind2nd(std::not_equal_to<double>(),0.))
00304 + count_if(paramUncertainties_.begin(), paramUncertainties_.end(),
00305 std::bind2nd(std::not_equal_to<double>(), 0.));
00306
00307 for (unsigned int iIOV = 0; iIOV < moduleGroupSelector_->numIovs(); ++iIOV) {
00308
00309 cond::Time_t firstRunOfIOV = moduleGroupSelector_->firstRunOfIOV(iIOV);
00310 SiPixelLorentzAngle *output = new SiPixelLorentzAngle;
00311
00312 for (auto iterIdValue = input->getLorentzAngles().begin();
00313 iterIdValue != input->getLorentzAngles().end(); ++iterIdValue) {
00314
00315 const unsigned int detId = iterIdValue->first;
00316 const double param = this->getParameterForDetId(detId, firstRunOfIOV);
00317
00318
00319 float value = iterIdValue->second + param;
00320 output->putLorentzAngle(detId, value);
00321 const int paramIndex = moduleGroupSelector_->getParameterIndexFromDetId(detId,firstRunOfIOV);
00322 treeInfo[detId] = TreeStruct(param, this->getParameterError(paramIndex), paramIndex);
00323 }
00324
00325 if (saveToDB_ || nonZeroParamsOrErrors != 0) {
00326 this->writeTree(output, treeInfo, (treeName + Form("result_%lld", firstRunOfIOV)).c_str());
00327 }
00328
00329 if (saveToDB_) {
00330 edm::Service<cond::service::PoolDBOutputService> dbService;
00331 if (dbService.isAvailable()) {
00332 dbService->writeOne(output, firstRunOfIOV, recordNameDBwrite_.c_str());
00333
00334 } else {
00335 delete output;
00336 edm::LogError("BadConfig") << "@SUB=SiPixelLorentzAngleCalibration::endOfJob"
00337 << "No PoolDBOutputService available, but saveToDB true!";
00338 }
00339 } else {
00340 delete output;
00341 }
00342 }
00343
00344 }
00345
00346
00347 bool SiPixelLorentzAngleCalibration::checkLorentzAngleInput(const edm::EventSetup &setup,
00348 const EventInfo &eventInfo)
00349 {
00350 edm::ESHandle<SiPixelLorentzAngle> lorentzAngleHandle;
00351 if (!siPixelLorentzAngleInput_) {
00352 setup.get<SiPixelLorentzAngleRcd>().get(lorentzAngleHandle);
00353 siPixelLorentzAngleInput_ = new SiPixelLorentzAngle(*lorentzAngleHandle);
00354 } else {
00355 if (watchLorentzAngleRcd_.check(setup)) {
00356 setup.get<SiPixelLorentzAngleRcd>().get(lorentzAngleHandle);
00357 if (lorentzAngleHandle->getLorentzAngles()
00358 != siPixelLorentzAngleInput_->getLorentzAngles()) {
00359
00360
00361 throw cms::Exception("BadInput")
00362 << "SiPixelLorentzAngleCalibration::checkLorentzAngleInput:\n"
00363 << "Content of SiPixelLorentzAngle changed at run " << eventInfo.eventId_.run()
00364 << ", but algorithm expects constant input!\n";
00365 return false;
00366 }
00367 }
00368 }
00369
00370 return true;
00371 }
00372
00373
00374 const SiPixelLorentzAngle* SiPixelLorentzAngleCalibration::getLorentzAnglesInput()
00375 {
00376
00377
00378
00379
00380 const std::string treeName(this->name() + "_input");
00381 for (auto iFile = mergeFileNames_.begin(); iFile != mergeFileNames_.end(); ++iFile) {
00382 SiPixelLorentzAngle* la = this->createFromTree(iFile->c_str(), treeName.c_str());
00383
00384
00385 if (!siPixelLorentzAngleInput_ || siPixelLorentzAngleInput_->getLorentzAngles().empty()) {
00386 delete siPixelLorentzAngleInput_;
00387 siPixelLorentzAngleInput_ = la;
00388 } else {
00389
00390 if (la && !la->getLorentzAngles().empty() &&
00391 la->getLorentzAngles() != siPixelLorentzAngleInput_->getLorentzAngles()) {
00392
00393 edm::LogError("NoInput") << "@SUB=SiPixelLorentzAngleCalibration::getLorentzAnglesInput"
00394 << "Different input values from tree " << treeName
00395 << " in file " << *iFile << ".";
00396
00397 }
00398 delete la;
00399 }
00400 }
00401
00402 if (!siPixelLorentzAngleInput_) {
00403 siPixelLorentzAngleInput_ = new SiPixelLorentzAngle;
00404 edm::LogError("NoInput") << "@SUB=SiPixelLorentzAngleCalibration::getLorentzAnglesInput"
00405 << "No input, create an empty one!";
00406 } else if (siPixelLorentzAngleInput_->getLorentzAngles().empty()) {
00407 edm::LogError("NoInput") << "@SUB=SiPixelLorentzAngleCalibration::getLorentzAnglesInput"
00408 << "Empty result!";
00409 }
00410
00411 return siPixelLorentzAngleInput_;
00412 }
00413
00414
00415 double SiPixelLorentzAngleCalibration::getParameterForDetId(unsigned int detId,
00416 edm::RunNumber_t run) const
00417 {
00418 const int index = moduleGroupSelector_->getParameterIndexFromDetId(detId, run);
00419 return (index < 0 ? 0. : parameters_[index]);
00420 }
00421
00422
00423 void SiPixelLorentzAngleCalibration::writeTree(const SiPixelLorentzAngle *lorentzAngle,
00424 const std::map<unsigned int,TreeStruct> &treeInfo,
00425 const char *treeName) const
00426 {
00427 if (!lorentzAngle) return;
00428
00429 TFile* file = TFile::Open(outFileName_.c_str(), "UPDATE");
00430 if (!file) {
00431 edm::LogError("BadConfig") << "@SUB=SiPixelLorentzAngleCalibration::writeTree"
00432 << "Could not open file '" << outFileName_ << "'.";
00433 return;
00434 }
00435
00436 TTree *tree = new TTree(treeName, treeName);
00437 unsigned int id = 0;
00438 float value = 0.;
00439 TreeStruct treeStruct;
00440 tree->Branch("detId", &id, "detId/i");
00441 tree->Branch("value", &value, "value/F");
00442 tree->Branch("treeStruct", &treeStruct, TreeStruct::LeafList());
00443
00444 for (auto iterIdValue = lorentzAngle->getLorentzAngles().begin();
00445 iterIdValue != lorentzAngle->getLorentzAngles().end(); ++iterIdValue) {
00446
00447 id = iterIdValue->first;
00448 value = iterIdValue->second;
00449
00450 auto treeStructIter = treeInfo.find(id);
00451 if (treeStructIter != treeInfo.end()) {
00452 treeStruct = treeStructIter->second;
00453 } else {
00454 const cond::Time_t run1of1stIov = moduleGroupSelector_->firstRunOfIOV(0);
00455 const int ind = moduleGroupSelector_->getParameterIndexFromDetId(id, run1of1stIov);
00456 treeStruct = TreeStruct(ind);
00457 }
00458 tree->Fill();
00459 }
00460 tree->Write();
00461 delete file;
00462
00463 }
00464
00465
00466 SiPixelLorentzAngle*
00467 SiPixelLorentzAngleCalibration::createFromTree(const char *fileName, const char *treeName) const
00468 {
00469
00470
00471 TFile* file = 0;
00472 FILE* testFile = fopen(fileName,"r");
00473 if (testFile) {
00474 fclose(testFile);
00475 file = TFile::Open(fileName, "READ");
00476 }
00477
00478 TTree *tree = 0;
00479 if (file) file->GetObject(treeName, tree);
00480
00481 SiPixelLorentzAngle *result = 0;
00482 if (tree) {
00483 unsigned int id = 0;
00484 float value = 0.;
00485 tree->SetBranchAddress("detId", &id);
00486 tree->SetBranchAddress("value", &value);
00487
00488 result = new SiPixelLorentzAngle;
00489 const Long64_t nEntries = tree->GetEntries();
00490 for (Long64_t iEntry = 0; iEntry < nEntries; ++iEntry) {
00491 tree->GetEntry(iEntry);
00492 result->putLorentzAngle(id, value);
00493 }
00494 } else {
00495 edm::LogWarning("Alignment") << "@SUB=SiPixelLorentzAngleCalibration::createFromTree"
00496 << "Could not get TTree '" << treeName << "' from file '"
00497 << fileName << (file ? "'." : "' (file does not exist).");
00498 }
00499
00500 delete file;
00501 return result;
00502 }
00503
00504
00505
00506
00507
00508
00509 #include "Alignment/CommonAlignmentAlgorithm/interface/IntegratedCalibrationPluginFactory.h"
00510
00511 DEFINE_EDM_PLUGIN(IntegratedCalibrationPluginFactory,
00512 SiPixelLorentzAngleCalibration, "SiPixelLorentzAngleCalibration");