CMS 3D CMS Logo

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

#include <ETLDeviceSim.h>

Public Member Functions

 ETLDeviceSim (const edm::ParameterSet &pset, edm::ConsumesCollector iC)
 
void getEvent (const edm::Event &evt)
 
void getEventSetup (const edm::EventSetup &evt)
 
void getHitsResponse (const std::vector< std::tuple< int, uint32_t, float > > &hitRefs, const edm::Handle< edm::PSimHitContainer > &hits, mtd_digitizer::MTDSimHitDataAccumulator *simHitAccumulator, CLHEP::HepRandomEngine *hre)
 

Private Attributes

const bool applyDegradation_
 
float bxTime_
 
const reco::FormulaEvaluator fluence_
 
const MTDGeometrygeom_
 
const edm::ESGetToken< MTDGeometry, MTDDigiGeometryRecordgeomToken_
 
const float integratedLum_
 
const reco::FormulaEvaluator lgadGain_
 
const reco::FormulaEvaluator lgadGainDegradation_
 
float MIPPerMeV_
 
const reco::FormulaEvaluator MPVElectron_
 
const reco::FormulaEvaluator MPVKaon_
 
const reco::FormulaEvaluator MPVMuon_
 
const reco::FormulaEvaluator MPVPion_
 
const reco::FormulaEvaluator MPVProton_
 
float tofDelay_
 

Detailed Description

Definition at line 23 of file ETLDeviceSim.h.

Constructor & Destructor Documentation

◆ ETLDeviceSim()

ETLDeviceSim::ETLDeviceSim ( const edm::ParameterSet pset,
edm::ConsumesCollector  iC 
)

Definition at line 15 of file ETLDeviceSim.cc.

16  : geomToken_(iC.esConsumes()),
17  geom_(nullptr),
18  MIPPerMeV_(1.0 / pset.getParameter<double>("meVPerMIP")),
19  integratedLum_(pset.getParameter<double>("IntegratedLuminosity")),
20  fluence_(pset.getParameter<std::string>("FluenceVsRadius")),
21  lgadGain_(pset.getParameter<std::string>("LGADGainVsFluence")),
22  lgadGainDegradation_(pset.getParameter<std::string>("LGADGainDegradation")),
23  applyDegradation_(pset.getParameter<bool>("applyDegradation")),
24  bxTime_(pset.getParameter<double>("bxTime")),
25  tofDelay_(pset.getParameter<double>("tofDelay")),
26  MPVMuon_(pset.getParameter<std::string>("MPVMuon")),
27  MPVPion_(pset.getParameter<std::string>("MPVPion")),
28  MPVKaon_(pset.getParameter<std::string>("MPVKaon")),
29  MPVElectron_(pset.getParameter<std::string>("MPVElectron")),
30  MPVProton_(pset.getParameter<std::string>("MPVProton")) {}
const MTDGeometry * geom_
Definition: ETLDeviceSim.h:38
const reco::FormulaEvaluator lgadGain_
Definition: ETLDeviceSim.h:42
const edm::ESGetToken< MTDGeometry, MTDDigiGeometryRecord > geomToken_
Definition: ETLDeviceSim.h:37
const reco::FormulaEvaluator MPVKaon_
Definition: ETLDeviceSim.h:49
const bool applyDegradation_
Definition: ETLDeviceSim.h:44
const reco::FormulaEvaluator MPVMuon_
Definition: ETLDeviceSim.h:47
const reco::FormulaEvaluator MPVElectron_
Definition: ETLDeviceSim.h:50
const reco::FormulaEvaluator fluence_
Definition: ETLDeviceSim.h:41
const float integratedLum_
Definition: ETLDeviceSim.h:40
const reco::FormulaEvaluator lgadGainDegradation_
Definition: ETLDeviceSim.h:43
const reco::FormulaEvaluator MPVProton_
Definition: ETLDeviceSim.h:51
float tofDelay_
Definition: ETLDeviceSim.h:46
float MIPPerMeV_
Definition: ETLDeviceSim.h:39
const reco::FormulaEvaluator MPVPion_
Definition: ETLDeviceSim.h:48

Member Function Documentation

◆ getEvent()

void ETLDeviceSim::getEvent ( const edm::Event evt)
inline

