LDMX Software
ldmx::HgcrocEmulator Class Reference

Emulate the digitization procedure performed by the HGCROC. More...

#include <HgcrocEmulator.h>

Public Member Functions

 HgcrocEmulator (const framework::config::Parameters &ps)
 Constructor.
 
 ~HgcrocEmulator ()
 Destructor.
 
bool hasSeed () const
 Check if emulator has been seeded.
 
void seedGenerator (uint64_t seed)
 Seed the emulator for random number generation.
 
void condition (const conditions::DoubleTableCondition &table)
 Set Conditions.
 
bool digitize (const int &channelID, std::vector< std::pair< double, double > > &arriving_pulses, std::vector< ldmx::HgcrocDigiCollection::Sample > &digiToAdd) const
 Digitize the signals from the simulated hits_.
 
std::vector< ldmx::HgcrocDigiCollection::SamplenoiseDigi (const int &channel, const double &soi_amplitude=0) const
 Generate a digi of pure noise.
 
double noise (const int &channelID) const
 Get random noise amplitdue for input channel [mV].
 
double gain (const int &channelID) const
 Gain for input channel.
 
double pedestal (const int &id) const
 Pedestal [ADC Counts] for input channel.
 
double readoutThreshold (const int &id) const
 Readout Threshold (ADC Counts)
 

Public Attributes

bool save_pulse_truth_info_ {false}
 
ldmx::HgcrocPulseTruthCollection * pulse_truth_coll_
 

Private Member Functions

double getCondition (int id, const std::string &name) const
 Get condition for input chip ID, condition name, and default value.
 

Private Attributes

bool noise_ {true}
 Put noise in channels, only configure to false if testing.
 
int n_ad_cs_
 Depth of ADC buffer.
 
int i_soi_
 Index for the Sample Of Interest in the list of digi samples.
 
double clock_cycle_
 Time interval for chip clock [ns].
 
double timing_jitter_
 Jitter of timing mechanism in the chip [ns].
 
double ns_
 Conversion from time [ns] to counts.
 
double rate_up_slope_
 Rate of Up Slope in Pulse Shape [1/ns].
 
double time_up_slope_
 Time of Up Slope relative to Pulse Shape Fit [ns].
 
double rate_dn_slope_
 Rate of Down Slope in Pulse Shape [1/ns].
 
double time_dn_slope_
 Time of Down Slope relative to Pulse Shape Fit [ns].
 
double time_peak_
 Time of Peak relative to pulse shape fit [ns].
 
double hit_merge_ns_
 Hit merging time [ns].
 
const conditions::DoubleTableConditionchip_conditions_ {nullptr}
 Handle to table of chip-dependent conditions.
 
std::map< std::string, int > condition_names_to_index_
 Map of condition names to column numbers.
 
std::unique_ptr< TRandom3 > noise_injector_
 Generates Gaussian noise on top of real hits_.
 
TF1 pulse_func_
 Functional shape of signal pulse in time.
 

Detailed Description

Emulate the digitization procedure performed by the HGCROC.

This object emulates how the chip converts the analog signal into DIGI samples. With that in mind, the digitize method converts a set of voltages and times into the DIGI.

This object does not do everything related to subsystem information. It does not simulate noise within the empty channels, and it does not convert simulated energy depositions into voltages. These tasks depend on the detector construction, so they are left to the individual subsystem producers.

@TODO time phase setting relative to target t=0ns using electronic IDs

@TODO accurately model recovering from saturation (TOT Mode). Currently, we just have all the samples after the sample triggering TOT mode rail.

Definition at line 42 of file HgcrocEmulator.h.

Constructor & Destructor Documentation

◆ HgcrocEmulator()

ldmx::HgcrocEmulator::HgcrocEmulator ( const framework::config::Parameters & ps)

Constructor.

Configures the chip emulator using the passed parameters.

Definition at line 6 of file HgcrocEmulator.cxx.

