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 "DTT0ChamberReferenceCorrection.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 #include "FWCore/Framework/interface/ConsumesCollector.h"
0012 
0013 #include "DataFormats/MuonDetId/interface/DTWireId.h"
0014 #include "CondFormats/DTObjects/interface/DTT0.h"
0015 #include "CondFormats/DataRecord/interface/DTT0Rcd.h"
0016 
0017 #include <string>
0018 #include <sstream>
0019 
0020 using namespace std;
0021 using namespace edm;
0022 
0023 namespace dtCalibration {
0024 
0025   DTT0ChamberReferenceCorrection::DTT0ChamberReferenceCorrection(const ParameterSet& pset, edm::ConsumesCollector cc)
0026       : calibChamber_(pset.getParameter<string>("calibChamber")), t0Token_(cc.esConsumes<edm::Transition::BeginRun>()) {
0027     //DTChamberId chosenChamberId;
0028     if (!calibChamber_.empty() && calibChamber_ != "None" && calibChamber_ != "All") {
0029       stringstream linestr;
0030       int selWheel, selStation, selSector;
0031       linestr << calibChamber_;
0032       linestr >> selWheel >> selStation >> selSector;
0033       chosenChamberId_ = DTChamberId(selWheel, selStation, selSector);
0034       LogVerbatim("Calibration") << "[DTT0ChamberReferenceCorrection] Chosen chamber: " << chosenChamberId_ << endl;
0035     }
0036     //FIXME: Check if chosen chamber is valid.
0037   }
0038 
0039   DTT0ChamberReferenceCorrection::~DTT0ChamberReferenceCorrection() {}
0040 
0041   void DTT0ChamberReferenceCorrection::setES(const EventSetup& setup) {
0042     // Get t0 record from DB
0043     ESHandle<DTT0> t0H;
0044     t0H = setup.getHandle(t0Token_);
0045     t0Map_ = &setup.getData(t0Token_);
0046     LogVerbatim("Calibration") << "[DTT0ChamberReferenceCorrection] T0 version: " << t0H->version();
0047   }
0048 
0049   DTT0Data DTT0ChamberReferenceCorrection::correction(const DTWireId& wireId) {
0050     // Compute for selected chamber (or All) correction using as reference chamber mean
0051 
0052     DTChamberId chamberId = wireId.layerId().superlayerId().chamberId();
0053 
0054     if (calibChamber_.empty() || calibChamber_ == "None")
0055       return defaultT0(wireId);
0056     if (calibChamber_ != "All" && chamberId != chosenChamberId_)
0057       return defaultT0(wireId);
0058 
0059     // Access DB
0060     float t0Mean, t0RMS;
0061     int status = t0Map_->get(wireId, t0Mean, t0RMS, DTTimeUnits::counts);
0062     if (status != 0)
0063       throw cms::Exception("[DTT0ChamberReferenceCorrection]") << "Could not find t0 entry in DB for" << wireId << endl;
0064     /*
0065     Leaving just the structure for future implementation
0066     ...
0067     ...
0068   */
0069     return DTT0Data(t0Mean, t0RMS);
0070   }
0071 
0072   DTT0Data DTT0ChamberReferenceCorrection::defaultT0(const DTWireId& wireId) {
0073     // Access default DB
0074     float t0Mean, t0RMS;
0075     int status = t0Map_->get(wireId, t0Mean, t0RMS, DTTimeUnits::counts);
0076     if (!status) {
0077       return DTT0Data(t0Mean, t0RMS);
0078     } else {
0079       //...
0080       throw cms::Exception("[DTT0ChamberReferenceCorrection]") << "Could not find t0 entry in DB for" << wireId << endl;
0081     }
0082   }
0083 
0084 }  // namespace dtCalibration