CMS 3D CMS Logo

Pixel3DDigitizerAlgorithm.cc
Go to the documentation of this file.
2 
3 // Framework infrastructure
7 
8 // Calibration & Conditions
10 
11 // Geometry
13 
14 //#include <iostream>
15 #include <cmath>
16 #include <vector>
17 #include <algorithm>
18 
19 using namespace sipixelobjects;
20 
21 namespace {
22  // Analogously to CMSUnits (no um defined)
23  constexpr double operator""_um(long double length) { return length * 1e-4; }
24  constexpr double operator""_um_inv(long double length) { return length * 1e4; }
25 } // namespace
26 
28  // XXX: Just copied from PixelDigitizer Algorithm
29  // CHECK if all these is needed
30 
31  if (use_ineff_from_db_) {
32  // load gain calibration service fromdb...
33  theSiPixelGainCalibrationService_->setESObjects(es);
34  }
35 
36  if (use_deadmodule_DB_) {
37  siPixelBadModule_ = &es.getData(siPixelBadModuleToken_);
38  }
39 
40  if (use_LorentzAngle_DB_) {
41  // Get Lorentz angle from DB record
42  siPixelLorentzAngle_ = &es.getData(siPixelLorentzAngleToken_);
43  }
44 
45  // gets the map and geometry from the DB (to kill ROCs)
46  fedCablingMap_ = &es.getData(fedCablingMapToken_);
47  geom_ = &es.getData(geomToken_);
48 }
49 
51  : Phase2TrackerDigitizerAlgorithm(conf.getParameter<edm::ParameterSet>("AlgorithmCommon"),
52  conf.getParameter<edm::ParameterSet>("Pixel3DDigitizerAlgorithm"),
53  iC),
54  np_column_radius_(
55  (conf.getParameter<edm::ParameterSet>("Pixel3DDigitizerAlgorithm").getParameter<double>("NPColumnRadius")) *
56  1.0_um),
57  ohm_column_radius_(
58  (conf.getParameter<edm::ParameterSet>("Pixel3DDigitizerAlgorithm").getParameter<double>("OhmicColumnRadius")) *
59  1.0_um),
60  np_column_gap_(
61  (conf.getParameter<edm::ParameterSet>("Pixel3DDigitizerAlgorithm").getParameter<double>("NPColumnGap")) *
62  1.0_um),
63  fedCablingMapToken_(iC.esConsumes()),
64  geomToken_(iC.esConsumes()) {
65  // XXX - NEEDED?
66  pixelFlag_ = true;
67 
72 
73  edm::LogInfo("Pixel3DDigitizerAlgorithm")
74  << "Algorithm constructed \n"
75  << "Configuration parameters:\n"
76  << "\n*** Threshold"
77  << "\n Endcap = " << theThresholdInE_Endcap_ << " electrons"
78  << "\n Barrel = " << theThresholdInE_Barrel_ << " electrons"
79  << "\n*** Gain"
80  << "\n Electrons per ADC:" << theElectronPerADC_ << "\n ADC Full Scale: " << theAdcFullScale_
81  << "\n*** The delta cut-off is set to " << tMax_ << "\n*** Pixel-inefficiency: " << addPixelInefficiency_;
82 }
83 
85 
86 //
87 // -- Select the Hit for Digitization
88 //
89 bool Pixel3DDigitizerAlgorithm::select_hit(const PSimHit& hit, double tCorr, double& sigScale) const {
90  double time = hit.tof() - tCorr;
91  return (time >= theTofLowerCut_ && time < theTofUpperCut_);
92 }
93 
94 const bool Pixel3DDigitizerAlgorithm::is_inside_n_column_(const LocalPoint& p, const float& sensor_thickness) const {
95  // The insensitive volume of the column: sensor thickness - column gap distance
96  return (p.perp() <= np_column_radius_ && p.z() <= (sensor_thickness - np_column_gap_));
97 }
98 
100  const std::pair<float, float>& half_pitch) const {
101  // The four corners of the cell
102  return ((p - LocalVector(half_pitch.first, half_pitch.second, 0)).perp() <= ohm_column_radius_) ||
103  ((p - LocalVector(-half_pitch.first, half_pitch.second, 0)).perp() <= ohm_column_radius_) ||
104  ((p - LocalVector(half_pitch.first, -half_pitch.second, 0)).perp() <= ohm_column_radius_) ||
105  ((p - LocalVector(-half_pitch.first, -half_pitch.second, 0)).perp() <= ohm_column_radius_);
106 }
107 
108 // Diffusion algorithm: Probably not needed,
109 // Assuming the position point is given in the reference system of the proxy
110 // cell, centered at the n-column.
111 // The algorithm assumes only 1-axis could produce the charge migration, this assumption
112 // could be enough given that the p-columns (5 um radius) are in the corners of the cell
113 // (no producing charge in there)
114 // The output is vector of newly created charge in the neighbour pixel i+1 or i-1,
115 // defined by its position higher than abs(half_pitch) and the the sign providing
116 // the addition or subtraction in the pixel (i+-1)
117 std::vector<DigitizerUtility::EnergyDepositUnit> Pixel3DDigitizerAlgorithm::diffusion(
118  const LocalPoint& pos,
119  const float& ncarriers,
120  const std::function<LocalVector(float, float)>& u_drift,
121  const std::pair<float, float> hpitches,
122  const float& thickness) const {
123  // FIXME -- DM : Note that with a 0.3 will be enough (if using current sigma formulae)
124  // With the current sigma, this value is dependent of the thickness,
125  // Note that this formulae is coming from planar sensors, a similar
126  // study with data will be needed to extract the sigma for 3D
127  const float max_migration_radius = 0.4_um;
128  // Need to know which axis is the relevant one
129  int displ_ind = -1;
130  float pitch = 0.0;
131 
132  // Check the group is near the edge of the pixel, so diffusion will
133  // be relevant in order to migrate between pixel cells
134  if (hpitches.first - std::abs(pos.x()) < max_migration_radius) {
135  displ_ind = 0;
136  pitch = hpitches.first;
137  } else if (hpitches.second - std::abs(pos.y()) < max_migration_radius) {
138  displ_ind = 1;
139  pitch = hpitches.second;
140  } else {
141  // Nothing to do, too far away
142  return std::vector<DigitizerUtility::EnergyDepositUnit>();
143  }
144 
145  // The new EnergyDeposits in the neighbour pixels
146  // (defined by +1 to the right (first axis) and +1 to the up (second axis)) <-- XXX
147  std::vector<DigitizerUtility::EnergyDepositUnit> migrated_charge;
148 
149  // FIXME -- DM
150  const float diffusion_step = 0.1_um;
151 
152  // The position while drifting
153  std::vector<float> pos_moving({pos.x(), pos.y(), pos.z()});
154  // The drifting: drift field and steps
155  std::function<std::vector<float>(int)> do_step =
156  [&pos_moving, &u_drift, diffusion_step](int i) -> std::vector<float> {
157  auto dd = u_drift(pos_moving[0], pos_moving[1]);
158  return std::vector<float>({i * diffusion_step * dd.x(), i * diffusion_step * dd.y(), i * diffusion_step * dd.z()});
159  };
160 
161  LogDebug("Pixel3DDigitizerAlgorithm::diffusion")
162  << "\nMax. radius from the pixel edge to migrate charge: " << max_migration_radius * 1.0_um_inv << " [um]"
163  << "\nMigration axis: " << displ_ind
164  << "\n(super-)Charge distance to the pixel edge: " << (pitch - pos_moving[displ_ind]) * 1.0_um_inv << " [um]";
165 
166  // How many sigmas (probably a configurable, to be decided not now)
167  const float N_SIGMA = 3.0;
168 
169  // Start the drift and check every step
170  // Some variables needed
171  float current_carriers = ncarriers;
172  std::vector<float> newpos({pos_moving[0], pos_moving[1], pos_moving[2]});
173  float distance_edge = 0.0_um;
174  // Current diffusion value
175  const float sigma = 0.4_um;
176  for (int i = 1;; ++i) {
177  std::transform(pos_moving.begin(), pos_moving.end(), do_step(i).begin(), pos_moving.begin(), std::plus<float>());
178  distance_edge = pitch - std::abs(pos_moving[displ_ind]);
179  // Get the amount of charge on the neighbor pixel: note the
180  // transformation to a Normal
181  float migrated_e = current_carriers * 0.5 * (1.0 - std::erf(distance_edge / (sigma * std::sqrt(2.0))));
182 
183  LogDebug("(super-)charge diffusion") << "step-" << i << ", Current carriers Ne= " << current_carriers << ","
184  << "r=(" << pos_moving[0] * 1.0_um_inv << ", " << pos_moving[1] * 1.0_um_inv
185  << ", " << pos_moving[2] * 1.0_um_inv << ") [um], "
186  << "Migrated charge: " << migrated_e;
187 
188  // Move the migrated charge
189  current_carriers -= migrated_e;
190 
191  // Either far away from the edge or almost half of the carriers already migrated
192  if (std::abs(distance_edge) >= max_migration_radius || current_carriers <= 0.5 * ncarriers) {
193  break;
194  }
195 
196  // Create the ionization point:
197  // First update the newpos vector: the new charge position at the neighbouring pixel
198  // is created in the same position as its "parent carriers"
199  // except the direction of migration
200  std::vector<float> newpos(pos_moving);
201  // Let's create the new charge carriers around 3 sigmas away
202  newpos[displ_ind] += std::copysign(N_SIGMA * sigma, newpos[displ_ind]);
203  migrated_charge.push_back(DigitizerUtility::EnergyDepositUnit(migrated_e, newpos[0], newpos[1], newpos[2]));
204  }
205  return migrated_charge;
206 }
207 
208 // ======================================================================
209 //
210 // Drift the charge segments to the column (collection surface)
211 // Include the effect of E-field and B-field
212 //
213 // =====================================================================
214 std::vector<DigitizerUtility::SignalPoint> Pixel3DDigitizerAlgorithm::drift(
215  const PSimHit& hit,
216  const Phase2TrackerGeomDetUnit* pixdet,
217  const GlobalVector& bfield,
218  const std::vector<DigitizerUtility::EnergyDepositUnit>& ionization_points) const {
219  return drift(hit, pixdet, bfield, ionization_points, true);
220 }
221 std::vector<DigitizerUtility::SignalPoint> Pixel3DDigitizerAlgorithm::drift(
222  const PSimHit& hit,
223  const Phase2TrackerGeomDetUnit* pixdet,
224  const GlobalVector& bfield,
225  const std::vector<DigitizerUtility::EnergyDepositUnit>& ionization_points,
226  bool diffusion_activated) const {
227  // -- Current reference system is placed in the center of the module
228  // -- The natural reference frame should be discribed taking advantatge of
229  // -- the cylindrical nature of the pixel geometry -->
230  // -- the new reference frame should be placed in the center of the n-column, and in the
231  // -- surface of the ROC using cylindrical coordinates
232 
233  // Get ROC pitch, half_pitch and sensor thickness to be used to create the
234  // proxy pixel cell reference frame
235  const auto pitch = pixdet->specificTopology().pitch();
236  const auto half_pitch = std::make_pair<float, float>(pitch.first * 0.5, pitch.second * 0.5);
237  const float thickness = pixdet->specificSurface().bounds().thickness();
238  const int nrows = pixdet->specificTopology().nrows();
239  const int ncolumns = pixdet->specificTopology().ncolumns();
240  const float pix_rounding = 0.99;
241 
242  // The maximum radial distance is going to be used to evaluate radiation damage XXX?
243  const float max_radial_distance =
244  std::sqrt(half_pitch.first * half_pitch.first + half_pitch.second * half_pitch.second);
245 
246  // All pixels are going to be translated to a proxy pixel cell (all pixels should behave
247  // equally no matter their position w.r.t. the module) and describe the movements there
248  // Define the center of the pixel local reference frame: the current cartesian local reference
249  // frame is centered at half_width_x,half_width_y,half_thickness
250  // XXX -- This info could be obtained at init/construction time?
251  LocalPoint center_proxy_cell(half_pitch.first, half_pitch.second, -0.5 * thickness);
252 
253  LogDebug("Pixel3DDigitizerAlgorithm::drift")
254  << "Pixel pitch:" << pitch.first * 1.0_um_inv << ", " << pitch.second * 1.0_um_inv << " [um]";
255 
256  // And the drift direction (assumed same for all the sensor)
257  // XXX call the function which will return a functional
258  std::function<LocalVector(float, float)> drift_direction = [](float x, float y) -> LocalVector {
259  const float theta = std::atan2(y, x);
260  return LocalVector(-std::cos(theta), -std::sin(theta), 0.0);
261  };
262  // The output
263  std::vector<DigitizerUtility::SignalPoint> collection_points;
264  //collection_points.resize(ionization_points.size());
265  collection_points.reserve(ionization_points.size());
266 
267  // Radiation damage limit of application
268  // (XXX: No sense for 3D, let this be until decided what algorithm to use)
269  const float RAD_DAMAGE = 0.001;
270 
271  for (const auto& super_charge : ionization_points) {
272  // Extract the pixel cell
273  auto current_pixel = pixdet->specificTopology().pixel(LocalPoint(super_charge.x(), super_charge.y()));
274  // `pixel` function does not check to be in the ROC bounds,
275  // so check it here and fix potential rounding problems.
276  // Careful, this is assuming a rounding problem (1 unit), more than 1 pixel
277  // away is probably showing some backward problem worth it to track.
278  // This is also correcting out of bounds migrated charge from diffusion.
279  // The charge will be moved to the edge of the row/column.
280  current_pixel.first = std::clamp(current_pixel.first, float(0.0), (nrows - 1) + pix_rounding);
281  current_pixel.second = std::clamp(current_pixel.second, float(0.0), (ncolumns - 1) + pix_rounding);
282 
283  const auto current_pixel_int = std::make_pair(std::floor(current_pixel.first), std::floor(current_pixel.second));
284 
285  // Convert to the 1x1 proxy pixel cell (pc), where all calculations are going to be
286  // performed. The pixel is scaled to the actual pitch
287  const auto relative_position_at_pc =
288  std::make_pair((current_pixel.first - current_pixel_int.first) * pitch.first,
289  (current_pixel.second - current_pixel_int.second) * pitch.second);
290 
291  // Changing the reference frame to the proxy pixel cell
292  LocalPoint position_at_pc(relative_position_at_pc.first - center_proxy_cell.x(),
293  relative_position_at_pc.second - center_proxy_cell.y(),
294  super_charge.z() - center_proxy_cell.z());
295 
296  LogDebug("Pixel3DDigitizerAlgorithm::drift")
297  << "(super-)Charge\nlocal position: (" << super_charge.x() * 1.0_um_inv << ", " << super_charge.y() * 1.0_um_inv
298  << ", " << super_charge.z() * 1.0_um_inv << ") [um]"
299  << "\nMeasurement Point (row,column) (" << current_pixel.first << ", " << current_pixel.second << ")"
300  << "\nProxy pixel-cell frame (centered at left-back corner): (" << relative_position_at_pc.first * 1.0_um_inv
301  << ", " << relative_position_at_pc.second * 1.0_um_inv << ") [um]"
302  << "\nProxy pixel-cell frame (centered at n-column): (" << position_at_pc.x() * 1.0_um_inv << ", "
303  << position_at_pc.y() * 1.0_um_inv << ") [um] "
304  << "\nNe=" << super_charge.energy() << " electrons";
305 
306  // Check if the point is inside any of the column --> no charge was actually created then
307  if (is_inside_n_column_(position_at_pc, thickness) || is_inside_ohmic_column_(position_at_pc, half_pitch)) {
308  LogDebug("Pixel3DDigitizerAlgorithm::drift") << "Remove charge, inside the n-column or p-column!!";
309  continue;
310  }
311 
312  float nelectrons = super_charge.energy();
313  // XXX -- Diffusion: using the center frame
314  if (diffusion_activated) {
315  auto migrated_charges = diffusion(position_at_pc, super_charge.energy(), drift_direction, half_pitch, thickness);
316  for (auto& mc : migrated_charges) {
317  // Remove the migrated charges
318  nelectrons -= mc.energy();
319  // and convert back to the pixel ref. system
320  // Low-left origin/pitch -> relative within the pixel (a)
321  // Adding the pixel
322  const float pixel_x = current_pixel_int.first + (mc.x() + center_proxy_cell.x()) / pitch.first;
323  const float pixel_y = current_pixel_int.second + (mc.y() + center_proxy_cell.y()) / pitch.second;
324  const auto lp = pixdet->specificTopology().localPosition(MeasurementPoint(pixel_x, pixel_y));
325  // Remember: the drift function will move the reference system to the top. We need to subtract
326  // (center_proxy_cell.z() is a constant negative value) what we previously added in order to
327  // avoid a double translation when calling the drift function below the drift function
328  // initially considers the reference system centered in the module at half thickness)
329  mc.migrate_position(LocalPoint(lp.x(), lp.y(), mc.z() + center_proxy_cell.z()));
330  }
331  if (!migrated_charges.empty()) {
332  LogDebug("Pixel3DDigitizerAlgorithm::drift") << "****************"
333  << "MIGRATING (super-)charges"
334  << "****************";
335  // Drift this charges on the other pixel
336  auto mig_colpoints = drift(hit, pixdet, bfield, migrated_charges, false);
337  collection_points.insert(std::end(collection_points), mig_colpoints.begin(), mig_colpoints.end());
338  LogDebug("Pixel3DDigitizerAlgorithm::drift") << "*****************"
339  << "DOME MIGRATION"
340  << "****************";
341  }
342  }
343 
344  // Perform the drift, and check a potential lost of carriers because
345  // they reach the pasivation region (-z < thickness/2)
346  // XXX: not doing nothing, the carriers reach the electrode surface,
347  const float drift_distance = position_at_pc.perp() - np_column_radius_;
348 
349  // Insert a charge loss due to Rad Damage here
350  // XXX: ??
351  float energyOnCollector = nelectrons;
352  // FIXME: is this correct? Not for 3D...
353 
354  if (pseudoRadDamage_ >= RAD_DAMAGE) {
355  const float module_radius = pixdet->surface().position().perp();
356  if (module_radius <= pseudoRadDamageRadius_) {
357  const float kValue = pseudoRadDamage_ / (module_radius * module_radius);
358  energyOnCollector = energyOnCollector * std::exp(-1.0 * kValue * drift_distance / max_radial_distance);
359  }
360  }
361  LogDebug("Pixel3DDigitizerAlgorithm::drift")
362  << "Drift distance = " << drift_distance * 1.0_um_inv << " [um], "
363  << "Initial electrons = " << super_charge.energy()
364  << " [electrons], Electrons after loss/diff= " << energyOnCollector << " [electrons] ";
365  // Load the Charge distribution parameters
366  // XXX -- probably makes no sense the SignalPoint anymore...
367  collection_points.push_back(DigitizerUtility::SignalPoint(
368  current_pixel_int.first, current_pixel_int.second, 0.0, 0.0, hit.tof(), energyOnCollector));
369  }
370 
371  return collection_points;
372 }
373 
374 // ====================================================================
375 // Signal is already "induced" (actually electrons transported to the
376 // n-column) at the electrode. Just collecting and adding-up all pixel
377 // signal and linking it to the simulated energy deposit (hit)
379  const size_t hitIndex,
380  const uint32_t tofBin,
381  const Phase2TrackerGeomDetUnit* pixdet,
382  const std::vector<DigitizerUtility::SignalPoint>& collection_points) {
383  // X - Rows, Left-Right
384  // Y - Columns, Down-Up
385  const uint32_t detId = pixdet->geographicalId().rawId();
386  // Accumulated signal at each channel of this detector
387  signal_map_type& the_signal = _signal[detId];
388 
389  // Choose the proper pixel-to-channel converter
390  std::function<int(int, int)> pixelToChannel =
392  : static_cast<std::function<int(int, int)> >(Phase2TrackerDigi::pixelToChannel);
393 
394  // Iterate over collection points on the collection plane
395  for (const auto& pt : collection_points) {
396  // Extract corresponding channel (position is already given in pixel indices)
397  const int channel = pixelToChannel(pt.position().x(), pt.position().y());
398 
399  float corr_time = hit.tof() - pixdet->surface().toGlobal(hit.localPosition()).mag() * c_inv;
400  if (makeDigiSimLinks_) {
401  the_signal[channel] +=
402  DigitizerUtility::Amplitude(pt.amplitude(), &hit, pt.amplitude(), corr_time, hitIndex, tofBin);
403  } else {
404  the_signal[channel] += DigitizerUtility::Amplitude(pt.amplitude(), nullptr, pt.amplitude());
405  }
406 
407  LogDebug("Pixel3DDigitizerAlgorithm::induce_signal")
408  << " Induce charge at row,col:" << pt.position() << " N_electrons:" << pt.amplitude() << " [Channel:" << channel
409  << "]\n [Accumulated signal in this channel:" << the_signal[channel].ampl() << "] "
410  << " Global index linked PSimHit:" << hitIndex;
411  }
412 }
Vector3DBase
Definition: Vector3DBase.h:8
Pixel3DDigitizerAlgorithm::select_hit
bool select_hit(const PSimHit &hit, double tCorr, double &sigScale) const override
Definition: Pixel3DDigitizerAlgorithm.cc:89
Phase2TrackerDigitizerAlgorithm::signal_map_type
std::map< int, DigitizerUtility::Amplitude, std::less< int > > signal_map_type
Definition: Phase2TrackerDigitizerAlgorithm.h:103
Phase2TrackerDigitizerAlgorithm::makeDigiSimLinks_
const bool makeDigiSimLinks_
Definition: Phase2TrackerDigitizerAlgorithm.h:111
Pixel3DDigitizerAlgorithm::is_inside_n_column_
const bool is_inside_n_column_(const LocalPoint &p, const float &sensor_thickness) const
Definition: Pixel3DDigitizerAlgorithm.cc:94
DDAxes::y
Phase2TrackerDigitizerAlgorithm::theThresholdInE_Barrel_
const float theThresholdInE_Barrel_
Definition: Phase2TrackerDigitizerAlgorithm.h:144
mps_fire.i
i
Definition: mps_fire.py:428
Phase2TrackerDigitizerAlgorithm::pseudoRadDamage_
const double pseudoRadDamage_
Definition: Phase2TrackerDigitizerAlgorithm.h:168
MessageLogger.h
Phase2TrackerDigitizerAlgorithm::theTofUpperCut_
const float theTofUpperCut_
Definition: Phase2TrackerDigitizerAlgorithm.h:153
Pixel3DDigitizerAlgorithm::ohm_column_radius_
const float ohm_column_radius_
Definition: Pixel3DDigitizerAlgorithm.h:66
PV3DBase::x
T x() const
Definition: PV3DBase.h:59
CaloTowersParam_cfi.mc
mc
Definition: CaloTowersParam_cfi.py:8
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
Pixel3DDigitizerAlgorithm.h
edm
HLT enums.
Definition: AlignableModifier.h:19
pos
Definition: PixelAliasList.h:18
PixelTopology::pitch
virtual std::pair< float, float > pitch() const =0
Pixel3DDigitizerAlgorithm::siPixelLorentzAngleToken_
edm::ESGetToken< SiPixelLorentzAngle, SiPixelLorentzAngleSimRcd > siPixelLorentzAngleToken_
Definition: Pixel3DDigitizerAlgorithm.h:71
Pixel3DDigitizerAlgorithm::~Pixel3DDigitizerAlgorithm
~Pixel3DDigitizerAlgorithm() override
Definition: Pixel3DDigitizerAlgorithm.cc:84
protons_cff.time
time
Definition: protons_cff.py:35
Phase2TrackerDigitizerAlgorithm::use_LorentzAngle_DB_
const bool use_LorentzAngle_DB_
Definition: Phase2TrackerDigitizerAlgorithm.h:116
Topology::localPosition
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
DDAxes::x
edm::LogInfo
Log< level::Info, false > LogInfo
Definition: MessageLogger.h:125
align::LocalPoint
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
SiPixelGainCalibrationOfflineSimService.h
B2GMonitoring_cff.nelectrons
nelectrons
Definition: B2GMonitoring_cff.py:145
edm::ConsumesCollector::esConsumes
auto esConsumes()
Definition: ConsumesCollector.h:97
Phase2TrackerDigi::pixelToChannel
static PackedDigiType pixelToChannel(unsigned int row, unsigned int col)
Definition: Phase2TrackerDigi.h:43
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Pixel3DDigitizerAlgorithm::drift
std::vector< DigitizerUtility::SignalPoint > drift(const PSimHit &hit, const Phase2TrackerGeomDetUnit *pixdet, const GlobalVector &bfield, const std::vector< DigitizerUtility::EnergyDepositUnit > &ionization_points) const override
Definition: Pixel3DDigitizerAlgorithm.cc:214
PV3DBase::z
T z() const
Definition: PV3DBase.h:61
DigitizerUtility::EnergyDepositUnit
Definition: DigitizerUtility.h:72
Phase2TrackerDigitizerAlgorithm::theTofLowerCut_
const float theTofLowerCut_
Definition: Phase2TrackerDigitizerAlgorithm.h:152
GeomDet::surface
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Phase2TrackerDigitizerAlgorithm::theAdcFullScale_
const int theAdcFullScale_
Definition: Phase2TrackerDigitizerAlgorithm.h:138
Calorimetry_cff.thickness
thickness
Definition: Calorimetry_cff.py:115
sipixelobjects
Definition: CablingPathToDetUnit.h:4
Phase2TrackerDigitizerAlgorithm::pseudoRadDamageRadius_
const double pseudoRadDamageRadius_
Definition: Phase2TrackerDigitizerAlgorithm.h:169
PixelGeomDetUnit
Definition: PixelGeomDetUnit.h:15
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
Surface::bounds
const Bounds & bounds() const
Definition: Surface.h:87
mps_fire.end
end
Definition: mps_fire.py:242
Surface::toGlobal
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
PixelTopology::ncolumns
virtual int ncolumns() const =0
Pixel3DDigitizerAlgorithm::diffusion
std::vector< DigitizerUtility::EnergyDepositUnit > diffusion(const LocalPoint &pos, const float &ncarriers, const std::function< LocalVector(float, float)> &u_drift, const std::pair< float, float > pitches, const float &thickness) const
Definition: Pixel3DDigitizerAlgorithm.cc:117
theta
Geom::Theta< T > theta() const
Definition: Basic3DVectorLD.h:150
HcalDetIdTransform::transform
unsigned transform(const HcalDetId &id, unsigned transformCode)
Definition: HcalDetIdTransform.cc:7
PixelDigi::pixelToChannel
static int pixelToChannel(int row, int col)
Definition: PixelDigi.h:75
Point3DBase< float, LocalTag >
createTree.dd
string dd
Definition: createTree.py:154
pixelgpudetails::pixelToChannel
constexpr uint32_t pixelToChannel(int row, int col)
Definition: SiPixelRawToClusterGPUKernel.h:113
MeasurementPoint
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
Definition: MeasurementPoint.h:12
GeomDet::geographicalId
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
Bounds::thickness
virtual float thickness() const =0
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:233
edm::ParameterSet
Definition: ParameterSet.h:47
AlCaHLTBitMon_ParallelJobs.p
def p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
ParameterSet
Definition: Functions.h:16
Phase2TrackerDigitizerAlgorithm
Definition: Phase2TrackerDigitizerAlgorithm.h:54
Phase2TrackerDigitizerAlgorithm::pixelFlag_
bool pixelFlag_
Definition: Phase2TrackerDigitizerAlgorithm.h:235
Phase2TrackerDigitizerAlgorithm::tMax_
const double tMax_
Definition: Phase2TrackerDigitizerAlgorithm.h:176
PixelGeomDetUnit::specificTopology
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
Definition: PixelGeomDetUnit.cc:17
Pixel3DDigitizerAlgorithm::induce_signal
void induce_signal(const PSimHit &hit, const size_t hitIndex, const uint32_t tofBin, const Phase2TrackerGeomDetUnit *pixdet, const std::vector< DigitizerUtility::SignalPoint > &collection_points) override
Definition: Pixel3DDigitizerAlgorithm.cc:378
PV3DBase::y
T y() const
Definition: PV3DBase.h:60
Pixel3DDigitizerAlgorithm::init
void init(const edm::EventSetup &es) override
Definition: Pixel3DDigitizerAlgorithm.cc:27
Pixel3DDigitizerAlgorithm::is_inside_ohmic_column_
const bool is_inside_ohmic_column_(const LocalPoint &p, const std::pair< float, float > &pitch) const
Definition: Pixel3DDigitizerAlgorithm.cc:99
LocalVector
Local3DVector LocalVector
Definition: LocalVector.h:12
createfilelist.int
int
Definition: createfilelist.py:10
Phase2TrackerDigitizerAlgorithm::theElectronPerADC_
const float theElectronPerADC_
Definition: Phase2TrackerDigitizerAlgorithm.h:137
edm::EventSetup
Definition: EventSetup.h:58
GeomDet::specificSurface
const Plane & specificSurface() const
Same as surface(), kept for backward compatibility.
Definition: GeomDet.h:40
Phase2TrackerDigitizerAlgorithm::theThresholdInE_Endcap_
const float theThresholdInE_Endcap_
Definition: Phase2TrackerDigitizerAlgorithm.h:143
edm::EventSetup::getData
bool getData(T &iHolder) const
Definition: EventSetup.h:127
Phase2TrackerDigitizerAlgorithm::use_deadmodule_DB_
const bool use_deadmodule_DB_
Definition: Phase2TrackerDigitizerAlgorithm.h:115
Phase2TrackerDigitizerAlgorithm::_signal
signalMaps _signal
Definition: Phase2TrackerDigitizerAlgorithm.h:109
GloballyPositioned::position
const PositionType & position() const
Definition: GloballyPositioned.h:36
mag
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Definition: Basic3DVectorLD.h:127
LaserClient_cfi.Amplitude
Amplitude
Definition: LaserClient_cfi.py:39
DetId::rawId
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
c_inv
constexpr double c_inv
Definition: Phase2TrackerDigitizerAlgorithm.h:52
DigitizerUtility::SignalPoint
Definition: DigitizerUtility.h:93
PixelTopology::pixel
virtual std::pair< float, float > pixel(const LocalPoint &p) const =0
PixelGeomDetUnit.h
HiBiasedCentrality_cfi.function
function
Definition: HiBiasedCentrality_cfi.py:4
Pixel3DDigitizerAlgorithm::siPixelBadModuleToken_
edm::ESGetToken< SiPixelQuality, SiPixelQualityRcd > siPixelBadModuleToken_
Definition: Pixel3DDigitizerAlgorithm.h:70
ConsumesCollector.h
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterSet.h
PSimHit
Definition: PSimHit.h:15
Phase2TrackerDigitizerAlgorithm::addPixelInefficiency_
const bool addPixelInefficiency_
Definition: Phase2TrackerDigitizerAlgorithm.h:163
JetChargeProducer_cfi.exp
exp
Definition: JetChargeProducer_cfi.py:6
Pixel3DDigitizerAlgorithm::Pixel3DDigitizerAlgorithm
Pixel3DDigitizerAlgorithm(const edm::ParameterSet &conf, edm::ConsumesCollector iC)
Definition: Pixel3DDigitizerAlgorithm.cc:50
DeDxTools::esConsumes
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
PixelTopology::nrows
virtual int nrows() const =0
hcaltb::N_SIGMA
static const double N_SIGMA
Definition: HcalTBTDCUnpacker.cc:263
PV3DBase::perp
T perp() const
Definition: PV3DBase.h:69
Pixel3DDigitizerAlgorithm::np_column_radius_
const float np_column_radius_
Definition: Pixel3DDigitizerAlgorithm.h:65
edm::ConsumesCollector
Definition: ConsumesCollector.h:45
hit
Definition: SiStripHitEffFromCalibTree.cc:88
vertexPlots.e4
e4
Definition: vertexPlots.py:64
Pixel3DDigitizerAlgorithm::np_column_gap_
const float np_column_gap_
Definition: Pixel3DDigitizerAlgorithm.h:68
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37