Belle II Software development
SVDDigitizerModule.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#include <svd/modules/svdSimulation/SVDDigitizerModule.h>
10#include <vxd/geometry/GeoCache.h>
11
12#include <framework/dataobjects/EventMetaData.h>
13#include <framework/logging/Logger.h>
14#include <framework/gearbox/Unit.h>
15#include <framework/gearbox/Const.h>
16#include <framework/datastore/DataStore.h>
17#include <framework/datastore/StoreArray.h>
18#include <framework/datastore/StoreObjPtr.h>
19#include <framework/datastore/RelationArray.h>
20#include <framework/datastore/RelationIndex.h>
21#include <mdst/dataobjects/MCParticle.h>
22#include <svd/dataobjects/SVDTrueHit.h>
23#include <svd/dataobjects/SVDShaperDigit.h>
24#include <svd/dataobjects/SVDEventInfo.h>
25#include <fstream>
26#include <sstream>
27#include <regex>
28#include <algorithm>
29#include <deque>
30#include <cmath>
31#include <root/TMath.h>
32#include <root/TRandom.h>
33
34using namespace std;
35using namespace std::placeholders;
36using namespace Belle2;
37using namespace Belle2::SVD;
38
39//-----------------------------------------------------------------
40// Auxiliaries for waveform storage
41//-----------------------------------------------------------------
42// FIXME: I believe this is wrong, these things must be created separately for
43// each instance of the digitizer, otherwise it can crash in parallel mode.
44int tree_vxdID = 0; // VXD ID of a sensor
45int tree_uv = 0; // U or V coordinate
46int tree_strip = 0; // number of strip
47double tree_signal[20]; // array for 20 samples of 10 ns
48
49//-----------------------------------------------------------------
50// Register the Module
51//-----------------------------------------------------------------
52REG_MODULE(SVDDigitizer);
53
54//-----------------------------------------------------------------
55// Implementation
56//-----------------------------------------------------------------
57
58std::string Belle2::SVD::SVDDigitizerModule::m_xmlFileName = std::string("SVDChannelMapping.xml");
60 m_currentSensorWaveforms(nullptr),
61 m_mapping(m_xmlFileName)
62{
63 //Set module properties
64 setDescription("Creates SVDShaperDigits from SVDSimHits");
66
67 // Define module parameters
68
69 // 1. Collections
70 addParam("MCParticles", m_storeMCParticlesName,
71 "MCParticle collection name", m_storeMCParticlesName);
72 addParam("SimHits", m_storeSimHitsName, "SimHit collection name",
74 addParam("TrueHits", m_storeTrueHitsName, "TrueHit collection name",
76 addParam("ShaperDigits", m_storeShaperDigitsName, "ShaperDigits collection name", m_storeShaperDigitsName);
77 addParam("SVDEventInfo", m_svdEventInfoName, "SVDEventInfo name", m_svdEventInfoName);
78
79 // 2. Physics
80 addParam("SegmentLength", m_segmentLength,
81 "Maximum segment length (in millimeters)", m_segmentLength);
82
83 // 3. Noise
84 addParam("PoissonSmearing", m_applyPoisson,
85 "Apply Poisson smearing on chargelets", m_applyPoisson);
86 addParam("ZeroSuppressionCut", m_SNAdjacent,
87 "Zero suppression cut in sigmas of strip noise", m_SNAdjacent);
88 addParam("FADCmode", m_roundZS,
89 "FADC mode: if True, ZS cut is rounded to nearest ADU ", m_roundZS);
90 addParam("numberOfSamples", m_nSamplesOverZS,
91 "Keep digit if numberOfSamples or more samples are over ZS threshold",
93
94 // 4. Timing
95 addParam("BetaPrimeDecayTimeU", m_betaPrimeDecayTimeU, "Decay time of betaprime waveform in ns, U-side",
97 addParam("BetaPrimeDecayTimeV", m_betaPrimeDecayTimeV, "Decay time of betaprime waveform in ns, V-side",
99 addParam("ADCSamplingTime", m_samplingTime,
100 "Interval between ADC samples in ns, if = -1 taken from HardwareClockSettings payload (default).", m_samplingTime);
101 addParam("StartSampling", m_startSampling,
102 "Start of the sampling window, in ns. Used to tune the SVD latency.", m_startSampling);
103 addParam("RandomizeEventTimes", m_randomizeEventTimes,
104 "Randomize event times over a frame interval", m_randomizeEventTimes);
105 addParam("TimeFrameLow", m_minTimeFrame,
106 "Left edge of event time randomization window, ns", m_minTimeFrame);
107 addParam("TimeFrameHigh", m_maxTimeFrame,
108 "Right edge of event time randomization window, ns", m_maxTimeFrame);
109
110
111 // 6. Reporting
112 addParam("statisticsFilename", m_rootFilename,
113 "ROOT Filename for statistics generation. If filename is empty, no statistics will be produced",
115 addParam("storeWaveforms", m_storeWaveforms,
116 "Store waveforms in a TTree in the statistics file.", m_storeWaveforms);
117 addParam("signalsList", m_signalsList,
118 "Store signals (time/charge/tau) in a tab-delimited file",
120}
121
123{
124
125 //Register all required collections
127
129
131 storeShaperDigits.registerInDataStore();
132 storeShaperDigits.registerRelationTo(storeMCParticles);
133 storeShaperDigits.registerRelationTo(storeTrueHits);
134
136 DataStore::arrayName<SVDShaperDigit>(m_storeShaperDigitsName),
137 DataStore::arrayName<MCParticle>(m_storeMCParticlesName));
139 DataStore::arrayName<SVDShaperDigit>(m_storeShaperDigitsName),
140 DataStore::arrayName<SVDTrueHit>(m_storeTrueHitsName));
141
143 StoreObjPtr<EventMetaData> storeEvents;
144
145 //Set names in case default was used. We need the names to initialize the RelationIndices.
147 DataStore::arrayName<MCParticle>(m_storeMCParticlesName),
148 DataStore::arrayName<SVDSimHit>(m_storeSimHitsName));
150 DataStore::arrayName<SVDTrueHit>(m_storeTrueHitsName),
151 DataStore::arrayName<SVDSimHit>(m_storeSimHitsName));
152
153 // Convert parameters to correct units
155 m_noiseFraction = TMath::Freq(m_SNAdjacent); // 0.9... !
161
162 B2DEBUG(29,
163 "SVDDigitizer parameters (in default system units, *=cannot be set directly):");
164 B2DEBUG(29, " DATASTORE COLLECTIONS:");
165 B2DEBUG(29,
166 " --> MCParticles: " << DataStore::arrayName<MCParticle>(m_storeMCParticlesName));
167 B2DEBUG(29,
168 " --> Digits: " << DataStore::arrayName<SVDShaperDigit>(m_storeShaperDigitsName));
169 B2DEBUG(29,
170 " --> SimHits: " << DataStore::arrayName<SVDSimHit>(m_storeSimHitsName));
171 B2DEBUG(29,
172 " --> TrueHits: " << DataStore::arrayName<SVDTrueHit>(m_storeTrueHitsName));
173 B2DEBUG(29, " --> MCSimHitRel: " << m_relMCParticleSimHitName);
174 B2DEBUG(29, " --> DigitMCRel: " << m_relShaperDigitMCParticleName);
175 B2DEBUG(29, " --> TrueSimRel: " << m_relTrueHitSimHitName);
176 B2DEBUG(29, " --> DigitTrueRel: " << m_relShaperDigitTrueHitName);
177 B2DEBUG(29, " PHYSICS: ");
178 B2DEBUG(29, " --> SegmentLength: " << m_segmentLength);
179 B2DEBUG(29, " NOISE: ");
180 B2DEBUG(29, " --> Add Poisson noise " << (m_applyPoisson ? "true" : "false"));
181 B2DEBUG(29, " --> Zero suppression cut" << m_SNAdjacent);
182 B2DEBUG(29, " --> Round ZS cut: " << (m_roundZS ? "true" : "false"));
183 B2DEBUG(29, " --> Samples over ZS cut:" << m_nSamplesOverZS);
184 B2DEBUG(29, " --> Noise fraction*: " << 1.0 - m_noiseFraction);
185 B2DEBUG(29, " TIMING: ");
186 B2DEBUG(29, " --> Sampling time: " << m_samplingTime);
187 B2DEBUG(29, " --> Start of int. wind.:" << m_startSampling);
188 B2DEBUG(29, " --> Random event times. " << (m_randomizeEventTimes ? "true" : "false"));
189 B2DEBUG(29, " REPORTING: ");
190 B2DEBUG(29, " --> statisticsFilename: " << m_rootFilename);
191 B2DEBUG(29,
192 " --> storeWaveforms: " << (m_storeWaveforms ? "true" : "false"));
193
194 if (!m_rootFilename.empty()) {
195 m_rootFile = new TFile(m_rootFilename.c_str(), "RECREATE");
196 m_rootFile->cd();
197 m_histChargeSharing_v = new TH1D("h_Diffusion_v", " 'Diffusion' distance, v",
198 200, -500, 500);
199 m_histChargeSharing_v->GetXaxis()->SetTitle(" distance v [um]");
200 m_histChargeSharing_u = new TH1D("h_Diffusion_u",
201 " 'Diffusion' distance, u", 100, -200, 200);
202 m_histChargeSharing_u->GetXaxis()->SetTitle("distance u [um]");
203 m_histLorentz_u = new TH1D("h_LorentzAngle_u", "Lorentz angle, holes",
204 100, -0.08, 0);
205 m_histLorentz_u->GetXaxis()->SetTitle("Lorentz angle");
206 m_histLorentz_v = new TH1D("h_LorentzAngle_v",
207 "Lorentz angle, electrons", 100, -0.002, 0.002);
208 m_histLorentz_v->GetXaxis()->SetTitle("Lorentz angle");
209 m_signalDist_u = new TH1D("h_signalDist_u",
210 "Strip signals vs. TrueHits, holes", 100, -400, 400);
211 m_signalDist_u->GetXaxis()->SetTitle("U strip position - TrueHit u [um]");
212 m_signalDist_v = new TH1D("h_signalDist_v",
213 "Strip signals vs. TrueHits, electrons", 100, -400, 400);
214 m_signalDist_v->GetXaxis()->SetTitle("V strip position - TrueHit v [um]");
215
216 m_histMobility_e = new TH1D("h_elecMobility", "electron Mobility",
217 30, 900, 1200);
218 m_histMobility_e->GetXaxis()->SetTitle("Electron Mobility");
219 m_histMobility_h = new TH1D("h_holeMobility", "hole Mobility",
220 30, 400, 500);
221 m_histMobility_h->GetXaxis()->SetTitle("Holes Mobility");
222
223 m_histDistanceToPlane_e = new TH1D("h_elecDistToPlane", "electron Distance to Plane",
224 50, -0.05, 0.05);
225 m_histDistanceToPlane_e->GetXaxis()->SetTitle("Electron Distance To Plane [cm]");
226 m_histDistanceToPlane_h = new TH1D("h_holeDistToPlane", "holes Distance to Plane",
227 50, -0.05, 0.05);
228 m_histDistanceToPlane_h->GetXaxis()->SetTitle("Holes Distance To Plane [cm]");
229
230 m_histVelocity_e = new TH1D("h_elecVelocity", "electrons Velocity (z)",
231 100, 0.001, 0.01);
232
233 m_histVelocity_e->GetXaxis()->SetTitle("Electron Velocity [cm/s]");
234 m_histVelocity_h = new TH1D("h_holeVelocity", "holes Velocity (z)",
235 30, -0.002, -0.0004);
236 m_histVelocity_h->GetXaxis()->SetTitle("holes Velocity [cm/s]");
237
238 m_histDriftTime_e = new TH1D("h_elecDriftTime", "electron Drift Time",
239 30, 0, 30);
240 m_histDriftTime_e->GetXaxis()->SetTitle("Electron Drift Time");
241 m_histDriftTime_h = new TH1D("h_holeDriftTime", "hole Drift Time",
242 30, 0, 30);
243 m_histDriftTime_h->GetXaxis()->SetTitle("Hole Drift Time");
244
245 m_histHitTime = new TH1D("h_startAPVTime", "start APV Time",
246 200, -100, 100);
247 m_histHitTime->GetXaxis()->SetTitle("time (ns)");
248 m_histHitTimeTB = new TH2F("h_startAPVTimeTB", "start APV Time vs TB",
249 200, -100, 100, 4, -0.5, 3.5);
250 m_histHitTimeTB->GetXaxis()->SetTitle("time (ns)");
251 m_histHitTimeTB->GetYaxis()->SetTitle("TB");
252
253 if (m_storeWaveforms) {
254 m_waveTree = new TTree("waveTree", "SVD waveforms");
255 m_waveTree->Branch("sensor", &tree_vxdID, "sensor/I");
256 m_waveTree->Branch("u_or_v", &tree_uv, "u_or_v/I");
257 m_waveTree->Branch("strip", &tree_strip, "strip/I");
258 m_waveTree->Branch("signal", tree_signal, "signal[20]/D");
259 }
260 } else {
261 // No waveforms can be stored if there is no statistics file.
262 m_storeWaveforms = false;
263 }
264}
265
267{
268
269 if (m_mapping.hasChanged()) { m_map = std::make_unique<SVDOnlineToOfflineMap>(m_mapping->getFileName()); }
270
271 //read sampling time from HardwareClockSettings
272 if (m_samplingTime == -1 && m_hwClock.isValid())
273 m_samplingTime = 1. / m_hwClock->getClockFrequency(Const::EDetector::SVD, "sampling");
274 else if (m_samplingTime == -1)
275 m_samplingTime = 16000. / 509;
276
277 //Fill map with all possible sensors This is too slow to be done every event so
278 //we fill it once and only clear the content of the sensors per event, not
279 //the whole map
280 m_waveforms.clear();
282 for (VxdID layer : geo.getLayers(SensorInfo::SVD)) {
283 for (VxdID ladder : geo.getLadders(layer)) {
284 for (VxdID sensor : geo.getSensors(ladder)) {
285 m_waveforms[sensor] = SensorWaveforms();
286 }
287 }
288 }
289
290 if (!m_MaskedStr.isValid())
291 B2WARNING("No valid SVDFADCMaskedStrip for the requested IoV -> no strips masked");
292 if (!m_map)
293 B2WARNING("No valid channel mapping -> all APVs will be enabled");
294
295
296}
297
299{
300
301 //get number of samples and relativeShift
303 SVDModeByte modeByte = storeSVDEvtInfo->getModeByte();
304 m_relativeShift = storeSVDEvtInfo->getRelativeShift();
305 m_nAPV25Samples = storeSVDEvtInfo->getNSamples();
306
307 //Compute time of the first sample, update latency
308 const double systemClockPeriod = 1. / m_hwClock->getGlobalClockFrequency();
309 int triggerBin = modeByte.getTriggerBin();
310
311 m_initTime = m_startSampling - systemClockPeriod * triggerBin;
312
313 m_is3sampleEvent = false;
314 if (m_nAPV25Samples == 3) {
315 m_is3sampleEvent = true;
317 B2DEBUG(25, "3-sample event, starting sample = " << m_startingSample);
318 } else m_startingSample = 0; //not used
319
320 // set APV mode for background overlay
322
323 // Generate current event time
327 // We have negative event times, so we have to encode!
328 storeEvent->setTime(static_cast<unsigned long>(1000 + m_currentEventTime));
329 } else
330 m_currentEventTime = 0.0;
331
332 // Clear sensors' waveforms and process SimHits
333 for (Waveforms::value_type& sensorWaveforms : m_waveforms) {
334 sensorWaveforms.second.first.clear(); // u-side channels
335 sensorWaveforms.second.second.clear(); // v-side channels
336 }
337 // m_currentSensorWaveforms = 0;
338 // m_currentSensorInfo = 0;
339
343
344 RelationArray mcParticlesToSimHits(storeMCParticles, storeSimHits, m_relMCParticleSimHitName);
345 RelationArray trueHitsToSimHits(storeTrueHits, storeSimHits, m_relTrueHitSimHitName);
346
347 RelationIndex<MCParticle, SVDSimHit> relMCParticleSimHit(storeMCParticles, storeSimHits, m_relMCParticleSimHitName);
348 RelationIndex<SVDTrueHit, SVDSimHit> relTrueHitSimHit(storeTrueHits, storeSimHits, m_relTrueHitSimHitName);
349
350 unsigned int nSimHits = storeSimHits.getEntries();
351 if (nSimHits == 0) {
352 return;
353 }
354
355 //Check sensor info and set pointers to current sensor
356 for (unsigned int i = 0; i < nSimHits; ++i) {
357 m_currentHit = storeSimHits[i];
359 relMCParticleSimHit.getFirstElementTo(m_currentHit);
360 if (mcRel) {
362 if (mcRel->weight < 0) {
363 //This simhit is from a particle which was not saved by the simulation
364 //so we do not take it into account for relations. Otherwise we might
365 //end up adding positive and negative weights together
367 }
368 } else {
369 // Don't bother with warnings for background SimHits
371 B2WARNING(
372 "Could not find MCParticle which produced SVDSimhit " << i);
374 }
376 relTrueHitSimHit.getFirstElementTo(m_currentHit);
377 //We only care about true hits from particles which have not been ignored
378 if (trueRel && trueRel->weight > 0) {
379 m_currentTrueHit = trueRel->indexFrom;
380 } else {
381 m_currentTrueHit = -1;
382 }
383
384 VxdID sensorID = m_currentHit->getSensorID();
385 if (!m_currentSensorInfo || sensorID != m_currentSensorInfo->getID()) {
387 dynamic_cast<const SensorInfo*>(&VXD::GeoCache::getInstance().getSensorInfo(sensorID));
389 B2FATAL(
390 "Sensor Information for Sensor " << sensorID << " not found, make sure that the geometry is set up correctly");
391
392 const SensorInfo& info = *m_currentSensorInfo;
393 // Publish some useful data
394 m_sensorThickness = info.getThickness();
396 B2DEBUG(29,
397 "Sensor Parameters for Sensor " << sensorID << ": " << endl
398 << " --> Width: " << m_currentSensorInfo->getWidth() << endl
399 << " --> Length: " << m_currentSensorInfo->getLength() << endl
400 << " --> uPitch: " << m_currentSensorInfo->getUPitch() << endl
401 << " --> vPitch: " << m_currentSensorInfo->getVPitch(-m_currentSensorInfo->getLength() / 2.0)
402 << ", " << m_currentSensorInfo->getVPitch(m_currentSensorInfo->getLength() / 2.0) << endl
403 << " --> Thickness: " << m_currentSensorInfo->getThickness() << endl
404 << " --> Deplet. voltage:" << m_currentSensorInfo->getDepletionVoltage() << endl
405 << " --> Bias voltage: " << m_currentSensorInfo->getBiasVoltage() << endl
406 );
407
408 }
409 B2DEBUG(28,
410 "Processing hit " << i << " in Sensor " << sensorID << ", related to MCParticle " << m_currentParticle);
411 processHit();
412 }
413 // If storage of waveforms is required, store them in the statistics file.
414 if (m_storeWaveforms) {
415 m_rootFile->cd();
417 }
418 if (m_signalsList != "")
419 saveSignals();
420
421
422 saveDigits();
423
424
425
426}
427
429{
430 // Set time of the hit
432
433 //Get Steplength and direction
434 const ROOT::Math::XYZVector& startPoint = m_currentHit->getPosIn();
435 const ROOT::Math::XYZVector& stopPoint = m_currentHit->getPosOut();
436 ROOT::Math::XYZVector direction = stopPoint - startPoint;
437 double trackLength = direction.R();
438
439 if (m_currentHit->getPDGcode() == Const::photon.getPDGCode() || trackLength < 0.1 * Unit::um) {
440 //Photons deposit energy at the end of their step
441 driftCharge(stopPoint, m_currentHit->getElectrons(), SVD::SensorInfo::electron);
443 } else {
444 //Otherwise, split into segments of (default) max. 5µm and
445 //drift the charges from the center of each segment
447 double lastFraction {0};
448 double lastElectrons {0};
449
450 for (auto& segment : segments) {
451 //Simhit returns step fraction and cumulative electrons. We want the
452 //center of these steps and electrons in this step
453 const double f = (segment.first + lastFraction) / 2;
454 const double e = segment.second - lastElectrons;
455 //Update last values
456 std::tie(lastFraction, lastElectrons) = segment;
457
458 //And drift charge from that position
459 const ROOT::Math::XYZVector position = startPoint + f * direction;
460 driftCharge(position, e, SVD::SensorInfo::electron);
462 }
463 }
464}
465
466
467void SVDDigitizerModule::driftCharge(const ROOT::Math::XYZVector& position, double carriers,
469{
470 bool have_electrons = (carrierType == SVD::SensorInfo::electron);
471
472 string carrierName = (have_electrons) ? "electron" : "hole";
473 B2DEBUG(29,
474 "Drifting " << carriers << " " << carrierName << "s at position (" << position.X() << ", " << position.Y() << ", " << position.Z()
475 << ").");
476 B2DEBUG(29, "@@@ driftCharge: drifting " << carriers << " " << carrierName << "s at position (" << position.X() << ", " <<
477 position.Y() << ", " << position.Z()
478 << ").");
479
480 // Get references to current sensor/info for ease of use
481 const SensorInfo& info = *m_currentSensorInfo;
482 StripWaveforms& waveforms = (!have_electrons) ? m_currentSensorWaveforms->first : m_currentSensorWaveforms->second;
483
484 double distanceToPlane = (have_electrons) ?
485 0.5 * m_sensorThickness - position.Z() :
486 -0.5 * m_sensorThickness - position.Z(); //cm
487
488 if (m_histDistanceToPlane_e && have_electrons) m_histDistanceToPlane_e->Fill(distanceToPlane);
489 if (m_histDistanceToPlane_h && !have_electrons) m_histDistanceToPlane_h->Fill(distanceToPlane);
490
491 // Approximation: calculate drift velocity at the point halfway towards
492 // the respective sensor surface.
493 ROOT::Math::XYZVector mean_pos(position.X(), position.Y(), position.Z() + 0.5 * distanceToPlane);
494
495 // Calculate drift times and widths of charge clouds.
496 ROOT::Math::XYZVector v = info.getVelocity(carrierType, mean_pos);
497 if (m_histVelocity_e && have_electrons) m_histVelocity_e->Fill(v.Z()); //Unit::cm/Unit::cm*Unit::eV/Unit::e*Unit::s);
498 if (m_histVelocity_h && !have_electrons) m_histVelocity_h->Fill(v.Z()); //Unit::cm/Unit::cm*Unit::eV/Unit::e*Unit::s);
499
500 double driftTime = distanceToPlane / v.Z(); //ns
501 if (m_histDriftTime_e && have_electrons) m_histDriftTime_e->Fill(driftTime); //ns
502 if (m_histDriftTime_h && !have_electrons) m_histDriftTime_h->Fill(driftTime); //ns
503
504 ROOT::Math::XYZVector center = position + driftTime * v; //cm
505 double mobility = (have_electrons) ?
506 info.getElectronMobility(info.getEField(mean_pos).R()) :
507 info.getHoleMobility(info.getEField(mean_pos).R());
508
509 if (m_histMobility_e && have_electrons) m_histMobility_e->Fill(mobility); //cm2/V/ns
510 if (m_histMobility_h && !have_electrons) m_histMobility_h->Fill(mobility); //cm2/V/ns
511
512 double D = Const::kBoltzmann * info.getTemperature() / Unit::e * mobility;
513 double sigma = std::max(1.0e-4, sqrt(2.0 * D * driftTime));
514 double tanLorentz = (!have_electrons) ? v.X() / v.Z() : v.Y() / v.Z();
515
516 B2DEBUG(29, "velocity (" << v.X() / Unit::um << ", " << v.Y() / Unit::um << ", " << v.Z() / Unit::um << ") um/ns");
517 B2DEBUG(29, "D = " << D << ", driftTime = " << driftTime / Unit::ns << " ns");
518 B2DEBUG(29, "sigma = " << sigma / Unit::um << " um");
519 B2DEBUG(29, "tan Lorentz = " << tanLorentz);
520
521 sigma *= sqrt(1.0 + tanLorentz * tanLorentz);
522 if (m_histLorentz_u && !have_electrons) m_histLorentz_u->Fill(tanLorentz);
523 if (m_histLorentz_v && have_electrons) m_histLorentz_v->Fill(tanLorentz);
524
525 //Distribute carrier cloud on strips
526 int vID = info.getVCellID(center.Y(), true);
527 int uID = info.getUCellID(center.X(), center.Y(), true);
528 int seedStrip = (!have_electrons) ? uID : vID;
529 double seedPos = (!have_electrons) ?
530 info.getUCellPosition(seedStrip, vID) :
531 info.getVCellPosition(seedStrip);
532 double geomPitch = (!have_electrons) ? 0.5 * info.getUPitch(center.Y()) : 0.5 * info.getVPitch();
533 int nCells = (!have_electrons) ? info.getUCells() : info.getVCells();
534 std::deque<double> stripCharges;
535 std::deque<double> strips; // intermediate strips will be half-integers, like 2.5.
536#define NORMAL_CDF(z) 0.5 * std::erfc( - (z) * 0.707107)
537 double current_pos = (!have_electrons) ? seedPos - center.X() : seedPos - center.Y();
538 double current_strip = seedStrip;
539 double cdf_low = NORMAL_CDF((current_pos - 0.5 * geomPitch) / sigma);
540 double cdf_high = NORMAL_CDF((current_pos + 0.5 * geomPitch) / sigma);
541 double charge = carriers * (cdf_high - cdf_low);
542
543 B2DEBUG(29, "geomPitch = " << geomPitch / Unit::um << " um");
544 B2DEBUG(29, "charge = " << charge << " = " << carriers << "(carriers) * (" << cdf_high << "(cdf_high) - " << cdf_low <<
545 "(cdf_low));");
546
547 stripCharges.push_back(charge);
548 strips.push_back(current_strip);
549 while (cdf_low > 1.0e-5) {
550 current_pos -= geomPitch;
551 current_strip -= 0.5;
552 double cdf_current = NORMAL_CDF((current_pos - 0.5 * geomPitch) / sigma);
553 charge = carriers * (cdf_low - cdf_current);
554 stripCharges.push_front(charge);
555 strips.push_front(current_strip);
556 cdf_low = cdf_current;
557 }
558 current_pos = (!have_electrons) ? seedPos - center.X() : seedPos - center.Y();
559 current_strip = seedStrip;
560 while (cdf_high < 1.0 - 1.0e-5) {
561 current_pos += geomPitch;
562 current_strip += 0.5;
563 double cdf_current = NORMAL_CDF((current_pos + 0.5 * geomPitch) / sigma);
564 charge = carriers * (cdf_current - cdf_high);
565 stripCharges.push_back(charge);
566 strips.push_back(current_strip);
567 cdf_high = cdf_current;
568 }
569#undef NORMAL_CDF
570
571 // Pad with zeros for smoothing
572 int npads = (strips.front() - floor(strips.front()) == 0) ? 5 : 4;
573 for (int i = 0; i < npads; ++i) {
574 strips.push_front(strips.front() - 0.5);
575 stripCharges.push_front(0);
576 }
577 npads = (strips.back() - floor(strips.back()) == 0) ? 5 : 4;
578 for (int i = 0; i < npads; ++i) {
579 strips.push_back(strips.back() + 0.5);
580 stripCharges.push_back(0);
581 }
582 // Charge sharing
583 B2DEBUG(29, " --> charge sharing simulation, # strips = " << strips.size());
584 std::deque<double> readoutCharges;
585 std::deque<int> readoutStrips;
586 VxdID currentSensorID = m_currentHit->getSensorID();
587 for (std::size_t index = 3; index < strips.size() - 3; index += 2) {
588 B2DEBUG(29, " index = " << index << ", strip = " << strips[index] << ", stripCharge = " << stripCharges[index]);
589 int currentStrip = static_cast<int>(strips[index]);
590
591 double c0 = m_ChargeSimCal.getCouplingConstant(currentSensorID, !have_electrons, "C0");
592 double c1 = m_ChargeSimCal.getCouplingConstant(currentSensorID, !have_electrons, "C1");
593 double c2 = m_ChargeSimCal.getCouplingConstant(currentSensorID, !have_electrons, "C2");
594 double c3 = m_ChargeSimCal.getCouplingConstant(currentSensorID, !have_electrons, "C3");
595
596 B2DEBUG(29, " current strip = " << currentStrip);
597 B2DEBUG(29, " index-3 = " << index - 3 << ", strip = " << strips[index - 3] << ", stripCharge = " << stripCharges[index - 3]);
598 B2DEBUG(29, " index-2 = " << index - 2 << ", strip = " << strips[index - 2] << ", stripCharge = " << stripCharges[index - 2]);
599 B2DEBUG(29, " index-1 = " << index - 1 << ", strip = " << strips[index - 1] << ", stripCharge = " << stripCharges[index - 1]);
600 B2DEBUG(29, " index = " << index << ", strip = " << strips[index] << ", stripCharge = " << stripCharges[index]);
601 B2DEBUG(29, " index+1 = " << index + 1 << ", strip = " << strips[index + 1] << ", stripCharge = " << stripCharges[index + 1]);
602 B2DEBUG(29, " index+2 = " << index + 2 << ", strip = " << strips[index + 2] << ", stripCharge = " << stripCharges[index + 2]);
603 B2DEBUG(29, " index+3 = " << index + 3 << ", strip = " << strips[index + 3] << ", stripCharge = " << stripCharges[index + 3]);
604
605 readoutCharges.push_back(c3 * stripCharges[index - 3]
606 + c2 * stripCharges[index - 2]
607 + c1 * stripCharges[index - 1]
608 + c0 * stripCharges[index]
609 + c1 * stripCharges[index + 1]
610 + c2 * stripCharges[index + 2]
611 + c3 * stripCharges[index + 3]
612 );
613 readoutStrips.push_back(currentStrip);
614 B2DEBUG(29, " post simulation: " << index << ", strip = " << currentStrip << ", readoutCharge = " <<
615 readoutCharges[readoutCharges.size() - 1]);
616 }
617 // Trim at sensor edges
618 double tail = 0;
619 while (readoutStrips.size() > 0 && readoutStrips.front() < 0) {
620 readoutStrips.pop_front();
621 tail += readoutCharges.front();
622 readoutCharges.pop_front();
623 }
624 readoutCharges.front() += tail;
625 tail = 0;
626 while (readoutStrips.size() > 0 && readoutStrips.back() > nCells - 1) {
627 readoutStrips.pop_back();
628 tail += readoutCharges.back();
629 readoutCharges.pop_back();
630 }
631 readoutCharges.back() += tail;
632 // Poisson smearing - Gaussian approximation
633 if (m_applyPoisson)
634 for (auto& c : readoutCharges)
635 c = (c <= 0) ? 0 : std::max(0.0, gRandom->Gaus(c, std::sqrt(info.c_fanoFactorSi * c)));
636
637 // Fill diagnostic charts
639 TH1D* histo = (!have_electrons) ? m_histChargeSharing_u : m_histChargeSharing_v;
640 double d = (!have_electrons) ? seedPos - center.X() : seedPos - center.Y();
641 for (std::size_t index = 0; index < readoutStrips.size(); ++ index) {
642 double dist = d + (readoutStrips[index] - seedStrip) * 2 * geomPitch;
643 histo->Fill(dist / Unit::um, readoutCharges[index]);
644 }
645 }
646 if (m_histHitTime) {
649 SVDModeByte modeByte = storeSVDEvtInfo->getModeByte();
651 }
652
653 // Store
654 B2DEBUG(29, "currentTime = " << m_currentTime << " + 0.5 driftTime = " << 0.5 * driftTime << " = " << m_currentTime + 0.5 *
655 driftTime);
656
657 // Specify beta prime decay time
658 double betaPrimeDecayTime = (!have_electrons) ? m_betaPrimeDecayTimeU : m_betaPrimeDecayTimeV;
659
660 // Specify coupling and adjacent-channel waveform shape
661 double apvCoupling = m_ChargeSimCal.getCouplingConstant(currentSensorID, !have_electrons, "APVCoupling");
662 WaveformShape w_adjacent = (!have_electrons) ? w_adjacentU : w_adjacentV;
663
664 double recoveredCharge = 0;
665 for (std::size_t index = 0; index < readoutStrips.size(); index ++) {
666 // NB> To first approximation, we assign to the signal 1/2*driftTime.
667 // This doesn't change the charge collection, only charge collection timing.
668 waveforms[readoutStrips[index]].add(m_currentTime + 0.5 * driftTime, readoutCharges[index],
669 betaPrimeDecayTime, m_currentParticle, m_currentTrueHit, w_betaprime);
670 // coupled signal left neighbour
671 if (index > 0)
672 waveforms[readoutStrips[index]].add(m_currentTime + 0.5 * driftTime, apvCoupling * readoutCharges[index - 1],
673 1, m_currentParticle, m_currentTrueHit, w_adjacent);
674 // coupled signal right neighbour
675 if (index < readoutStrips.size() - 1)
676 waveforms[readoutStrips[index]].add(m_currentTime + 0.5 * driftTime, apvCoupling * readoutCharges[index + 1],
677 1, m_currentParticle, m_currentTrueHit, w_adjacent);
678 recoveredCharge += readoutCharges[index];
679 B2DEBUG(29, "strip: " << readoutStrips[index] << " charge: " << readoutCharges[index]);
680 }
681 B2DEBUG(29, "Digitized " << recoveredCharge << " of " << carriers << " original carriers.");
682}
683
684double SVDDigitizerModule::addNoise(double charge, double noise)
685{
686 charge += gRandom->Gaus(0., noise);
687 return charge;
688}
689
691{
692
696 RelationArray relShaperDigitMCParticle(storeShaperDigits, storeMCParticles,
698 RelationArray relShaperDigitTrueHit(storeShaperDigits, storeTrueHits,
700
701 //Get SVD config from SVDEventInfo
702 // int runType = (int) modeByte.getRunType();
703 // int eventType = (int) modeByte.getEventType();
704
705
706 // ... to store digit-digit relations
707 vector<pair<unsigned int, float> > digit_weights;
708
709 // Take samples at the desired times, add noise, zero-suppress and save digits.
710 for (Waveforms::value_type& sensorWaveforms : m_waveforms) {
711 int sensorID = sensorWaveforms.first;
712 // u-side digits:
713
714 // Cycle through signals and generate samples
715 for (StripWaveforms::value_type& stripWaveform : sensorWaveforms.second.first) {
716 short int iStrip = stripWaveform.first;
717 SVDWaveform& s = stripWaveform.second;
718 // Now generate samples in time and save as digits.
719 vector<double> samples;
720 // ... to store digit-digit relations
721 digit_weights.clear();
722 digit_weights.reserve(SVDShaperDigit::c_nAPVSamples);
723
724 double elNoise = m_NoiseCal.getNoiseInElectrons(sensorID, true, iStrip);
725 double gain = 1 / m_PulseShapeCal.getChargeFromADC(sensorID, true, iStrip, 1);
726 double electronWeight = m_ChargeSimCal.getElectronWeight(sensorID, true);
727
728 double t = m_initTime;
729 B2DEBUG(25, "start sampling at " << m_initTime);
730 for (int iSample = 0; iSample < (int) SVDShaperDigit::c_nAPVSamples; iSample ++) {
731 samples.push_back(addNoise(electronWeight * s(t), elNoise));
732 t += m_samplingTime;
733 }
734
735 SVDWaveform::relations_map particles = s.getMCParticleRelations();
736 SVDWaveform::relations_map truehits = s.getTrueHitRelations();
737
738 // Save SVDShaperDigits
739
740 // 1. Convert to ADU
742 std::transform(samples.begin(), samples.end(), rawSamples.begin(),
743 [&](double x)->SVDShaperDigit::APVRawSampleType {
744 return SVDShaperDigit::trimToSampleRange(x * gain);
745 });
746
747 // 2.a Check if over threshold
748 auto rawThreshold = m_SNAdjacent * elNoise * gain;
749 if (m_roundZS) rawThreshold = round(rawThreshold);
750 auto n_over = std::count_if(rawSamples.begin(), rawSamples.end(),
751 std::bind(std::greater<double>(), _1, rawThreshold)
752 );
753 if (n_over < m_nSamplesOverZS) continue;
754
755 // 2.b check if the strip is masked
756 if (m_MaskedStr.isMasked(sensorID, true, iStrip)) continue;
757
758 // 2.c check if the APV is disabled
759 if (!m_map->isAPVinMap(sensorID, true, iStrip)) continue;
760
761 // 2.d.1 check if it's a 3-sample event
762 if (m_is3sampleEvent) {
763 rawSamples[0] = rawSamples[m_startingSample];
764 rawSamples[1] = rawSamples[m_startingSample + 1];
765 rawSamples[2] = rawSamples[m_startingSample + 2];
766 rawSamples[3] = 0.;
767 rawSamples[4] = 0.;
768 rawSamples[5] = 0.;
769 //2.d.2 check if still over threshold
770 n_over = std::count_if(rawSamples.begin(), rawSamples.end(),
771 std::bind(std::greater<double>(), _1, rawThreshold)
772 );
773 if (n_over < m_nSamplesOverZS) continue;
774
775 }
776
777 // 3. Save as a new digit
778 int digIndex = storeShaperDigits.getEntries();
779 storeShaperDigits.appendNew(sensorID, true, iStrip, rawSamples, 0);
780
781 //If the digit has any relations to MCParticles, add the Relation
782 if (particles.size() > 0) {
783 relShaperDigitMCParticle.add(digIndex, particles.begin(), particles.end());
784 }
785 //If the digit has any relations to truehits, add the Relations.
786 if (truehits.size() > 0) {
787 relShaperDigitTrueHit.add(digIndex, truehits.begin(), truehits.end());
788 }
789 // generate SVDShaperDigits
790 } // for stripSignals
791
792 // v-side digits:
793
794 // Cycle through signals and generate samples
795 for (StripWaveforms::value_type& stripWaveform : sensorWaveforms.second.second) {
796 short int iStrip = stripWaveform.first;
797 SVDWaveform& s = stripWaveform.second;
798 // Now generate samples in time and save as digits.
799 vector<double> samples;
800 // ... to store digit-digit relations
801 digit_weights.clear();
802 digit_weights.reserve(SVDShaperDigit::c_nAPVSamples);
803
804 double elNoise = m_NoiseCal.getNoiseInElectrons(sensorID, false, iStrip);
805 double gain = 1 / m_PulseShapeCal.getChargeFromADC(sensorID, false, iStrip, 1);
806 double electronWeight = m_ChargeSimCal.getElectronWeight(sensorID, false);
807
808 double t = m_initTime;
809 for (int iSample = 0; iSample < (int)SVDShaperDigit::c_nAPVSamples; iSample ++) {
810 samples.push_back(addNoise(electronWeight * s(t), elNoise));
811 t += m_samplingTime;
812 }
813
814 SVDWaveform::relations_map particles = s.getMCParticleRelations();
815 SVDWaveform::relations_map truehits = s.getTrueHitRelations();
816
817 // Save SVDShaperDigits
818 // 1. Convert to ADU
820 std::transform(samples.begin(), samples.end(), rawSamples.begin(),
821 [&](double x)->SVDShaperDigit::APVRawSampleType {
822 return SVDShaperDigit::trimToSampleRange(x * gain);
823 });
824
825 // 2.a Check if over threshold
826 auto rawThreshold = m_SNAdjacent * elNoise * gain;
827 if (m_roundZS) rawThreshold = round(rawThreshold);
828 auto n_over = std::count_if(rawSamples.begin(), rawSamples.end(),
829 std::bind(std::greater<double>(), _1, rawThreshold)
830 );
831 if (n_over < m_nSamplesOverZS) continue;
832
833 // 2.b check if the strip is masked
834 if (m_MaskedStr.isMasked(sensorID, false, iStrip)) continue;
835
836 // 2.c check if the APV is disabled
837 if (!m_map->isAPVinMap(sensorID, false, iStrip)) continue;
838
839 // 2.d.1 check if it's a 3-sample event
840 if (m_is3sampleEvent) {
841 rawSamples[0] = rawSamples[m_startingSample];
842 rawSamples[1] = rawSamples[m_startingSample + 1];
843 rawSamples[2] = rawSamples[m_startingSample + 2];
844 rawSamples[3] = 0.;
845 rawSamples[4] = 0.;
846 rawSamples[5] = 0.;
847 //2.d.2 check if still over threshold
848 n_over = std::count_if(rawSamples.begin(), rawSamples.end(),
849 std::bind(std::greater<double>(), _1, rawThreshold)
850 );
851 if (n_over < m_nSamplesOverZS) continue;
852 }
853
854 // 3. Save as a new digit
855 int digIndex = storeShaperDigits.getEntries();
856 storeShaperDigits.appendNew(sensorID, false, iStrip, rawSamples, 0);
857
858 //If the digit has any relations to MCParticles, add the Relation
859 if (particles.size() > 0) {
860 relShaperDigitMCParticle.add(digIndex, particles.begin(), particles.end());
861 }
862 //If the digit has any relations to truehits, add the Relations.
863 if (truehits.size() > 0) {
864 relShaperDigitTrueHit.add(digIndex, truehits.begin(), truehits.end());
865 }
866 } // for stripSignals
867 } // FOREACH sensor
868}
869
871{
872 for (Waveforms::value_type& sensorWaveforms : m_waveforms) {
873 tree_vxdID = sensorWaveforms.first;
874 const SensorInfo& info =
875 dynamic_cast<const SensorInfo&>(VXD::GeoCache::getInstance().getSensorInfo(sensorWaveforms.first));
876 // u-side digits:
877 tree_uv = 1;
878 double thresholdU = 3.0 * info.getElectronicNoiseU();
879 for (StripWaveforms::value_type& stripWaveform : sensorWaveforms.second.first) {
880 tree_strip = stripWaveform.first;
881 SVDWaveform& s = stripWaveform.second;
882 // Read the value only if the signal is large enough.
883 if (s.getCharge() < thresholdU)
884 continue;
885 for (int iTime = 0; iTime < 20; ++iTime) {
886 tree_signal[iTime] = s(10 * iTime);
887 }
888 m_waveTree->Fill();
889 } // FOREACH stripSignal
890 // v-side digits:
891 tree_uv = 0;
892 double thresholdV = 3.0 * info.getElectronicNoiseV();
893 for (StripWaveforms::value_type& stripWaveform : sensorWaveforms.second.second) {
894 tree_strip = stripWaveform.first;
895 SVDWaveform& s = stripWaveform.second;
896 // Read the values only if the signal is large enough
897 if (s.getCharge() < thresholdV)
898 continue;
899 for (int iTime = 0; iTime < 20; ++iTime) {
900 tree_signal[iTime] = s(10. * iTime);
901 }
902 m_waveTree->Fill();
903 } // FOREACH stripSignal
904 } // FOREACH sensor
905 m_rootFile->Flush();
906}
907
909{
910 static size_t recordNo = 0;
911 static const string header("Event\tSensor\tSide\tStrip\tContrib\tTime\tCharge\tTau");
912 regex startLine("^|\n"); // for inserting event/sensor/etc info
913 ofstream outfile(m_signalsList, ios::out | ios::app);
914 if (recordNo == 0) outfile << header << endl;
915 for (Waveforms::value_type& sensorWaveforms : m_waveforms) {
916 VxdID sensorID(sensorWaveforms.first);
917 const SensorInfo& info =
918 dynamic_cast<const SensorInfo&>(VXD::GeoCache::getInstance().getSensorInfo(sensorWaveforms.first));
919 // u-side digits:
920 size_t isU = 1;
921 double thresholdU = 3.0 * info.getElectronicNoiseU();
922 for (StripWaveforms::value_type& stripWaveform : sensorWaveforms.second.first) {
923 size_t strip = stripWaveform.first;
924 SVDWaveform& s = stripWaveform.second;
925 // Read the value only if the signal is large enough.
926 if (s.getCharge() < thresholdU)
927 continue;
928 // Else print to a string
929 ostringstream preamble;
930 // We don't have eventNo, but we don't care about event boundaries.
931 preamble << "$&" << recordNo << '\t' << sensorID << '\t' << isU << '\t' << strip << '\t';
932 string signalString = s.toString();
933 signalString.pop_back(); // remove the last newline!!!
934 string tableString = regex_replace(signalString, startLine, preamble.str());
935 outfile << tableString << endl; // now we have to add the newline back.
936 } // FOREACH stripSignal
937 // x-side digits:
938 isU = 0;
939 double thresholdV = 3.0 * info.getElectronicNoiseV();
940 for (StripWaveforms::value_type& stripWaveform : sensorWaveforms.second.second) {
941 size_t strip = stripWaveform.first;
942 SVDWaveform& s = stripWaveform.second;
943 // Read the value only if the signal is large enough.
944 if (s.getCharge() < thresholdV)
945 continue;
946 // Else print to a string
947 ostringstream preamble;
948 // We don't have eventNo, but we don't care about event boundaries.
949 preamble << "$&" << recordNo << '\t' << sensorID << '\t' << isU << '\t' << strip << '\t';
950 string signalString = s.toString();
951 signalString.pop_back(); // remove the last newline!!!
952 string tableString = regex_replace(signalString, startLine, preamble.str());
953 outfile << tableString << endl; // now we have to add the newline back.
954 } // FOREACH stripSignal
955 } // for sensors
956 outfile.close();
957 recordNo++;
958}
959
961{
962 if (m_rootFile) {
963 m_rootFile->Write();
964 m_rootFile->Close();
965 }
966}
967
968int SVDDigitizerModule::getFirstSample(int triggerBin, int relativeShift)
969{
970 int nTriggerClocks = triggerBin + relativeShift;
971 return floor(nTriggerClocks / 4);
972}
int getPDGCode() const
PDG code.
Definition: Const.h:473
static const double kBoltzmann
Boltzmann constant in GeV/K.
Definition: Const.h:696
static const ParticleType photon
photon particle
Definition: Const.h:673
bool hasChanged()
Check whether the object has changed since the last call to hasChanged of the accessor).
static std::string relationName(const std::string &fromName, const std::string &toName, std::string const &namedRelation="")
Return storage name for a relation between two arrays of the given names.
Definition: DataStore.h:180
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
std::string getFileName() const
Get the name of the downloaded payload file.
Definition: PayloadFile.h:35
Low-level class to create/modify relations between StoreArrays.
Definition: RelationArray.h:62
void add(index_type from, index_type to, weight_type weight=1.0)
Add a new element to the relation.
Provides access to fast ( O(log n) ) bi-directional lookups on a specified relation.
Definition: RelationIndex.h:76
const Element * getFirstElementTo(const TO &to) const
Return a pointer to the first relation Element of the given object.
float getCouplingConstant(const VxdID &sensorID, const bool &isU, const std::string &couplingName) const
Return coupling constant.
float getElectronWeight(const VxdID &sensorID, const bool &isU) const
Return Geant4 electron weight.
float isMasked(const VxdID &sensorID, const bool &isU, const unsigned short &strip) const
This is the method for getting the comprehensive list of masked strips at FADC level.
bool isValid()
returns true if the m_aDBObtPtr is valid in the requested IoV
Class to store SVD mode information.
Definition: SVDModeByte.h:69
baseType getTriggerBin() const
Get the triggerBin id.
Definition: SVDModeByte.h:140
float getNoiseInElectrons(const VxdID &sensorID, const bool &isU, const unsigned short &strip) const
This method provides the correct noise conversion into electrons, taking into account that the noise ...
double getChargeFromADC(const Belle2::VxdID &sensorID, const bool &isU, const unsigned short &strip, const double &pulseADC) const
Return the charge (number of electrons/holes) collected on a specific strip, given the number of ADC ...
static void setAPVMode(size_t mode, size_t firstSample)
set APV mode for the event
static const std::size_t c_nAPVSamples
Number of APV samples stored.
std::array< APVRawSampleType, c_nAPVSamples > APVRawSamples
array of APVRawSamplesType objects
double m_startSampling
Time window start, excluding trigger bin effect.
TH1D * m_histHitTime
Histogram showing the hit time.
TH1D * m_histDistanceToPlane_e
Histogram showing the distance to plane for e.
TH1D * m_histChargeSharing_u
Histogram showing the charge sharing + diffusion in u (r-phi).
double m_sensorThickness
Thickness of current sensor (read from m_currentSensorInfo)
void processHit()
Process one SVDSimHit by dividing the step in smaller steps and drifting the charge.
TH1D * m_histVelocity_h
Histogram showing the velocity of h.
SVDFADCMaskedStrips m_MaskedStr
FADC masked strip payload.
TH1D * m_histVelocity_e
Histogram showing the velocity of e-.
Waveforms m_waveforms
Structure containing waveforms in all existing sensors.
double m_samplingTime
Interval between two waveform samples, by default taken from HardwareClockSettings.
std::string m_relTrueHitSimHitName
Name of the relation between SVDTrueHits and SVDSimHits.
TH1D * m_histLorentz_v
Histogram showing the Lorentz angles in v (z).
void saveDigits()
Save digits to the DataStore Saves samples of generated waveforms.
bool m_randomizeEventTimes
Randomize event times? If set to true, event times will be randomized uniformly from m_minTimeFrame t...
int getFirstSample(int triggerBin, int relativShift)
return the starting sample
int m_startingSample
Starting sample for the selection of 3 samples in 3-mixed-6.
int m_currentParticle
Index of the particle which caused the current hit.
std::string m_relShaperDigitMCParticleName
Name of the relation between SVDShaperDigits and MCParticles.
int m_nSamplesOverZS
Keep digit if at least m_nSamplesOverZS are over threshold.
TH1D * m_histChargeSharing_v
Histogram showing the charge sharing + diffusion in v (z).
TH1D * m_histMobility_h
Histogram showing the mobility of h.
virtual void initialize() override
Initialize the module and check module parameters.
bool m_roundZS
Round ZS cut to nearest ADU.
std::string m_storeShaperDigitsName
Name of the collection for the SVDShaperDigits.
virtual void event() override
Digitize one event.
bool m_is3sampleEvent
True if the event should be simulated with 3 sample.
double m_SNAdjacent
Zero-suppression cut.
double m_currentTime
Time of the current SimHit.
TFile * m_rootFile
Pointer to the ROOT filename for statistics.
SVDNoiseCalibrations m_NoiseCal
SVDNoise calibrations db object.
SensorWaveforms * m_currentSensorWaveforms
Pointer to the sensor in which the current hit occurred.
std::string m_storeTrueHitsName
Name of the collection for the SVDTrueHits.
TH1D * m_signalDist_u
Histogram showing the distribution of digit signals in u (r-phi).
DBObjPtr< PayloadFile > m_mapping
channel mapping payload
virtual void terminate() override
Terminate the module.
void saveSignals()
Save signals to a root-delimited file (to be analyzed in Python).
std::string m_storeMCParticlesName
Name of the collection for the MCParticles.
TH1D * m_signalDist_v
Histogram showing the distribution of digit signals in v (z).
double m_noiseFraction
(derived from SNAdjacent) Fraction of noisy strips per sensor.
static std::string m_xmlFileName
< channel mapping xml filename
double m_betaPrimeDecayTimeU
Decay time of betaprime waveform U-side.
SVDPulseShapeCalibrations m_PulseShapeCal
SVDPulseShapeCalibrations calibrations db object.
TH1D * m_histMobility_e
Histogram showing the mobility of e-.
std::string m_relShaperDigitTrueHitName
Name of the relation between SVDShaperDigits and SVDTrueHits.
SVDChargeSimulationCalibrations m_ChargeSimCal
SVDChargeSimulationCalibrations calibrations db object.
int m_nAPV25Samples
number of digitized samples read from SVDEventInfo
TTree * m_waveTree
Tree for waveform storage.
bool m_applyPoisson
Whether or not to apply poisson fluctuation of charge (Fano factor)
virtual void beginRun() override
Initialize the list of existing SVD Sensors.
TH1D * m_histDistanceToPlane_h
Histogram showing the distance to plane for h.
TH1D * m_histDriftTime_e
Histogram showing the drift time of e.
DBObjPtr< HardwareClockSettings > m_hwClock
Hardware Clocks.
std::string m_svdEventInfoName
Name of the SVDEventInfo object.
std::string m_signalsList
Name of the tab-delimited listing of waveforms.
double addNoise(double charge, double noise)
Calculate the noise contribution to one strip with given charge.
const SensorInfo * m_currentSensorInfo
Pointer to the SensorInfo of the current sensor.
TH1D * m_histLorentz_u
Histogram showing the Lorentz angles in u (r-phi).
int m_relativeShift
relative shift in SVDEventInfo obj
TH2F * m_histHitTimeTB
Histogram showing the hit time vs TB.
bool m_storeWaveforms
Store waveform data in the reporting file?
void driftCharge(const ROOT::Math::XYZVector &position, double carriers, SVD::SensorInfo::CarrierType carrierType)
Drift the charge inside the silicon.
void saveWaveforms()
Save waveforms to the statistics file.
std::unique_ptr< SVDOnlineToOfflineMap > m_map
channel mapping map
float m_maxTimeFrame
High edge of randomization time frame.
TH1D * m_histDriftTime_h
Histogram showing the drift time of h.
double m_initTime
Time window start, including the triggerBin effect.
int m_currentTrueHit
Index of the TrueHit the current hit belongs to.
float m_minTimeFrame
Low edge of randomization time frame.
const SVDSimHit * m_currentHit
Pointer to the SVDSimhit currently digitized.
std::string m_relMCParticleSimHitName
Name of the relation between MCParticles and SVDSimHits.
double m_betaPrimeDecayTimeV
Decay time of betaprime waveform V-side.
std::string m_rootFilename
Name of the ROOT filename to output statistics.
std::string m_storeSimHitsName
Name of the collection for the SVDSimhits.
float m_currentEventTime
Current event time.
The SVD waveform class.
Definition: SVDWaveform.h:38
std::map< RelationElement::index_type, RelationElement::weight_type > relations_map
Type to store contributions to strip signal by different particles on output of SVDWaveform.
Definition: SVDWaveform.h:80
Specific implementation of SensorInfo for SVD Sensors which provides additional sensor specific infor...
Definition: SensorInfo.h:25
CarrierType
Enum to flag charge carriers.
Definition: SensorInfo.h:38
double getBiasVoltage() const
Return the bias voltage on the sensor.
Definition: SensorInfo.h:163
double getDepletionVoltage() const
Return the depletion voltage of the sensor.
Definition: SensorInfo.h:161
virtual unsigned short getBackgroundTag() const
Get background tag.
Definition: SimHitBase.h:46
bool registerInDataStore(DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut)
Register the object/array in the DataStore.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
T * appendNew()
Construct a new T object at the end of the array.
Definition: StoreArray.h:246
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
bool registerRelationTo(const StoreArray< TO > &toArray, DataStore::EDurability durability=DataStore::c_Event, DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut, const std::string &namedRelation="") const
Register a relation to the given StoreArray.
Definition: StoreArray.h:140
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
static const double mm
[millimeters]
Definition: Unit.h:70
static const double e
Standard of [electric charge].
Definition: Unit.h:53
static const double um
[micrometers]
Definition: Unit.h:71
static const double ns
Standard of [time].
Definition: Unit.h:48
float getGlobalTime() const override
Return the time of the electron deposition.
Definition: VXDSimHit.h:78
int getPDGcode() const
Return the PDG code of the particle causing the electron deposition.
Definition: VXDSimHit.h:68
ROOT::Math::XYZVector getPosOut() const
Return the end point of the electron deposition in local coordinates.
Definition: VXDSimHit.h:72
std::vector< std::pair< float, float > > getElectronsConstantDistance(double length) const
Get the electron deposition along constant stepsize.
Definition: VXDSimhit.cc:32
float getElectrons() const
Return the number of created electrons.
Definition: VXDSimhit.cc:15
ROOT::Math::XYZVector getPosIn() const
Return the start point of the electron deposition in local coordinates.
Definition: VXDSimHit.h:70
VxdID getSensorID() const
Return the sensorID of the sensor the electron was deposited in.
Definition: VXDSimHit.h:66
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:39
const std::set< Belle2::VxdID > getLayers(SensorInfoBase::SensorType sensortype=SensorInfoBase::VXD)
Return a set of all known Layers.
Definition: GeoCache.cc:176
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:67
const std::set< Belle2::VxdID > & getSensors(Belle2::VxdID ladder) const
Return a set of all sensor IDs belonging to a given ladder.
Definition: GeoCache.cc:204
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
const std::set< Belle2::VxdID > & getLadders(Belle2::VxdID layer) const
Return a set of all ladder IDs belonging to a given layer.
Definition: GeoCache.cc:193
double getUPitch(double v=0) const
Return the pitch of the sensor.
double getWidth(double v=0) const
Return the width of the sensor.
VxdID getID() const
Return the ID of the Sensor.
double getVPitch(double v=0) const
Return the pitch of the sensor.
double getThickness() const
Return the thickness of the sensor.
double getLength() const
Return the length of the sensor.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Namespace to encapsulate code needed for simulation and reconstrucion of the SVD.
Definition: GeoSVDCreator.h:23
double w_adjacentU(double t)
Adjacent-channel waveform U-side.
double w_betaprime(double t)
Beta-prime waveform shape, x^alpha/(1+x)^beta.
std::pair< StripWaveforms, StripWaveforms > SensorWaveforms
Waveforms of u- and v- channels in one sensor.
std::function< double(double)> WaveformShape
WaveformShape type.
double w_adjacentV(double t)
Adjacent-channel waveform V-side.
std::map< short int, SVDWaveform > StripWaveforms
Map of all channels' waveforms in one sensor side.
Abstract base class for different kinds of events.
STL namespace.
RelationElement::index_type indexFrom
index of the element from which the relation points.
RelationElement::weight_type weight
weight of the relation.