Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:58:27

0001 /*
0002  *  See header file for a description of this class.
0003  *
0004  */
0005 
0006 #include "DTTTrigConstantShift.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include "FWCore/Framework/interface/EventSetup.h"
0009 #include "FWCore/Framework/interface/ESHandle.h"
0010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0011 
0012 #include "DataFormats/MuonDetId/interface/DTSuperLayerId.h"
0013 #include "CondFormats/DTObjects/interface/DTTtrig.h"
0014 #include "CondFormats/DataRecord/interface/DTTtrigRcd.h"
0015 #include "FWCore/Framework/interface/ConsumesCollector.h"
0016 
0017 #include <cmath>
0018 
0019 using namespace std;
0020 using namespace edm;
0021 
0022 namespace dtCalibration {
0023 
0024   DTTTrigConstantShift::DTTTrigConstantShift(const ParameterSet& pset, edm::ConsumesCollector cc)
0025       : calibChamber_(pset.getParameter<string>("calibChamber")), value_(pset.getParameter<double>("value")) {
0026     ttrigToken_ =
0027         cc.esConsumes<edm::Transition::BeginRun>(edm::ESInputTag("", pset.getUntrackedParameter<string>("dbLabel")));
0028     LogVerbatim("Calibration") << "[DTTTrigConstantShift] Applying constant correction value: " << value_ << endl;
0029 
0030     if (!calibChamber_.empty() && calibChamber_ != "None" && calibChamber_ != "All") {
0031       stringstream linestr;
0032       int selWheel, selStation, selSector;
0033       linestr << calibChamber_;
0034       linestr >> selWheel >> selStation >> selSector;
0035       chosenChamberId_ = DTChamberId(selWheel, selStation, selSector);
0036       LogVerbatim("Calibration") << "[DTTTrigConstantShift] Chosen chamber: " << chosenChamberId_ << endl;
0037     }
0038     //FIXME: Check if chosen chamber is valid.
0039   }
0040 
0041   DTTTrigConstantShift::~DTTTrigConstantShift() {}
0042 
0043   void DTTTrigConstantShift::setES(const EventSetup& setup) {
0044     // Get tTrig record from DB
0045     ESHandle<DTTtrig> tTrig = setup.getHandle(ttrigToken_);
0046     tTrigMap_ = &*tTrig;
0047   }
0048 
0049   DTTTrigData DTTTrigConstantShift::correction(const DTSuperLayerId& slId) {
0050     float tTrigMean, tTrigSigma, kFactor;
0051     int status = tTrigMap_->get(slId, tTrigMean, tTrigSigma, kFactor, DTTimeUnits::ns);
0052     if (status != 0)
0053       throw cms::Exception("[DTTTrigConstantShift]") << "Could not find tTrig entry in DB for" << slId << endl;
0054 
0055     float tTrigMeanNew = tTrigMean;
0056     if (!calibChamber_.empty() && calibChamber_ != "None") {
0057       if ((calibChamber_ == "All") || (calibChamber_ != "All" && slId.chamberId() == chosenChamberId_)) {
0058         tTrigMeanNew = tTrigMean + value_;
0059       }
0060     }
0061 
0062     return DTTTrigData(tTrigMeanNew, tTrigSigma, kFactor);
0063   }
0064 
0065 }  // namespace dtCalibration