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