CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes | Static Private Attributes
CTPPSDiamondRecHitProducerAlgorithm Class Reference

#include <CTPPSDiamondRecHitProducerAlgorithm.h>

Public Member Functions

void build (const CTPPSGeometry &, const edm::DetSetVector< CTPPSDiamondDigi > &, edm::DetSetVector< CTPPSDiamondRecHit > &)
 
 CTPPSDiamondRecHitProducerAlgorithm (const edm::ParameterSet &conf)
 
void setCalibration (const PPSTimingCalibration &)
 

Private Attributes

bool apply_calib_
 Switch on/off the timing calibration. More...
 
PPSTimingCalibration calib_
 
std::unique_ptr< reco::FormulaEvaluatorcalib_fct_
 
double ts_to_ns_
 Conversion constant between HPTDC time slice and absolute time (in ns) More...
 

Static Private Attributes

static constexpr unsigned short MAX_CHANNEL = 20
 

Detailed Description

Definition at line 25 of file CTPPSDiamondRecHitProducerAlgorithm.h.

Constructor & Destructor Documentation

◆ CTPPSDiamondRecHitProducerAlgorithm()

CTPPSDiamondRecHitProducerAlgorithm::CTPPSDiamondRecHitProducerAlgorithm ( const edm::ParameterSet conf)

Definition at line 13 of file CTPPSDiamondRecHitProducerAlgorithm.cc.

15  : ts_to_ns_(iConfig.getParameter<double>("timeSliceNs")),

Member Function Documentation

◆ build()

void CTPPSDiamondRecHitProducerAlgorithm::build ( const CTPPSGeometry geom,
const edm::DetSetVector< CTPPSDiamondDigi > &  input,
edm::DetSetVector< CTPPSDiamondRecHit > &  output 
)

Definition at line 22 of file CTPPSDiamondRecHitProducerAlgorithm.cc.