Definition at line 27 of file ETLDeviceSim.h.

27 {}

◆ getEventSetup()

void ETLDeviceSim::getEventSetup ( const edm::EventSetup evt)

Definition at line 32 of file ETLDeviceSim.cc.

References geom_, geomToken_, and edm::EventSetup::getData().

32 { geom_ = &evs.getData(geomToken_); }
const MTDGeometry * geom_
Definition: ETLDeviceSim.h:38
const edm::ESGetToken< MTDGeometry, MTDDigiGeometryRecord > geomToken_
Definition: ETLDeviceSim.h:37

◆ getHitsResponse()

void ETLDeviceSim::getHitsResponse ( const std::vector< std::tuple< int, uint32_t, float > > &  hitRefs,
const edm::Handle< edm::PSimHitContainer > &  hits,
mtd_digitizer::MTDSimHitDataAccumulator simHitAccumulator,
CLHEP::HepRandomEngine *  hre 
)

Definition at line 34 of file ETLDeviceSim.cc.

References funct::abs(), applyDegradation_, bxTime_, ALCARECOTkAlJpsiMuMu_cff::charge, cuy::col, angle_units::operators::convertGeVToMeV(), angle_units::operators::convertMmToCm(), TauDecayModes::dec, hcalRecHitTable_cff::detId, reco::FormulaEvaluator::evaluate(), Exception, fluence_, DetId::Forward, PedestalClient_cfi::gain, geom_, hfClusterShapes_cfi::hits, mps_fire::i, MTDGeometry::idToDet(), integratedLum_, RectangularMTDTopology::isInPixel(), lgadGain_, lgadGainDegradation_, LogTrace, METSkim_cff::Max, MIPPerMeV_, MPVElectron_, MPVKaon_, MPVMuon_, MPVPion_, MPVProton_, PbPb_ZMuSkimMuonDPG_cff::particleType, RectangularMTDTopology::pixelIndex(), position, CosmicsPD_Skims::radius, DetId::rawId(), ProxyMTDTopology::specificTopology(), tofDelay_, GeomDet::toGlobal(), GeomDet::topology(), PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