6 {
7 // settings of readout chip that are the same for all chips
8 // used in actual digitization
9 noise_ = ps.get<bool>("noise");
10 timing_jitter_ = ps.get<double>("timingJitter");
11 rate_up_slope_ = ps.get<double>("rateUpSlope");
12 time_up_slope_ = ps.get<double>("timeUpSlope");
13 rate_dn_slope_ = ps.get<double>("rateDnSlope");
14 time_dn_slope_ = ps.get<double>("timeDnSlope");
15 time_peak_ = ps.get<double>("timePeak");
16 clock_cycle_ = ps.get<double>("clockCycle");
17 n_ad_cs_ = ps.get<int>("nADCs");
18 i_soi_ = ps.get<int>("iSOI");
19
20 // Time -> clock counts conversion
21 // time [ns] * ( 2^10 / max time in ns ) = clock counts
22 ns_ = 1024. / clock_cycle_;
23
24 hit_merge_ns_ = 0.05; // combine at 50 ps level
25
26 // Configure the pulse shape function
28 TF1("pulseFunc",
29 "[0]*((1.0+exp([1]*(-[2]+[3])))*(1.0+exp([5]*(-[6]+[3]))))/"
30 "((1.0+exp([1]*(x-[2]+[3]-[4])))*(1.0+exp([5]*(x-[6]+[3]-[4]))))",
31 0.0, (double)n_ad_cs_ * clock_cycle_);
32 pulse_func_.FixParameter(0, 1.0); // amplitude is set externally
33 pulse_func_.FixParameter(1, rate_up_slope_);
34 pulse_func_.FixParameter(2, time_up_slope_);
35 pulse_func_.FixParameter(3, time_peak_);
36 pulse_func_.FixParameter(4, 0); // not using time offset in this way
37 pulse_func_.FixParameter(5, rate_dn_slope_);
38 pulse_func_.FixParameter(6, time_dn_slope_);
39}
const T & get(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:78
double clock_cycle_
Time interval for chip clock [ns].
int n_ad_cs_
Depth of ADC buffer.
double timing_jitter_
Jitter of timing mechanism in the chip [ns].
double time_up_slope_
Time of Up Slope relative to Pulse Shape Fit [ns].
double hit_merge_ns_
Hit merging time [ns].
double rate_dn_slope_
Rate of Down Slope in Pulse Shape [1/ns].
int i_soi_
Index for the Sample Of Interest in the list of digi samples.
double time_dn_slope_
Time of Down Slope relative to Pulse Shape Fit [ns].
TF1 pulse_func_
Functional shape of signal pulse in time.
bool noise_
Put noise in channels, only configure to false if testing.
double time_peak_
Time of Peak relative to pulse shape fit [ns].
double ns_
Conversion from time [ns] to counts.
double rate_up_slope_
Rate of Up Slope in Pulse Shape [1/ns].

References clock_cycle_, framework::config::Parameters::get(), hit_merge_ns_, i_soi_, n_ad_cs_, noise_, ns_, pulse_func_, rate_dn_slope_, rate_up_slope_, time_dn_slope_, time_peak_, time_up_slope_, and timing_jitter_.

◆ ~HgcrocEmulator()

ldmx::HgcrocEmulator::~HgcrocEmulator ( )
inline

Destructor.

Definition at line 52 of file HgcrocEmulator.h.

52{}

Member Function Documentation

◆ condition()

void ldmx::HgcrocEmulator::condition ( const conditions::DoubleTableCondition & table)
inline

Set Conditions.

Passes the chips conditions to be cached here and used later in digitization.

Parameters
tableconditions::DoubleTableConditions to be used for chip parameters

Definition at line 75 of file HgcrocEmulator.h.

75 {
76 // reset cache of column numbers if table changes
77 if (&table != chip_conditions_) condition_names_to_index_.clear();
78 chip_conditions_ = &table;
79 }
std::map< std::string, int > condition_names_to_index_
Map of condition names to column numbers.
const conditions::DoubleTableCondition * chip_conditions_
Handle to table of chip-dependent conditions.

References chip_conditions_, and condition_names_to_index_.

◆ digitize()

bool ldmx::HgcrocEmulator::digitize ( const int & channelID,
std::vector< std::pair< double, double > > & arriving_pulses,
std::vector< ldmx::HgcrocDigiCollection::Sample > & digiToAdd ) const

Digitize the signals from the simulated hits_.

This is where the hefty amount of work is done.

  1. Prepare for Emulation
    • Clear input list of digi samples
    • Get conditions for the current chip
    • Sort the input sim voltage hits by amplitude
  2. Combine input simulated hits into one CompositePulse to digitize.
    • This composite pulse decides whether to merge two simulated hits_ into one larger pulse depending on how close they are in time.
  3. Add a timing jitter TODO
  4. Go through sampling baskets one-by-one
    • If enter pulse goes above tot threshold, then enter TOT readout mode, digitize, and return true.
    • If pulse never goes above tot threshold, then only return true if the voltage sample taken in the SOI is above the readout threshold.

ADC Mode

Here, we measure the height of the pulse once per clock cycle. This leaves us with nADCs_ samples for each digitized hit. The voltage measurements are converted to ADC counts using the parameter gain_.

The time of arrival (TOA) is zero unless the amplitude is greater than toaThreshold_. Then the TOA is set to the point the pulse crosses the toaThreshold_ with respect to the current clock window. The time measurements are converted to clock counts using 2^10=1024 and clock_cycle_.

Both the tot_complete_ and tot_progress_ flags are set to false for all the samples.

TOT Mode

Here, we measure how long the chip is in saturation. This is calculated using the drain rate and assuming a linear drain of the charge off of the chip.

The charge deposited on the chip is converted from the pulse triggering TOT without including the pedestal. The pedestal is included in the pulse when determining if the pulse should enter TOT mode.

Thus, TOT = charge deposited / drain rate and TOA is calculated as before, seeing when the pulse crossed the TOA threshold.

Pulse Measurement

All "measurements" of the pulse use the object CompositePulse. This function incorporates the pedestal_ and linearly adds all the pulses at a variety of times. The pulses at different times are "merged" upon addition to the composite pulse depending on how close they are. This is done in a single pass, so we might end up with pulses closer than the provided separation time.

Note
For more realism, some chip parameters should change depending on the chip they are coming from. This should be modified here, with a package of "HgcrocConditions" passed to this function to configure the emulator before digitizing.
Parameters
[in]channelIDraw integer ID for this readout channel
[in]arriving_pulsespairs of (voltage,time) of hits arriving at the chip
[out]digiToAdddigi that will be filled with the samples from the chip
Returns
true if digis were constructed (false if hit was below readout)

the time here is nominal (zero gives peak if hit.second is zero)

Definition at line 45 of file HgcrocEmulator.cxx.

48 {
49 // step 0: prepare ourselves for emulation
50 digiToAdd.clear(); // make sure it is clean
51
52 // Configure chip settings based off of table (that may have been passed)
53 double tot_max = getCondition(channelID, "TOT_MAX");
54 double pad_capacitance = getCondition(channelID, "PAD_CAPACITANCE");
55 double gain = this->gain(channelID);
56 double pedestal = this->pedestal(channelID);
57 double toa_threshold = getCondition(channelID, "TOA_THRESHOLD");
58 double tot_threshold = getCondition(channelID, "TOT_THRESHOLD");
59 // measTime defines the point in the BX where an in-time
60 // (time=0 in times vector) hit would arrive.
61 // Used to determine BX boundaries and TOA behavior.
62 double meas_time = getCondition(channelID, "MEAS_TIME");
63 double drain_rate = getCondition(channelID, "DRAIN_RATE");
64 double readout_threshold_float = this->readoutThreshold(channelID);
65 int readout_threshold = int(readout_threshold_float);
66
67 // sort by amplitude
68 // ==> makes sure that puleses are merged towards higher ones
69 std::sort(
70 arriving_pulses.begin(), arriving_pulses.end(),
71 [](const std::pair<double, double> &a,
72 const std::pair<double, double> &b) { return a.first > b.first; });
73
74 // step 1: gather voltages into groups separated by (programmable) ns, single
75 // pass
77
78 for (auto hit : arriving_pulses) pulse.addOrMerge(hit, hit_merge_ns_);
79
80 // TODO step 2: add timing jitter
81 // if (noise_) pulse.jitter();
82
84
85 // step 3: go through each BX sample one by one
86 bool was_toa = false;
87 for (int i_adc = 0; i_adc < n_ad_cs_; i_adc++) {
88 double start_bx = (i_adc - i_soi_) * clock_cycle_ - meas_time;
89 ldmx_log(trace) << " iADC = " << i_adc << " at startBX = " << start_bx;
90
91 // step 3b: check each merged hit to see if it peaks in this BX. If so,
92 // check its peak time to see if it's over TOT or TOA.
93 bool start_tot = false;
94 bool over_toa = false;
95 double tover_toa = -1;
96 double tover_tot = -1;
97 for (auto hit : pulse.hits()) {
98 int hit_bx = int((hit.second + meas_time) / clock_cycle_ + i_soi_);
99 // if this hit wasn't in the current BX, continue...
100 if (hit_bx != i_adc) {
101 continue;
102 }
103
104 double vpeak = pulse(hit.second);
105
106 if (vpeak > tot_threshold) {
107 start_tot = true;
108 // use the latest time in the window
109 if (tover_tot < hit.second) {
110 tover_tot = hit.second;
111 }
112 }
113
114 if (vpeak > toa_threshold) {
115 if (!over_toa || hit.second < tover_toa) tover_toa = hit.second;
116 over_toa = true;
117 }
118
119 } // loop over sim hits_
120
121 // check for the case of a TOA even though the peak is in the next BX
122 if (!over_toa && pulse(start_bx + clock_cycle_) > toa_threshold) {
123 if (pulse(start_bx) < toa_threshold) {
124 // pulse crossed TOA threshold somewhere between the start of this
125 // basket and the end
126 over_toa = true;
127 tover_toa = start_bx + clock_cycle_;
128 }
129 }
130
131 if (start_tot) {
132 // above TOT threshold -> do TOT readout mode
133
134 // @TODO NO NOISE
135 // CompositePulse includes pedestal, we need to remove it
136 // when calculating the charge deposited.
137 double charge_deposited =
138 (pulse(tover_tot) - gain * pedestal) * pad_capacitance;
139
140 // Measure Time Over Threshold (TOT) by using the drain rate.
141 // 1. Use drain rate to see how long it takes for the charge to drain off
142 // 2. Translate this into DIGI samples
143
144 // Assume linear drain with slope drain rate:
145 // y_-intercept = pulse amplitude
146 // slope = drain rate
147 // ==> x-intercept = amplitude / rate
148 // actual time over threshold using the real signal voltage amplitude
149 double tot = charge_deposited / drain_rate;
150 ldmx_log(trace) << " we are in TOT read-out mode, TOT = " << tot;
151
152 // calculate the TDC counts for this tot measurement
153 // internally, the chip uses 12 bits (2^12 = 4096)
154 // to measure a maximum of tot Max [ns]
155 int tdc_counts = int(tot * 4096 / tot_max) + pedestal;
156
157 // were we already over TOA? TOT is reported in BX where TOA went over
158 // threshold...
159 int toa{0};
160 if (was_toa) {
161 // TOA was in the past
162 toa = digiToAdd.back().toa();
163 } else {
164 // TOA is here and we need to find it
165 double timecross =
166 pulse.findCrossing(start_bx, tover_tot, toa_threshold);
167 toa = int((timecross - start_bx) * ns_);
168 // keep inside valid limits
169 if (toa == 0) toa = 1;
170 if (toa > 1023) toa = 1023;
171 }
172 ldmx_log(trace) << " Adding TOT hit with toa = " << toa
173 << ", tdc_counts = " << tdc_counts
174 << " adcT at prev iADC = "
175 << digiToAdd.at(i_adc - 1).adcT();
176 // ADC at t-1
177 auto adc_at_tminus1 =
178 (i_adc > 0) ? digiToAdd.at(i_adc - 1).adcT() : pedestal;
179 auto i_tot_sample = digiToAdd.size();
180 // mark as a TOT measurement with 2nd boolean as true
181 digiToAdd.emplace_back(false, true, adc_at_tminus1, tdc_counts, toa);
182
183 // TODO: properly handle saturation and recovery, eventually.
184 // Now just kill everything...
185 ldmx_log(trace)
186 << " Adding further hits_ with ADC [t-1] = 0x3FF, toa = "
187 "0x3FF, until digiToAdd.size() = "
188 << digiToAdd.size() << " < n_ad_cs_(" << n_ad_cs_ << ")";
189 while (digiToAdd.size() < n_ad_cs_) {
190 // flags to mark type of sample
191 digiToAdd.emplace_back(true, false, 0x3FF, 0x3FF, 0);
192 }
193 // Read out if the toa is within one Bx after nominal
194 return (i_tot_sample <= i_soi_ + 1);
195 } else {
196 // determine the voltage at the sampling time
197 double bxvolts = pulse((i_adc - i_soi_) * clock_cycle_);
198 // add noise if requested
199 if (noise_) bxvolts += noise(channelID);
200 // convert to integer and keep in range (handle low and high saturation)
201 int adc = bxvolts / gain;
202 ldmx_log(trace) << " we are in ADC read-out mode, adc = " << adc;
203 if (adc < 0) adc = 0;
204 if (adc > 1023) adc = 1023;
205
206 // check for TOA
207 int toa(0);
208 if (pulse(start_bx) < toa_threshold && over_toa) {
209 double timecross =
210 pulse.findCrossing(start_bx, tover_toa, toa_threshold);
211 toa = int((timecross - start_bx) * ns_);
212 // keep inside valid limits
213 if (toa == 0) toa = 1;
214 if (toa > 1023) toa = 1023;
215 was_toa = true;
216 } else {
217 was_toa = false;
218 }
219 // ADC at t-1
220 auto adc_t_minus1 =
221 (i_adc > 0) ? digiToAdd.at(i_adc - 1).adcT() : pedestal;
222
223 digiToAdd.emplace_back(false, false, adc_t_minus1, adc, toa);
224 } // TOT or ADC Mode
225 } // sampling baskets
226
227 if (save_pulse_truth_info_)
228 pulse_truth_coll_->push_back(ldmx::HgcrocPulseTruth(channelID, pulse));
229
230 // we only get here if we never went into TOT mode
231 // check the SOI to see if we should read out
232 ldmx_log(trace) << " we are adding the hit IFF iSOI= " << i_soi_
233 << "'s adc_t = " << digiToAdd.at(i_soi_).adcT()
234 << " >= thresh (" << readout_threshold << ")";
235 return digiToAdd.at(i_soi_).adcT() >= readout_threshold;
236} // HgcrocEmulator::digitize
CompositePulse.
double readoutThreshold(const int &id) const
Readout Threshold (ADC Counts)
double noise(const int &channelID) const
Get random noise amplitdue for input channel [mV].
double pedestal(const int &id) const
Pedestal [ADC Counts] for input channel.
double getCondition(int id, const std::string &name) const
Get condition for input chip ID, condition name, and default value.
double gain(const int &channelID) const
Gain for input channel.

References ldmx::CompositePulse::addOrMerge(), clock_cycle_, ldmx::CompositePulse::findCrossing(), gain(), getCondition(), hit_merge_ns_, ldmx::CompositePulse::hits(), i_soi_, n_ad_cs_, noise(), noise_, ns_, pedestal(), pulse_func_, and readoutThreshold().

◆ gain()

double ldmx::HgcrocEmulator::gain ( const int & channelID) const
inline

Gain for input channel.

Definition at line 194 of file HgcrocEmulator.h.

194 {
195 return getCondition(channelID, "GAIN");
196 }

References getCondition().

Referenced by digitize(), noise(), and noiseDigi().

◆ getCondition()

double ldmx::HgcrocEmulator::getCondition ( int id,
const std::string & name ) const
inlineprivate

Get condition for input chip ID, condition name, and default value.

Parameters
[in]idchip global integer ID used in condition table
[in]namestd::string name of chip parameter in table
[in]defdefault value for parameter if not found in table (or table not set)
Returns
value of chip parameter

Definition at line 219 of file HgcrocEmulator.h.

219 {
220 // check if emulator has been passed a table of conditions
221 if (!chip_conditions_) {
222 EXCEPTION_RAISE("HgcrocCond",
223 "HGC ROC Emulator was not given a conditions table.");
224 }
225
226 // cache column index for the input name
227 if (condition_names_to_index_.count(name) == 0)
229
230 // get condition
231 return chip_conditions_->get(id, condition_names_to_index_.at(name));
232 }
unsigned int getColumnNumber(const std::string &colname) const
Get a the column number for the given column name.
T get(unsigned int id, unsigned int col) const
Get an entry by DetectorId and number.

References chip_conditions_, condition_names_to_index_, conditions::HomogenousTableCondition< T >::get(), and conditions::BaseTableCondition::getColumnNumber().

Referenced by digitize(), gain(), noise(), pedestal(), and readoutThreshold().

◆ hasSeed()

bool ldmx::HgcrocEmulator::hasSeed ( ) const
inline

Check if emulator has been seeded.

Returns
true if random generator has been seeded

Definition at line 58 of file HgcrocEmulator.h.

58{ return noise_injector_.get() != nullptr; }
std::unique_ptr< TRandom3 > noise_injector_
Generates Gaussian noise on top of real hits_.

References noise_injector_.

◆ noise()

double ldmx::HgcrocEmulator::noise ( const int & channelID) const
inline

Get random noise amplitdue for input channel [mV].

Parameters
[in]channelID
Returns
electronic noise amplitude [mV] above pedestal

Definition at line 188 of file HgcrocEmulator.h.

188 {
189 return noise_injector_->Gaus(
190 0, getCondition(channelID, "NOISE") * gain(channelID));
191 };

References gain(), getCondition(), and noise_injector_.

Referenced by digitize(), and noiseDigi().

◆ noiseDigi()

std::vector< ldmx::HgcrocDigiCollection::Sample > ldmx::HgcrocEmulator::noiseDigi ( const int & channel,
const double & soi_amplitude = 0 ) const

Generate a digi of pure noise.

The (optional) input soi amplitude is added on-top of pedestal in the SOI after being divided by the gain. This is designed to be close to the output amplitudes generated by the noise generator.

Note
We don't go into TOT mode. We simply fill a DIGI with ADC samples with noise on top of the pedestal, including the SOI amplitude if provided.

We set the TOA for each of the samples created to 0 (or "not found"), so this also indirectly assumes we are below the TOA threshold as well as below the TOT threshold.

Parameters
[in]channelraw integer ID for this readout channel
[in]soi_amplitudeamplitude of noise "pulse" in mV
Returns
DIGI of pure noise

Definition at line 238 of file HgcrocEmulator.cxx.

239 {
240 // get chip conditions from emulator
241 double pedestal{this->pedestal(channel)};
242 double gain{this->gain(channel)};
243 // fill a digi with noise samples
244 std::vector<ldmx::HgcrocDigiCollection::Sample> noise_digi;
245 for (int i_adc{0}; i_adc < n_ad_cs_; i_adc++) {
246 // gen noise for ADC samples
247 // ADC at t-1
248 int adc_tm1{static_cast<int>(pedestal)};
249 if (i_adc > 0) {
250 adc_tm1 = noise_digi.at(i_adc - 1).adcT();
251 } else {
252 adc_tm1 += noise(channel) / gain;
253 }
254 // ADC at t
255 int adc_t{static_cast<int>(pedestal + noise(channel) / gain)};
256
257 if (i_adc == i_soi_) adc_t += soi_amplitude / gain;
258
259 // set toa to 0 (not determined)
260 // put new sample into noise digi
261 noise_digi.emplace_back(false, false, adc_tm1, adc_t, 0);
262 } // samples in noise digi
263 return noise_digi;
264}

References gain(), i_soi_, n_ad_cs_, noise(), and pedestal().

◆ pedestal()

double ldmx::HgcrocEmulator::pedestal ( const int & id) const
inline

Pedestal [ADC Counts] for input channel.

Definition at line 199 of file HgcrocEmulator.h.

199{ return getCondition(id, "PEDESTAL"); }

References getCondition().

Referenced by digitize(), and noiseDigi().

◆ readoutThreshold()

double ldmx::HgcrocEmulator::readoutThreshold ( const int & id) const
inline

Readout Threshold (ADC Counts)

Definition at line 202 of file HgcrocEmulator.h.

202 {
203 return getCondition(id, "READOUT_THRESHOLD");
204 }

References getCondition().

Referenced by digitize().

◆ seedGenerator()

void ldmx::HgcrocEmulator::seedGenerator ( uint64_t seed)

Seed the emulator for random number generation.

Parameters
[in]seedinteger to use as random seed

Definition at line 41 of file HgcrocEmulator.cxx.

41 {
42 noise_injector_ = std::make_unique<TRandom3>(seed);
43}

References noise_injector_.

Member Data Documentation

◆ chip_conditions_

const conditions::DoubleTableCondition* ldmx::HgcrocEmulator::chip_conditions_ {nullptr}
private

Handle to table of chip-dependent conditions.

The defaults are listed below and are separate parameters passed through the python configuration.

Definition at line 285 of file HgcrocEmulator.h.

285{nullptr};

Referenced by condition(), and getCondition().

◆ clock_cycle_

double ldmx::HgcrocEmulator::clock_cycle_
private

Time interval for chip clock [ns].

Definition at line 249 of file HgcrocEmulator.h.

Referenced by digitize(), and HgcrocEmulator().

◆ condition_names_to_index_

std::map<std::string, int> ldmx::HgcrocEmulator::condition_names_to_index_
mutableprivate

Map of condition names to column numbers.

mutable so that we can update the cached column values in getCondition during processing.

Definition at line 293 of file HgcrocEmulator.h.

Referenced by condition(), and getCondition().

◆ hit_merge_ns_

double ldmx::HgcrocEmulator::hit_merge_ns_
private

Hit merging time [ns].

Definition at line 273 of file HgcrocEmulator.h.

Referenced by digitize(), and HgcrocEmulator().

◆ i_soi_

int ldmx::HgcrocEmulator::i_soi_
private

Index for the Sample Of Interest in the list of digi samples.

Definition at line 246 of file HgcrocEmulator.h.

Referenced by digitize(), HgcrocEmulator(), and noiseDigi().

◆ n_ad_cs_

int ldmx::HgcrocEmulator::n_ad_cs_
private

Depth of ADC buffer.

Definition at line 243 of file HgcrocEmulator.h.

Referenced by digitize(), HgcrocEmulator(), and noiseDigi().

◆ noise_

bool ldmx::HgcrocEmulator::noise_ {true}
private

Put noise in channels, only configure to false if testing.

Definition at line 240 of file HgcrocEmulator.h.

240{true};

Referenced by digitize(), and HgcrocEmulator().

◆ noise_injector_

std::unique_ptr<TRandom3> ldmx::HgcrocEmulator::noise_injector_
private

Generates Gaussian noise on top of real hits_.

Definition at line 300 of file HgcrocEmulator.h.

Referenced by hasSeed(), noise(), and seedGenerator().

◆ ns_

double ldmx::HgcrocEmulator::ns_
private

Conversion from time [ns] to counts.

Definition at line 255 of file HgcrocEmulator.h.

Referenced by digitize(), and HgcrocEmulator().

◆ pulse_func_

TF1 ldmx::HgcrocEmulator::pulse_func_
mutableprivate

Functional shape of signal pulse in time.

Shape parameters are hardcoded into the function currently. Pulse Shape: [0]*((1.0+exp([1]*(-[2]+[3])))*(1.0+exp([5]*(-[6]+[3]))))/((1.0+exp([1]*(x-[2]+[3]-[4])))*(1.0+exp([5]*(x-[6]+[3]-[4])))) p[0] = amplitude (height of peak in mV) p[1] = rate of up slope - rateUpSlope_ p[2] = time of up slope relative to shape fit - timeUpSlope_ p[3] = time of peak relative to shape fit - timePeak_ p[4] = peak time (related to time of hit [ns]) p[5] = rate of down slope - rateDnSlope_ p[6] = time of down slope relative to shape fit - timeDnSlope_

\[
 V(t) =
 p_0\frac{(1+\exp(p_1(-p_2+p_3)))(1+\exp(p_5*(-p_6+p_3)))}
         {(1+\exp(p_1(t-p_2+p_3-p_4)))(1+\exp(p_5*(t-p_6+p_3-p_4)))}
\]

Definition at line 322 of file HgcrocEmulator.h.

Referenced by digitize(), and HgcrocEmulator().

◆ pulse_truth_coll_

ldmx::HgcrocPulseTruthCollection* ldmx::HgcrocEmulator::pulse_truth_coll_

Definition at line 207 of file HgcrocEmulator.h.

◆ rate_dn_slope_

double ldmx::HgcrocEmulator::rate_dn_slope_
private

Rate of Down Slope in Pulse Shape [1/ns].

Definition at line 264 of file HgcrocEmulator.h.

Referenced by HgcrocEmulator().

◆ rate_up_slope_

double ldmx::HgcrocEmulator::rate_up_slope_
private

Rate of Up Slope in Pulse Shape [1/ns].

Definition at line 258 of file HgcrocEmulator.h.

Referenced by HgcrocEmulator().

◆ save_pulse_truth_info_

bool ldmx::HgcrocEmulator::save_pulse_truth_info_ {false}

Definition at line 206 of file HgcrocEmulator.h.

206{false};

◆ time_dn_slope_

double ldmx::HgcrocEmulator::time_dn_slope_
private

Time of Down Slope relative to Pulse Shape Fit [ns].

Definition at line 267 of file HgcrocEmulator.h.

Referenced by HgcrocEmulator().

◆ time_peak_

double ldmx::HgcrocEmulator::time_peak_
private

Time of Peak relative to pulse shape fit [ns].

Definition at line 270 of file HgcrocEmulator.h.

Referenced by HgcrocEmulator().

◆ time_up_slope_

double ldmx::HgcrocEmulator::time_up_slope_
private

Time of Up Slope relative to Pulse Shape Fit [ns].

Definition at line 261 of file HgcrocEmulator.h.

Referenced by HgcrocEmulator().

◆ timing_jitter_

double ldmx::HgcrocEmulator::timing_jitter_
private

Jitter of timing mechanism in the chip [ns].

Definition at line 252 of file HgcrocEmulator.h.

Referenced by HgcrocEmulator().


The documentation for this class was generated from the following files: