CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiStripDigitizerAlgorithm.cc
Go to the documentation of this file.
1 // File: SiStripDigitizerAlgorithm.cc
2 // Description: Steering class for digitization.
3 
4 // Modified 15/May/2013 mark.grimes@bristol.ac.uk - Modified so that the digi-sim link has the correct
5 // index for the sim hits stored. It was previously always set to zero (I won't mention that it was
6 // me who originally wrote that).
7 
8 #include <vector>
9 #include <algorithm>
10 #include <iostream>
27 #include "CLHEP/Random/RandFlat.h"
28 
30  lorentzAngleName(conf.getParameter<std::string>("LorentzAngle")),
31  theThreshold(conf.getParameter<double>("NoiseSigmaThreshold")),
32  cmnRMStib(conf.getParameter<double>("cmnRMStib")),
33  cmnRMStob(conf.getParameter<double>("cmnRMStob")),
34  cmnRMStid(conf.getParameter<double>("cmnRMStid")),
35  cmnRMStec(conf.getParameter<double>("cmnRMStec")),
36  APVSaturationProb(conf.getParameter<double>("APVSaturationProb")),
37  makeDigiSimLinks_(conf.getUntrackedParameter<bool>("makeDigiSimLinks", false)),
38  peakMode(conf.getParameter<bool>("APVpeakmode")),
39  noise(conf.getParameter<bool>("Noise")),
40  RealPedestals(conf.getParameter<bool>("RealPedestals")),
41  SingleStripNoise(conf.getParameter<bool>("SingleStripNoise")),
42  CommonModeNoise(conf.getParameter<bool>("CommonModeNoise")),
43  BaselineShift(conf.getParameter<bool>("BaselineShift")),
44  APVSaturationFromHIP(conf.getParameter<bool>("APVSaturationFromHIP")),
45  theFedAlgo(conf.getParameter<int>("FedAlgorithm")),
46  zeroSuppression(conf.getParameter<bool>("ZeroSuppression")),
47  theElectronPerADC(conf.getParameter<double>( peakMode ? "electronPerAdcPeak" : "electronPerAdcDec" )),
48  theTOFCutForPeak(conf.getParameter<double>("TOFCutForPeak")),
49  theTOFCutForDeconvolution(conf.getParameter<double>("TOFCutForDeconvolution")),
50  tofCut(peakMode ? theTOFCutForPeak : theTOFCutForDeconvolution),
51  cosmicShift(conf.getUntrackedParameter<double>("CosmicDelayShift")),
52  inefficiency(conf.getParameter<double>("Inefficiency")),
53  pedOffset((unsigned int)conf.getParameter<double>("PedestalsOffset")),
54  theSiHitDigitizer(new SiHitDigitizer(conf)),
55  theSiPileUpSignals(new SiPileUpSignals()),
56  theSiNoiseAdder(new SiGaussianTailNoiseAdder(theThreshold)),
57  theSiDigitalConverter(new SiTrivialDigitalConverter(theElectronPerADC)),
58  theSiZeroSuppress(new SiStripFedZeroSuppression(theFedAlgo)) {
59 
60  if (peakMode) {
61  LogDebug("StripDigiInfo")<<"APVs running in peak mode (poor time resolution)";
62  } else {
63  LogDebug("StripDigiInfo")<<"APVs running in deconvolution mode (good time resolution)";
64  };
65 }
66 
68 }
69 
70 void
72  edm::ESHandle<SiStripBadStrip> deadChannelHandle;
73  iSetup.get<SiStripBadChannelRcd>().get(deadChannelHandle);
74 
75  unsigned int detId = det->geographicalId().rawId();
76  int numStrips = (det->specificTopology()).nstrips();
77 
78  SiStripBadStrip::Range detBadStripRange = deadChannelHandle->getRange(detId);
79  //storing the bad strip of the the module. the module is not removed but just signal put to 0
80  std::vector<bool>& badChannels = allBadChannels[detId];
81  badChannels.clear();
82  badChannels.insert(badChannels.begin(), numStrips, false);
83  for(SiStripBadStrip::ContainerIterator it = detBadStripRange.first; it != detBadStripRange.second; ++it) {
84  SiStripBadStrip::data fs = deadChannelHandle->decode(*it);
85  for(int strip = fs.firstStrip; strip < fs.firstStrip + fs.range; ++strip) badChannels[strip] = true;
86  }
87  firstChannelsWithSignal[detId] = numStrips;
88  lastChannelsWithSignal[detId]= 0;
89 }
90 
91 void
93  theSiPileUpSignals->reset();
94  // This should be clear by after all calls to digitize(), but I might as well make sure
96 
97  //get gain noise pedestal lorentzAngle from ES handle
99  iSetup.getData(pdt);
100  setParticleDataTable(&*pdt);
102 }
103 
104 // Run the algorithm for a given module
105 // ------------------------------------
106 
107 void
108 SiStripDigitizerAlgorithm::accumulateSimHits(std::vector<PSimHit>::const_iterator inputBegin,
109  std::vector<PSimHit>::const_iterator inputEnd,
110  size_t inputBeginGlobalIndex,
111  unsigned int tofBin,
112  const StripGeomDetUnit* det,
113  const GlobalVector& bfield,
114  const TrackerTopology *tTopo,
115  CLHEP::HepRandomEngine* engine) {
116  // produce SignalPoints for all SimHits in detector
117  unsigned int detID = det->geographicalId().rawId();
118  int numStrips = (det->specificTopology()).nstrips();
119 
120  std::vector<bool>& badChannels = allBadChannels[detID];
121  size_t thisFirstChannelWithSignal = numStrips;
122  size_t thisLastChannelWithSignal = 0;
123 
124  float langle = (lorentzAngleHandle.isValid()) ? lorentzAngleHandle->getLorentzAngle(detID) : 0.;
125 
126  std::vector<float> locAmpl(numStrips, 0.);
127 
128  // Loop over hits
129 
130  uint32_t detId = det->geographicalId().rawId();
131  // First: loop on the SimHits
132  if(CLHEP::RandFlat::shoot(engine) > inefficiency) {
133  AssociationInfoForChannel* pDetIDAssociationInfo; // I only need this if makeDigiSimLinks_ is true...
134  if( makeDigiSimLinks_ ) pDetIDAssociationInfo=&(associationInfoForDetId_[detId]); // ...so only search the map if that is the case
135  std::vector<float> previousLocalAmplitude; // Only used if makeDigiSimLinks_ is true. Needed to work out the change in amplitude.
136 
137  size_t simHitGlobalIndex=inputBeginGlobalIndex; // This needs to stored to create the digi-sim link later
138  for (std::vector<PSimHit>::const_iterator simHitIter = inputBegin; simHitIter != inputEnd; ++simHitIter, ++simHitGlobalIndex ) {
139  // skip hits not in this detector.
140  if((*simHitIter).detUnitId() != detId) {
141  continue;
142  }
143  // check TOF
144  if (std::fabs(simHitIter->tof() - cosmicShift - det->surface().toGlobal(simHitIter->localPosition()).mag()/30.) < tofCut && simHitIter->energyLoss()>0) {
145  if( makeDigiSimLinks_ ) previousLocalAmplitude=locAmpl; // Not needed except to make the sim link association.
146  size_t localFirstChannel = numStrips;
147  size_t localLastChannel = 0;
148  // process the hit
149  theSiHitDigitizer->processHit(&*simHitIter, *det, bfield, langle, locAmpl, localFirstChannel, localLastChannel, tTopo, engine);
150 
151  //APV Killer to simulate HIP effect
152  //------------------------------------------------------
153 
155  int pdg_id = simHitIter->particleType();
156  particle = pdt->particle(pdg_id);
157  if(particle != NULL){
158  float charge = particle->charge();
159  bool isHadron = particle->isHadron();
160  if(charge!=0 && isHadron){
161  if(CLHEP::RandFlat::shoot(engine) < APVSaturationProb){
162  int FirstAPV = localFirstChannel/128;
163  int LastAPV = localLastChannel/128;
164  std::cout << "-------------------HIP--------------" << std::endl;
165  std::cout << "Killing APVs " << FirstAPV << " - " <<LastAPV << " " << detID <<std::endl;
166  for(int strip = FirstAPV*128; strip < LastAPV*128 +128; ++strip) {
167  badChannels[strip] = true;
168  }
169  //doing like that I remove the signal information only after the
170  //stip that got the HIP but it remains the signal of the previous
171  //one. I'll make a further loop to remove all signal
172  }
173  }
174  }
175  }
176 
177 
178  if(thisFirstChannelWithSignal > localFirstChannel) thisFirstChannelWithSignal = localFirstChannel;
179  if(thisLastChannelWithSignal < localLastChannel) thisLastChannelWithSignal = localLastChannel;
180 
181  if( makeDigiSimLinks_ ) { // No need to do any of this if truth association was turned off in the configuration
182  for( size_t stripIndex=0; stripIndex<locAmpl.size(); ++stripIndex ) {
183  // Work out the amplitude from this SimHit from the difference of what it was before and what it is now
184  float signalFromThisSimHit=locAmpl[stripIndex]-previousLocalAmplitude[stripIndex];
185  if( signalFromThisSimHit!=0 ) { // If this SimHit had any contribution I need to record it.
186  auto& associationVector=(*pDetIDAssociationInfo)[stripIndex];
187  bool addNewEntry=true;
188  // Make sure the hit isn't in already. I've seen this a few times, it always seems to happen in pairs so I think
189  // it's something to do with the stereo strips.
190  for( auto& associationInfo : associationVector ) {
191  if( associationInfo.trackID==simHitIter->trackId() && associationInfo.eventID==simHitIter->eventId() ) {
192  // The hit is already in, so add this second contribution and move on
193  associationInfo.contributionToADC+=signalFromThisSimHit;
194  addNewEntry=false;
195  break;
196  }
197  } // end of loop over associationVector
198  // If the hit wasn't already in create a new association info structure.
199  if( addNewEntry ) associationVector.push_back( AssociationInfo{ simHitIter->trackId(), simHitIter->eventId(), signalFromThisSimHit, simHitGlobalIndex, tofBin } );
200  } // end of "if( signalFromThisSimHit!=0 )"
201  } // end of loop over locAmpl strips
202  } // end of "if( makeDigiSimLinks_ )"
203  } // end of TOF check
204  } // end for
205  }
206  theSiPileUpSignals->add(detID, locAmpl, thisFirstChannelWithSignal, thisLastChannelWithSignal);
207 
208  if(firstChannelsWithSignal[detID] > thisFirstChannelWithSignal) firstChannelsWithSignal[detID] = thisFirstChannelWithSignal;
209  if(lastChannelsWithSignal[detID] < thisLastChannelWithSignal) lastChannelsWithSignal[detID] = thisLastChannelWithSignal;
210 }
211 
212 void
214  edm::DetSet<SiStripDigi>& outdigi,
215  edm::DetSet<SiStripRawDigi>& outrawdigi,
217  const StripGeomDetUnit *det,
218  edm::ESHandle<SiStripGain> & gainHandle,
219  edm::ESHandle<SiStripThreshold> & thresholdHandle,
220  edm::ESHandle<SiStripNoises> & noiseHandle,
221  edm::ESHandle<SiStripPedestals> & pedestalHandle,
222  CLHEP::HepRandomEngine* engine) {
223  unsigned int detID = det->geographicalId().rawId();
224  int numStrips = (det->specificTopology()).nstrips();
225 
226  const SiPileUpSignals::SignalMapType* theSignal(theSiPileUpSignals->getSignal(detID));
227 
228  std::vector<float> detAmpl(numStrips, 0.);
229  if(theSignal) {
230  for(const auto& amp : *theSignal) {
231  detAmpl[amp.first] = amp.second;
232  }
233  }
234 
235  //removing signal from the dead (and HIP effected) strips
236  std::vector<bool>& badChannels = allBadChannels[detID];
237  for(int strip =0; strip < numStrips; ++strip) if(badChannels[strip]) detAmpl[strip] = 0.;
238 
239  SiStripNoises::Range detNoiseRange = noiseHandle->getRange(detID);
240  SiStripApvGain::Range detGainRange = gainHandle->getRange(detID);
241  SiStripPedestals::Range detPedestalRange = pedestalHandle->getRange(detID);
242 
243 // -----------------------------------------------------------
244 
245  auto& firstChannelWithSignal = firstChannelsWithSignal[detID];
246  auto& lastChannelWithSignal = lastChannelsWithSignal[detID];
247  auto iAssociationInfoByChannel=associationInfoForDetId_.find(detID); // Use an iterator so that I can easily remove it once finished
248 
249  if(zeroSuppression){
250  if(noise){
251  int RefStrip = int(numStrips/2.);
252  while(RefStrip<numStrips&&badChannels[RefStrip]){ //if the refstrip is bad, I move up to when I don't find it
253  RefStrip++;
254  }
255  if(RefStrip<numStrips){
256  float noiseRMS = noiseHandle->getNoise(RefStrip,detNoiseRange);
257  float gainValue = gainHandle->getStripGain(RefStrip, detGainRange);
258  theSiNoiseAdder->addNoise(detAmpl,firstChannelWithSignal,lastChannelWithSignal,numStrips,noiseRMS*theElectronPerADC/gainValue, engine);
259  }
260  }
261  DigitalVecType digis;
262  theSiZeroSuppress->suppress(theSiDigitalConverter->convert(detAmpl, gainHandle, detID), digis, detID,noiseHandle,thresholdHandle);
263  // Now do the association to truth. Note that if truth association was turned off in the configuration this map
264  // will be empty and the iterator will always equal associationInfoForDetId_.end().
265  if( iAssociationInfoByChannel!=associationInfoForDetId_.end() ) { // make sure the readings for this DetID aren't completely from noise
266  for( const auto& iDigi : digis ) {
267  auto& associationInfoByChannel=iAssociationInfoByChannel->second;
268  const std::vector<AssociationInfo>& associationInfo=associationInfoByChannel[iDigi.channel()];
269 
270  // Need to find the total from all sim hits, because this might not be the same as the total
271  // digitised due to noise or whatever.
272  float totalSimADC=0;
273  for( const auto& iAssociationInfo : associationInfo ) totalSimADC+=iAssociationInfo.contributionToADC;
274  // Now I know that I can loop again and create the links
275  for( const auto& iAssociationInfo : associationInfo ) {
276  // Note simHitGlobalIndex used to have +1 because TrackerHitAssociator (the only place I can find this value being used)
277  // expected counting to start at 1, not 0. Now changed.
278  outLink.push_back( StripDigiSimLink( iDigi.channel(), iAssociationInfo.trackID, iAssociationInfo.simHitGlobalIndex, iAssociationInfo.tofBin, iAssociationInfo.eventID, iAssociationInfo.contributionToADC/totalSimADC ) );
279  } // end of loop over associationInfo
280  } // end of loop over the digis
281  } // end of check that iAssociationInfoByChannel is a valid iterator
282  outdigi.data = digis;
283  }
284 
285 
286  if(!zeroSuppression){
287  //if(noise){
288  // the constant pedestal offset is needed because
289  // negative adc counts are not allowed in case
290  // Pedestal and CMN subtraction is performed.
291  // The pedestal value read from the conditions
292  // is pedValue and after the pedestal subtraction
293  // the baseline is zero. The Common Mode Noise
294  // is not subtracted from the negative adc counts
295  // channels. Adding pedOffset the baseline is set
296  // to pedOffset after pedestal subtraction and CMN
297  // is subtracted to all the channels since none of
298  // them has negative adc value. The pedOffset is
299  // treated as a constant component in the CMN
300  // estimation and subtracted as CMN.
301 
302 
303  //calculating the charge deposited on each APV and subtracting the shift
304  //------------------------------------------------------
305  if(BaselineShift){
306  theSiNoiseAdder->addBaselineShift(detAmpl, badChannels);
307  }
308 
309  //Adding the strip noise
310  //------------------------------------------------------
311  if(noise){
312  std::vector<float> noiseRMSv;
313  noiseRMSv.clear();
314  noiseRMSv.insert(noiseRMSv.begin(),numStrips,0.);
315 
316  if(SingleStripNoise){
317  for(int strip=0; strip< numStrips; ++strip){
318  if(!badChannels[strip]) noiseRMSv[strip] = (noiseHandle->getNoise(strip,detNoiseRange))* theElectronPerADC;
319  }
320 
321  } else {
322  int RefStrip = 0; //int(numStrips/2.);
323  while(RefStrip<numStrips&&badChannels[RefStrip]){ //if the refstrip is bad, I move up to when I don't find it
324  RefStrip++;
325  }
326  if(RefStrip<numStrips){
327  float noiseRMS = noiseHandle->getNoise(RefStrip,detNoiseRange) *theElectronPerADC;
328  for(int strip=0; strip< numStrips; ++strip){
329  if(!badChannels[strip]) noiseRMSv[strip] = noiseRMS;
330  }
331  }
332  }
333 
334  theSiNoiseAdder->addNoiseVR(detAmpl, noiseRMSv, engine);
335  }
336 
337  //adding the CMN
338  //------------------------------------------------------
339  if(CommonModeNoise){
340  float cmnRMS = 0.;
341  DetId detId(detID);
342  uint32_t SubDet = detId.subdetId();
343  if(SubDet==3){
344  cmnRMS = cmnRMStib;
345  }else if(SubDet==4){
346  cmnRMS = cmnRMStid;
347  }else if(SubDet==5){
348  cmnRMS = cmnRMStob;
349  }else if(SubDet==6){
350  cmnRMS = cmnRMStec;
351  }
352  cmnRMS *= theElectronPerADC;
353  theSiNoiseAdder->addCMNoise(detAmpl, cmnRMS, badChannels, engine);
354  }
355 
356 
357  //Adding the pedestals
358  //------------------------------------------------------
359 
360  std::vector<float> vPeds;
361  vPeds.clear();
362  vPeds.insert(vPeds.begin(),numStrips,0.);
363 
364  if(RealPedestals){
365  for(int strip=0; strip< numStrips; ++strip){
366  if(!badChannels[strip]) vPeds[strip] = (pedestalHandle->getPed(strip,detPedestalRange)+pedOffset)* theElectronPerADC;
367  }
368  } else {
369  for(int strip=0; strip< numStrips; ++strip){
370  if(!badChannels[strip]) vPeds[strip] = pedOffset* theElectronPerADC;
371  }
372  }
373 
374  theSiNoiseAdder->addPedestals(detAmpl, vPeds);
375 
376 
377  //if(!RealPedestals&&!CommonModeNoise&&!noise&&!BaselineShift&&!APVSaturationFromHIP){
378  // edm::LogWarning("SiStripDigitizer")<<"You are running the digitizer without Noise generation and without applying Zero Suppression. ARE YOU SURE???";
379  //}else{
380 
381  DigitalRawVecType rawdigis = theSiDigitalConverter->convertRaw(detAmpl, gainHandle, detID);
382 
383  // Now do the association to truth. Note that if truth association was turned off in the configuration this map
384  // will be empty and the iterator will always equal associationInfoForDetId_.end().
385  if( iAssociationInfoByChannel!=associationInfoForDetId_.end() ) { // make sure the readings for this DetID aren't completely from noise
386  // N.B. For the raw digis the channel is inferred from the position in the vector.
387  // I'VE NOT TESTED THIS YET!!!!!
388  // ToDo Test this properly.
389  for( size_t channel=0; channel<rawdigis.size(); ++channel ) {
390  auto& associationInfoByChannel=iAssociationInfoByChannel->second;
391  const auto iAssociationInfo=associationInfoByChannel.find(channel);
392  if( iAssociationInfo==associationInfoByChannel.end() ) continue; // Skip if there is no sim information for this channel (i.e. it's all noise)
393  const std::vector<AssociationInfo>& associationInfo=iAssociationInfo->second;
394 
395  // Need to find the total from all sim hits, because this might not be the same as the total
396  // digitised due to noise or whatever.
397  float totalSimADC=0;
398  for( const auto& iAssociationInfo : associationInfo ) totalSimADC+=iAssociationInfo.contributionToADC;
399  // Now I know that I can loop again and create the links
400  for( const auto& iAssociationInfo : associationInfo ) {
401  // Note simHitGlobalIndex used to have +1 because TrackerHitAssociator (the only place I can find this value being used)
402  // expected counting to start at 1, not 0. Now changed.
403  outLink.push_back( StripDigiSimLink( channel, iAssociationInfo.trackID, iAssociationInfo.simHitGlobalIndex, iAssociationInfo.tofBin, iAssociationInfo.eventID, iAssociationInfo.contributionToADC/totalSimADC ) );
404  } // end of loop over associationInfo
405  } // end of loop over the digis
406  } // end of check that iAssociationInfoByChannel is a valid iterator
407 
408  outrawdigi.data = rawdigis;
409 
410  //}
411  }
412 
413  // Now that I've finished with this entry in the map of associations, I can remove it.
414  // Note that there might not be an association if the ADC reading is from noise in which
415  // case associationIsValid will be false.
416  if( iAssociationInfoByChannel!=associationInfoForDetId_.end() ) associationInfoForDetId_.erase(iAssociationInfoByChannel);
417 }
#define LogDebug(id)
unsigned short range
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:114
const std::unique_ptr< SiPileUpSignals > theSiPileUpSignals
edm::ESHandle< SiStripLorentzAngle > lorentzAngleHandle
void push_back(const T &t)
Definition: DetSet.h:68
AssociationInfoForDetId associationInfoForDetId_
Structure that holds the information on the SimTrack contributions. Only filled if makeDigiSimLinks_ ...
SiDigitalConverter::DigitalRawVecType DigitalRawVecType
SiDigitalConverter::DigitalVecType DigitalVecType
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
std::vector< unsigned int >::const_iterator ContainerIterator
void initializeEvent(const edm::EventSetup &iSetup)
#define NULL
Definition: scimark2.h:8
std::pair< ContainerIterator, ContainerIterator > Range
void initializeDetUnit(StripGeomDetUnit const *det, const edm::EventSetup &iSetup)
void digitize(edm::DetSet< SiStripDigi > &outDigis, edm::DetSet< SiStripRawDigi > &outRawDigis, edm::DetSet< StripDigiSimLink > &outLink, const StripGeomDetUnit *stripdet, edm::ESHandle< SiStripGain > &, edm::ESHandle< SiStripThreshold > &, edm::ESHandle< SiStripNoises > &, edm::ESHandle< SiStripPedestals > &, CLHEP::HepRandomEngine *)
virtual const StripTopology & specificTopology() const
Returns a reference to the strip proxy topology.
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:40
const std::unique_ptr< const SiGaussianTailNoiseAdder > theSiNoiseAdder
const std::unique_ptr< SiHitDigitizer > theSiHitDigitizer
void getData(T &iHolder) const
Definition: EventSetup.h:78
std::map< unsigned int, size_t > firstChannelsWithSignal
std::map< int, Amplitude > SignalMapType
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
const std::unique_ptr< SiTrivialDigitalConverter > theSiDigitalConverter
const std::unique_ptr< SiStripFedZeroSuppression > theSiZeroSuppress
std::pair< ContainerIterator, ContainerIterator > Range
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:77
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
std::map< unsigned int, std::vector< bool > > allBadChannels
tuple conf
Definition: dbtoconf.py:185
void setParticleDataTable(const ParticleDataTable *pardt)
Definition: DetId.h:18
unsigned short firstStrip
std::map< int, std::vector< AssociationInfo > > AssociationInfoForChannel
const T & get() const
Definition: EventSetup.h:55
SiStripDigitizerAlgorithm(const edm::ParameterSet &conf)
collection_type data
Definition: DetSet.h:78
std::pair< ContainerIterator, ContainerIterator > Range
std::map< unsigned int, size_t > lastChannelsWithSignal
void accumulateSimHits(const std::vector< PSimHit >::const_iterator inputBegin, const std::vector< PSimHit >::const_iterator inputEnd, size_t inputBeginGlobalIndex, unsigned int tofBin, const StripGeomDetUnit *stripdet, const GlobalVector &bfield, const TrackerTopology *tTopo, CLHEP::HepRandomEngine *)
tuple cout
Definition: gather_cfg.py:121
volatile std::atomic< bool > shutdown_flag false
std::pair< ContainerIterator, ContainerIterator > Range
Definition: SiStripNoises.h:48
bool isValid() const
Definition: ESHandle.h:47
const ParticleDataTable * pdt