37  {
38  using namespace geant_units::operators;
39 
40  std::vector<double> emptyV;
41  std::vector<double> radius(1);
42  std::vector<double> fluence(1);
43  std::vector<double> gain(1);
44  std::vector<double> param(2);
45  std::vector<double> momentum(1);
46 
47  //loop over sorted hits
48  const int nchits = hitRefs.size();
49  LogTrace("ETLDeviceSim") << "Processing " << nchits << " SIM hits";
50 
51  for (int i = 0; i < nchits; ++i) {
52  const int hitidx = std::get<0>(hitRefs[i]);
53  const uint32_t id = std::get<1>(hitRefs[i]);
54  const MTDDetId detId(id);
55 
56  // Safety check (this should never happen, it should be an exception
57  if (detId.det() != DetId::Forward || detId.mtdSubDetector() != 2) {
58  throw cms::Exception("ETLDeviceSim")
59  << "got a DetId that was not ETL: Det = " << detId.det() << " subDet = " << detId.mtdSubDetector();
60  }
61 
62  if (id == 0)
63  continue; // to be ignored at RECO level
64 
65  ETLDetId etlid(detId);
66  DetId geoId = ETLDetId(etlid);
67  const MTDGeomDet* thedet = geom_->idToDet(geoId);
68  if (thedet == nullptr) {
69  throw cms::Exception("ETLDeviceSim") << "GeographicalID: " << std::hex << geoId.rawId() << " (" << detId.rawId()
70  << ") is invalid!" << std::dec << std::endl;
71  }
72  const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(thedet->topology());
73  const RectangularMTDTopology& topo = static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
74 
75  const float toa = std::get<2>(hitRefs[i]) + tofDelay_;
76  const PSimHit& hit = hits->at(hitidx);
77  float charge = convertGeVToMeV(hit.energyLoss()) * MIPPerMeV_;
78 
79  momentum[0] = hit.pabs();
80 
81  // particle type
82  int particleType = abs(hit.particleType());
83  float MPV_ = 0;
84  if (particleType == 11) {
85  MPV_ = MPVElectron_.evaluate(momentum, emptyV);
86  } else if (particleType == 13) {
87  MPV_ = MPVMuon_.evaluate(momentum, emptyV);
88  } else if (particleType == 211) {
89  MPV_ = MPVPion_.evaluate(momentum, emptyV);
90  } else if (particleType == 321) {
91  MPV_ = MPVKaon_.evaluate(momentum, emptyV);
92  } else {
93  MPV_ = MPVProton_.evaluate(momentum, emptyV);
94  }
95  float MPV_charge = convertGeVToMeV(MPV_) * MIPPerMeV_;
96 
97  // calculate the simhit row and column
98  const auto& position = hit.localPosition();
99 
100  // ETL is already in module-local coordinates so just scale to cm from mm
102  const auto& global_point = thedet->toGlobal(simscaled);
103 
104  radius[0] = global_point.perp();
105  fluence[0] = integratedLum_ * fluence_.evaluate(radius, emptyV);
106  gain[0] = lgadGain_.evaluate(fluence, emptyV);
107 
108  //The following lines check whether the pixel point is actually out of the active area.
109  if (topo.isInPixel(simscaled)) {
110  charge *= gain[0];
111  MPV_charge *= gain[0];
112  } else {
113  if (applyDegradation_) {
114  double dGapCenter = TMath::Max(TMath::Abs(simscaled.x()), TMath::Abs(simscaled.y()));
115  param[0] = gain[0];
116  param[1] = dGapCenter;
117  gain[0] = lgadGainDegradation_.evaluate(param, emptyV);
118  charge *= gain[0];
119  MPV_charge *= gain[0];
120  }
121  }
122 
123  const auto& thepixel = topo.pixelIndex(simscaled);
124  const uint8_t row = static_cast<uint8_t>(thepixel.first);
125  const uint8_t col = static_cast<uint8_t>(thepixel.second);
126  LogTrace("ETLDeviceSim") << "Processing hit in pixel # " << hitidx << " DetId " << etlid.rawId() << " row/col "
127  << (uint32_t)row << " " << (uint32_t)col << " inPixel " << topo.isInPixel(simscaled)
128  << " tof " << toa << " ene " << hit.energyLoss() << " MIP " << MIPPerMeV_ << " gain "
129  << gain[0] << " charge " << charge << " MPV " << MPV_charge;
130 
131  auto simHitIt =
132  simHitAccumulator->emplace(mtd_digitizer::MTDCellId(id, row, col), mtd_digitizer::MTDCellInfo()).first;
133 
134  // Accumulate in 15 buckets of 25ns (9 pre-samples, 1 in-time, 5 post-samples)
135  const int itime = std::floor(toa / bxTime_) + 9;
136  if (itime < 0 || itime > 14)
137  continue;
138 
139  // Check if time index is ok and store energy
140  if (itime >= (int)simHitIt->second.hit_info[0].size())
141  continue;
142  (simHitIt->second).hit_info[0][itime] += charge;
143  // Store the time of the first SimHit in the right DataFrame bucket
144  const float tof = toa - (itime - 9) * bxTime_;
145 
146  if ((simHitIt->second).hit_info[1][itime] == 0. || tof < (simHitIt->second).hit_info[1][itime]) {
147  (simHitIt->second).hit_info[1][itime] = tof;
148  }
149  (simHitIt->second).hit_info[2][itime] += MPV_charge;
150  }
151 }
const MTDGeometry * geom_
Definition: ETLDeviceSim.h:38
const reco::FormulaEvaluator lgadGain_
Definition: ETLDeviceSim.h:42
bool isInPixel(const LocalPoint &p) const
double evaluate(V const &iVariables, P const &iParameters) const
virtual const Topology & topology() const
Definition: GeomDet.cc:67
virtual const PixelTopology & specificTopology() const
const reco::FormulaEvaluator MPVKaon_
Definition: ETLDeviceSim.h:49
const bool applyDegradation_
Definition: ETLDeviceSim.h:44
#define LogTrace(id)
Detector identifier base class for the MIP Timing Layer.
Definition: MTDDetId.h:21
const reco::FormulaEvaluator MPVMuon_
Definition: ETLDeviceSim.h:47
const reco::FormulaEvaluator MPVElectron_
Definition: ETLDeviceSim.h:50
constexpr NumType convertGeVToMeV(NumType gev)
Definition: angle_units.h:74
const reco::FormulaEvaluator fluence_
Definition: ETLDeviceSim.h:41
const MTDGeomDet * idToDet(DetId) const override
Definition: MTDGeometry.cc:171
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const float integratedLum_
Definition: ETLDeviceSim.h:40
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
Definition: DetId.h:17
const reco::FormulaEvaluator lgadGainDegradation_
Definition: ETLDeviceSim.h:43
constexpr NumType convertMmToCm(NumType millimeters)
Definition: angle_units.h:44
const reco::FormulaEvaluator MPVProton_
Definition: ETLDeviceSim.h:51
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
float tofDelay_
Definition: ETLDeviceSim.h:46
Detector identifier class for the Endcap Timing Layer.
Definition: ETLDetId.h:16
static int position[264][3]
Definition: ReadPGInfo.cc:289
col
Definition: cuy.py:1009
std::pair< int, int > pixelIndex(const LocalPoint &p) const
float MIPPerMeV_
Definition: ETLDeviceSim.h:39
const reco::FormulaEvaluator MPVPion_
Definition: ETLDeviceSim.h:48