25  {
26  for (const auto& vec : input) {
27  const CTPPSDiamondDetId detid(vec.detId());
28 
29  if (detid.channel() > MAX_CHANNEL) // VFAT-like information, to be ignored
30  continue;
31 
32  // retrieve the geometry element associated to this DetID
33  const DetGeomDesc* det = geom.sensor(detid);
34 
35  const float x_pos = det->translation().x(), y_pos = det->translation().y();
36  float z_pos = 0.;
37  z_pos = det->parentZPosition(); // retrieve the plane position;
38 
39  // parameters stand for half the size
40  const float x_width = 2.0 * det->params().at(0);
41  const float y_width = 2.0 * det->params().at(1);
42  const float z_width = 2.0 * det->params().at(2);
43 
44  // retrieve the timing calibration part for this channel
45  const int sector = detid.arm(), station = detid.station(), plane = detid.plane(), channel = detid.channel();
46  const auto& ch_params = (apply_calib_) ? calib_.parameters(sector, station, plane, channel) : std::vector<double>{};
47  // default values for offset + time precision if calibration object not found
48  const double ch_t_offset = (apply_calib_) ? calib_.timeOffset(sector, station, plane, channel) : 0.;
49  const double ch_t_precis = (apply_calib_) ? calib_.timePrecision(sector, station, plane, channel) : 0.;
50 
51  edm::DetSet<CTPPSDiamondRecHit>& rec_hits = output.find_or_insert(detid);
52 
53  for (const auto& digi : vec) {
54  const int t_lead = digi.leadingEdge(), t_trail = digi.trailingEdge();
55  // skip invalid digis
56  if (t_lead == 0 && t_trail == 0)
57  continue;
58 
59  double tot = -1., ch_t_twc = 0.;
60  if (t_lead != 0 && t_trail != 0) {
61  tot = (t_trail - t_lead) * ts_to_ns_; // in ns
62  if (calib_fct_ && apply_calib_) {
63  // compute the time-walk correction
64  ch_t_twc = calib_fct_->evaluate(std::vector<double>{tot}, ch_params);
65  if (edm::isNotFinite(ch_t_twc))
66  ch_t_twc = 0.;
67  }
68  }
69 
70  const int time_slice =
71  (t_lead != 0) ? (t_lead - ch_t_offset / ts_to_ns_) / 1024 : CTPPSDiamondRecHit::TIMESLICE_WITHOUT_LEADING;
72 
73  // calibrated time of arrival
74  const double t0 = (t_lead % 1024) * ts_to_ns_ - ch_t_twc;
75 
76  rec_hits.emplace_back(
77  // spatial information
78  x_pos,
79  x_width,
80  y_pos,
81  y_width,
82  z_pos,
83  z_width,
84  // timing information
85  t0,
86  tot,
87  ch_t_precis,
88  time_slice,
89  // readout information
90  digi.hptdcErrorFlags(),
91  digi.multipleHit());
92  }
93  }

References apply_calib_, calib_, calib_fct_, edm::DetSet< T >::emplace_back(), relativeConstraints::geom, input, edm::isNotFinite(), MAX_CHANNEL, convertSQLitetoXML_cfg::output, PPSTimingCalibration::parameters(), DetGeomDesc::params(), DetGeomDesc::parentZPosition(), relativeConstraints::station, FrontierCondition_GT_autoExpress_cfi::t0, PPSTimingCalibration::timeOffset(), PPSTimingCalibration::timePrecision(), CTPPSDiamondRecHit::TIMESLICE_WITHOUT_LEADING, DetGeomDesc::translation(), and ts_to_ns_.

Referenced by CTPPSDiamondRecHitProducer::produce().

◆ setCalibration()

void CTPPSDiamondRecHitProducerAlgorithm::setCalibration ( const PPSTimingCalibration calib)

Member Data Documentation

◆ apply_calib_

bool CTPPSDiamondRecHitProducerAlgorithm::apply_calib_
private

Switch on/off the timing calibration.

Definition at line 42 of file CTPPSDiamondRecHitProducerAlgorithm.h.

Referenced by build().

◆ calib_

PPSTimingCalibration CTPPSDiamondRecHitProducerAlgorithm::calib_
private

Definition at line 43 of file CTPPSDiamondRecHitProducerAlgorithm.h.

Referenced by build(), and setCalibration().

◆ calib_fct_

std::unique_ptr<reco::FormulaEvaluator> CTPPSDiamondRecHitProducerAlgorithm::calib_fct_
private

Definition at line 44 of file CTPPSDiamondRecHitProducerAlgorithm.h.

Referenced by build(), and setCalibration().

◆ MAX_CHANNEL

constexpr unsigned short CTPPSDiamondRecHitProducerAlgorithm::MAX_CHANNEL = 20
staticconstexprprivate

Definition at line 38 of file CTPPSDiamondRecHitProducerAlgorithm.h.

Referenced by build().

◆ ts_to_ns_

double CTPPSDiamondRecHitProducerAlgorithm::ts_to_ns_
private

Conversion constant between HPTDC time slice and absolute time (in ns)

Definition at line 40 of file CTPPSDiamondRecHitProducerAlgorithm.h.

Referenced by build().

CTPPSDiamondRecHitProducerAlgorithm::calib_
PPSTimingCalibration calib_
Definition: CTPPSDiamondRecHitProducerAlgorithm.h:43
input
static const std::string input
Definition: EdmProvDump.cc:48
edm::isNotFinite
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
edm::DetSet
Definition: DetSet.h:23
convertSQLitetoXML_cfg.output
output
Definition: convertSQLitetoXML_cfg.py:32
relativeConstraints.station
station
Definition: relativeConstraints.py:67
DetGeomDesc::parentZPosition
float parentZPosition() const
Definition: DetGeomDesc.h:56
PPSTimingCalibration::formula
const std::string & formula() const
Definition: PPSTimingCalibration.h:46
CTPPSDiamondRecHitProducerAlgorithm::ts_to_ns_
double ts_to_ns_
Conversion constant between HPTDC time slice and absolute time (in ns)
Definition: CTPPSDiamondRecHitProducerAlgorithm.h:40
CTPPSDiamondRecHitProducerAlgorithm::apply_calib_
bool apply_calib_
Switch on/off the timing calibration.
Definition: CTPPSDiamondRecHitProducerAlgorithm.h:42
CTPPSDiamondRecHitProducerAlgorithm::MAX_CHANNEL
static constexpr unsigned short MAX_CHANNEL
Definition: CTPPSDiamondRecHitProducerAlgorithm.h:38
PPSTimingCalibration::timeOffset
double timeOffset(int key1, int key2, int key3, int key4=-1) const
Definition: PPSTimingCalibration.cc:39
reco::FormulaEvaluator
Definition: FormulaEvaluator.h:67
FrontierCondition_GT_autoExpress_cfi.t0
t0
Definition: FrontierCondition_GT_autoExpress_cfi.py:148
PPSTimingCalibration::timePrecision
double timePrecision(int key1, int key2, int key3, int key4=-1) const
Definition: PPSTimingCalibration.cc:47
relativeConstraints.geom
geom
Definition: relativeConstraints.py:72
PPSTimingCalibration::parameters
std::vector< double > parameters(int key1, int key2, int key3, int key4) const
Definition: PPSTimingCalibration.cc:31
CTPPSDiamondDetId
Detector ID class for CTPPS Timing Diamond detectors. Bits [19:31] : Assigend in CTPPSDetId Calss Bit...
Definition: CTPPSDiamondDetId.h:24
edm::DetSet::emplace_back
decltype(auto) emplace_back(Args &&... args)
Definition: DetSet.h:68
CTPPSDiamondRecHit::TIMESLICE_WITHOUT_LEADING
static constexpr int TIMESLICE_WITHOUT_LEADING
Definition: CTPPSDiamondRecHit.h:44
DetGeomDesc
Geometrical description of a sensor.
Definition: DetGeomDesc.h:34
L1EGammaCrystalsEmulatorProducer_cfi.calib
calib
Definition: L1EGammaCrystalsEmulatorProducer_cfi.py:6
CTPPSDiamondRecHitProducerAlgorithm::calib_fct_
std::unique_ptr< reco::FormulaEvaluator > calib_fct_
Definition: CTPPSDiamondRecHitProducerAlgorithm.h:44
DetGeomDesc::params
std::vector< double > params() const
Definition: DetGeomDesc.h:66
DetGeomDesc::translation
Translation translation() const
Definition: DetGeomDesc.h:64