Public Member Functions | Private Attributes | List of all members
pdsp::EMCNNCheck Class Reference
Inheritance diagram for pdsp::EMCNNCheck:
art::EDAnalyzer art::detail::Analyzer art::detail::LegacyModule art::Observer art::ModuleBase

Public Member Functions

 EMCNNCheck (fhicl::ParameterSet const &p)
 
 EMCNNCheck (EMCNNCheck const &)=delete
 
 EMCNNCheck (EMCNNCheck &&)=delete
 
EMCNNCheckoperator= (EMCNNCheck const &)=delete
 
EMCNNCheckoperator= (EMCNNCheck &&)=delete
 
void analyze (art::Event const &e) override
 
void beginJob () override
 
- Public Member Functions inherited from art::EDAnalyzer
 EDAnalyzer (fhicl::ParameterSet const &pset)
 
template<typename Config >
 EDAnalyzer (Table< Config > const &config)
 
std::string workerType () const
 
- Public Member Functions inherited from art::detail::Analyzer
virtual ~Analyzer () noexcept
 
 Analyzer (fhicl::ParameterSet const &pset)
 
template<typename Config >
 Analyzer (Table< Config > const &config)
 
void doBeginJob (SharedResources const &resources)
 
void doEndJob ()
 
void doRespondToOpenInputFile (FileBlock const &fb)
 
void doRespondToCloseInputFile (FileBlock const &fb)
 
void doRespondToOpenOutputFiles (FileBlock const &fb)
 
void doRespondToCloseOutputFiles (FileBlock const &fb)
 
bool doBeginRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doEndRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doBeginSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEndSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEvent (EventPrincipal &ep, ModuleContext const &mc, std::atomic< std::size_t > &counts_run, std::atomic< std::size_t > &counts_passed, std::atomic< std::size_t > &counts_failed)
 
- Public Member Functions inherited from art::Observer
 ~Observer () noexcept
 
 Observer (Observer const &)=delete
 
 Observer (Observer &&)=delete
 
Observeroperator= (Observer const &)=delete
 
Observeroperator= (Observer &&)=delete
 
void registerProducts (ProductDescriptions &, ModuleDescription const &)
 
void fillDescriptions (ModuleDescription const &)
 
fhicl::ParameterSetID selectorConfig () const
 
- Public Member Functions inherited from art::ModuleBase
virtual ~ModuleBase () noexcept
 
 ModuleBase ()
 
ModuleDescription const & moduleDescription () const
 
void setModuleDescription (ModuleDescription const &)
 
std::array< std::vector< ProductInfo >, NumBranchTypes > const & getConsumables () const
 
void sortConsumables (std::string const &current_process_name)
 
template<typename T , BranchType BT>
ViewToken< T > consumesView (InputTag const &tag)
 
template<typename T , BranchType BT>
ViewToken< T > mayConsumeView (InputTag const &tag)
 

Private Attributes

std::string fGeneratorTag
 
std::string fCNNTag
 
protoana::ProtoDUNEBeamCuts beam_cuts
 
protoana::ProtoDUNEBeamlineUtils fBeamlineUtils
 
TTree * ftree
 
int run
 
int subrun
 
int event
 
int beampdg
 
double average_score_em
 
double average_score_trk
 
double average_score_mic
 
double track_endz
 
int ndaughterhits
 
double average_daughter_score_mic
 
double vtxx
 
double vtxy
 
double vtxz
 
double endx
 
double endy
 
double endz
 
double dirx
 
double diry
 
double dirz
 
std::vector< short > channel
 
std::vector< short > tpc
 
std::vector< short > plane
 
std::vector< short > wire
 
std::vector< double > charge
 
std::vector< double > peakt
 
std::vector< double > score_em
 
std::vector< double > score_trk
 
std::vector< double > score_mic
 
std::vector< short > daughter_channel
 
std::vector< short > daughter_tpc
 
std::vector< short > daughter_plane
 
std::vector< short > daughter_wire
 
std::vector< double > daughter_charge
 
std::vector< double > daughter_peakt
 
std::vector< double > daughter_score_em
 
std::vector< double > daughter_score_trk
 
std::vector< double > daughter_score_mic
 
std::vector< int > pdg
 
std::vector< int > origin
 
std::vector< std::stringprocess
 

Additional Inherited Members

- Public Types inherited from art::EDAnalyzer
using WorkerType = WorkerT< EDAnalyzer >
 
using ModuleType = EDAnalyzer
 
- Protected Member Functions inherited from art::Observer
std::string const & processName () const
 
bool wantAllEvents () const noexcept
 
bool wantEvent (ScheduleID id, Event const &e) const
 
Handle< TriggerResultsgetTriggerResults (Event const &e) const
 
 Observer (fhicl::ParameterSet const &config)
 
 Observer (std::vector< std::string > const &select_paths, std::vector< std::string > const &reject_paths, fhicl::ParameterSet const &config)
 
