OptDetDigitizer_module.cc
Go to the documentation of this file.
1 // OptDetDigitizer_module.cc
2 // Kazuhiro Terao <kazuhiro@nevis.columbia.edu>, Jul 2013
3 // based on code by Ben Jones and Christie Chiu, MIT, Sept 2012
4 // bjpjones@mit.edu, cschiu@mit.edu
5 //
6 // This module starts from MC truth sim::OnePhoton objects
7 // and produces a digitized waveform.
8 
9 // LArSoft includes
17 
18 // ART includes
22 #include "fhiclcpp/ParameterSet.h"
24 
25 // nurandom
26 #include "nurandom/RandomUtils/NuRandomService.h"
27 
28 // CLHEP includes
29 #include "CLHEP/Random/RandFlat.h"
30 #include "CLHEP/Random/RandPoisson.h"
31 
32 // C++ language includes
33 #include <cstring>
34 
35 namespace opdet {
36 
38  public:
39  explicit OptDetDigitizer(const fhicl::ParameterSet&);
40 
41  private:
42  void produce(art::Event&) override;
43 
44  // The parameters we'll read from the .fcl file.
45  std::string fInputModule; // Input tag for OpDet collection
46  float fSampleFreq; // in MHz
47  float fTimeBegin; // in us
48  float fTimeEnd; // in us
49  float fQE; // quantum efficiency of opdet
50  optdata::ADC_Count_t fSaturationScale; // adc count w/ saturation occurs
51  std::vector<optdata::ADC_Count_t> fPedMeanArray; // Array of pedestal baseline (per ch)
52  float fDarkRate; // Noise rate in Hz
53  optdata::ADC_Count_t fPedFlucAmp; // Pedestal fluctuation amplitude
54  float fPedFlucRate; // Pedestal fluctuation rate
55  // float fWFRandTimeOffsetLow; // The lower bound of WF's T=0 offset from Trigger
56  // float fWFRandTimeOffsetHigh; // The upper bound of WF's T=0 offset from Trigger
57  std::vector<double> fSinglePEWaveform;
58 
60 
61  CLHEP::HepRandomEngine& fEngine;
62  CLHEP::RandFlat fFlatRandom;
63  CLHEP::RandPoisson fPoissonRandom;
64  void AddDarkNoise (std::vector<double> &RawWF,double gain);
66  std::vector<double>& OldPulse,
67  std::vector<double>& NewPulse,
68  double factor,
69  bool extend=false);
70  optdata::ChannelData ApplyDigitization (std::vector<double> const RawWF,
71  optdata::Channel_t const ch) const;
74 
75  };
76 } // namespace opdet
77 
78 namespace opdet{
79 
81 
82 } //end namespace opdet
83 
84 namespace opdet {
85 
87  : EDProducer{pset}
91  {
92  // Infrastructure piece
93  produces<std::vector< optdata::ChannelDataGroup> >();
94 
96  // Input Module and histogram parameters come from .fcl
97  fInputModule = pset.get<std::string>("InputModule");
98  fSimGainSpread = pset.get<bool >("SimGainSpread");
99  fTimeBegin = fOpDigiProperties->TimeBegin();
100  fTimeEnd = fOpDigiProperties->TimeEnd();
101  fSampleFreq = fOpDigiProperties->SampleFreq();
102  fQE = fOpDigiProperties->QE();
103  fDarkRate = fOpDigiProperties->DarkRate();
104  fPedFlucAmp = fOpDigiProperties->PedFlucAmp();
105  fPedFlucRate= fOpDigiProperties->PedFlucRate();
106  fSaturationScale = fOpDigiProperties->SaturationScale();
107  fPedMeanArray = fOpDigiProperties->PedMeanArray();
108 
109  fSinglePEWaveform = fOpDigiProperties->SinglePEWaveform();
110  }
111 
112  //-------------------------------------------------
113 
115  std::vector<double> &OldPulse,
116  std::vector<double> &NewPulse,
117  double const factor,
118  bool const extend)
119  {
120  if( (time+NewPulse.size()) > OldPulse.size() && extend )
121  OldPulse.resize(time + NewPulse.size());
122  for(size_t i = 0; i<NewPulse.size() && (time+i)<OldPulse.size(); ++i)
123  OldPulse[time+i] += NewPulse[i] * factor;
124  }
125 
126  //-------------------------------------------------
127 
128  void OptDetDigitizer::AddDarkNoise(std::vector<double> &RawWF, double gain){
129  // Add dark noise
130  double MeanDarkPulses = fDarkRate * (fTimeEnd-fTimeBegin) / 1000000;
131 
132  unsigned int NumberOfPulses = fPoissonRandom.fire(MeanDarkPulses);
133  for(size_t i=0; i!=NumberOfPulses; ++i)
134  {
135  double PulseTime_ns = fTimeBegin*1000 + (fTimeEnd-fTimeBegin)*1000*(fFlatRandom.fire(1.0)); // Should be in ns
136  optdata::TimeSlice_t PulseTime_ts = fOpDigiProperties->GetTimeSlice(PulseTime_ns);
137  AddWaveform( PulseTime_ts,
138  RawWF,
140  gain);
141  }
142 
143  }
144 
146  optdata::Channel_t const ch ) const
147  {
148  //
149  // Digitization includes...
150  // (a) amplitude digitization
151  // (b) saturation
152  // (c) pedestal fluctuation
153  //
154 
155  // prepare return data container
156  optdata::ChannelData chData(ch);
157  chData.reserve(rawWF.size());
158  optdata::ADC_Count_t baseMean(fPedMeanArray.at(ch));
159  for(optdata::TimeSlice_t time=0; time<rawWF.size(); ++time)
160  {
161  double thisSample = rawWF[time];
162 
163  optdata::ADC_Count_t thisCount = (optdata::ADC_Count_t)(thisSample)+baseMean;
164 
165  // (a) amplitude digitization
166  if(CLHEP::RandFlat::shoot(1.0) < (thisSample - int(thisSample)))
167  thisCount+=1;
168 
169  // (b) saturation
170  if(thisCount > fSaturationScale) thisCount = fSaturationScale;
171 
172  chData.push_back(thisCount);
173  }
174 
175  // (c) pedestal fluctuation
176  double timeSpan = chData.size() * 1.e-6/(fOpDigiProperties->SampleFreq());
177  unsigned int nFluc = CLHEP::RandPoisson::shoot(fPedFlucRate * timeSpan);
178  for(size_t i=0; i<nFluc; ++i)
179  {
180  optdata::TimeSlice_t pulseTime(CLHEP::RandFlat::shoot(0.0,(double)(chData.size())));
181  optdata::ADC_Count_t amp = chData[pulseTime];
182  if( CLHEP::RandFlat::shoot(0.,1.) > 0.5)
183  {
184  amp += fPedFlucAmp;
185  if(amp > fSaturationScale) amp=fSaturationScale;
186  }
187  else
188  amp -= fPedFlucAmp;
189  chData[pulseTime] = amp;
190  }
191 
192  return chData;
193  }
194 
195  //-------------------------------------------------
196 
198  {
199 
200  //
201  // Event-wise initialization
202  //
203 
204  // Infrastructure piece
205  std::unique_ptr< std::vector<optdata::ChannelDataGroup > > StoragePtr (new std::vector<optdata::ChannelDataGroup>);
206 
207  // Read in the Sim Photons
209 
210  // Convert units into ns from us/MHz
211  double timeBegin_ns = fTimeBegin * 1000;
212  double timeEnd_ns = fTimeEnd * 1000;
213  double sampleFreq_ns = fSampleFreq / 1000;
214 
215  // Compute # of timeslices to be stored in the output. This is defined by a user input (fcl file)
216  optdata::TimeSlice_t timeSliceWindow(fOpDigiProperties->GetTimeSlice(timeEnd_ns));
217 
218  /*
219  Create output data product, optdata::ChannelDataGroup for each gain channel.
220  Note : Although the frame + sample number in DATA should have a reference of T=0 @ DAQ start time, this is
221  not handled in MC. Hence we do not set them here (use constructor default)
222  */
223  optdata::ChannelDataGroup rawWFGroup_HighGain(optdata::kHighGain);
224  optdata::ChannelDataGroup rawWFGroup_LowGain(optdata::kLowGain);
225  // Reserve entries equal to # of channels
226  rawWFGroup_HighGain.reserve(fGeom->NOpChannels());
227  rawWFGroup_LowGain.reserve(fGeom->NOpChannels());
228 
229  /*
230  Define "raw" waveform container which will be filled based on G4 photon timing + SPE waveform information.
231  Note this is not completely an analog waveform because it is digitized in terms of time (as it is using std::vector).
232  */
233  std::vector<std::vector<double> > rawWF_HighGain(fGeom->NOpChannels(),std::vector<double>(timeSliceWindow,0.0));
234  std::vector<std::vector<double> > rawWF_LowGain(fGeom->NOpChannels(),std::vector<double>(timeSliceWindow,0.0));
235 
236  /*
237  Start data processing ... see following steps
238  (1) Loop over input array of optical photons & fill "raw" waveform container w/ corresponding SPE waveform
239  (2) Loop over filled "raw" waveform and process (digitization, adding noise, baseline spread, etc)
240  */
241 
242  //
243  // Step (1) ... loop over G4 optical photons
244  //
245 
246  // For every OpDet, convert PE into waveform and combine all together
247  for(sim::SimPhotonsCollection::const_iterator itOpDet=ThePhotCollection.begin(); itOpDet!=ThePhotCollection.end(); itOpDet++)
248  {
249  const sim::SimPhotons& ThePhot=itOpDet->second;
250 
251  int ch = ThePhot.OpChannel();
252  // For every photon in the hit:
253  for(const sim::OnePhoton& Phot: ThePhot)
254  {
255  // Sample a random subset according to QE
256  if(fFlatRandom.fire(1.0)<=fQE)
257  {
258  optdata::TimeSlice_t PhotonTime(fOpDigiProperties->GetTimeSlice(Phot.Time));
259  if( Phot.Time > timeBegin_ns && Phot.Time < timeEnd_ns )
260  {
261  if(fSimGainSpread)
262  {
263  AddWaveform( PhotonTime, rawWF_HighGain[ch], fSinglePEWaveform, fOpDigiProperties->HighGain(ch));
264  AddWaveform( PhotonTime, rawWF_LowGain[ch], fSinglePEWaveform, fOpDigiProperties->LowGain(ch));
265  }
266  else
267  {
268  AddWaveform( PhotonTime, rawWF_HighGain[ch], fSinglePEWaveform, fOpDigiProperties->HighGainMean(ch));
269  AddWaveform( PhotonTime, rawWF_LowGain[ch], fSinglePEWaveform, fOpDigiProperties->LowGainMean(ch));
270  }
271  }
272  } // random QE cut
273  } // for each Photon in SimPhotons
274  }
275 
276  //
277  // Loop over "raw" waveform (channel-wise)
278  //
279  for(unsigned short iCh = 0; iCh < rawWF_LowGain.size(); ++iCh){
280  rawWF_LowGain[iCh].resize((timeEnd_ns - timeBegin_ns) * sampleFreq_ns);
281  rawWF_HighGain[iCh].resize((timeEnd_ns - timeBegin_ns) * sampleFreq_ns);
282 
283  // Add dark noise
284  if(fSimGainSpread){
285  AddDarkNoise(rawWF_LowGain[iCh],fOpDigiProperties->LowGain(iCh));
286  AddDarkNoise(rawWF_HighGain[iCh],fOpDigiProperties->HighGain(iCh));
287  }else{
288  AddDarkNoise(rawWF_LowGain[iCh],fOpDigiProperties->LowGainMean(iCh));
289  AddDarkNoise(rawWF_HighGain[iCh],fOpDigiProperties->HighGainMean(iCh));
290  }
291 
292  // Apply digitization and make channel data
293  optdata::ChannelData chData_HighGain(ApplyDigitization(rawWF_HighGain[iCh],iCh));
294  optdata::ChannelData chData_LowGain(ApplyDigitization(rawWF_LowGain[iCh],iCh));
295 
296  rawWFGroup_HighGain.push_back(chData_HighGain);
297  rawWFGroup_LowGain.push_back(chData_LowGain);
298  } // for each OpDet in SimPhotonsCollection
299 
300  StoragePtr->push_back(rawWFGroup_HighGain);
301  StoragePtr->push_back(rawWFGroup_LowGain);
302 
303  evt.put(std::move(StoragePtr));
304  }
305 }
art::ServiceHandle< geo::Geometry const > fGeom
optdata::ChannelData ApplyDigitization(std::vector< double > const RawWF, optdata::Channel_t const ch) const
base_engine_t & createEngine(seed_t seed)
optdata::ADC_Count_t fSaturationScale
int OpChannel() const
Returns the optical channel number this object is associated to.
Definition: SimPhotons.h:254
std::string string
Definition: nybbler.cc:12
EDProducer(fhicl::ParameterSet const &pset)
Definition: EDProducer.h:20
All information of a photon entering the sensitive optical detector volume.
Definition: SimPhotons.h:64
OptDetDigitizer(const fhicl::ParameterSet &)
unsigned int NOpChannels() const
Number of electronics channels for all the optical detectors.
void produce(art::Event &) override
art framework interface to geometry description
uint16_t ADC_Count_t
Definition: OpticalTypes.h:16
Simulation objects for optical detectors.
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:67
def move(depos, offset)
Definition: depos.py:107
optdata::ADC_Count_t fPedFlucAmp
ProductID put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: DataViewImpl.h:686
void AddWaveform(optdata::TimeSlice_t time, std::vector< double > &OldPulse, std::vector< double > &NewPulse, double factor, bool extend=false)
CLHEP::RandPoisson fPoissonRandom
art::ServiceHandle< OpDigiProperties > fOpDigiProperties
Collection of photons which recorded on one channel.
Definition: SimPhotons.h:136
list_type::const_iterator const_iterator
Definition: SimPhotons.h:206
unsigned int TimeSlice_t
Definition: OpticalTypes.h:20
std::vector< double > fSinglePEWaveform
std::vector< optdata::ADC_Count_t > fPedMeanArray
TCEvent evt
Definition: DataStructs.cxx:7
unsigned int Channel_t
Definition: OpticalTypes.h:19
CLHEP::HepRandomEngine & fEngine
Collection of sim::SimPhotons, indexed by channel number.
Definition: SimPhotons.h:192
void AddDarkNoise(std::vector< double > &RawWF, double gain)
static sim::SimPhotonsCollection GetSimPhotonsCollection(const art::Event &evt, std::string moduleLabel)