CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
dtCalibration::DTTTrigConstantShift Class Reference

#include <DTTTrigConstantShift.h>

Inheritance diagram for dtCalibration::DTTTrigConstantShift:
dtCalibration::DTTTrigBaseCorrection

Public Member Functions

DTTTrigData correction (const DTSuperLayerId &) override
 
 DTTTrigConstantShift (const edm::ParameterSet &, edm::ConsumesCollector)
 
void setES (const edm::EventSetup &setup) override
 
 ~DTTTrigConstantShift () override
 
- Public Member Functions inherited from dtCalibration::DTTTrigBaseCorrection
 DTTTrigBaseCorrection ()
 
virtual ~DTTTrigBaseCorrection ()
 

Private Attributes

std::string calibChamber_
 
DTChamberId chosenChamberId_
 
const DTTtrigtTrigMap_
 
edm::ESGetToken< DTTtrig, DTTtrigRcdttrigToken_
 
double value_
 

Detailed Description

Definition at line 27 of file DTTTrigConstantShift.h.

Constructor & Destructor Documentation

◆ DTTTrigConstantShift()

DTTTrigConstantShift::DTTTrigConstantShift ( const edm::ParameterSet pset,
edm::ConsumesCollector  cc 
)

Definition at line 24 of file DTTTrigConstantShift.cc.

References edm::BeginRun, calibChamber_, gpuPixelDoublets::cc, chosenChamberId_, muonDTDigis_cfi::pset, ttrigToken_, and value_.

25  : calibChamber_(pset.getParameter<string>("calibChamber")), value_(pset.getParameter<double>("value")) {
26  ttrigToken_ =
27  cc.esConsumes<edm::Transition::BeginRun>(edm::ESInputTag("", pset.getUntrackedParameter<string>("dbLabel")));
28  LogVerbatim("Calibration") << "[DTTTrigConstantShift] Applying constant correction value: " << value_ << endl;
29 
30  if (!calibChamber_.empty() && calibChamber_ != "None" && calibChamber_ != "All") {
31  stringstream linestr;
32  int selWheel, selStation, selSector;
33  linestr << calibChamber_;
34  linestr >> selWheel >> selStation >> selSector;
35  chosenChamberId_ = DTChamberId(selWheel, selStation, selSector);
36  LogVerbatim("Calibration") << "[DTTTrigConstantShift] Chosen chamber: " << chosenChamberId_ << endl;
37  }
38  //FIXME: Check if chosen chamber is valid.
39  }
Log< level::Info, true > LogVerbatim
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
edm::ESGetToken< DTTtrig, DTTtrigRcd > ttrigToken_

◆ ~DTTTrigConstantShift()

DTTTrigConstantShift::~DTTTrigConstantShift ( )
override

Definition at line 41 of file DTTTrigConstantShift.cc.

41 {}

Member Function Documentation

◆ correction()

DTTTrigData DTTTrigConstantShift::correction ( const DTSuperLayerId slId)
overridevirtual

Implements dtCalibration::DTTTrigBaseCorrection.

Definition at line 49 of file DTTTrigConstantShift.cc.

References calibChamber_, DTSuperLayerId::chamberId(), chosenChamberId_, Exception, DTTtrig::get(), dtTTrigCalibration_cfi::kFactor, DTTimeUnits::ns, mps_update::status, tTrigMap_, and value_.

49  {
50  float tTrigMean, tTrigSigma, kFactor;
51  int status = tTrigMap_->get(slId, tTrigMean, tTrigSigma, kFactor, DTTimeUnits::ns);
52  if (status != 0)
53  throw cms::Exception("[DTTTrigConstantShift]") << "Could not find tTrig entry in DB for" << slId << endl;
54 
55  float tTrigMeanNew = tTrigMean;
56  if (!calibChamber_.empty() && calibChamber_ != "None") {
57  if ((calibChamber_ == "All") || (calibChamber_ != "All" && slId.chamberId() == chosenChamberId_)) {
58  tTrigMeanNew = tTrigMean + value_;
59  }
60  }
61 
62  return DTTTrigData(tTrigMeanNew, tTrigSigma, kFactor);
63  }
DTChamberId chamberId() const
Return the corresponding ChamberId.
int get(int wheelId, int stationId, int sectorId, int slId, float &tTrig, float &tTrms, float &kFact, DTTimeUnits::type unit) const
get content
Definition: DTTtrig.cc:59

◆ setES()

void DTTTrigConstantShift::setES ( const edm::EventSetup setup)
overridevirtual

Member Data Documentation

◆ calibChamber_

std::string dtCalibration::DTTTrigConstantShift::calibChamber_
private

Definition at line 39 of file DTTTrigConstantShift.h.

Referenced by correction(), and DTTTrigConstantShift().

◆ chosenChamberId_

DTChamberId dtCalibration::DTTTrigConstantShift::chosenChamberId_
private

Definition at line 43 of file DTTTrigConstantShift.h.

Referenced by correction(), and DTTTrigConstantShift().

◆ tTrigMap_

const DTTtrig* dtCalibration::DTTTrigConstantShift::tTrigMap_
private

Definition at line 42 of file DTTTrigConstantShift.h.

Referenced by correction(), and setES().

◆ ttrigToken_

edm::ESGetToken<DTTtrig, DTTtrigRcd> dtCalibration::DTTTrigConstantShift::ttrigToken_
private

Definition at line 44 of file DTTTrigConstantShift.h.

Referenced by DTTTrigConstantShift(), and setES().

◆ value_

double dtCalibration::DTTTrigConstantShift::value_
private

Definition at line 40 of file DTTTrigConstantShift.h.

Referenced by correction(), and DTTTrigConstantShift().