- Protected Member Functions inherited from art::ModuleBase
ConsumesCollectorconsumesCollector ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > consumes (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > consumesView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void consumesMany ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > mayConsume (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > mayConsumeView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void mayConsumeMany ()
 

Detailed Description

Definition at line 46 of file EMCNNCheck_module.cc.

Constructor & Destructor Documentation

pdsp::EMCNNCheck::EMCNNCheck ( fhicl::ParameterSet const &  p)
explicit

Definition at line 113 of file EMCNNCheck_module.cc.

114  : EDAnalyzer{p},
115 //fselectpdg(p.get<int>("selectpdg")),
116  fGeneratorTag(p.get<std::string>("GeneratorTag")),
117  fCNNTag(p.get<std::string>("CNNTag","emtrkmichelid:emtrkmichel")),
118  beam_cuts(p.get<fhicl::ParameterSet>("BeamCuts")),
119  fBeamlineUtils(p.get<fhicl::ParameterSet>("BeamlineUtils"))
120  //BeamCuts(p.get<fhicl::ParameterSet>("BeamCuts"))
121  {
122  //beam_cuts = protoana::ProtoDUNEBeamCuts(BeamCuts);
123  }
std::string string
Definition: nybbler.cc:12
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.h:25
protoana::ProtoDUNEBeamlineUtils fBeamlineUtils
p
Definition: test.py:223
std::string fGeneratorTag
protoana::ProtoDUNEBeamCuts beam_cuts
pdsp::EMCNNCheck::EMCNNCheck ( EMCNNCheck const &  )
delete
pdsp::EMCNNCheck::EMCNNCheck ( EMCNNCheck &&  )
delete

Member Function Documentation

void pdsp::EMCNNCheck::analyze ( art::Event const &  e)
overridevirtual

Implements art::EDAnalyzer.

Definition at line 125 of file EMCNNCheck_module.cc.

126 {
127 
128  run = e.run();
129  subrun = e.subRun();
130  event = e.id().event();
131  beampdg = 0;
132  average_score_em = 0.;
133  average_score_trk = 0.;
134  average_score_mic = 0.;
135  track_endz = -1;
136  ndaughterhits = 0;
138  vtxx = -9999;
139  vtxy = -9999;
140  vtxz = -9999;
141  endx = -9999;
142  endy = -9999;
143  endz = -9999;
144  dirx = -9999;
145  diry = -9999;
146  dirz = -9999;
147  channel.clear();
148  tpc.clear();
149  plane.clear();
150  wire.clear();
151  charge.clear();
152  peakt.clear();
153  score_em.clear();
154  score_trk.clear();
155  score_mic.clear();
156 
157  daughter_channel.clear();
158  daughter_tpc.clear();
159  daughter_plane.clear();
160  daughter_wire.clear();
161  daughter_charge.clear();
162  daughter_peakt.clear();
163  daughter_score_em.clear();
164  daughter_score_trk.clear();
165  daughter_score_mic.clear();
166 
167  pdg.clear();
168  origin.clear();
169  process.clear();
170 
171  //Services
174  auto const* SCE = lar::providerFrom<spacecharge::SpaceChargeService>();
175  auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService>()->DataFor(e);
176  auto const detProp = art::ServiceHandle<detinfo::DetectorPropertiesService>()->DataFor(e, clockData);
178 
179  std::vector < art::Ptr < recob::Slice > > sliceList;
180  auto sliceListHandle = e.getHandle < std::vector < recob::Slice > > ("pandora");
181  if (sliceListHandle) {
182  art::fill_ptr_vector(sliceList, sliceListHandle);
183  }
184  else return;
185 
186  // Get all pfparticles
187  std::vector < art::Ptr < recob::PFParticle > > pfpList;
188  auto pfpListHandle = e.getHandle < std::vector < recob::PFParticle > >("pandora");
189  if (pfpListHandle) {
190  art::fill_ptr_vector(pfpList, pfpListHandle);
191  }
192 
193  // Get all clusters
194  std::vector < art::Ptr < recob::Cluster > > cluList;
195  auto cluListHandle = e.getHandle < std::vector < recob::Cluster > >("pandora");
196  if (cluListHandle) {
197  art::fill_ptr_vector(cluList, cluListHandle);
198  }
199 
200  // Get all tracks
201  std::vector < art::Ptr < recob::Track > > trkList;
202  auto trkListHandle = e.getHandle < std::vector < recob::Track > >("pandoraTrack");
203  if (trkListHandle) {
204  art::fill_ptr_vector(trkList, trkListHandle);
205  }
206 
207  // Get all hits
208  std::vector < art::Ptr < recob::Hit > > hitList;
209  auto hitListHandle = e.getHandle < std::vector < recob::Hit > >("hitpdune");
210  if (hitListHandle) {
211  art::fill_ptr_vector(hitList, hitListHandle);
212  }
213 
214  // Get cluster-PFParticle association
215  art::FindManyP<recob::Cluster> fmcpfp(pfpListHandle, e, "pandora");
216 
217  // Get vertex-PFParticle association
218  art::FindManyP<recob::Vertex> fmvpfp(pfpListHandle, e, "pandora");
219 
220  // Get hit-cluster association
221  art::FindManyP<recob::Hit> fmhc(cluListHandle, e, "pandora");
222 
223  art::FindManyP <recob::Hit> hitsFromSlice(sliceListHandle, e, "pandora");
224 
225  // Get track-hit association
226  art::FindManyP<recob::Hit, recob::TrackHitMeta> fmthm(trkListHandle, e,"pandoraTrack"); // to associate tracks and hits
227 
228  // Get hit-track association
229  art::FindManyP<recob::Track> thass(hitListHandle, e, "pandoraTrack"); //to associate hit just trying
230 
232 
233  if (!e.isRealData()){
234  // Get the truth utility to help us out
236  // Firstly we need to get the list of MCTruth objects from the generator. The standard protoDUNE
237  // simulation has fGeneratorTag = "generator"
238  auto mcTruths = e.getValidHandle<std::vector<simb::MCTruth>>(fGeneratorTag);
239  // mcTruths is basically a pointer to an std::vector of simb::MCTruth objects. There should only be one
240  // of these, so we pass the first element into the function to get the good particle
241  const simb::MCParticle* geantGoodParticle = truthUtil.GetGeantGoodParticle((*mcTruths)[0],e);
242  if(geantGoodParticle != 0x0){
243  //std::cout << "Found GEANT particle corresponding to the good particle with pdg = " << geantGoodParticle->PdgCode() << std::endl;
244 // if (fselectpdg==211){
245 // if (geantGoodParticle->PdgCode()!=211 &&
246 // geantGoodParticle->PdgCode()!=13){
247 // return;
248 // }
249 // }
250 // else{
251 // if (geantGoodParticle->PdgCode()!=fselectpdg) return;
252 // }
253  beampdg = geantGoodParticle->PdgCode();
254  }
255  }
256  else{
257  //Access the Beam Event
258  auto beamHandle = e.getValidHandle<std::vector<beam::ProtoDUNEBeamEvent>>("beamevent");
259 
260  std::vector<art::Ptr<beam::ProtoDUNEBeamEvent>> beamVec;
261  if( beamHandle.isValid()){
262  art::fill_ptr_vector(beamVec, beamHandle);
263  }
264 
265  if (beamVec.empty()) return;
266 
267  const beam::ProtoDUNEBeamEvent & beamEvent = *(beamVec.at(0)); //Should just have one
268  /////////////////////////////////////////////////////////////
269 
270 
271  //Check the quality of the event
272  std::cout << "Timing Trigger: " << beamEvent.GetTimingTrigger() << std::endl;
273  std::cout << "Is Matched: " << beamEvent.CheckIsMatched() << std::endl << std::endl;
274 
276  std::cout << "Failed quality check" << std::endl;
277  return;
278  }
279  //Access PID
280  std::vector< int > pids = fBeamlineUtils.GetPID( beamEvent, 1. );
281 // bool foundparticle = false;
282 // for( size_t i = 0; i < pids.size(); ++i ){
283 // //std::cout << pids[i] << std::endl;
284 // if (pids[i] == fselectpdg) foundparticle = true;
285 // }
286 // if (!foundparticle) return;
287  if (pids.empty()) return;
288  beampdg = pids[0];
289  }
290 
291  if (!beampdg) return;
292 
293  //std::cout<<"Found pion"<<std::endl;
294 
295  // Get the PFParticle utility
297 
298  std::vector<const recob::PFParticle*> beamParticles = pfpUtil.GetPFParticlesFromBeamSlice(e,"pandora");
299 
300  if(beamParticles.size() == 0){
301  std::cerr << "We found no beam particles for this event... moving on" << std::endl;
302  return;
303  }
304 
305  // We can now look at these particles
306  int trackid = -1;
307  int endwire = -1;
308  int endtpc = -1;
309  double endpeakt = -1;
310  std::vector<int> wirekeys;
311  for(const recob::PFParticle* particle : beamParticles){
312 
313  const recob::Track* thisTrack = pfpUtil.GetPFParticleTrack(*particle,e,"pandora","pandoraTrack");
314  const recob::Shower* thisShower = pfpUtil.GetPFParticleShower(*particle,e,"pandora","pandoraShower");
315  if (thisTrack){
316  //if (!beam_cuts.IsBeamlike(*thisTrack, e, "1")) return;
317  // Track ID
318  trackid = thisTrack->ID();
319  // Track end point z
320  track_endz = thisTrack->End().Z();
321  vtxx = thisTrack->Vertex().X();
322  vtxy = thisTrack->Vertex().Y();
323  vtxz = thisTrack->Vertex().Z();
324  if (!geom->FindTPCAtPosition(geo::Point_t(vtxx, vtxy, vtxz)).isValid) return;
325  geo::Vector_t offset = {0., 0., 0.};
326  if (SCE->EnableCalSpatialSCE()){
327  offset = SCE->GetCalPosOffsets(geo::Point_t(vtxx, vtxy, vtxz), (geom->FindTPCAtPosition(geo::Point_t(vtxx, vtxy, vtxz))).TPC);
328  }
329  //std::cout<<"track "<<offset.X()<<" "<<offset.Y()<<" "<<offset.Z()<<std::endl;
330  vtxx -= offset.X();
331  vtxy += offset.Y();
332  vtxz += offset.Z();
333  endx = thisTrack->End().X();
334  endy = thisTrack->End().Y();
335  endz = thisTrack->End().Z();
336  if (!geom->FindTPCAtPosition(geo::Point_t(endx, endy, endz)).isValid) return;
337  offset = {0., 0., 0.};
338  if (SCE->EnableCalSpatialSCE()){
339  offset = SCE->GetCalPosOffsets(geo::Point_t(endx, endy, endz), (geom->FindTPCAtPosition(geo::Point_t(endx, endy, endz))).TPC);
340  }
341  endx -= offset.X();
342  endy += offset.Y();
343  endz += offset.Z();
344  TVector3 dir(endx-vtxx, endy-vtxy, endz-vtxz);
345  dir = dir.Unit();
346  dirx = dir.X();
347  diry = dir.Y();
348  dirz = dir.Z();
349  // Find the last wire number and peak time on the track
350  if (fmthm.isValid()){
351  float zlast0=-99999;
352  auto vhit=fmthm.at(trackid);
353  auto vmeta=fmthm.data(trackid);
354  for (size_t ii = 0; ii<vhit.size(); ++ii){ //loop over all meta data hit
355  bool fBadhit = false;
356  if (vmeta[ii]->Index() == static_cast<unsigned int>(std::numeric_limits<int>::max())){
357  fBadhit = true;
358  continue;
359  }
360  if (vmeta[ii]->Index()>=thisTrack->NumberTrajectoryPoints()){
361  throw cet::exception("Calorimetry_module.cc") << "Requested track trajectory index "<<vmeta[ii]->Index()<<" exceeds the total number of trajectory points "<<thisTrack->NumberTrajectoryPoints()<<" for track index "<<trackid<<". Something is wrong with the track reconstruction. Please contact tjyang@fnal.gov!!";
362  }
363  if (!thisTrack->HasValidPoint(vmeta[ii]->Index())){
364  fBadhit = true;
365  continue;
366  }
367  auto loc = thisTrack->LocationAtPoint(vmeta[ii]->Index());
368  if (fBadhit) continue; //HY::If BAD hit, skip this hit and go next
369  if (loc.Z()<-100) continue; //hit not on track
370  if(vhit[ii]->WireID().Plane==2){
371  wirekeys.push_back(vhit[ii].key());
372  float zlast=loc.Z();
373  if(zlast>zlast0){
374  zlast0=zlast;
375  endwire=vhit[ii]->WireID().Wire;
376  endpeakt=vhit[ii]->PeakTime();
377  endtpc=vhit[ii]->WireID().TPC;
378  }
379  }
380  }
381  }
382  }
383  if (thisShower){
384  //if (!beam_cuts.IsBeamlike(*thisShower, e, "1")) return;
385  vtxx = thisShower->ShowerStart().X();
386  vtxy = thisShower->ShowerStart().Y();
387  vtxz = thisShower->ShowerStart().Z();
388  endx = vtxx + thisShower->Direction().X();
389  endy = vtxy + thisShower->Direction().Y();
390  endz = vtxz + thisShower->Direction().Z();
391  if (!geom->FindTPCAtPosition(geo::Point_t(vtxx, vtxy, vtxz)).isValid) return;
392  geo::Vector_t offset = {0., 0., 0.};
393  if (SCE->EnableCalSpatialSCE()){
394  offset = SCE->GetCalPosOffsets(geo::Point_t(vtxx, vtxy, vtxz), (geom->FindTPCAtPosition(geo::Point_t(vtxx, vtxy, vtxz))).TPC);
395  }
396  //std::cout<<"shower "<<offset.X()<<" "<<offset.Y()<<" "<<offset.Z()<<std::endl;
397  vtxx -= offset.X();
398  vtxy += offset.Y();
399  vtxz += offset.Z();
400  if (!geom->FindTPCAtPosition(geo::Point_t(endx, endy, endz)).isValid) return;
401  offset = {0., 0., 0.};
402  if (SCE->EnableCalSpatialSCE()){
403  offset = SCE->GetCalPosOffsets(geo::Point_t(endx, endy, endz), (geom->FindTPCAtPosition(geo::Point_t(endx, endy, endz))).TPC);
404  }
405  endx -= offset.X();
406  endy += offset.Y();
407  endz += offset.Z();
408  TVector3 dir(endx-vtxx, endy-vtxy, endz-vtxz);
409  dir = dir.Unit();
410  dirx = dir.X();
411  diry = dir.Y();
412  dirz = dir.Z();
413  }
414  }
415 
416  //int sliceid = pfpUtil.GetBeamSlice(e, "pandora");
417 
418  //if (sliceid!=9999){
419  if (fmcpfp.isValid()){
420  // Get clusters associated with pfparticle
421  auto const& clusters = fmcpfp.at(beamParticles[0]->Self());
422  for (auto const & cluster : clusters){
423  if (fmhc.isValid()){
424  // Get hits associated with cluster
425  auto const& hits = fmhc.at(cluster.key());
426  //auto const& hits = hitsFromSlice.at(sliceid);
427  //std::cout<<hits.size()<<std::endl;
428  for (auto & hit : hits){
429  std::array<float,4> cnn_out = hitResults.getOutput(hit);
430  if (hit->WireID().Plane == 2){
431  channel.push_back(hit->Channel());
432  tpc.push_back(hit->WireID().TPC);
433  plane.push_back(hit->WireID().Plane);
434  wire.push_back(hit->WireID().Wire);
435  charge.push_back(hit->Integral());
436  peakt.push_back(hit->PeakTime());
437  score_em.push_back(cnn_out[hitResults.getIndex("em")]);
438  score_trk.push_back(cnn_out[hitResults.getIndex("track")]);
439  score_mic.push_back(cnn_out[hitResults.getIndex("michel")]);
440  int this_pdg = 0;
441  int this_origin = -1;
442  std::string this_process = "null";
443  if (!e.isRealData()){
444  int TrackID = 0;
445  std::map<int,double> trkide;
446  std::vector<sim::TrackIDE> TrackIDs = bt_serv->HitToEveTrackIDEs(clockData, hit);
447  for(size_t e = 0; e < TrackIDs.size(); ++e){
448  trkide[TrackIDs[e].trackID] += TrackIDs[e].energy;
449  }
450 
451  // Work out which IDE despoited the most charge in the hit if there was more than one.
452  double maxe = -1;
453  double tote = 0;
454  for (std::map<int,double>::iterator ii = trkide.begin(); ii!=trkide.end(); ++ii){
455  tote += ii->second;
456  if ((ii->second)>maxe){
457  maxe = ii->second;
458  TrackID = ii->first;
459  }
460  }
461  // Now have trackID, so get PdG code and T0 etc.
462  const simb::MCParticle *particle = pi_serv->TrackIdToParticle_P(TrackID);
463  if (particle){
464  this_pdg = particle->PdgCode();
465  this_process = particle->Process();
466  this_origin = pi_serv->ParticleToMCTruth_P(particle)->Origin();
467  }
468  }
469  pdg.push_back(this_pdg);
470  origin.push_back(this_origin);
471  process.push_back(this_process);
472  }
473  }
474  }
475  }
476  }
477 
478  // Get the average of the collection plane scores
479  unsigned int nCollectionHits = 0;
480  for(unsigned int h = 0; h < plane.size(); ++h){
481  if(plane.at(h) == 2){
482  ++nCollectionHits;
483  average_score_em += score_em.at(h);
486  }
487  }
488  if(nCollectionHits > 0){
489  average_score_em /= static_cast<double>(nCollectionHits);
490  average_score_trk /= static_cast<double>(nCollectionHits);
491  average_score_mic /= static_cast<double>(nCollectionHits);
492  }
493 
494  // Get the hits near the track end (Michel candidate)
495  //std::cout<<trackid<<std::endl;
496  if (trackid!=-1){
497  for(size_t hitl=0;hitl<hitList.size();hitl++){
498  std::array<float,4> cnn_out=hitResults.getOutput(hitList[hitl]);
499  auto & tracks = thass.at(hitList[hitl].key());
500  // hit not on the track
501  if (std::find(wirekeys.begin(), wirekeys.end(), hitl) != wirekeys.end()) continue;
502  // hit not on a long track
503  if (!tracks.empty() && int(tracks[0].key()) != trackid && trkList[tracks[0].key()]->Length()>25) continue;
504  int planeid=hitList[hitl]->WireID().Plane;
505  if (planeid!=2) continue;
506  int tpcid=hitList[hitl]->WireID().TPC;
507  if (tpcid!=endtpc) continue;
508  float peakth1=hitList[hitl]->PeakTime();
509  int wireh1=hitList[hitl]->WireID().Wire;
510  if(std::abs(wireh1-endwire)<15 && std::abs(peakth1-endpeakt)<100 && tpcid==endtpc){
511  daughter_channel.push_back(hitList[hitl]->Channel());
512  daughter_tpc.push_back(hitList[hitl]->WireID().TPC);
513  daughter_plane.push_back(hitList[hitl]->WireID().Plane);
514  daughter_wire.push_back(hitList[hitl]->WireID().Wire);
515  daughter_charge.push_back(hitList[hitl]->Integral());
516  daughter_peakt.push_back(hitList[hitl]->PeakTime());
517  daughter_score_em.push_back(cnn_out[hitResults.getIndex("em")]);
518  daughter_score_trk.push_back(cnn_out[hitResults.getIndex("track")]);
519  daughter_score_mic.push_back(cnn_out[hitResults.getIndex("michel")]);
520 
521  ++ndaughterhits;
522  average_daughter_score_mic += cnn_out[hitResults.getIndex("michel")];
523  //std::cout<<hitList[hitl]->WireID().Wire<<" "<<hitList[hitl]->PeakTime()<<" "<<hitList[hitl]->Integral()<<" "<<cnn_out[hitResults.getIndex("michel")]<<std::endl;
524  }
525  }
526  }
528 
529  if (!channel.empty()) ftree->Fill();
530 }
std::vector< double > daughter_score_trk
intermediate_table::iterator iterator
std::vector< int > origin
const TVector3 & ShowerStart() const
Definition: Shower.h:192
int PdgCode() const
Definition: MCParticle.h:212
const int & GetTimingTrigger() const
std::vector< int > pdg
std::vector< double > daughter_score_em
const simb::MCParticle * GetGeantGoodParticle(const simb::MCTruth &genTruth, const art::Event &evt) const
std::vector< short > wire
const recob::Shower * GetPFParticleShower(const recob::PFParticle &particle, art::Event const &evt, const std::string particleLabel, const std::string showerLabel) const
Get the shower associated to this particle. Returns a null pointer if not found.
std::string string
Definition: nybbler.cc:12
std::vector< TrackID > TrackIDs
const simb::MCParticle * TrackIdToParticle_P(int id) const
Point_t const & LocationAtPoint(size_t i) const
Definition: Track.h:126
simb::Origin_t Origin() const
Definition: MCTruth.h:74
bool isValid
Whether this ID points to a valid element.
Definition: geo_types.h:211
bool HasValidPoint(size_t i) const
Definition: Track.h:111
size_t NumberTrajectoryPoints() const
Various functions related to the presence and the number of (valid) points.
Definition: Track.h:102
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Vector_t
Type for representation of momenta in 3D space.
Definition: geo_vectors.h:164
string dir
unsigned int Index
Cluster finding and building.
std::string Process() const
Definition: MCParticle.h:215
std::vector< short > daughter_plane
const art::Ptr< simb::MCTruth > & ParticleToMCTruth_P(const simb::MCParticle *p) const
std::vector< double > daughter_charge
geo::TPCID FindTPCAtPosition(double const worldLoc[3]) const
Returns the ID of the TPC at specified location.
std::vector< double > score_em
T abs(T value)
const double e
IDparameter< geo::WireID > WireID
Member type of validated geo::WireID parameter.
protoana::ProtoDUNEBeamlineUtils fBeamlineUtils
def key(type, name=None)
Definition: graph.py:13
std::vector< double > peakt
const recob::Track * GetPFParticleTrack(const recob::PFParticle &particle, art::Event const &evt, const std::string particleLabel, const std::string trackLabel) const
Get the track associated to this particle. Returns a null pointer if not found.
Point_t const & Vertex() const
Definition: Track.h:124
double average_daughter_score_mic
const std::vector< const recob::PFParticle * > GetPFParticlesFromBeamSlice(art::Event const &evt, const std::string particleLabel) const
Return the pointers for the PFParticles in the beam slice. Returns an empty vector is no beam slice w...
const TVector3 & Direction() const
Definition: Shower.h:189
bool IsGoodBeamlineTrigger(art::Event const &evt) const
static int max(int a, int b)
std::string fGeneratorTag
Detector simulation of raw signals on wires.
std::vector< double > score_mic
std::vector< short > plane
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Point_t
Type for representation of position in physical 3D space.
Definition: geo_vectors.h:184
std::vector< short > channel
Definition: tracks.py:1
int ID() const
Definition: Track.h:198
std::vector< short > tpc
std::vector< int > GetPID(beam::ProtoDUNEBeamEvent const &beamevt, double nominal_momentum)
Hierarchical representation of particle flow.
Definition: PFParticle.h:44
std::vector< short > daughter_wire
std::vector< sim::TrackIDE > HitToEveTrackIDEs(detinfo::DetectorClocksData const &clockData, recob::Hit const &hit) const
Point_t const & End() const
Definition: Track.h:125
std::vector< double > daughter_peakt
std::vector< double > daughter_score_mic
std::vector< short > daughter_channel
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:297
recob::tracking::Plane Plane
Definition: TrackState.h:17
std::vector< double > score_trk
const bool & CheckIsMatched() const
std::vector< double > charge
std::vector< std::string > process
Track from a non-cascading particle.A recob::Track consists of a recob::TrackTrajectory, plus additional members relevant for a "fitted" track:
Definition: Track.h:49
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
QTextStream & endl(QTextStream &s)
std::vector< short > daughter_tpc
void pdsp::EMCNNCheck::beginJob ( )
overridevirtual

Reimplemented from art::EDAnalyzer.

Definition at line 532 of file EMCNNCheck_module.cc.

532  {
533 
534  art::ServiceHandle<art::TFileService> fileServiceHandle;
535  ftree = fileServiceHandle->make<TTree>("ftree", "hit info");
536  ftree->Branch("run", &run, "run/I");
537  ftree->Branch("event", &event, "event/I");
538  ftree->Branch("beampdg", &beampdg, "beampdg/I");
539  ftree->Branch("average_score_em" , &average_score_em , "average_score_em/D");
540  ftree->Branch("average_score_trk", &average_score_trk, "average_score_trk/D");
541  ftree->Branch("average_score_mic", &average_score_mic, "average_score_mic/D");
542  ftree->Branch("track_endz", &track_endz, "track_endz/D");
543  ftree->Branch("ndaughterhits", &ndaughterhits, "ndaughterhits/I");
544  ftree->Branch("average_daughter_score_mic", &average_daughter_score_mic, "average_daughter_score_mic/D");
545  ftree->Branch("vtxx", &vtxx, "vtxx/D");
546  ftree->Branch("vtxy", &vtxy, "vtxy/D");
547  ftree->Branch("vtxz", &vtxz, "vtxz/D");
548  ftree->Branch("endx", &endx, "endx/D");
549  ftree->Branch("endy", &endy, "endy/D");
550  ftree->Branch("endz", &endz, "endz/D");
551  ftree->Branch("dirx", &dirx, "dirx/D");
552  ftree->Branch("diry", &diry, "diry/D");
553  ftree->Branch("dirz", &dirz, "dirz/D");
554  ftree->Branch("channel", &channel);
555  ftree->Branch("tpc", &tpc);
556  ftree->Branch("plane", &plane);
557  ftree->Branch("wire", &wire);
558  ftree->Branch("charge", &charge);
559  ftree->Branch("peakt", &peakt);
560  ftree->Branch("score_em", &score_em);
561  ftree->Branch("score_trk", &score_trk);
562  ftree->Branch("score_mic", &score_mic);
563 
564  ftree->Branch("daughter_channel", &daughter_channel);
565  ftree->Branch("daughter_tpc", &daughter_tpc);
566  ftree->Branch("daughter_plane", &daughter_plane);
567  ftree->Branch("daughter_wire", &daughter_wire);
568  ftree->Branch("daughter_charge", &daughter_charge);
569  ftree->Branch("daughter_peakt", &daughter_peakt);
570  ftree->Branch("daughter_score_em", &daughter_score_em);
571  ftree->Branch("daughter_score_trk", &daughter_score_trk);
572  ftree->Branch("daughter_score_mic", &daughter_score_mic);
573 
574  ftree->Branch("pdg", &pdg);
575  ftree->Branch("origin", &origin);
576  ftree->Branch("process", &process);
577 
578 }
std::vector< double > daughter_score_trk
std::vector< int > origin
std::vector< int > pdg
std::vector< double > daughter_score_em
std::vector< short > wire
std::vector< short > daughter_plane
std::vector< double > daughter_charge
std::vector< double > score_em
std::vector< double > peakt
double average_daughter_score_mic
std::vector< double > score_mic
std::vector< short > plane
std::vector< short > channel
std::vector< short > tpc
std::vector< short > daughter_wire
std::vector< double > daughter_peakt
std::vector< double > daughter_score_mic
std::vector< short > daughter_channel
std::vector< double > score_trk
std::vector< double > charge
std::vector< std::string > process
Event finding and building.
std::vector< short > daughter_tpc
EMCNNCheck& pdsp::EMCNNCheck::operator= ( EMCNNCheck const &  )
delete
EMCNNCheck& pdsp::EMCNNCheck::operator= ( EMCNNCheck &&  )
delete

Member Data Documentation

double pdsp::EMCNNCheck::average_daughter_score_mic
private

Definition at line 82 of file EMCNNCheck_module.cc.

double pdsp::EMCNNCheck::average_score_em
private

Definition at line 77 of file EMCNNCheck_module.cc.

double pdsp::EMCNNCheck::average_score_mic
private

Definition at line 79 of file EMCNNCheck_module.cc.

double pdsp::EMCNNCheck::average_score_trk
private

Definition at line 78 of file EMCNNCheck_module.cc.

protoana::ProtoDUNEBeamCuts pdsp::EMCNNCheck::beam_cuts
private

Definition at line 69 of file EMCNNCheck_module.cc.

int pdsp::EMCNNCheck::beampdg
private

Definition at line 76 of file EMCNNCheck_module.cc.

std::vector<short> pdsp::EMCNNCheck::channel
private

Definition at line 86 of file EMCNNCheck_module.cc.

std::vector<double> pdsp::EMCNNCheck::charge
private

Definition at line 90 of file EMCNNCheck_module.cc.

std::vector<short> pdsp::EMCNNCheck::daughter_channel
private

Definition at line 96 of file EMCNNCheck_module.cc.

std::vector<double> pdsp::EMCNNCheck::daughter_charge
private

Definition at line 100 of file EMCNNCheck_module.cc.

std::vector<double> pdsp::EMCNNCheck::daughter_peakt
private

Definition at line 101 of file EMCNNCheck_module.cc.

std::vector<short> pdsp::EMCNNCheck::daughter_plane
private

Definition at line 98 of file EMCNNCheck_module.cc.

std::vector<double> pdsp::EMCNNCheck::daughter_score_em
private

Definition at line 102 of file EMCNNCheck_module.cc.

std::vector<double> pdsp::EMCNNCheck::daughter_score_mic
private

Definition at line 104 of file EMCNNCheck_module.cc.

std::vector<double> pdsp::EMCNNCheck::daughter_score_trk
private

Definition at line 103 of file EMCNNCheck_module.cc.

std::vector<short> pdsp::EMCNNCheck::daughter_tpc
private

Definition at line 97 of file EMCNNCheck_module.cc.

std::vector<short> pdsp::EMCNNCheck::daughter_wire
private

Definition at line 99 of file EMCNNCheck_module.cc.

double pdsp::EMCNNCheck::dirx
private

Definition at line 85 of file EMCNNCheck_module.cc.

double pdsp::EMCNNCheck::diry
private

Definition at line 85 of file EMCNNCheck_module.cc.

double pdsp::EMCNNCheck::dirz
private

Definition at line 85 of file EMCNNCheck_module.cc.

double pdsp::EMCNNCheck::endx
private

Definition at line 84 of file EMCNNCheck_module.cc.

double pdsp::EMCNNCheck::endy
private

Definition at line 84 of file EMCNNCheck_module.cc.

double pdsp::EMCNNCheck::endz
private

Definition at line 84 of file EMCNNCheck_module.cc.

int pdsp::EMCNNCheck::event
private

Definition at line 75 of file EMCNNCheck_module.cc.

protoana::ProtoDUNEBeamlineUtils pdsp::EMCNNCheck::fBeamlineUtils
private

Definition at line 70 of file EMCNNCheck_module.cc.

std::string pdsp::EMCNNCheck::fCNNTag
private

Definition at line 67 of file EMCNNCheck_module.cc.

std::string pdsp::EMCNNCheck::fGeneratorTag
private

Definition at line 66 of file EMCNNCheck_module.cc.

TTree* pdsp::EMCNNCheck::ftree
private

Definition at line 72 of file EMCNNCheck_module.cc.

int pdsp::EMCNNCheck::ndaughterhits
private

Definition at line 81 of file EMCNNCheck_module.cc.

std::vector<int> pdsp::EMCNNCheck::origin
private

Definition at line 107 of file EMCNNCheck_module.cc.

std::vector<int> pdsp::EMCNNCheck::pdg
private

Definition at line 106 of file EMCNNCheck_module.cc.

std::vector<double> pdsp::EMCNNCheck::peakt
private

Definition at line 91 of file EMCNNCheck_module.cc.

std::vector<short> pdsp::EMCNNCheck::plane
private

Definition at line 88 of file EMCNNCheck_module.cc.

std::vector<std::string> pdsp::EMCNNCheck::process
private

Definition at line 108 of file EMCNNCheck_module.cc.

int pdsp::EMCNNCheck::run
private

Definition at line 73 of file EMCNNCheck_module.cc.

std::vector<double> pdsp::EMCNNCheck::score_em
private

Definition at line 92 of file EMCNNCheck_module.cc.

std::vector<double> pdsp::EMCNNCheck::score_mic
private

Definition at line 94 of file EMCNNCheck_module.cc.

std::vector<double> pdsp::EMCNNCheck::score_trk
private

Definition at line 93 of file EMCNNCheck_module.cc.

int pdsp::EMCNNCheck::subrun
private

Definition at line 74 of file EMCNNCheck_module.cc.

std::vector<short> pdsp::EMCNNCheck::tpc
private

Definition at line 87 of file EMCNNCheck_module.cc.

double pdsp::EMCNNCheck::track_endz
private

Definition at line 80 of file EMCNNCheck_module.cc.

double pdsp::EMCNNCheck::vtxx
private

Definition at line 83 of file EMCNNCheck_module.cc.

double pdsp::EMCNNCheck::vtxy
private

Definition at line 83 of file EMCNNCheck_module.cc.

double pdsp::EMCNNCheck::vtxz
private

Definition at line 83 of file EMCNNCheck_module.cc.

std::vector<short> pdsp::EMCNNCheck::wire
private

Definition at line 89 of file EMCNNCheck_module.cc.


The documentation for this class was generated from the following file: