LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
CellTree_module.cc
Go to the documentation of this file.
1 // Read data from MC raw files and convert it into ROOT tree
2 // Chao Zhang (chao@bnl.gov) 5/13/2014
3 // Added optical info --- Brooke Russell (brussell@yale.edu) 1/31/2017
4 
5 #ifndef CELLTREE_MODULE
6 #define CELLTREE_MODULE
7 
8 // LArSoft includes
10 
23 // #include "MCCheater/BackTracker.h"
24 #include "lardataobj/RawData/raw.h"
33 
34 // Framework includes
45 #include "fhiclcpp/ParameterSet.h"
46 
47 // ROOT includes.
48 #include "TFile.h"
49 #include "TTree.h"
50 #include "TDirectory.h"
51 #include "TMath.h"
52 #include "TString.h"
53 #include "TClonesArray.h"
54 #include "TH1F.h"
55 #include "TH2F.h"
56 #include "TLorentzVector.h"
57 #include "TUnixSystem.h"
58 #include "TDatabasePDG.h"
59 #include "TObjArray.h"
60 #include "TTimeStamp.h"
61 
62 // C++ Includes
63 #include <map>
64 #include <vector>
65 #include <fstream>
66 #include <iostream>
67 #include <cstdio>
68 
69 #define MAX_TRACKS 30000
70 
71 using namespace std;
72 
73 namespace wc {
74 
75 class CellTree : public art::EDAnalyzer {
76 public:
77 
78  explicit CellTree(fhicl::ParameterSet const& pset);
79  virtual ~CellTree();
80 
81  void beginJob();
82  void endJob();
83  void beginRun(const art::Run& run);
84  void analyze(const art::Event& evt);
85 
86  void reconfigure(fhicl::ParameterSet const& pset);
87  void initOutput();
88  void printEvent();
89  void print_vector(ostream& out, vector<double>& v, TString desc, bool end=false);
90 
91  void processRaw(const art::Event& evt);
92  void processCalib(const art::Event& evt);
93  void processOpHit(const art::Event& evt);
94  void processOpFlash(const art::Event& evt);
95  void processSpacePoint( const art::Event& event, TString option, ostream& out=cout);
96  void processSimChannel(const art::Event& evt);
97  void processMC(const art::Event& evt);
98  void processMCTracks();
99  void processTrigger(const art::Event& evt);
100 
101  void reset();
102  void InitProcessMap();
103 
104  bool IsPrimary(int i) { return mc_mother[i] == 0 ; }
105  bool KeepMC(int i);
106  double KE(float* momentum); // KE
107  TString PDGName(int pdg);
108  bool DumpMCJSON(int id, ostream& out);
109  void DumpMCJSON(ostream& out=cout);
110 
111 
112 private:
113 
114  // the parameters we'll read from the .fcl
115  std::string fRawDigitLabel;
116  std::string fCalibLabel;
117  std::string fOpHitLabel;
118  std::string fOpFlashLabel;
119  std::string fTriggerLabel;
120  std::vector<std::string> fSpacePointLabels;
121  std::string fOutFileName;
122  std::string mcOption;
127  bool fSaveRaw;
131  bool fSaveMC;
133  bool fSaveJSON;
134  art::ServiceHandle<geo::Geometry> fGeometry; // pointer to Geometry service
135 
136  // art::ServiceHandle<geo::Geometry> fGeom;
137  // // auto const* larp = lar::providerFrom<detinfo::LArPropertiesService>();
138 
139  TFile *fOutFile;
140  TTree *fEventTree;
141  std::map<std::string, int> processMap;
142  std::map<int, int> savedMCTrackIdMap; // key: id; value: pdg;
143 
144  int entryNo;
145 
146  // Event Tree Leafs
147  int fEvent;
148  int fRun;
149  int fSubRun;
150  double fEventTime;
151 
152  unsigned int fTriggernumber; //trigger counter
153  double fTriggertime; //trigger time w.r.t. electronics clock T0
154  double fBeamgatetime; //beamgate time w.r.t. electronics clock T0
155  unsigned int fTriggerbits; //trigger bits
156 
158  // int fCalib_channelId[MAX_CHANNEL]; // hit channel id; size == fCalib_Nhit
159  // // FIXEME:: cannot save e.g std::vector<std::vector<float> > in ttree
160  std::vector<int> fCalib_channelId;
161  // std::vector<std::vector<float> > fCalib_wf;
162  TClonesArray *fCalib_wf;
163  // std::vector<std::vector<int> > fCalib_wfTDC;
164 
165  int oh_nHits;
166  vector<int> oh_channel;
167  vector<double> oh_bgtime;
168  vector<double> oh_trigtime;
169  vector<double> oh_pe;
170 
172  vector<float> of_t;
173  vector<float> of_peTotal;
174  vector<int> of_multiplicity;
175  TClonesArray *fPEperOpDet;
176 
178  std::vector<int> fRaw_channelId;
179  TClonesArray *fRaw_wf;
180 
182  vector<int> fSIMIDE_channelIdY;
183  vector<int> fSIMIDE_trackId;
184  vector<unsigned short> fSIMIDE_tdc;
185  vector<float> fSIMIDE_x;
186  vector<float> fSIMIDE_y;
187  vector<float> fSIMIDE_z;
188  vector<float> fSIMIDE_numElectrons;
189 
190  int mc_Ntrack; // number of tracks in MC
191  int mc_id[MAX_TRACKS]; // track id; size == mc_Ntrack
192  int mc_pdg[MAX_TRACKS]; // track particle pdg; size == mc_Ntrack
193  int mc_process[MAX_TRACKS]; // track generation process code; size == mc_Ntrack
194  int mc_mother[MAX_TRACKS]; // mother id of this track; size == mc_Ntrack
195  float mc_startXYZT[MAX_TRACKS][4]; // start position of this track; size == mc_Ntrack
196  float mc_endXYZT[MAX_TRACKS][4]; // end position of this track; size == mc_Ntrack
197  float mc_startMomentum[MAX_TRACKS][4]; // start momentum of this track; size == mc_Ntrack
198  float mc_endMomentum[MAX_TRACKS][4]; // end momentum of this track; size == mc_Ntrack
199  std::vector<std::vector<int> > mc_daughters; // daughters id of this track; vector
200  TObjArray *fMC_trackPosition;
201 
202  int mc_isnu; // is neutrino interaction
203  int mc_nGeniePrimaries; // number of Genie primaries
204  int mc_nu_pdg; // pdg code of neutrino
205  int mc_nu_ccnc; // cc or nc
206  int mc_nu_mode; // mode: http://nusoft.fnal.gov/larsoft/doxsvn/html/MCNeutrino_8h_source.html
207  int mc_nu_intType; // interaction type
208  int mc_nu_target; // target interaction
209  int mc_hitnuc; // hit nucleon
210  int mc_hitquark; // hit quark
211  double mc_nu_Q2; // Q^2
212  double mc_nu_W; // W
213  double mc_nu_X; // X
214  double mc_nu_Y; // Y
215  double mc_nu_Pt; // Pt
216  double mc_nu_Theta; // angle relative to lepton
217  float mc_nu_pos[4]; // interaction position of nu
218  float mc_nu_mom[4]; // interaction momentum of nu
219 
220  // ----- derived ---
221  std::map<int, int> trackIndex;
222  std::vector<std::vector<int> > trackParents;
223  std::vector<std::vector<int> > trackChildren;
224  std::vector<std::vector<int> > trackSiblings;
225  TDatabasePDG *dbPDG;
226 
227 }; // class CellTree
228 
229 
230 //-----------------------------------------------------------------------
231 CellTree::CellTree(fhicl::ParameterSet const& parameterSet)
232  : EDAnalyzer(parameterSet)
233 {
234  dbPDG = new TDatabasePDG();
235 
236  reconfigure(parameterSet);
237  InitProcessMap();
238  initOutput();
239  entryNo = 0;
240 }
241 
242 //-----------------------------------------------------------------------
244 {
245 }
246 
247 //-----------------------------------------------------------------------
249  fRawDigitLabel = p.get<std::string>("RawDigitLabel");
250  fCalibLabel = p.get<std::string>("CalibLabel");
251  fOpHitLabel = p.get<std::string>("OpHitLabel");
252  fOpFlashLabel = p.get<std::string>("OpFlashLabel");
253  fTriggerLabel = p.get<std::string>("TriggerLabel");
254  fSpacePointLabels= p.get<std::vector<std::string> >("SpacePointLabels");
255  fOutFileName = p.get<std::string>("outFile");
256  mcOption = p.get<std::string>("mcOption");
257  fSaveMCTrackPoints = p.get<bool>("saveMCTrackPoints");
258  fSaveRaw = p.get<bool>("saveRaw");
259  fSaveCalib = p.get<bool>("saveCalib");
260  fSaveOpHit = p.get<bool>("saveOpHit");
261  fSaveOpFlash = p.get<bool>("saveOpFlash");
262  fSaveMC = p.get<bool>("saveMC");
263  fSaveSimChannel = p.get<bool>("saveSimChannel");
264  fSaveTrigger = p.get<bool>("saveTrigger");
265  fSaveJSON = p.get<bool>("saveJSON");
266  opMultPEThresh = p.get<float>("opMultPEThresh");
267  nRawSamples = p.get<int>("nRawSamples");
268 }
269 
270 //-----------------------------------------------------------------------
272 {
273  TDirectory* tmpDir = gDirectory;
274 
275  fOutFile = new TFile(fOutFileName.c_str(), "recreate");
276 
277  // 3.1: add mc_trackPosition
278  TNamed version("version", "4.0");
279  version.Write();
280 
281  // init Event TTree
282  TDirectory* subDir = fOutFile->mkdir("Event");
283  subDir->cd();
284  fEventTree = new TTree("Sim", "Event Tree from Simulation");
285  fEventTree->Branch("eventNo", &fEvent);
286  fEventTree->Branch("runNo", &fRun);
287  fEventTree->Branch("subRunNo", &fSubRun);
288  fEventTree->Branch("eventTime", &fEventTime); // timestamp
289 
290  fEventTree->Branch("triggerNo", &fTriggernumber); // timestamp
291  fEventTree->Branch("triggerTime", &fTriggertime); // timestamp
292  fEventTree->Branch("beamgateTime", &fBeamgatetime); // timestamp
293  fEventTree->Branch("triggerBits", &fTriggerbits); // timestamp
294 
295  fEventTree->Branch("raw_nChannel", &fRaw_nChannel); // number of hit channels above threshold
296  fEventTree->Branch("raw_channelId" , &fRaw_channelId); // hit channel id; size == raw_nChannel
297  fRaw_wf = new TClonesArray("TH1F");
298  fEventTree->Branch("raw_wf", &fRaw_wf, 256000, 0); // raw waveform adc of each channel
299 
300 
301  fEventTree->Branch("calib_nChannel", &fCalib_nChannel); // number of hit channels above threshold
302  fEventTree->Branch("calib_channelId" , &fCalib_channelId); // hit channel id; size == calib_Nhit
303  fCalib_wf = new TClonesArray("TH1F");
304  fEventTree->Branch("calib_wf", &fCalib_wf, 256000, 0); // calib waveform adc of each channel
305  // fCalib_wf->BypassStreamer();
306  // fEventTree->Branch("calib_wfTDC", &fCalib_wfTDC); // calib waveform tdc of each channel
307 
308  fEventTree->Branch("oh_nHits", &oh_nHits); // number of op hits
309  fEventTree->Branch("oh_channel", &oh_channel); //opchannel id; size == no ophits
310  fEventTree->Branch("oh_bgtime", &oh_bgtime); // optical pulse peak time w.r.t. start of beam gate
311  fEventTree->Branch("oh_trigtime", &oh_trigtime); // optical pulse peak time w.r.t. trigger
312  fEventTree->Branch("oh_pe", &oh_pe); // pulse PE
313 
314  fEventTree->Branch("of_nFlash", &of_nFlash);
315  fEventTree->Branch("of_t", &of_t); // time in us w.r.t. the trigger for each flash
316  fEventTree->Branch("of_peTotal", &of_peTotal); // total PE (sum of all PMTs) for each flash
317  fEventTree->Branch("of_multiplicity", &of_multiplicity); // total number of PMTs above threshold for each flash
318  fPEperOpDet = new TClonesArray("TH1F");
319  fEventTree->Branch("pe_opdet", &fPEperOpDet, 256000, 0);
320 
321  fEventTree->Branch("simide_size", &fSIMIDE_size); // size of stored sim:IDE
322  fEventTree->Branch("simide_channelIdY", &fSIMIDE_channelIdY);
323  fEventTree->Branch("simide_trackId", &fSIMIDE_trackId);
324  fEventTree->Branch("simide_tdc", &fSIMIDE_tdc);
325  fEventTree->Branch("simide_x", &fSIMIDE_x);
326  fEventTree->Branch("simide_y", &fSIMIDE_y);
327  fEventTree->Branch("simide_z", &fSIMIDE_z);
328  fEventTree->Branch("simide_numElectrons", &fSIMIDE_numElectrons);
329 
330  fEventTree->Branch("mc_Ntrack", &mc_Ntrack); // number of tracks in MC
331  fEventTree->Branch("mc_id", &mc_id, "mc_id[mc_Ntrack]/I"); // track id; size == mc_Ntrack
332  fEventTree->Branch("mc_pdg", &mc_pdg, "mc_pdg[mc_Ntrack]/I"); // track particle pdg; size == mc_Ntrack
333  fEventTree->Branch("mc_process", &mc_process, "mc_process[mc_Ntrack]/I"); // track generation process code; size == mc_Ntrack
334  fEventTree->Branch("mc_mother", &mc_mother, "mc_mother[mc_Ntrack]/I"); // mother id of this track; size == mc_Ntrack
335  fEventTree->Branch("mc_daughters", &mc_daughters); // daughters id of this track; vector
336  fEventTree->Branch("mc_startXYZT", &mc_startXYZT, "mc_startXYZT[mc_Ntrack][4]/F"); // start position of this track; size == mc_Ntrack
337  fEventTree->Branch("mc_endXYZT", &mc_endXYZT, "mc_endXYZT[mc_Ntrack][4]/F"); // start position of this track; size == mc_Ntrack
338  fEventTree->Branch("mc_startMomentum", &mc_startMomentum, "mc_startMomentum[mc_Ntrack][4]/F"); // start momentum of this track; size == mc_Ntrack
339  fEventTree->Branch("mc_endMomentum", &mc_endMomentum, "mc_endMomentum[mc_Ntrack][4]/F"); // start momentum of this track; size == mc_Ntrack
340  fMC_trackPosition = new TObjArray();
341  fMC_trackPosition->SetOwner(kTRUE);
342  fEventTree->Branch("mc_trackPosition", &fMC_trackPosition);
343 
344  fEventTree->Branch("mc_isnu", &mc_isnu);
345  fEventTree->Branch("mc_nGeniePrimaries", &mc_nGeniePrimaries);
346  fEventTree->Branch("mc_nu_pdg", &mc_nu_pdg);
347  fEventTree->Branch("mc_nu_ccnc", &mc_nu_ccnc);
348  fEventTree->Branch("mc_nu_mode", &mc_nu_mode);
349  fEventTree->Branch("mc_nu_intType", &mc_nu_intType);
350  fEventTree->Branch("mc_nu_target", &mc_nu_target);
351  fEventTree->Branch("mc_hitnuc", &mc_hitnuc);
352  fEventTree->Branch("mc_hitquark", &mc_hitquark);
353  fEventTree->Branch("mc_nu_Q2", &mc_nu_Q2);
354  fEventTree->Branch("mc_nu_W", &mc_nu_W);
355  fEventTree->Branch("mc_nu_X", &mc_nu_X);
356  fEventTree->Branch("mc_nu_Y", &mc_nu_Y);
357  fEventTree->Branch("mc_nu_Pt", &mc_nu_Pt);
358  fEventTree->Branch("mc_nu_Theta", &mc_nu_Theta);
359  fEventTree->Branch("mc_nu_pos", &mc_nu_pos, "mc_nu_pos[4]/F");
360  fEventTree->Branch("mc_nu_mom", &mc_nu_mom, "mc_nu_mom[4]/F");
361 
362  gDirectory = tmpDir;
363 
364  if (fSaveJSON) {
365  system("rm -rf bee");
366  gSystem->MakeDirectory("bee");
367  gSystem->ChangeDirectory("bee");
368  gSystem->MakeDirectory("data");
369  }
370 
371 }
372 
373 //-----------------------------------------------------------------------
375 {
376 
377 
378 }
379 
380 
381 //-----------------------------------------------------------------------
383 {
384  // Write fEventTree to file
385  TDirectory* tmpDir = gDirectory;
386  fOutFile->cd("/Event");
387 
388  fEventTree->Write();
389 
390  gDirectory = tmpDir;
391 
392  fOutFile->Close();
393 
394  if (fSaveJSON) {
395  system("zip -r bee_upload data");
396  gSystem->ChangeDirectory("..");
397  }
398 }
399 
400 //-----------------------------------------------------------------------
401 void CellTree::beginRun(const art::Run& /*run*/)
402 {
403  mf::LogInfo("CellTree") << "begin run";
404 }
405 
406 //-----------------------------------------------------------------------
408 {
409  reset();
410  fEvent = event.id().event();
411  fRun = event.run();
412  fSubRun = event.subRun();
413  art::Timestamp ts = event.time();
414  TTimeStamp tts(ts.timeHigh(), ts.timeLow());
415  fEventTime = tts.AsDouble();
416 
417  if (fSaveRaw) processRaw(event);
418  if (fSaveCalib) processCalib(event);
419  if (fSaveOpHit) processOpHit(event);
420  if (fSaveOpFlash) processOpFlash(event);
422  if (fSaveMC) processMC(event);
423  if (fSaveTrigger) processTrigger(event);
424 
425  if (fSaveJSON) {
426  gSystem->MakeDirectory(TString::Format("data/%i", entryNo).Data());
427  int nSp = fSpacePointLabels.size();
428  for (int i=0; i<nSp; i++) {
429  TString jsonfile;
430  jsonfile.Form("data/%i/%i-%s.json", entryNo, entryNo, fSpacePointLabels[i].c_str());
431  std::ofstream out(jsonfile.Data());
432  processSpacePoint(event, fSpacePointLabels[i], out);
433  out.close();
434  }
435 
436  if(fSaveMC) {
437  processMCTracks();
438  TString jsonfile;
439  jsonfile.Form("data/%i/%i-mc.json", entryNo, entryNo);
440  std::ofstream out(jsonfile.Data());
441  DumpMCJSON(out);
442  out.close();
443  }
444  }
445 
446  // printEvent();
447  fEventTree->Fill();
448 
449  entryNo++;
450 }
451 
452 //-----------------------------------------------------------------------
454 {
455 
456  fRaw_channelId.clear();
457  // fRaw_wf->Clear();
458  fRaw_wf->Delete();
459 
460  fCalib_channelId.clear();
461  fCalib_wf->Clear();
462 
463  oh_channel.clear();
464  oh_bgtime.clear();
465  oh_trigtime.clear();
466  oh_pe.clear();
467 
468  of_t.clear();
469  of_peTotal.clear();
470  of_multiplicity.clear();
471  fPEperOpDet->Delete();
472 
473  fSIMIDE_channelIdY.clear();
474  fSIMIDE_trackId.clear();
475  fSIMIDE_tdc.clear();
476  fSIMIDE_x.clear();
477  fSIMIDE_y.clear();
478  fSIMIDE_z.clear();
479  fSIMIDE_numElectrons.clear();
480 
481  mc_Ntrack = 0;
482  for (int i=0; i<MAX_TRACKS; i++) {
483  mc_id[i] = 0;
484  mc_pdg[i] = 0;
485  mc_mother[i] = 0;
486  for (int j=0; j<4; j++) {
487  mc_startXYZT[i][j] = 0;
488  mc_endXYZT[i][j] = 0;
489  mc_startMomentum[i][j] = 0;
490  mc_endMomentum[i][j] = 0;
491  }
492  }
493  mc_daughters.clear();
494  savedMCTrackIdMap.clear();
495  fMC_trackPosition->Clear();
496 
497  mc_isnu = 0;
498  mc_nGeniePrimaries = -1;
499  mc_nu_pdg = -1;
500  mc_nu_ccnc = -1;
501  mc_nu_mode = -1;
502  mc_nu_intType = -1;
503  mc_nu_target = -1;
504  mc_hitnuc = -1;
505  mc_hitquark = -1;
506  mc_nu_Q2 = -1;
507  mc_nu_W = -1;
508  mc_nu_X = -1;
509  mc_nu_Y = -1;
510  mc_nu_Pt = -1;
511  mc_nu_Theta = -1;
512  for (int i=0; i<4; i++) {
513  mc_nu_pos[i] = 0;
514  mc_nu_mom[i] = 0;
515  }
516 
517  trackIndex.clear();
518  trackParents.clear();
519  trackChildren.clear();
520  trackSiblings.clear();
521 }
522 
523 //-----------------------------------------------------------------------
525 {
527  if (! event.getByLabel(fRawDigitLabel, rawdigit)) {
528  cout << "WARNING: no label " << fRawDigitLabel << endl;
529  return;
530  }
531  std::vector< art::Ptr<raw::RawDigit> > wires;
532  art::fill_ptr_vector(wires, rawdigit);
533 
534  fRaw_nChannel = wires.size();
535 
536  int i=0;
537  for (auto const& wire: wires) {
538  int chanId = wire->Channel();
539  fRaw_channelId.push_back(chanId);
540 
541  int nSamples = wire->Samples();
542  std::vector<short> uncompressed(nSamples);
543  raw::Uncompress(wire->ADCs(), uncompressed, wire->Compression());
544 
545  TH1F *h = new((*fRaw_wf)[i]) TH1F("", "", nRawSamples, 0, nRawSamples);
546  for (int j=1; j<=nSamples; j++) {
547  h->SetBinContent(j, uncompressed[j-1]);
548  }
549  i++;
550  if (i==1) {
551  cout << nSamples << " samples expanding to " << nRawSamples << endl;
552  }
553  }
554 
555 }
556 
557 //-----------------------------------------------------------------------
559 {
560 
562  if (! event.getByLabel(fCalibLabel, wires_handle)) {
563  cout << "WARNING: no label " << fCalibLabel << endl;
564  return;
565  }
566  std::vector< art::Ptr<recob::Wire> > wires;
567  art::fill_ptr_vector(wires, wires_handle);
568 
569  // wires size should == Nchannels == 1992; (no hit channel has a flat 0-waveform)
570  // cout << "\n wires size: " << wires.size() << endl;
571  fCalib_nChannel = wires.size();
572 
573  int i=0;
574  for (auto const& wire: wires) {
575  std::vector<float> calibwf = wire->Signal();
576  int chanId = wire->Channel();
577  fCalib_channelId.push_back(chanId);
578  TH1F *h = new((*fCalib_wf)[i]) TH1F("", "", nRawSamples, 0, nRawSamples);
579  for (int j=1; j<=nRawSamples; j++) {
580  h->SetBinContent(j, calibwf[j]);
581  }
582  // fCalib_wf.push_back(calibwf);
583  // cout << chanId << ", " << nSamples << endl;
584  i++;
585  }
586 
587 }
588 
589  //----------------------------------------------------------------------
591 {
593  if(! event.getByLabel(fOpHitLabel, ophit_handle)){
594  cout << "WARNING: no label " << fOpHitLabel << endl;
595  return;
596  }
597  std::vector<art::Ptr<recob::OpHit> > ophits;
598  art::fill_ptr_vector(ophits, ophit_handle);
599  oh_nHits = (int)ophits.size();
600 
601  for(auto const& oh : ophits){
602  oh_channel.push_back(oh->OpChannel());
603  oh_bgtime.push_back(oh->PeakTime());
604  oh_trigtime.push_back(oh->PeakTimeAbs());
605  oh_pe.push_back(oh->PE());
606  }
607 
608 }
609 
610  //----------------------------------------------------------------------
612  {
614  if(! event.getByLabel(fOpFlashLabel, flash_handle)){
615  cout << "WARNING: no label " << fOpFlashLabel << endl;
616  return;
617  }
618  std::vector<art::Ptr<recob::OpFlash> > flashes;
619  art::fill_ptr_vector(flashes, flash_handle);
620  of_nFlash = (int)flashes.size();
621 
622  int a=0;
623  int nOpDet = fGeometry->NOpDets();
624 
625  for(auto const& flash: flashes){
626  of_t.push_back(flash->Time());
627  of_peTotal.push_back(flash->TotalPE());
628  TH1F *h = new ((*fPEperOpDet)[a]) TH1F("","",nOpDet,0,nOpDet);
629 
630  int mult = 0;
631  for(int i=0; i<nOpDet; ++i){
632  if(flash->PE(i) >= opMultPEThresh){
633  mult++;
634  }
635  h->SetBinContent(i, flash->PE(i));
636  }
637  of_multiplicity.push_back(mult);
638  a++;
639  }
640  }
641 
642 //-----------------------------------------------------------------------
644 {
646  event.getByLabel("largeant", simChannelHandle);
647  // cout << "total simChannel: " << (*simChannelHandle).size() << endl;
648  fSIMIDE_size = 0;
649  for ( auto const& channel : (*simChannelHandle) ) {
650  auto channelNumber = channel.Channel();
651  // cout << channelNumber << endl;
652  // if (! (fGeometry->SignalType( channelNumber ) == geo::kCollection) ) {
653  // continue;
654  // }
655  auto const& timeSlices = channel.TDCIDEMap();
656  for ( auto const& timeSlice : timeSlices ) {
657  auto const& energyDeposits = timeSlice.second;
658  for ( auto const& energyDeposit : energyDeposits ) {
659  fSIMIDE_size++;
660  fSIMIDE_channelIdY.push_back(channelNumber);
661  fSIMIDE_tdc.push_back(timeSlice.first);
662  fSIMIDE_trackId.push_back(energyDeposit.trackID);
663  fSIMIDE_x.push_back(energyDeposit.x);
664  fSIMIDE_y.push_back(energyDeposit.y);
665  fSIMIDE_z.push_back(energyDeposit.z);
666  fSIMIDE_numElectrons.push_back(energyDeposit.numElectrons);
667  // cout << channelNumber << ": " << energyDeposit.trackID << ": " << timeSlice.first << ": "
668  // << energyDeposit.x << ", " << energyDeposit.y << ", " << energyDeposit.z << ", "
669  // << energyDeposit.numElectrons << endl;
670  }
671  }
672 
673  }
674  cout << "total IDEs: " << fSIMIDE_size << endl;
675 
676 
677 }
678 
679 //-----------------------------------------------------------------------
681 {
683  if (! event.getByLabel("largeant", particleHandle)) return;
684  std::vector< art::Ptr<simb::MCParticle> > particles;
685  art::fill_ptr_vector(particles, particleHandle);
686 
688  event.getByLabel("largeant", simChannelHandle);
689 
690  // art::ServiceHandle<cheat::BackTracker> bt;
691  art::FindOneP<simb::MCTruth> fo(particleHandle, event, "largeant");
692 
693  int i=0; // track index in saved MCParticles;
694  int i_all=0; // track index in all MCParticles;
695  for (auto const& particle: particles ) {
696  art::Ptr<simb::MCTruth> mctruth = fo.at(i_all);
697  i_all++;
698 
699  if (mcOption == "nuOnly") {
700  if ( !(mctruth->Origin() == 1 && particle->Mother() == 0) ) {
701  continue;
702  }
703  }
704 
705  // if (mctruth->Origin() == 1 || mc_mother[i] == 0) {
706  // cout << "process: " << particle->Process()
707  // << ", id: " << mc_id[i]
708  // << ", pdg: " << mc_pdg[i]
709  // << ", mother: " << mc_mother[i]
710  // << ", nDaughter: " << (particle->NumberDaughters())
711  // << ", truth: " << mctruth->Origin()
712  // << endl;
713  // }
714  // const art::Ptr<simb::MCTruth> mctruth = bt_serv->TrackIDToMCTruth(mc_id[i]);
715 
716  mc_process[i] = processMap[particle->Process()];
717  if (mc_process[i] == 0) cout << "unknown process: " << particle->Process() << endl;
718  mc_id[i] = particle->TrackId();
719  mc_pdg[i] = particle->PdgCode();
720  mc_mother[i] = particle->Mother();
721  savedMCTrackIdMap[mc_id[i]] = mc_pdg[i];
722 
723  int Ndaughters = particle->NumberDaughters();
724  vector<int> daughters;
725  for (int i=0; i<Ndaughters; i++) {
726  daughters.push_back(particle->Daughter(i));
727  }
728  mc_daughters.push_back(daughters);
729  size_t numberTrajectoryPoints = particle->NumberTrajectoryPoints();
730  int last = numberTrajectoryPoints - 1;
731  const TLorentzVector& positionStart = particle->Position(0);
732  const TLorentzVector& positionEnd = particle->Position(last);
733  const TLorentzVector& momentumStart = particle->Momentum(0);
734  const TLorentzVector& momentumEnd = particle->Momentum(last);
735  positionStart.GetXYZT(mc_startXYZT[i]);
736  positionEnd.GetXYZT(mc_endXYZT[i]);
737  momentumStart.GetXYZT(mc_startMomentum[i]);
738  momentumEnd.GetXYZT(mc_endMomentum[i]);
739 
740  if (fSaveMCTrackPoints) {
741  TClonesArray *Lposition = new TClonesArray("TLorentzVector", numberTrajectoryPoints);
742  // Read the position and momentum along this particle track
743  for(unsigned int j=0; j<numberTrajectoryPoints; j++) {
744  new ((*Lposition)[j]) TLorentzVector(particle->Position(j));
745  }
746  fMC_trackPosition->Add(Lposition);
747  }
748 
749  i++;
750  if (i==MAX_TRACKS) {
751  cout << "WARNING:: # tracks exceeded MAX_TRACKS " << MAX_TRACKS << endl;
752  break;
753  }
754  } // particle loop done
755  mc_Ntrack = i;
756  // cout << "MC_Ntracks:" << mc_Ntrack << endl;
757 
758  // Generator Info
759  art::Handle< std::vector<simb::MCTruth> > mctruthListHandle;
760  event.getByLabel("generator",mctruthListHandle);
761  std::vector<art::Ptr<simb::MCTruth> > mclist;
762  art::fill_ptr_vector(mclist, mctruthListHandle);
763  art::Ptr<simb::MCTruth> mctruth;
764 
765  if (mclist.size()>0) {
766  mctruth = mclist.at(0);
767  if (mctruth->NeutrinoSet()) {
768  simb::MCNeutrino nu = mctruth->GetNeutrino();
769  mc_isnu = 1;
770  mc_nGeniePrimaries = mctruth->NParticles();
771  mc_nu_pdg = nu.Nu().PdgCode();
772  mc_nu_ccnc = nu.CCNC();
773  mc_nu_mode = nu.Mode();
775  mc_nu_target = nu.Target();
776  mc_hitnuc = nu.HitNuc();
777  mc_hitquark = nu.HitQuark();
778  mc_nu_Q2 = nu.QSqr();
779  mc_nu_W = nu.W();
780  mc_nu_X = nu.X();
781  mc_nu_Y = nu.Y();
782  mc_nu_Pt = nu.Pt();
783  mc_nu_Theta = nu.Theta();
784 
785  const TLorentzVector& position = nu.Nu().Position(0);
786  const TLorentzVector& momentum = nu.Nu().Momentum(0);
787  position.GetXYZT(mc_nu_pos);
788  momentum.GetXYZT(mc_nu_mom);
789  // cout << "nu: " << mc_nu_pdg << ", nPrim: " << mc_nGeniePrimaries
790  // << ", ccnc: " << mc_nu_ccnc << endl;
791  // for (int i=0; i<mc_nGeniePrimaries; i++) {
792  // simb::MCParticle particle = mctruth->GetParticle(i);
793  // cout << "id: " << particle.TrackId()
794  // << ", pdg: " << particle.PdgCode()
795  // << endl;
796  // }
797  }
798  }
799 
800 }
801 
802 
803 //-----------------------------------------------------------------------
804 void CellTree::processSpacePoint( const art::Event& event, TString option, ostream& out)
805 {
806 
808  if (! event.getByLabel(option.Data(), handle)) {
809  cout << "WARNING: no label " << option << endl;
810  return;
811  }
812  std::vector< art::Ptr<recob::SpacePoint> > sps;
813  art::fill_ptr_vector(sps, handle);
814 
815  double x=0, y=0, z=0, q=0, nq=1;
816  vector<double> vx, vy, vz, vq, vnq;
817 
818  for (auto const& sp: sps ) {
819  // cout << sp->XYZ()[0] << ", " << sp->XYZ()[1] << ", " << sp->XYZ()[2] << endl;
820  x = sp->XYZ()[0];
821  y = sp->XYZ()[1];
822  z = sp->XYZ()[2];
823  vx.push_back(x);
824  vy.push_back(y);
825  vz.push_back(z);
826  vq.push_back(q);
827  vnq.push_back(nq);
828  }
829 
830  out << fixed << setprecision(1);
831  out << "{" << endl;
832 
833  out << '"' << "runNo" << '"' << ":" << '"' << fRun << '"' << "," << endl;
834  out << '"' << "subRunNo" << '"' << ":" << '"' << fSubRun << '"' << "," << endl;
835  out << '"' << "eventNo" << '"' << ":" << '"' << fEvent << '"' << "," << endl;
836 
837  TString geomName(fGeometry->DetectorName().c_str());
838  if (geomName.Contains("35t")) { geomName = "dune35t"; }
839  else if (geomName.Contains("protodune")) { geomName = "protodune"; }
840  else if (geomName.Contains("workspace")) { geomName = "dune10kt_workspace"; }
841  else { geomName = "uboone"; } // use uboone as default
842  out << '"' << "geom" << '"' << ":" << '"' << geomName << '"' << "," << endl;
843 
844 
845  print_vector(out, vx, "x");
846  print_vector(out, vy, "y");
847  print_vector(out, vz, "z");
848 
849  out << fixed << setprecision(0);
850  print_vector(out, vq, "q");
851  print_vector(out, vnq, "nq");
852 
853  out << '"' << "type" << '"' << ":" << '"' << option << '"' << endl;
854  out << "}" << endl;
855 }
856 
857 //-----------------------------------------------------------------------
858 void CellTree::print_vector(ostream& out, vector<double>& v, TString desc, bool end)
859 {
860  int N = v.size();
861 
862  out << '"' << desc << '"' << ":[";
863  for (int i=0; i<N; i++) {
864  out << v[i];
865  if (i!=N-1) {
866  out << ",";
867  }
868  }
869  out << "]";
870  if (!end) out << ",";
871  out << endl;
872 }
873 
874 //-----------------------------------------------------------------------
876 {
877  // map track id to track index in the array
878  for (int i=0; i<mc_Ntrack; i++) {
879  trackIndex[mc_id[i]] = i;
880  }
881 
882  // in trackParents, trackChildren, trackSiblings vectors, store track index (not track id)
883  for (int i=0; i<mc_Ntrack; i++) {
884  // currently, parent size == 1;
885  // for primary particle, parent id = 0;
886  vector<int> parents;
887  if ( !IsPrimary(i) ) {
888  parents.push_back(trackIndex[mc_mother[i]]);
889  }
890  trackParents.push_back(parents); // primary track will have 0 parents
891 
892  vector<int> children;
893  int nChildren = mc_daughters.at(i).size();
894  for (int j=0; j<nChildren; j++) {
895  children.push_back(trackIndex[mc_daughters.at(i).at(j)]);
896  }
897  trackChildren.push_back(children);
898 
899  }
900 
901  // siblings
902  for (int i=0; i<mc_Ntrack; i++) {
903  vector<int> siblings;
904  if ( IsPrimary(i) ) {
905  for (int j=0; j<mc_Ntrack; j++) {
906  if( IsPrimary(j) ) {
907  siblings.push_back(j);
908  }
909  }
910  }
911  else {
912  // siblings are simply children of the mother
913  int mother = trackIndex[mc_mother[i]];
914  int nSiblings = trackChildren.at(mother).size();
915  for (int j=0; j<nSiblings; j++) {
916  siblings.push_back(trackChildren.at(mother).at(j));
917  }
918  }
919  trackSiblings.push_back(siblings);
920  }
921 
922 }
923 
924 //-----------------------------------------------------------------------
926 {
927  art::Handle< std::vector<raw::Trigger>> triggerListHandle;
928  std::vector<art::Ptr<raw::Trigger>> triggerlist;
929  if (event.getByLabel(fTriggerLabel, triggerListHandle)) {
930  art::fill_ptr_vector(triggerlist, triggerListHandle);
931  }
932  else {
933  cout << "WARNING: no label " << fTriggerLabel << endl;
934  }
935  if (triggerlist.size()){
936  fTriggernumber = triggerlist[0]->TriggerNumber();
937  fTriggertime = triggerlist[0]->TriggerTime();
938  fBeamgatetime = triggerlist[0]->BeamGateTime();
939  fTriggerbits = triggerlist[0]->TriggerBits();
940  }
941  else {
942  fTriggernumber = 0;
943  fTriggertime = 0;
944  fBeamgatetime = 0;
945  fTriggerbits = 0;
946  }
947 }
948 
949 //-----------------------------------------------------------------------
950 bool CellTree::DumpMCJSON(int id, ostream& out)
951 {
952  int i = trackIndex[id];
953  if (!KeepMC(i)) return false;
954 
955  int e = KE(mc_startMomentum[i])*1000;
956 
957  int nDaughter = mc_daughters.at(i).size();
958  vector<int> saved_daughters;
959  for (int j=0; j<nDaughter; j++) {
960  int daughter_id = mc_daughters.at(i).at(j);
961  // int e_daughter = KE(mc_startMomentum[ trackIndex[daughter_id] ])*1000;
962  // if (e_daughter >= thresh_KE) {
963  if ( KeepMC(trackIndex[daughter_id]) ) {
964  saved_daughters.push_back(daughter_id);
965  }
966  }
967 
968  out << fixed << setprecision(1);
969  out << "{";
970 
971  out << "\"id\":" << id << ",";
972  out << "\"text\":" << "\"" << PDGName(mc_pdg[i]) << " " << e << " MeV\",";
973  out << "\"data\":{";
974  out << "\"start\":[" << mc_startXYZT[i][0] << ", " << mc_startXYZT[i][1] << ", " << mc_startXYZT[i][2] << "],";
975  out << "\"end\":[" << mc_endXYZT[i][0] << ", " << mc_endXYZT[i][1] << ", " << mc_endXYZT[i][2] << "]";
976  out << "},";
977  out << "\"children\":[";
978  int nSavedDaughter = saved_daughters.size();
979  if (nSavedDaughter == 0) {
980  out << "],";
981  out << "\"icon\":" << "\"jstree-file\"";
982  out << "}";
983  return true;
984  }
985  else {
986  for (int j=0; j<nSavedDaughter; j++) {
987  DumpMCJSON(saved_daughters.at(j), out);
988  if (j!=nSavedDaughter-1) {
989  out << ",";
990  }
991  }
992  out << "]";
993  out << "}";
994  return true;
995  }
996 }
997 
998 //-----------------------------------------------------------------------
999 void CellTree::DumpMCJSON(ostream& out)
1000 {
1001  out << "[";
1002  vector<int> primaries;
1003  for (int i=0; i<mc_Ntrack; i++) {
1004  if (IsPrimary(i)) {
1005  // int e = KE(mc_startMomentum[i])*1000;
1006  // if (e<thresh_KE) continue;
1007  if (KeepMC(i)) {
1008  primaries.push_back(i);
1009  }
1010  }
1011  }
1012  int size = primaries.size();
1013  // cout << size << endl;
1014  for (int i=0; i<size; i++) {
1015  if (DumpMCJSON(mc_id[primaries[i]], out) && i!=size-1) {
1016  out << ", ";
1017  }
1018  }
1019 
1020  out << "]";
1021 }
1022 
1023 //-----------------------------------------------------------------------
1024 double CellTree::KE(float* momentum)
1025 {
1026  TLorentzVector particle(momentum);
1027  return particle.E()-particle.M();
1028 }
1029 
1030 //-----------------------------------------------------------------------
1031 bool CellTree::KeepMC(int i)
1032 {
1033  double e = KE(mc_startMomentum[i])*1000;
1034  double thresh_KE_em = 5.; // MeV
1035  double thresh_KE_np = 10; // MeV
1036  if (mc_pdg[i]==22 || mc_pdg[i]==11 || mc_pdg[i]==-11) {
1037  if (e>=thresh_KE_em) return true;
1038  else return false;
1039  }
1040  else if (mc_pdg[i]==2112 || mc_pdg[i]==2212 || mc_pdg[i]>1e9) {
1041  if (e>=thresh_KE_np) return true;
1042  else return false;
1043  }
1044  return true;
1045 }
1046 
1047 //-----------------------------------------------------------------------
1048 TString CellTree::PDGName(int pdg)
1049 {
1050  TParticlePDG *p = dbPDG->GetParticle(pdg);
1051  if (p == 0) {
1052  if (pdg>1e9) {
1053  int z = (pdg - 1e9) / 10000;
1054  int a = (pdg - 1e9 - z*1e4) / 10;
1055  TString name;
1056  if (z == 18) name = "Ar";
1057 
1058  else if (z == 17) name = "Cl";
1059  else if (z == 19) name = "Ca";
1060  else if (z == 16) name = "S";
1061  else if (z == 15) name = "P";
1062  else if (z == 14) name = "Si";
1063  else if (z == 1) name = "H";
1064  else if (z == 2) name = "He";
1065 
1066  else return Form("%i", pdg);
1067  return Form("%s-%i", name.Data(), a);
1068  }
1069  return Form("%i", pdg);
1070  }
1071  else {
1072  return p->GetName();
1073  }
1074 }
1075 
1076 //-----------------------------------------------------------------------
1078 {
1079  cout << " Run/SubRun/Event: " << fRun << "/" << fSubRun << "/" << fEvent << endl;
1080  cout << " Ntracks:" << mc_Ntrack << endl;
1081 
1082  for (int i=0; i<mc_Ntrack; i++) {
1083  cout << "\n id: " << mc_id[i];
1084  cout << "\n pdg: " << mc_pdg[i];
1085  cout << "\n mother: " << mc_mother[i];
1086  cout << "\n Ndaughters: " << mc_daughters.at(i).size();
1087  cout << "\n start XYZT: (" << mc_startXYZT[i][0] << ", " << mc_startXYZT[i][1] << ", " << mc_startXYZT[i][2] << ", " << mc_startXYZT[i][3] << ")";
1088  cout << "\n end XYZT: (" << mc_endXYZT[i][0] << ", " << mc_endXYZT[i][1] << ", " << mc_endXYZT[i][2] << ", " << mc_endXYZT[i][3] << ")";
1089  cout << "\n start momentum: (" << mc_startMomentum[i][0] << ", " << mc_startMomentum[i][1] << ", " << mc_startMomentum[i][2] << ", " << mc_startMomentum[i][3] << ")";
1090  cout << "\n end momentum: (" << mc_endMomentum[i][0] << ", " << mc_endMomentum[i][1] << ", " << mc_endMomentum[i][2] << ", " << mc_endMomentum[i][3] << ")";
1091 
1092  cout << endl;
1093  }
1094 }
1095 
1096 //-----------------------------------------------------------------------
1098 {
1099  processMap["unknown"] = 0;
1100  processMap["primary"] = 1;
1101  processMap["compt"] = 2;
1102  processMap["phot"] = 3;
1103  processMap["annihil"] = 4;
1104  processMap["eIoni"] = 5;
1105  processMap["eBrem"] = 6;
1106  processMap["conv"] = 7;
1107  processMap["muIoni"] = 8;
1108  processMap["muMinusCaptureAtRest"] = 9;
1109  processMap["NeutronInelastic"] = 10;
1110  processMap["nCapture"] = 11;
1111  processMap["hadElastic"] = 12;
1112  processMap["Decay"] = 13;
1113  processMap["CoulombScat"] = 14;
1114  processMap["muPairProd"] = 15;
1115  processMap["muBrems"] = 16;
1116  processMap["muPairProd"] = 17;
1117  processMap["PhotonInelastic"] = 18;
1118  processMap["hIoni"] = 19;
1119  processMap["ProtonInelastic"] = 20;
1120  processMap["PionPlusInelastic"] = 21;
1121  processMap["CHIPSNuclearCaptureAtRest"] = 22;
1122  processMap["PionMinusInelastic"] = 23;
1123 }
1124 
1125 //-----------------------------------------------------------------------
1127 } // namespace microboone
1128 
1129 
1130 #endif
Float_t x
Definition: compare.C:6
float mc_startMomentum[MAX_TRACKS][4]
Store parameters for running LArG4.
std::vector< std::vector< int > > trackChildren
vector< double > oh_pe
vector< int > fSIMIDE_channelIdY
const TLorentzVector & Position(const int i=0) const
Definition: MCParticle.h:223
int PdgCode() const
Definition: MCParticle.h:216
int CCNC() const
Definition: MCNeutrino.h:152
constexpr std::uint32_t timeLow() const
Definition: Timestamp.h:30
double QSqr() const
Definition: MCNeutrino.h:161
double Theta() const
angle between incoming and outgoing leptons, in radians
Definition: MCNeutrino.cxx:63
const simb::MCNeutrino & GetNeutrino() const
Definition: MCTruth.h:74
system("rm -rf microbeam.root")
std::map< int, int > trackIndex
vector< int > fSIMIDE_trackId
int HitQuark() const
Definition: MCNeutrino.h:157
void reconfigure(fhicl::ParameterSet const &pset)
float mc_endXYZT[MAX_TRACKS][4]
void processMCTracks()
void processRaw(const art::Event &evt)
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
int mc_id[MAX_TRACKS]
TString PDGName(int pdg)
Declaration of signal hit object.
Float_t y
Definition: compare.C:6
const simb::MCParticle & Nu() const
Definition: MCNeutrino.h:150
constexpr std::uint32_t timeHigh() const
Definition: Timestamp.h:35
Double_t z
Definition: plot.C:279
simb::Origin_t Origin() const
Definition: MCTruth.h:71
double Pt() const
transverse momentum of interaction, in GeV/c
Definition: MCNeutrino.cxx:74
int HitNuc() const
Definition: MCNeutrino.h:156
int mc_pdg[MAX_TRACKS]
std::vector< std::vector< int > > trackSiblings
vector< int > of_multiplicity
STL namespace.
Definition of basic raw digits.
vector< float > fSIMIDE_y
int NParticles() const
Definition: MCTruth.h:72
void processSimChannel(const art::Event &evt)
Particle class.
std::string fOutFileName
std::string fRawDigitLabel
Definition: Run.h:30
vector< float > fSIMIDE_x
bool DumpMCJSON(int id, ostream &out)
int mc_mother[MAX_TRACKS]
vector< unsigned short > fSIMIDE_tdc
std::map< int, int > savedMCTrackIdMap
float mc_endMomentum[MAX_TRACKS][4]
vector< int > oh_channel
std::string fOpHitLabel
std::string fOpFlashLabel
int InteractionType() const
Definition: MCNeutrino.h:154
int mc_process[MAX_TRACKS]
void analyze(const art::Event &evt)
void print_vector(ostream &out, vector< double > &v, TString desc, bool end=false)
double W() const
Definition: MCNeutrino.h:158
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:42
void beginJob()
Definition: Breakpoints.cc:14
vector< double > oh_bgtime
bool KeepMC(int i)
std::string DetectorName() const
Returns a string with the name of the detector, as configured.
double Y() const
Definition: MCNeutrino.h:160
TDatabasePDG * dbPDG
void processOpFlash(const art::Event &evt)
vector< double > oh_trigtime
TClonesArray * fCalib_wf
Collect all the RawData header files together.
T get(std::string const &key) const
Definition: ParameterSet.h:231
TObjArray * fMC_trackPosition
std::vector< int > fRaw_channelId
double X() const
Definition: MCNeutrino.h:159
vector< float > fSIMIDE_z
Declaration of cluster object.
unsigned int fTriggerbits
Provides recob::Track data product.
Definition of data types for geometry description.
#define MAX_TRACKS
void processOpHit(const art::Event &evt)
std::vector< std::vector< int > > mc_daughters
TClonesArray * fRaw_wf
void beginRun(const art::Run &run)
unsigned int NOpDets() const
Number of OpDets in the whole detector.
unsigned int fTriggernumber
double KE(float *momentum)
std::string mcOption
std::vector< int > fCalib_channelId
virtual ~CellTree()
int Target() const
Definition: MCNeutrino.h:155
art::ServiceHandle< geo::Geometry > fGeometry
void processCalib(const art::Event &evt)
Encapsulate the construction of a single detector plane.
std::string fCalibLabel
bool getByLabel(std::string const &label, std::string const &productInstanceName, Handle< PROD > &result) const
Definition: DataViewImpl.h:344
float mc_nu_mom[4]
void processMC(const art::Event &evt)
std::vector< std::string > fSpacePointLabels
const TLorentzVector & Momentum(const int i=0) const
Definition: MCParticle.h:224
object containing MC truth information necessary for making RawDigits and doing back tracking ...
TClonesArray * fPEperOpDet
vector< float > fSIMIDE_numElectrons
std::string fTriggerLabel
std::map< std::string, int > processMap
Declaration of basic channel signal object.
float mc_startXYZT[MAX_TRACKS][4]
std::vector< evd::details::RawDigitInfo_t >::const_iterator end(RawDigitCacheDataClass const &cache)
bool NeutrinoSet() const
Definition: MCTruth.h:75
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:464
void Uncompress(const std::vector< short > &adc, std::vector< short > &uncompressed, raw::Compress_t compress)
Uncompresses a raw data buffer.
Definition: raw.cxx:756
Float_t e
Definition: plot.C:34
void processSpacePoint(const art::Event &event, TString option, ostream &out=cout)
bool IsPrimary(int i)
Event generator information.
Definition: MCNeutrino.h:18
vector< float > of_peTotal
int Mode() const
Definition: MCNeutrino.h:153
float mc_nu_pos[4]
art framework interface to geometry description
Event finding and building.
void processTrigger(const art::Event &evt)
vector< float > of_t
std::vector< std::vector< int > > trackParents