CMS 3D CMS Logo

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