RawWaveformClnSigDump_module.cc
Go to the documentation of this file.
1 //////////////////////////////////////////////
2 //
3 // Module to dump raw data for ROI training
4 //
5 // mwang@fnal.gov
6 // tjyang@fnal.gov
7 // wwu@fnal.gov
8 //
9 // This version:
10 // a) uses only RawDigits, not WireProducer
11 // b) reads in two separate RawDigits:
12 // - signal+noise
13 // - same signal as above but no noise
14 // c) saves full waveform or short waveform
15 // - for short waveform, picks signal
16 // associated with largest energy
17 // deposit & saves only tdcmin/tdcmaxes
18 // that are within window extents
19 // d) selects only signal waveforms that have
20 // minimum ADC count associated with the
21 // pure signal
22 // e) produces a 2nd npy file for pure signal
23 // waveform
24 //
25 //////////////////////////////////////////////
26 
27 #include <random>
28 
29 // Framework includes
37 #include "canvas/Persistency/Common/FindManyP.h"
39 #include "fhiclcpp/ParameterSet.h"
41 
42 // LArSoft libraries
45 #include "larcoreobj/SimpleTypesAndConstants/RawTypes.h" // raw::ChannelID_t
49 #include "lardataobj/RawData/raw.h"
56 
59 
60 #include "nurandom/RandomUtils/NuRandomService.h"
61 #include "CLHEP/Random/RandFlat.h"
62 #include "c2numpy.h"
63 
64 using std::cout;
65 using std::endl;
66 using std::ofstream;
67 using std::string;
68 
69 struct WireSigInfo {
70  int pdgcode;
73  unsigned int tdcmin;
74  unsigned int tdcmax;
75  int tdcpeak;
76  int adcpeak;
77  int numel;
78  double edep;
79 };
80 
81 namespace nnet {
82  class RawWaveformClnSigDump;
83 }
84 
86 
87 public:
89 
90  // Plugins should not be copied or assigned.
93  RawWaveformClnSigDump& operator=(RawWaveformClnSigDump const&) = delete;
94  RawWaveformClnSigDump& operator=(RawWaveformClnSigDump&&) = delete;
95 
96  // Required functions.
97  void analyze(art::Event const& e) override;
98 
99  //void reconfigure(fhicl::ParameterSet const & p);
100 
101  void beginJob() override;
102  void endJob() override;
103 
104 private:
107 
108  std::string fSimulationProducerLabel; ///< producer that tracked simulated part. through detector
109  std::string fSimChannelLabel; ///< module that made simchannels
110  std::string fDigitModuleLabel; ///< module that made digits
111  std::string fCleanSignalDigitModuleLabel; ///< module that made the signal-only digits
113  unsigned int fShortWaveformSize;
116 
131 
132  CLHEP::RandFlat fRandFlat;
133 
136 };
137 
138 //-----------------------------------------------------------------------
139 struct genFinder {
140 private:
141  typedef std::pair<int, std::string> track_id_to_string;
142  std::vector<track_id_to_string> track_id_map;
143  std::set<std::string> generator_names;
144  bool isSorted = false;
145 
146 public:
147  void
149  {
150  std::sort(this->track_id_map.begin(),
151  this->track_id_map.end(),
152  [](const auto& a, const auto& b) { return (a.first < b.first); });
153  isSorted = true;
154  }
155  void
156  add(const int& track_id, const std::string& gname)
157  {
158  this->track_id_map.push_back(std::make_pair(track_id, gname));
159  generator_names.emplace(gname);
160  isSorted = false;
161  }
162  bool
164  {
165  return static_cast<bool>(generator_names.count(gname));
166  };
168  get_gen(int tid)
169  {
170  if (!isSorted) { this->sort_now(); }
171  return std::lower_bound(track_id_map.begin(),
172  track_id_map.end(),
173  tid,
174  [](const auto& a, const auto& b) { return (a.first < b); })
175  ->second;
176  };
177 };
178 
179 // Create the random number generator
180 namespace {
181  std::string const instanceName = "RawWaveformClnSigDump";
182 }
183 
184 //-----------------------------------------------------------------------
186  : EDAnalyzer{p}
187  , fDumpWaveformsFileName(p.get<std::string>("DumpWaveformsFileName", "dumpwaveforms"))
188  , fDumpCleanSignalFileName(p.get<std::string>("CleanSignalFileName", "dumpcleansignal"))
189  , fSimulationProducerLabel(p.get<std::string>("SimulationProducerLabel", "larg4Main"))
190  , fSimChannelLabel(p.get<std::string>("SimChannelLabel", "elecDrift"))
191  , fDigitModuleLabel(p.get<std::string>("DigitModuleLabel", "simWire"))
192  , fCleanSignalDigitModuleLabel(p.get<std::string>("CleanSignalDigitModuleLabel", "simWire:signal"))
193  , fUseFullWaveform(p.get<bool>("UseFullWaveform", true))
194  , fShortWaveformSize(p.get<unsigned int>("ShortWaveformSize"))
195  , fEstIndFWForOffset(p.get<int>("EstIndFWForOffset"))
196  , fEstColFWForOffset(p.get<int>("EstIndFWForOffset"))
197  , fSelectGenLabel(p.get<std::string>("SelectGenLabel", "ANY"))
198  , fSelectProcID(p.get<std::string>("SelectProcID", "ANY"))
199  , fSelectPDGCode(p.get<int>("SelectPDGCode", 0))
200  , fPlaneToDump(p.get<std::string>("PlaneToDump"))
201  , fMinParticleEnergyGeV(p.get<double>("MinParticleEnergyGeV", 0.))
202  , fMinEnergyDepositedMeV(p.get<double>("MinEnergyDepositedMeV", 0.))
203  , fMinNumberOfElectrons(p.get<int>("MinNumberOfElectrons", 600))
204  , fMaxNumberOfElectrons(p.get<int>("MaxNumberOfElectrons", -1))
205  , fMinPureSignalADCs(p.get<int>("MinPureSignalADCs", 0))
206  , fSaveSignal(p.get<bool>("SaveSignal", true))
207  , fMaxNoiseChannelsPerEvent(p.get<int>("MaxNoiseChannelsPerEvent"))
208  , fCollectionPlaneLabel(p.get<std::string>("CollectionPlaneLabel"))
210  instanceName, p, "SeedForRawWaveformDump"),
211  "HepJamesRandom", instanceName)}
212 {
213  if (std::getenv("CLUSTER") && std::getenv("PROCESS")) {
214  fDumpWaveformsFileName += string(std::getenv("CLUSTER")) + "-" + string(std::getenv("PROCESS")) + "-";
215  fDumpCleanSignalFileName += string(std::getenv("CLUSTER")) + "-" + string(std::getenv("PROCESS")) + "-";
216  }
217 
218  if (fDigitModuleLabel.empty() && fCleanSignalDigitModuleLabel.empty()) {
219  throw cet::exception("RawWaveformClnSigDump")
220  << "Both DigitModuleLabel and CleanSignalModuleLabel are empty";
221  }
222 }
223 
224 //-----------------------------------------------------------------------
225 void
227 {
228  auto const detProp = art::ServiceHandle<detinfo::DetectorPropertiesService const>()->DataForJob();
229 
235 
236  for (unsigned int i = 0; i < 5; i++) {
237  std::ostringstream name;
238 
239  name.str("");
240  name << "tid" << i;
241  c2numpy_addcolumn(&npywriter, name.str().c_str(), C2NUMPY_INT32);
242 
243  name.str("");
244  name << "pdg" << i;
245  c2numpy_addcolumn(&npywriter, name.str().c_str(), C2NUMPY_INT32);
246 
247  name.str("");
248  name << "gen" << i;
249  c2numpy_addcolumn(&npywriter, name.str().c_str(), (c2numpy_type)((int)C2NUMPY_STRING + 6));
250 
251  name.str("");
252  name << "pid" << i;
253  c2numpy_addcolumn(&npywriter, name.str().c_str(), (c2numpy_type)((int)C2NUMPY_STRING + 7));
254 
255  name.str("");
256  name << "edp" << i;
257  c2numpy_addcolumn(&npywriter, name.str().c_str(), C2NUMPY_FLOAT32);
258 
259  name.str("");
260  name << "nel" << i;
261  c2numpy_addcolumn(&npywriter, name.str().c_str(), C2NUMPY_UINT32);
262 
263  name.str("");
264  name << "sti" << i;
265  c2numpy_addcolumn(&npywriter, name.str().c_str(), C2NUMPY_UINT16);
266 
267  name.str("");
268  name << "stf" << i;
269  c2numpy_addcolumn(&npywriter, name.str().c_str(), C2NUMPY_UINT16);
270 
271  name.str("");
272  name << "stp" << i;
273  c2numpy_addcolumn(&npywriter, name.str().c_str(), C2NUMPY_INT32);
274 
275  name.str("");
276  name << "adc" << i;
277  c2numpy_addcolumn(&npywriter, name.str().c_str(), C2NUMPY_INT32);
278  }
279 
280  for (unsigned int i = 0;
281  i < (fUseFullWaveform ? detProp.ReadOutWindowSize() : fShortWaveformSize);
282  i++) {
283  std::ostringstream name;
284  name << "tck_" << i;
285  c2numpy_addcolumn(&npywriter, name.str().c_str(), C2NUMPY_INT16);
286  }
287 
288  // ... this is for storing the clean signal (no noise) waveform
290 
291  for (unsigned int i = 0;
292  i < (fUseFullWaveform ? detProp.ReadOutWindowSize() : fShortWaveformSize);
293  i++) {
294  std::ostringstream name;
295  name << "tck_" << i;
296  c2numpy_addcolumn(&npywriter2, name.str().c_str(), C2NUMPY_INT16);
297  }
298 }
299 
300 //-----------------------------------------------------------------------
301 void
303 {
306 }
307 
308 //-----------------------------------------------------------------------
309 void
311 {
312  cout << "Event "
313  << " " << evt.id().run() << " " << evt.id().subRun() << " " << evt.id().event() << endl;
314 
315  std::unique_ptr<genFinder> gf(new genFinder());
316 
317  // ... Read in the digit List object(s).
319  std::vector<art::Ptr<raw::RawDigit>> rawdigitlist;
320  if (evt.getByLabel(fDigitModuleLabel, digitVecHandle)) {
321  //std::cout << " !!!! RawWaveformClnSigDump: fDigitModuleLabel -> " << fDigitModuleLabel << std::endl;
322  art::fill_ptr_vector(rawdigitlist, digitVecHandle);
323  }
324 
325  // ... Read in the signal-only digit List object(s).
327  std::vector<art::Ptr<raw::RawDigit>> rawdigitlist2;
328  if (evt.getByLabel(fCleanSignalDigitModuleLabel, digitVecHandle2)) {
329  //std::cout << " !!!! RawWaveformClnSigDump: fCleanSignalDigitModuleLabel -> " << fCleanSignalDigitModuleLabel << std::endl;
330  art::fill_ptr_vector(rawdigitlist2, digitVecHandle2);
331  }
332 
333  if (rawdigitlist.empty() && rawdigitlist2.empty()) return;
334 
335  auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(evt);
336  auto const detProp =
338 
339  // ... Use the handle to get a particular (0th) element of collection.
340  unsigned int dataSize;
341  art::Ptr<raw::RawDigit> digitVec0(digitVecHandle, 0);
342  dataSize = digitVec0->Samples(); //size of raw data vectors
343  if (dataSize != detProp.ReadOutWindowSize()) {
344  throw cet::exception("RawWaveformClnSigDump") << "Bad dataSize: " << dataSize;
345  }
346  art::Ptr<raw::RawDigit> digitVec20(digitVecHandle2, 0);
347  unsigned int dataSize2 = digitVec20->Samples();
348  if (dataSize != dataSize2) {
349  throw cet::exception("RawWaveformClnSigDump")
350  << "RawDigits from the 2 data products have different dataSizes: " << dataSize << "not eq to" << dataSize2;
351  }
352 
353  // ... Build a map from channel number -> rawdigitVec
354  std::map<raw::ChannelID_t, art::Ptr<raw::RawDigit>> rawdigitMap;
355  raw::ChannelID_t chnum = raw::InvalidChannelID; // channel number
356  if (rawdigitlist.size()) {
357  for (size_t rdIter = 0; rdIter < digitVecHandle->size(); ++rdIter) {
358  art::Ptr<raw::RawDigit> digitVec(digitVecHandle, rdIter);
359  chnum = digitVec->Channel();
360  if (chnum == raw::InvalidChannelID) continue;
361  rawdigitMap[chnum] = digitVec;
362  }
363  }
364  std::map<raw::ChannelID_t, art::Ptr<raw::RawDigit>> rawdigitMap2;
365  raw::ChannelID_t chnum2 = raw::InvalidChannelID; // channel number
366  if (rawdigitlist2.size()) {
367  for (size_t rdIter = 0; rdIter < digitVecHandle2->size(); ++rdIter) {
368  art::Ptr<raw::RawDigit> digitVec2(digitVecHandle2, rdIter);
369  chnum2 = digitVec2->Channel();
370  if (chnum2 == raw::InvalidChannelID) continue;
371  rawdigitMap2[chnum2] = digitVec2;
372  }
373  }
374 
375  // ... Read in MC particle list
377  if (!evt.getByLabel(fSimulationProducerLabel, particleHandle)) {
378  throw cet::exception("AnalysisExample")
379  << " No simb::MCParticle objects in this event - "
380  << " Line " << __LINE__ << " in file " << __FILE__ << std::endl;
381  }
382 
383  // ... Read in sim channel list
384  auto simChannelHandle =
385  evt.getValidHandle<std::vector<sim::SimChannel>>(fSimChannelLabel);
386 
387  if (!simChannelHandle->size()) return;
388 
389  // ... Create a map of track IDs to generator labels
390  //Get a list of generator names.
391  //std::vector<art::Handle<std::vector<simb::MCTruth>>> mcHandles;
392  //evt.getManyByType(mcHandles);
393  auto mcHandles = evt.getMany<std::vector<simb::MCTruth>>();
394  std::vector<std::pair<int, std::string>> track_id_to_label;
395 
396  for (auto const& mcHandle : mcHandles) {
397  const std::string& sModuleLabel = mcHandle.provenance()->moduleLabel();
398  art::FindManyP<simb::MCParticle> findMCParts(mcHandle, evt, fSimulationProducerLabel);
399  std::vector<art::Ptr<simb::MCParticle>> mcParts = findMCParts.at(0);
400  for (const art::Ptr<simb::MCParticle> ptr : mcParts) {
401  int track_id = ptr->TrackId();
402  gf->add(track_id, sModuleLabel);
403  }
404  }
405 
406  std::string dummystr6 = "none ";
407  std::string dummystr7 = "none ";
408 
409  if (fSaveSignal) {
410  // .. create a channel number to trackid-wire signal info map
411  std::map<raw::ChannelID_t, std::map<int, WireSigInfo>> Ch2TrkWSInfoMap;
412 
413  // .. create a track ID to vector of channel numbers (in w/c this track deposited energy) map
414  std::map<int, std::vector<raw::ChannelID_t>> Trk2ChVecMap;
415 
416  // ... Loop over simChannels
417  for (auto const& channel : (*simChannelHandle)) {
418 
419  // .. get simChannel channel number
420  const raw::ChannelID_t ch1 = channel.Channel();
421  if (ch1 == raw::InvalidChannelID) continue;
422  if (geo::PlaneGeo::ViewName(fgeom->View(ch1)) != fPlaneToDump[0]) continue;
423 
424  bool selectThisChannel = false;
425 
426  // .. create a track ID to wire signal info map
427  std::map<int, WireSigInfo> Trk2WSInfoMap;
428 
429  // ... Loop over all ticks with ionization energy deposited
430  auto const& timeSlices = channel.TDCIDEMap();
431  for (auto const& timeSlice : timeSlices) {
432 
433  auto const& energyDeposits = timeSlice.second;
434  auto const tpctime = timeSlice.first;
435  unsigned int tdctick = static_cast<unsigned int>(clockData.TPCTDC2Tick(double(tpctime)));
436  if (tdctick < 0 || tdctick > (dataSize - 1)) continue;
437 
438  // ... Loop over all energy depositions in this tick
439  for (auto const& energyDeposit : energyDeposits) {
440 
441  if (!energyDeposit.trackID) continue;
442  int trkid = energyDeposit.trackID;
443  simb::MCParticle particle = PIS->TrackIdToMotherParticle(trkid);
444  //std::cout << energyDeposit.trackID << " " << trkid << " " << particle.TrackId() << std::endl;
445 
446  // .. ignore this energy deposition if incident particle energy below some threshold
447  if (particle.E() < fMinParticleEnergyGeV) continue;
448 
449  int eve_id = PIS->TrackIdToEveTrackId(trkid);
450  if (!eve_id) continue;
451  std::string genlab = gf->get_gen(eve_id);
452 
453  if (Trk2WSInfoMap.find(trkid) == Trk2WSInfoMap.end()) {
454  WireSigInfo wsinf;
455  wsinf.pdgcode = particle.PdgCode();
456  wsinf.genlab = genlab;
457  wsinf.procid = particle.Process();
458  wsinf.tdcmin = dataSize - 1;
459  wsinf.tdcmax = 0;
460  wsinf.tdcpeak = -1;
461  wsinf.adcpeak = 0;
462  wsinf.edep = 0.;
463  wsinf.numel = 0;
464  Trk2WSInfoMap.insert(std::pair<int, WireSigInfo>(trkid, wsinf));
465  }
466  if (tdctick < Trk2WSInfoMap.at(trkid).tdcmin) Trk2WSInfoMap.at(trkid).tdcmin = tdctick;
467  if (tdctick > Trk2WSInfoMap.at(trkid).tdcmax) Trk2WSInfoMap.at(trkid).tdcmax = tdctick;
468  Trk2WSInfoMap.at(trkid).edep += energyDeposit.energy;
469  Trk2WSInfoMap.at(trkid).numel += energyDeposit.numElectrons;
470  }
471  } // loop over timeSlices
472 
473  auto search2 = rawdigitMap2.find(ch1);
474  if (search2 == rawdigitMap2.end()) continue;
475  art::Ptr<raw::RawDigit> rawdig2 = (*search2).second;
476  std::vector<short> rawadc(dataSize);
477  raw::Uncompress(rawdig2->ADCs(), rawadc, rawdig2->GetPedestal(), rawdig2->Compression());
478  std::vector<short> adcvec2(dataSize);
479  for (size_t j = 0; j < rawadc.size(); ++j) {
480  adcvec2[j] = rawadc[j] - rawdig2->GetPedestal();
481  }
482 
483  if (!Trk2WSInfoMap.empty()) {
484  for (std::pair<int, WireSigInfo> itmap : Trk2WSInfoMap) {
485  // find the peak adc value in the signal-only raw digits within the range tdcmin->tdcmax
486  int pkadc = 0;
487  int pktdc = -1;
488  for (size_t i = itmap.second.tdcmin; i <= itmap.second.tdcmax; i++) {
489  if (abs(adcvec2[i]) > abs(pkadc)) {
490  pkadc = adcvec2[i];
491  pktdc = i;
492  }
493  }
494  Trk2WSInfoMap.at(itmap.first).tdcpeak = pktdc;
495  Trk2WSInfoMap.at(itmap.first).adcpeak = pkadc;
496 
497  if (fSelectGenLabel != "ANY") {
498  if (itmap.second.genlab != fSelectGenLabel) continue;
499  }
500  if (fSelectProcID != "ANY") {
501  if (itmap.second.procid != fSelectProcID) continue;
502  }
503  if (fSelectPDGCode != 0) {
504  if (itmap.second.pdgcode != fSelectPDGCode) continue;
505  }
506  itmap.second.genlab.resize(6, ' ');
507  itmap.second.procid.resize(7, ' ');
508  if (itmap.second.numel >= fMinNumberOfElectrons &&
509  itmap.second.edep >= fMinEnergyDepositedMeV && abs(pkadc) >= fMinPureSignalADCs) {
510  if (fMaxNumberOfElectrons >= 0 && itmap.second.numel >= fMaxNumberOfElectrons) {
511  continue;
512  }
513  else {
514  int trkid = itmap.first;
515  if (Trk2ChVecMap.find(trkid) == Trk2ChVecMap.end()) {
516  std::vector<raw::ChannelID_t> chvec;
517  Trk2ChVecMap.insert(std::pair<int, std::vector<raw::ChannelID_t>>(trkid, chvec));
518  }
519  Trk2ChVecMap.at(trkid).push_back(ch1);
520  selectThisChannel = true;
521  }
522  }
523  } // loop over Trk2WSinfoMap
524  if (selectThisChannel) {
525  Ch2TrkWSInfoMap.insert(
526  std::pair<raw::ChannelID_t, std::map<int, WireSigInfo>>(ch1, Trk2WSInfoMap));
527  }
528  } // if Trk2WSInfoMap not empty
529 
530  } // loop over SimChannels
531 
532  std::set<raw::ChannelID_t> selected_channels;
533 
534  // ... Now write out the signal waveforms for each track
535  if (!Trk2ChVecMap.empty()) {
536  for (auto const& ittrk : Trk2ChVecMap) {
537  int i = fRandFlat.fireInt(ittrk.second.size()); // randomly select one channel with a signal from this particle
538  chnum = ittrk.second[i];
539 
540  if (not selected_channels.insert(chnum).second) {
541  continue;
542  }
543 
544  std::map<raw::ChannelID_t, std::map<int, WireSigInfo>>::iterator itchn;
545  itchn = Ch2TrkWSInfoMap.find(chnum);
546  if (itchn != Ch2TrkWSInfoMap.end()) {
547 
548  std::vector<short> adcvec(dataSize); // vector to hold zero-padded full waveform
549 
550  auto search = rawdigitMap.find(chnum);
551  if (search == rawdigitMap.end()) continue;
552  art::Ptr<raw::RawDigit> rawdig = (*search).second;
553  std::vector<short> rawadc(dataSize); // vector to hold uncompressed adc values later
554  raw::Uncompress(rawdig->ADCs(), rawadc, rawdig->GetPedestal(), rawdig->Compression());
555  for (size_t j = 0; j < rawadc.size(); ++j) {
556  adcvec[j] = rawadc[j] - rawdig->GetPedestal();
557  }
558 
559  std::vector<short> adcvec2(dataSize); // vector to hold zero-padded full signal-only waveform
560 
561  auto search2 = rawdigitMap2.find(chnum);
562  if (search2 == rawdigitMap2.end()) continue;
563  art::Ptr<raw::RawDigit> rawdig2 = (*search2).second;
564  raw::Uncompress(rawdig2->ADCs(), rawadc, rawdig2->GetPedestal(), rawdig2->Compression());
565  for (size_t j = 0; j < rawadc.size(); ++j) {
566  adcvec2[j] = rawadc[j] - rawdig2->GetPedestal();
567  }
568 
569  // .. write out info for each peak
570  // a full waveform has at least one peak; the output will save up to 5 peaks (if there is
571  // only 1 peak, will fill the other 4 with 0);
572  // for fShortWaveformSize: only use the first peak's start_tick
573 
574  if (fUseFullWaveform) {
575 
576  c2numpy_uint32(&npywriter, evt.id().event());
577  c2numpy_uint32(&npywriter, chnum);
579  c2numpy_uint16(&npywriter, itchn->second.size()); // size of Trk2WSInfoMap, or #peaks
580  unsigned int icnt = 0;
581  for (auto & it : itchn->second) {
582  c2numpy_int32(&npywriter, it.first); // trackid
583  c2numpy_int32(&npywriter, it.second.pdgcode); // pdgcode
584  c2numpy_string(&npywriter, it.second.genlab.c_str()); // genlab
585  c2numpy_string(&npywriter, it.second.procid.c_str()); // procid
586  c2numpy_float32(&npywriter, it.second.edep); // edepo
587  c2numpy_uint32(&npywriter, it.second.numel); // numelec
588 
589  c2numpy_uint16(&npywriter, it.second.tdcmin); // stck1
590  c2numpy_uint16(&npywriter, it.second.tdcmax); // stck2
591  c2numpy_int32(&npywriter, it.second.tdcpeak); // pktdc
592  c2numpy_int32(&npywriter, it.second.adcpeak); // pkadc
593 
594  icnt++;
595  if (icnt == 5) break;
596  }
597 
598  // .. pad with 0's if number of peaks less than 5
599  for (unsigned int i = icnt; i < 5; ++i) {
602  c2numpy_string(&npywriter, dummystr6.c_str());
603  c2numpy_string(&npywriter, dummystr7.c_str());
610  }
611 
612  for (unsigned int itck = 0; itck < dataSize; ++itck) {
613  c2numpy_int16(&npywriter, adcvec[itck]);
614  }
615  for (unsigned int itck = 0; itck < dataSize; ++itck) {
616  c2numpy_int16(&npywriter2, adcvec2[itck]);
617  }
618 
619  } else {
620 
621  // .. first loop to find largest signal
622  double EDep = 0.;
623  unsigned int TDCMin, TDCMax;
624  bool foundmaxsig = false;
625  for (auto & it : itchn->second) {
626  if (it.second.edep > EDep && it.second.adcpeak != 0 && it.second.numel > 0){
627  EDep = it.second.edep;
628  TDCMin = it.second.tdcmin;
629  TDCMax = it.second.tdcmax;
630  foundmaxsig = true;
631  }
632  }
633  if (foundmaxsig) {
634  int sigtdc1, sigtdc2, sighwid, sigfwid, sigtdcm;
636  sigtdc1 = TDCMin - fEstIndFWForOffset/2;
637  sigtdc2 = TDCMax + 3*fEstIndFWForOffset/2;
638  } else {
639  sigtdc1 = TDCMin - fEstColFWForOffset/2;
640  sigtdc2 = TDCMax + fEstColFWForOffset/2;
641  }
642  sigfwid=sigtdc2 - sigtdc1;
643  sighwid=sigfwid/2;
644  sigtdcm=sigtdc1+sighwid;
645 
646  int start_tick = -1;
647  int end_tick = -1;
648  // .. set window edges to contain the largest signal
649  if (sigfwid < (int)fShortWaveformSize) {
650  // --> case 1: signal range fits within window
651  int dt = fShortWaveformSize - sigfwid;
652  start_tick = sigtdc1 - dt * fRandFlat.fire(0,1);
653  }
654  else {
655  // --> case 2: signal range larger than window
656  int mrgn = fShortWaveformSize/20;
657  int dt = fShortWaveformSize - 2*mrgn;
658  start_tick = sigtdcm - mrgn - dt * fRandFlat.fire(0,1);
659  }
660  if (start_tick < 0) start_tick = 0;
661  end_tick = start_tick + fShortWaveformSize - 1;
662  if (end_tick > int (dataSize - 1)) {
663  end_tick = dataSize - 1;
664  start_tick = end_tick - fShortWaveformSize + 1;
665  }
666 
667  c2numpy_uint32(&npywriter, evt.id().event());
668  c2numpy_uint32(&npywriter, chnum);
670 
671  // .. second loop to select only signals that are within the window
672 
673  int it_trk[5],it_pdg[5],it_nel[5],pk_tdc[5],pk_adc[5];
674  unsigned int stck_1[5],stck_2[5];
675  std::string it_glb[5], it_prc[5];
676  double it_edp[5];
677 
678  unsigned int icnt = 0;
679 
680  for (auto & it : itchn->second) {
681  if (abs(it.second.adcpeak) < fMinPureSignalADCs) continue;
682  if (( it.second.tdcmin >= (unsigned int)start_tick && it.second.tdcmin < (unsigned int)end_tick) ||
683  ( it.second.tdcmax > (unsigned int)start_tick && it.second.tdcmax <= (unsigned int)end_tick)) {
684 
685  it_trk[icnt] = it.first;
686  it_pdg[icnt] = it.second.pdgcode;
687  it_glb[icnt] = it.second.genlab;
688  it_prc[icnt] = it.second.procid;
689  it_edp[icnt] = it.second.edep;
690  it_nel[icnt] = it.second.numel;
691 
692  unsigned int mintdc = it.second.tdcmin;
693  unsigned int maxtdc = it.second.tdcmax;
694  if (mintdc < (unsigned int)start_tick)mintdc = start_tick;
695  if (maxtdc > (unsigned int)end_tick)maxtdc = end_tick;
696 
697  stck_1[icnt] = mintdc - start_tick;
698  stck_2[icnt] = maxtdc - start_tick;
699  pk_tdc[icnt] = it.second.tdcpeak - start_tick;
700  pk_adc[icnt] = it.second.adcpeak;
701 
702  icnt++;
703  if (icnt == 5) break;
704  }
705  }
706 
707  c2numpy_uint16(&npywriter, icnt); // number of peaks
708 
709  for (unsigned int i = 0; i < icnt; ++i) {
710  c2numpy_int32(&npywriter, it_trk[i]); // trackid
711  c2numpy_int32(&npywriter, it_pdg[i]); // pdgcode
712  c2numpy_string(&npywriter, it_glb[i].c_str()); // genlab
713  c2numpy_string(&npywriter, it_prc[i].c_str()); // procid
714  c2numpy_float32(&npywriter,it_edp[i]); // edepo
715  c2numpy_uint32(&npywriter, it_nel[i]); // numelec
716  c2numpy_uint16(&npywriter, stck_1[i]); // stck1
717  c2numpy_uint16(&npywriter, stck_2[i]); // stck2
718  c2numpy_int32(&npywriter, pk_tdc[i]); // pktdc
719  c2numpy_int32(&npywriter, pk_adc[i]); // pkadc
720  }
721 
722  // .. pad with 0's if number of peaks less than 5
723  for (unsigned int i = icnt; i < 5; ++i) {
726  c2numpy_string(&npywriter, dummystr6.c_str());
727  c2numpy_string(&npywriter, dummystr7.c_str());
734  }
735 
736  for (unsigned int itck = start_tick; itck < (start_tick + fShortWaveformSize); ++itck) {
737  c2numpy_int16(&npywriter, adcvec[itck]);
738  }
739  for (unsigned int itck = start_tick; itck < (start_tick + fShortWaveformSize); ++itck) {
740  c2numpy_int16(&npywriter2, adcvec2[itck]);
741  }
742 
743  } // foundmaxsig
744  }
745  }
746  }
747  }
748  }
749  else {
750  //save noise
751  int noisechancount = 0;
752  std::map<raw::ChannelID_t, bool> signalMap;
753  for (auto const& channel : (*simChannelHandle)) {
754  signalMap[channel.Channel()] = true;
755  }
756  // .. create a vector for shuffling the wire channel indices
757  auto seed = std::chrono::high_resolution_clock::now().time_since_epoch().count();
758  std::vector<size_t> randigitmap;
759  for (size_t i=0; i<rawdigitlist.size(); ++i) randigitmap.push_back(i);
760  std::shuffle ( randigitmap.begin(), randigitmap.end(), std::mt19937(seed) );
761 
762  for (size_t rdIter = 0; rdIter < rawdigitlist.size(); ++rdIter) {
763 
764  if (noisechancount==fMaxNoiseChannelsPerEvent)break;
765 
766  std::vector<float> adcvec(dataSize); // vector to wire adc values
767 
768  size_t ranIdx=randigitmap[rdIter];
769  art::Ptr<raw::RawDigit> digitVec(digitVecHandle, ranIdx);
770  if (signalMap[digitVec->Channel()]) continue;
771 
772  std::vector<short> rawadc(dataSize); // vector to hold uncompressed adc values later
773  if (geo::PlaneGeo::ViewName(fgeom->View(digitVec->Channel())) != fPlaneToDump[0]) continue;
774  raw::Uncompress(digitVec->ADCs(), rawadc, digitVec->GetPedestal(), digitVec->Compression());
775  for (size_t j = 0; j < rawadc.size(); ++j) {
776  adcvec[j] = rawadc[j] - digitVec->GetPedestal();
777  }
778  c2numpy_uint32(&npywriter, evt.id().event());
779  c2numpy_uint32(&npywriter, digitVec->Channel());
781  geo::PlaneGeo::ViewName(fgeom->View(digitVec->Channel())).c_str());
782 
783  c2numpy_uint16(&npywriter, 0); //number of peaks
784  for (unsigned int i = 0; i < 5; ++i) {
787  c2numpy_string(&npywriter, dummystr6.c_str());
788  c2numpy_string(&npywriter, dummystr7.c_str());
795  }
796 
797  if (fUseFullWaveform) {
798  for (unsigned int itck = 0; itck < dataSize; ++itck) {
799  c2numpy_int16(&npywriter, short(adcvec[itck]));
800  }
801  }
802  else {
803  int start_tick = int((dataSize - fShortWaveformSize) * fRandFlat.fire(0, 1));
804  for (unsigned int itck = start_tick; itck < (start_tick + fShortWaveformSize); ++itck) {
805  c2numpy_int16(&npywriter, short(adcvec[itck]));
806  }
807  }
808 
809  ++noisechancount;
810  }
811  std::cout << "Total number of noise channels " << noisechancount << std::endl;
812  }
813 }
def analyze(root, level, gtrees, gbranches, doprint)
Definition: rootstat.py:69
static QCString name
Definition: declinfo.cpp:673
double E(const int i=0) const
Definition: MCParticle.h:233
float GetPedestal() const
Definition: RawDigit.h:214
base_engine_t & createEngine(seed_t seed)
bool has_gen(std::string gname)
int c2numpy_int32(c2numpy_writer *writer, int32_t data)
Definition: c2numpy.h:290
const ADCvector_t & ADCs() const
Reference to the compressed ADC count vector.
Definition: RawDigit.h:210
int c2numpy_init(c2numpy_writer *writer, const std::string outputFilePrefix, int32_t numRowsPerFile)
Definition: c2numpy.h:142
int PdgCode() const
Definition: MCParticle.h:212
static std::string ViewName(geo::View_t view)
Returns the name of the specified view.
Definition: PlaneGeo.cxx:763
ULong64_t Samples() const
Number of samples in the uncompressed ADC data.
Definition: RawDigit.h:213
int c2numpy_int16(c2numpy_writer *writer, int16_t data)
Definition: c2numpy.h:282
void add(const int &track_id, const std::string &gname)
std::string string
Definition: nybbler.cc:12
ChannelID_t Channel() const
DAQ channel this raw data was read from.
Definition: RawDigit.h:212
uint8_t channel
Definition: CRTFragment.hh:201
std::string Process() const
Definition: MCParticle.h:215
int c2numpy_close(c2numpy_writer *writer)
Definition: c2numpy.h:407
Particle class.
int c2numpy_uint16(c2numpy_writer *writer, uint16_t data)
Definition: c2numpy.h:314
RunNumber_t run() const
Definition: EventID.h:98
art framework interface to geometry description
std::pair< int, std::string > track_id_to_string
std::string fDigitModuleLabel
module that made digits
int c2numpy_addcolumn(c2numpy_writer *writer, const std::string name, c2numpy_type type)
Definition: c2numpy.h:158
T abs(T value)
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
Definition: DataViewImpl.h:633
constexpr ChannelID_t InvalidChannelID
ID of an invalid channel.
Definition: RawTypes.h:32
const double e
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:67
void beginJob()
Definition: Breakpoints.cc:14
const double a
std::string getenv(std::string const &name)
Definition: getenv.cc:15
std::vector< Handle< PROD > > getMany(SelectorBase const &selector=MatchAllSelector{}) const
Definition: DataViewImpl.h:479
std::vector< track_id_to_string > track_id_map
int c2numpy_uint32(c2numpy_writer *writer, uint32_t data)
Definition: c2numpy.h:322
RawWaveformClnSigDump(fhicl::ParameterSet const &p)
ValidHandle< PROD > getValidHandle(InputTag const &tag) const
Definition: DataViewImpl.h:441
p
Definition: test.py:223
Definition: search.py:1
View_t View(geo::PlaneID const &pid) const
Returns the view (wire orientation) on the channels of specified TPC plane.
std::string fSimChannelLabel
module that made simchannels
std::string fSimulationProducerLabel
producer that tracked simulated part. through detector
int c2numpy_string(c2numpy_writer *writer, const char *data)
Definition: c2numpy.h:394
art::ServiceHandle< geo::Geometry > fgeom
raw::Compress_t Compression() const
Compression algorithm used to store the ADC counts.
Definition: RawDigit.h:216
c2numpy_type
Definition: c2numpy.h:29
std::set< std::string > generator_names
simb::MCParticle TrackIdToMotherParticle(int const id) const
Declaration of signal hit object.
Encapsulate the construction of a single detector plane.
int c2numpy_float32(c2numpy_writer *writer, float data)
Definition: c2numpy.h:354
void analyze(art::Event const &e) override
genFinder * gf
Interface for experiment-specific channel quality info provider.
static bool * b
Definition: config.cpp:1043
EventNumber_t event() const
Definition: EventID.h:116
Declaration of basic channel signal object.
TCEvent evt
Definition: DataStructs.cxx:7
unsigned int ChannelID_t
Type representing the ID of a readout channel.
Definition: RawTypes.h:28
std::string fCleanSignalDigitModuleLabel
module that made the signal-only digits
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:297
Interface for experiment-specific service for channel quality info.
void Uncompress(const std::vector< short > &adc, std::vector< short > &uncompressed, raw::Compress_t compress)
Uncompresses a raw data buffer.
Definition: raw.cxx:776
second_as<> second
Type of time stored in seconds, in double precision.
Definition: spacetime.h:85
std::string get_gen(int tid)
SubRunNumber_t subRun() const
Definition: EventID.h:110
EventID id() const
Definition: Event.cc:34
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
QTextStream & endl(QTextStream &s)
art::ServiceHandle< cheat::ParticleInventoryService > PIS