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