Member Data Documentation

◆ applyDegradation_

const bool ETLDeviceSim::applyDegradation_
private

Definition at line 44 of file ETLDeviceSim.h.

Referenced by getHitsResponse().

◆ bxTime_

float ETLDeviceSim::bxTime_
private

Definition at line 45 of file ETLDeviceSim.h.

Referenced by getHitsResponse().

◆ fluence_

const reco::FormulaEvaluator ETLDeviceSim::fluence_
private

Definition at line 41 of file ETLDeviceSim.h.

Referenced by getHitsResponse().

◆ geom_

const MTDGeometry* ETLDeviceSim::geom_
private

Definition at line 38 of file ETLDeviceSim.h.

Referenced by getEventSetup(), and getHitsResponse().

◆ geomToken_

const edm::ESGetToken<MTDGeometry, MTDDigiGeometryRecord> ETLDeviceSim::geomToken_
private

Definition at line 37 of file ETLDeviceSim.h.

Referenced by getEventSetup().

◆ integratedLum_

const float ETLDeviceSim::integratedLum_
private

Definition at line 40 of file ETLDeviceSim.h.

Referenced by getHitsResponse().

◆ lgadGain_

const reco::FormulaEvaluator ETLDeviceSim::lgadGain_
private

Definition at line 42 of file ETLDeviceSim.h.

Referenced by getHitsResponse().

◆ lgadGainDegradation_

const reco::FormulaEvaluator ETLDeviceSim::lgadGainDegradation_
private

Definition at line 43 of file ETLDeviceSim.h.

Referenced by getHitsResponse().

◆ MIPPerMeV_

float ETLDeviceSim::MIPPerMeV_
private

Definition at line 39 of file ETLDeviceSim.h.

Referenced by getHitsResponse().

◆ MPVElectron_

const reco::FormulaEvaluator ETLDeviceSim::MPVElectron_
private

Definition at line 50 of file ETLDeviceSim.h.

Referenced by getHitsResponse().

◆ MPVKaon_

const reco::FormulaEvaluator ETLDeviceSim::MPVKaon_
private

Definition at line 49 of file ETLDeviceSim.h.

Referenced by getHitsResponse().

◆ MPVMuon_

const reco::FormulaEvaluator ETLDeviceSim::MPVMuon_
private

Definition at line 47 of file ETLDeviceSim.h.

Referenced by getHitsResponse().

◆ MPVPion_

const reco::FormulaEvaluator ETLDeviceSim::MPVPion_
private

Definition at line 48 of file ETLDeviceSim.h.

Referenced by getHitsResponse().

◆ MPVProton_

const reco::FormulaEvaluator ETLDeviceSim::MPVProton_
private

Definition at line 51 of file ETLDeviceSim.h.

Referenced by getHitsResponse().

◆ tofDelay_

float ETLDeviceSim::tofDelay_
private

Definition at line 46 of file ETLDeviceSim.h.

Referenced by getHitsResponse().