19 using namespace sipixelobjects;
23 constexpr
double operator""_um(
long double length) {
return length * 1
e-4; }
24 constexpr
double operator""_um_inv(
long double length) {
return length * 1e4; }
31 if (use_ineff_from_db_) {
33 theSiPixelGainCalibrationService_->setESObjects(es);
36 if (use_deadmodule_DB_) {
37 siPixelBadModule_ = &es.
getData(siPixelBadModuleToken_);
40 if (use_LorentzAngle_DB_) {
42 siPixelLorentzAngle_ = &es.
getData(siPixelLorentzAngleToken_);
46 fedCablingMap_ = &es.
getData(fedCablingMapToken_);
47 geom_ = &es.
getData(geomToken_);
52 conf.getParameter<edm::
ParameterSet>(
"Pixel3DDigitizerAlgorithm"),
55 (conf.getParameter<edm::
ParameterSet>(
"Pixel3DDigitizerAlgorithm").getParameter<double>(
"NPColumnRadius")) *
58 (conf.getParameter<edm::
ParameterSet>(
"Pixel3DDigitizerAlgorithm").getParameter<double>(
"OhmicColumnRadius")) *
61 (conf.getParameter<edm::
ParameterSet>(
"Pixel3DDigitizerAlgorithm").getParameter<double>(
"NPColumnGap")) *
74 <<
"Algorithm constructed \n"
75 <<
"Configuration parameters:\n"
90 double time = hit.
tof() - tCorr;
100 const std::pair<float, float>& half_pitch)
const {
119 const float& ncarriers,
121 const std::pair<float, float> hpitches,
127 const float max_migration_radius = 0.4_um;
134 if (hpitches.first -
std::abs(pos.
x()) < max_migration_radius) {
136 pitch = hpitches.first;
137 }
else if (hpitches.second -
std::abs(pos.
y()) < max_migration_radius) {
139 pitch = hpitches.second;
142 return std::vector<DigitizerUtility::EnergyDepositUnit>();
147 std::vector<DigitizerUtility::EnergyDepositUnit> migrated_charge;
150 const float diffusion_step = 0.1_um;
153 std::vector<float> pos_moving({pos.
x(), pos.
y(), pos.
z()});
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()});
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]";
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;
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]);
181 float migrated_e = current_carriers * 0.5 * (1.0 - std::erf(distance_edge / (sigma *
std::sqrt(2.0))));
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;
189 current_carriers -= migrated_e;
192 if (std::abs(distance_edge) >= max_migration_radius || current_carriers <= 0.5 * ncarriers) {
200 std::vector<float> newpos(pos_moving);
202 newpos[displ_ind] += std::copysign(N_SIGMA * sigma, newpos[displ_ind]);
205 return migrated_charge;
218 const std::vector<DigitizerUtility::EnergyDepositUnit>& ionization_points)
const {
219 return drift(hit, pixdet, bfield, ionization_points,
true);
225 const std::vector<DigitizerUtility::EnergyDepositUnit>& ionization_points,
226 bool diffusion_activated)
const {
236 const auto half_pitch = std::make_pair<float, float>(pitch.first * 0.5, pitch.second * 0.5);
240 const float pix_rounding = 0.99;
243 const float max_radial_distance =
244 std::sqrt(half_pitch.first * half_pitch.first + half_pitch.second * half_pitch.second);
251 LocalPoint center_proxy_cell(half_pitch.first, half_pitch.second, -0.5 * thickness);
253 LogDebug(
"Pixel3DDigitizerAlgorithm::drift")
254 <<
"Pixel pitch:" << pitch.first * 1.0_um_inv <<
", " << pitch.second * 1.0_um_inv <<
" [um]";
258 std::function<LocalVector(float, float)> drift_direction = [](
float x,
float y) ->
LocalVector {
259 const float theta = std::atan2(
y, x);
263 std::vector<DigitizerUtility::SignalPoint> collection_points;
265 collection_points.reserve(ionization_points.size());
269 const float RAD_DAMAGE = 0.001;
271 for (
const auto& super_charge : ionization_points) {
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);
283 const auto current_pixel_int = std::make_pair(std::floor(current_pixel.first), std::floor(current_pixel.second));
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);
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());
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";
308 LogDebug(
"Pixel3DDigitizerAlgorithm::drift") <<
"Remove charge, inside the n-column or p-column!!";
312 float nelectrons = super_charge.energy();
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) {
318 nelectrons -= mc.energy();
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;
329 mc.migrate_position(
LocalPoint(lp.x(), lp.y(), mc.z() + center_proxy_cell.
z()));
331 if (!migrated_charges.empty()) {
332 LogDebug(
"Pixel3DDigitizerAlgorithm::drift") <<
"****************"
333 <<
"MIGRATING (super-)charges"
334 <<
"****************";
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") <<
"*****************"
340 <<
"****************";
351 float energyOnCollector = nelectrons;
358 energyOnCollector = energyOnCollector *
std::exp(-1.0 * kValue * drift_distance / max_radial_distance);
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] ";
368 current_pixel_int.first, current_pixel_int.second, 0.0, 0.0, hit.
tof(), energyOnCollector));
371 return collection_points;
379 const size_t hitIndex,
380 const uint32_t tofBin,
382 const std::vector<DigitizerUtility::SignalPoint>& collection_points) {
395 for (
const auto&
pt : collection_points) {
401 the_signal[channel] +=
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;
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Pixel3DDigitizerAlgorithm(const edm::ParameterSet &conf, edm::ConsumesCollector iC)
Local3DVector LocalVector
float tof() const
deprecated name for timeOfFlight()
const float theElectronPerADC_
Point3DBase< Scalar, LocalTag > LocalPoint
virtual std::pair< float, float > pixel(const LocalPoint &p) const =0
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
virtual int ncolumns() const =0
const double pseudoRadDamage_
const double pseudoRadDamageRadius_
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Sin< T >::type sin(const T &t)
void init(const edm::EventSetup &es) override
constexpr uint32_t rawId() const
get the raw id
Geom::Theta< T > theta() const
virtual int nrows() const =0
const Bounds & bounds() const
const float theThresholdInE_Endcap_
Exp< T >::type exp(const T &t)
const Plane & surface() const
The nominal surface of the GeomDet.
edm::ESGetToken< SiPixelLorentzAngle, SiPixelLorentzAngleSimRcd > siPixelLorentzAngleToken_
const bool use_LorentzAngle_DB_
static int pixelToChannel(int row, int col)
virtual float thickness() const =0
bool getData(T &iHolder) const
const bool use_deadmodule_DB_
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
Local3DPoint localPosition() const
const bool is_inside_ohmic_column_(const LocalPoint &p, const std::pair< float, float > &pitch) const
static PackedDigiType pixelToChannel(unsigned int row, unsigned int col)
const bool addPixelInefficiency_
Cos< T >::type cos(const T &t)
Abs< T >::type abs(const T &t)
const float theThresholdInE_Barrel_
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
DetId geographicalId() const
The label of this GeomDet.
std::map< int, DigitizerUtility::Amplitude, std::less< int > > signal_map_type
bool select_hit(const PSimHit &hit, double tCorr, double &sigScale) const override
const float theTofLowerCut_
Log< level::Info, false > LogInfo
const float ohm_column_radius_
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
edm::ESGetToken< SiPixelQuality, SiPixelQualityRcd > siPixelBadModuleToken_
static const double N_SIGMA
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
const int theAdcFullScale_
~Pixel3DDigitizerAlgorithm() override
virtual std::pair< float, float > pitch() const =0
const bool makeDigiSimLinks_
const float np_column_radius_
constexpr uint32_t pixelToChannel(int row, int col)
const PositionType & position() const
const float np_column_gap_
const bool is_inside_n_column_(const LocalPoint &p, const float &sensor_thickness) const
const float theTofUpperCut_
const Plane & specificSurface() const
Same as surface(), kept for backward compatibility.
std::vector< DigitizerUtility::SignalPoint > drift(const PSimHit &hit, const Phase2TrackerGeomDetUnit *pixdet, const GlobalVector &bfield, const std::vector< DigitizerUtility::EnergyDepositUnit > &ionization_points) const override