Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-10-25 09:34:07

0001 
0002 /*
0003  *  See header file for a description of this class.
0004  *
0005  *  Author of original version: M. Giunta
0006  *  \author A. Vilela Pereira
0007  */
0008 
0009 #include "DTVDriftWriter.h"
0010 
0011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0012 #include "FWCore/Framework/interface/ESHandle.h"
0013 #include "FWCore/Framework/interface/EventSetup.h"
0014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0015 
0016 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0017 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0018 #include "Geometry/DTGeometry/interface/DTSuperLayer.h"
0019 
0020 #include "CondFormats/DTObjects/interface/DTMtime.h"
0021 #include "CondFormats/DataRecord/interface/DTMtimeRcd.h"
0022 #include "CondFormats/DTObjects/interface/DTRecoConditions.h"
0023 #include "CondFormats/DataRecord/interface/DTRecoConditionsVdriftRcd.h"
0024 
0025 #include "CalibMuon/DTCalibration/interface/DTCalibDBUtils.h"
0026 #include "CalibMuon/DTCalibration/interface/DTVDriftPluginFactory.h"
0027 #include "CalibMuon/DTCalibration/interface/DTVDriftBaseAlgo.h"
0028 
0029 #include <string>
0030 #include <vector>
0031 
0032 using namespace std;
0033 using namespace edm;
0034 
0035 DTVDriftWriter::DTVDriftWriter(const ParameterSet& pset)
0036     : mTimeMapToken_(esConsumes<edm::Transition::BeginRun>()),
0037       vDriftMapToken_(esConsumes<edm::Transition::BeginRun>()),
0038       dtGeomToken_(esConsumes<edm::Transition::BeginRun>()),
0039       granularity_(pset.getUntrackedParameter<string>("calibGranularity", "bySL")),
0040       mTimeMap_(nullptr),
0041       vDriftMap_(nullptr),
0042       vDriftAlgo_{DTVDriftPluginFactory::get()->create(pset.getParameter<string>("vDriftAlgo"),
0043                                                        pset.getParameter<ParameterSet>("vDriftAlgoConfig"),
0044                                                        consumesCollector())} {
0045   LogVerbatim("Calibration") << "[DTVDriftWriter]Constructor called!";
0046 
0047   if (granularity_ != "bySL")
0048     throw cms::Exception("Configuration")
0049         << "[DTVDriftWriter] Check parameter calibGranularity: " << granularity_ << " option not available.";
0050 
0051   readLegacyVDriftDB = pset.getParameter<bool>("readLegacyVDriftDB");
0052   writeLegacyVDriftDB = pset.getParameter<bool>("writeLegacyVDriftDB");
0053 }
0054 
0055 DTVDriftWriter::~DTVDriftWriter() { LogVerbatim("Calibration") << "[DTVDriftWriter]Destructor called!"; }
0056 
0057 void DTVDriftWriter::beginRun(const edm::Run& run, const edm::EventSetup& setup) {
0058   // Get the map of vdrift from the Setup
0059   if (readLegacyVDriftDB) {
0060     mTimeMap_ = &setup.getData(mTimeMapToken_);
0061   } else {
0062     vDriftMap_ = &setup.getData(vDriftMapToken_);
0063     // Consistency check: no parametrization is implemented for the time being
0064     int version = vDriftMap_->version();
0065     if (version != 1) {
0066       throw cms::Exception("Configuration") << "only version 1 is presently supported for VDriftDB";
0067     }
0068   }
0069 
0070   // Get geometry from Event Setup
0071   dtGeom_ = setup.getHandle(dtGeomToken_);
0072 
0073   // Pass EventSetup to concrete implementation
0074   vDriftAlgo_->setES(setup);
0075 }
0076 
0077 void DTVDriftWriter::endJob() {
0078   // Create the object to be written to DB
0079   std::unique_ptr<DTMtime> mTimeNewMap;
0080   std::unique_ptr<DTRecoConditions> vDriftNewMap;
0081   if (writeLegacyVDriftDB) {
0082     mTimeNewMap = std::make_unique<DTMtime>();
0083   } else {
0084     vDriftNewMap = std::make_unique<DTRecoConditions>();
0085     vDriftNewMap->setFormulaExpr("[0]");
0086     //vDriftNewMap->setFormulaExpr("[0]*(1-[1]*x)"); // add parametrization for dependency along Y
0087     vDriftNewMap->setVersion(1);
0088   }
0089 
0090   if (granularity_ == "bySL") {
0091     // Get all the sls from the geometry
0092     const vector<const DTSuperLayer*>& superLayers = dtGeom_->superLayers();
0093     auto sl = superLayers.begin();
0094     auto sl_end = superLayers.end();
0095     for (; sl != sl_end; ++sl) {
0096       DTSuperLayerId slId = (*sl)->id();
0097 
0098       // Compute vDrift
0099       float vDriftNew = -1.;
0100       float resolutionNew = -1;
0101       try {
0102         dtCalibration::DTVDriftData vDriftData = vDriftAlgo_->compute(slId);
0103         vDriftNew = vDriftData.vdrift;
0104         resolutionNew = vDriftData.resolution;
0105         LogVerbatim("Calibration") << "vDrift for: " << slId << " Mean " << vDriftNew << " Resolution "
0106                                    << resolutionNew;
0107       } catch (cms::Exception& e) {  // Failure to compute new value, fall back to old table
0108         LogError("Calibration") << e.explainSelf();
0109         if (readLegacyVDriftDB) {  //...reading old db format...
0110           int status = mTimeMap_->get(slId, vDriftNew, resolutionNew, DTVelocityUnits::cm_per_ns);
0111           if (status == 0) {  // not found; silently skip this SL
0112             continue;
0113           }
0114         } else {  //...reading new db format
0115           try {
0116             vDriftNew = vDriftMap_->get(DTWireId(slId.rawId()));
0117           } catch (cms::Exception& e2) {
0118             // not found; silently skip this SL
0119             continue;
0120           }
0121         }
0122         LogVerbatim("Calibration") << "Keep original vDrift for: " << slId << " Mean " << vDriftNew << " Resolution "
0123                                    << resolutionNew;
0124       }
0125 
0126       // Add value to the vdrift table
0127       if (writeLegacyVDriftDB) {
0128         mTimeNewMap->set(slId, vDriftNew, resolutionNew, DTVelocityUnits::cm_per_ns);
0129       } else {
0130         vector<double> params = {vDriftNew};
0131         vDriftNewMap->set(DTWireId(slId.rawId()), params);
0132       }
0133     }  // End of loop on superlayers
0134   }
0135 
0136   // Write the vDrift object to DB
0137   LogVerbatim("Calibration") << "[DTVDriftWriter]Writing vdrift object to DB!";
0138   if (writeLegacyVDriftDB) {
0139     string record = "DTMtimeRcd";
0140     DTCalibDBUtils::writeToDB<DTMtime>(record, *mTimeNewMap);
0141   } else {
0142     DTCalibDBUtils::writeToDB<DTRecoConditions>("DTRecoConditionsVdriftRcd", *vDriftNewMap);
0143   }
0144 }