All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Namespaces | Classes | Typedefs | Enumerations | Functions
numu Namespace Reference

Namespaces

 flat
 

Classes

struct  CRTMatch
 
struct  CRTHit
 
struct  FlashTriggerPrimitive
 
struct  FlashMatch
 
struct  RecoSlice
 
struct  RecoInteraction
 
struct  RecoEvent
 
struct  RecoParticle
 
struct  MCSFitResult
 
struct  RecoTrack
 
struct  TrueParticle
 
struct  TrackTruthMatch
 
struct  TruthMatch
 
class  MuonS2NStudy
 
class  MuPVertexStudy
 
class  NuMuEfficiencyStudy
 

Typedefs

using TrackFunction = std::function< uscript::Value(const numu::RecoTrack *, const numu::TrueParticle *, const unsigned *)>
 
using TrackSelector = std::function< bool(const numu::RecoTrack &, const numu::TrueParticle &, const unsigned &)>
 

Enumerations

enum  MCType { fUnknown =0, fOverlay =1, fIntimeCosmic =2, fCosmic =3 }
 
enum  InteractionMode {
  mCC = 0, mCCNonPrimary = 1, mNC = 2, mNCNonPrimary = 3,
  mCosmic = 4, mIntimeCosmic = 5, mOther = 6, mAll = 7
}
 
enum  TrackMode { tmOther = -1, tmCosmic = 1, tmNeutrino = 2 }
 
enum  Wall {
  wNone =0, wTop =1, wBottom =2, wLeft =3,
  wRight =4, wFront =5, wBack =6
}
 

Functions

float dist2Match (const event::Interaction &truth, const std::vector< numu::RecoInteraction > &candidates)
 
float trackMatchCompletion (unsigned truth_index, const numu::RecoEvent &event)
 
std::vector< std::string > MultiplyNames (const std::vector< std::vector< std::string >> &strings, char join='_')
 
std::vector< TrackSelectorMultiplyTrackSelectors (const std::vector< std::vector< std::string >> &track_function_strings)
 
int SelectLongestTrack (const std::map< size_t, RecoTrack > &tracks, const RecoSlice &slice)
 
int SelectLongestIDdMuon (const std::map< size_t, RecoTrack > &tracks, const RecoSlice &slice)
 
float MeanTruncateddQdx (const anab::Calorimetry &calo)
 
float TrackMomentum (const numu::RecoTrack &track)
 
float RangeMomentum (const numu::RecoTrack &track)
 
float MCSMomentum (const numu::RecoTrack &track)
 
TruthMatch InteractionTruthMatch (const std::vector< event::Interaction > &truth, const std::map< size_t, numu::RecoTrack > &reco_tracks, const numu::RecoInteraction &reco)
 
void CorrectMultiMatches (RecoEvent &event, std::vector< RecoInteraction > &recos)
 
InteractionMode GetMode (const event::Interaction &truth)
 
std::vector
< FlashTriggerPrimitive
TriggerPrimitives (const std::vector< raw::OpDetWaveform > &waveforms, double tick_period, std::pair< double, double > &window, int thresh, bool is_sbnd)
 
bool HasTrigger (const std::vector< FlashTriggerPrimitive > &primitives, int threshold, unsigned n_above_threshold)
 
std::vector< int > TriggerThresholds (const std::vector< numu::FlashTriggerPrimitive > &primitives, unsigned size)
 

Typedef Documentation

using numu::TrackFunction = typedef std::function<uscript::Value (const numu::RecoTrack *, const numu::TrueParticle *, const unsigned *)>

Definition at line 12 of file DynamicSelector.h.

using numu::TrackSelector = typedef std::function<bool (const numu::RecoTrack &, const numu::TrueParticle &, const unsigned &)>

Definition at line 13 of file DynamicSelector.h.

Enumeration Type Documentation

Enum to hold each different typoe of reconstructed event

Enumerator
mCC 
mCCNonPrimary 
mNC 
mNCNonPrimary 
mCosmic 
mIntimeCosmic 
mOther 
mAll 

Definition at line 8 of file Mode.h.

8  {
9  mCC = 0,
10  mCCNonPrimary = 1,
11  mNC = 2,
12  mNCNonPrimary = 3,
13  mCosmic = 4,
14  mIntimeCosmic = 5,
15  mOther = 6,
16  mAll = 7
17 };
Definition: Mode.h:11
Definition: Mode.h:9
Enumerator
fUnknown 
fOverlay 
fIntimeCosmic 
fCosmic 

Definition at line 5 of file MCType.h.

5  {
6  fUnknown=0,
7  fOverlay=1,
9  fCosmic=3
10 };

Enum to hold different types of tracks.

Enumerator
tmOther 
tmCosmic 
tmNeutrino 

Definition at line 22 of file Mode.h.

22  {
23  tmOther = -1,
24  tmCosmic = 1,
25  tmNeutrino = 2
26 };
enum numu::Wall
Enumerator
wNone 
wTop 
wBottom 
wLeft 
wRight 
wFront 
wBack 

Definition at line 10 of file TrueParticle.h.

10  {
11  wNone=0,
12  wTop=1,
13  wBottom=2,
14  wLeft=3,
15  wRight=4,
16  wFront=5,
17  wBack=6
18 };

Function Documentation

void numu::CorrectMultiMatches ( numu::RecoEvent event,
std::vector< RecoInteraction > &  recos 
)

Corrects a list of reco interaction objects when some are matched to the same truth neutrino interaction by making some of them non-primary

Parameters
eventThe reconstructed event
recosThe list of reconstructed vertices. The "match" field of each RecoInteraction may be edited by this function

Definition at line 110 of file TruthMatch.cc.

110  {
111  // check for no two reco matched to same truth
112  std::map<int, unsigned> matched_vertices;
113  for (unsigned reco_i = 0; reco_i < recos.size(); reco_i++) {
114  numu::RecoInteraction &reco = recos[reco_i];
115  unsigned set_i = reco_i;
116 
117  numu::RecoTrack &primary_track = event.tracks.at(reco.slice.primary_track_index);
118 
119  if (primary_track.match.has_match && primary_track.match.is_primary && reco.match.mctruth_track_id >= 0) {
120  if (matched_vertices.count(reco.match.mctruth_track_id)) {
121  std::cout << "BAAAAAAAAD: two reco matched to same truth." << std::endl;
122  std::cout << "This id: " << reco.match.mctruth_track_id << std::endl;
123  for (const numu::RecoInteraction &reco: recos) {
124  std::cout << "Interaction x: " << reco.position.X() << " match: " << reco.match.mctruth_track_id << " primary track pdg: " << primary_track.match.match_pdg << std::endl;
125  }
126 
127  numu::RecoInteraction &other = recos[matched_vertices[reco.match.mctruth_track_id]];
128  numu::RecoTrack &other_primary_track = event.tracks.at(other.slice.primary_track_index);
129 
130  float dist_reco = std::min( (primary_track.start - event.particles.at(primary_track.match.mcparticle_id).start).Mag(),
131  (primary_track.end - event.particles.at(primary_track.match.mcparticle_id).start).Mag());
132 
133  float dist_other = std::min( (other_primary_track.start - event.particles.at(other_primary_track.match.mcparticle_id).start).Mag(),
134  (other_primary_track.end - event.particles.at(other_primary_track.match.mcparticle_id).start).Mag());
135 
136  // Handle if one of the matched particles in not a muon
137  if (abs(primary_track.match.match_pdg) != 13 && abs(other_primary_track.match.match_pdg) == 13) {
138  set_i = matched_vertices[reco.match.mctruth_track_id];
139  primary_track.match.is_primary = false;
140  std::cout << "Corrected -- matched track is non-muon\n";
141  }
142  else if (abs(primary_track.match.match_pdg) == 13 && abs(other_primary_track.match.match_pdg) != 13) {
143  other_primary_track.match.is_primary = false;
144  std::cout << "Corrected -- matched track is non-muon\n";
145  }
146  // Handle both two different particles neither of which are muons
147  else if (primary_track.match.mcparticle_id != other_primary_track.match.mcparticle_id) {
148  std::cout << "Reco track 1 length: " << primary_track.length << " pdg: " << primary_track.match.match_pdg << std::endl;
149  std::cout << "Reco track 2 length: " << other_primary_track.length << " pdg: " << other_primary_track.match.match_pdg << std::endl;
150  if (primary_track.length > other_primary_track.length) {
151  other_primary_track.match.is_primary = false;
152  std::cout << "Corrected -- used longer track. Track 1 is primary.\n";
153  }
154  else {
155  primary_track.match.is_primary = false;
156  std::cout << "Corrected -- used longer track. Track 2 is primary.\n";
157  }
158  }
159  // Handle if the priamry track was split in two
160  else if (primary_track.match.mcparticle_id == other_primary_track.match.mcparticle_id) {
161  std::cout << "Reco track 1 pos: " << primary_track.start.X() << " " << primary_track.start.Y() << " " << primary_track.start.Z() << " to: " << primary_track.end.X() << " " << primary_track.end.Y() << " " << primary_track.end.Z() << std::endl;
162  std::cout << "Reco track 2 pos: " << other_primary_track.start.X() << " " << other_primary_track.start.Y() << " " << other_primary_track.start.Z() << " to: " << other_primary_track.end.X() << " " << other_primary_track.end.Y() << " " << other_primary_track.end.Z() << std::endl;
163  TVector3 match_start = event.particles.at(other_primary_track.match.mcparticle_id).start;
164  std::cout << "Match start pos: " << match_start.X() << " " << match_start.Y() << " " << match_start.Z() << std::endl;
165  if (dist_reco <= dist_other) {
166  std::cout << "Corrected -- split muon. Track 1 is primary.\n";
167  other_primary_track.match.is_primary = false;
168  }
169  else {
170  std::cout << "Corrected -- split muon. Track 2 is primary.\n";
171  set_i = matched_vertices[reco.match.mctruth_track_id];
172  primary_track.match.is_primary = false;
173  }
174  }
175  else {
176  std::cout << "Unable to correct!\n";
177  std::cerr << "Exiting on failure.\n";
178  assert(false);
179  }
180  }
181  matched_vertices[reco.match.mctruth_track_id] = set_i;
182  }
183  }
184 }
BEGIN_PROLOG could also be cerr
TruthMatch match
Info for mathing to truth.
Definition: RecoEvent.h:42
int mcparticle_id
MCParticle ID of the particle this track matches to (same as the ID of the RecoTrack of that particle...
int primary_track_index
Index of the primary track.
Definition: RecoEvent.h:27
std::map< size_t, TrueParticle > particles
Map of indices to True particle information.
Definition: RecoEvent.h:49
T abs(T value)
TVector3 position
location of the vertex
Definition: RecoEvent.h:40
process_name standard_reco_uboone reco
TVector3 start
start position of track
Definition: RecoTrack.h:54
float length
Length of track.
Definition: RecoTrack.h:48
bool has_match
Whether a track match exists.
bool is_primary
Whether this track was produced as the &quot;primary&quot; process.
TrackTruthMatch match
Truth matching information.
Definition: RecoTrack.h:57
RecoSlice slice
Particle content of the interaction.
Definition: RecoEvent.h:39
TVector3 end
end position of track
Definition: RecoTrack.h:55
int mctruth_track_id
index of the primary track into the list of MCTruths. -1 if no match.
int match_pdg
PDG of the MCParticle this track matches to.
BEGIN_PROLOG could also be cout
float numu::dist2Match ( const event::Interaction truth,
const std::vector< numu::RecoInteraction > &  candidates 
)

Distance between one interation vertex and a list of candidate matching vertices

Parameters
truthThe true neutrino interaction to match
candidatesThe list of candidate vertices matching to this one
Returns
The minimum distance between the list of candidates and the vertex

Definition at line 3 of file Derived.cc.

3  {
4  float dist = -1;
5  for (const numu::RecoInteraction &candidate: candidates) {
6  float this_dist = (candidate.position - truth.neutrino.position).Mag();
7  if (dist < 0 || this_dist < dist) dist = this_dist;
8  }
9  return dist;
10 
11 }
Neutrino neutrino
The neutrino.
Definition: Event.hh:152
constexpr double dist(const TReal *x, const TReal *y, const unsigned int dimension)
TVector3 position
Neutrino interaction position.
Definition: Event.hh:92
numu::InteractionMode numu::GetMode ( const event::Interaction truth)

Definition at line 5 of file TruthMatch.cc.

5  {
6  if (truth.neutrino.iscc) return numu::mCC;
7  else return numu::mNC;
8 }
Definition: Mode.h:11
Neutrino neutrino
The neutrino.
Definition: Event.hh:152
bool iscc
CC (true) or NC/interference (false)
Definition: Event.hh:75
Definition: Mode.h:9
bool numu::HasTrigger ( const std::vector< FlashTriggerPrimitive > &  primitives,
int  threshold,
unsigned  n_above_threshold 
)

Definition at line 65 of file PMTTrigger.cc.

65  {
66  if (n_above_threshold == 0) return true;
67 
68  std::map<int, std::vector<unsigned>> above_threshold;
69 
70  for (const numu::FlashTriggerPrimitive &primitive: primitives) {
71  for (const numu::FlashTriggerPrimitive::Trig &trig: primitive.triggers) {
72  if (trig.adc <= threshold) {
73  above_threshold[trig.tdc].push_back(primitive.channel);
74  }
75  }
76  }
77 
78  for (auto const &pair: above_threshold) {
79  std::cout << "At time: " << pair.first << " above thresh: ";
80  for (unsigned ch: pair.second) std::cout << ch << " ";
81  std::cout << std::endl;
82  if (pair.second.size() >= n_above_threshold) {
83  return true;
84  }
85  }
86 
87  return false;
88 }
BEGIN_PROLOG could also be cout
numu::TruthMatch numu::InteractionTruthMatch ( const std::vector< event::Interaction > &  truth,
const std::map< size_t, numu::RecoTrack > &  reco_tracks,
const numu::RecoInteraction reco 
)

Returns the truth match interaction in a RecoInteraction

Parameters
truthThe list of true interactions in the event
reco_tracksThe set of reconstructed tracks in the event
recoThe reconstructed interaction in the event. The existing content of the reco.match instance does not affect the algorithm
Returns
The TruthMatch object for the RecoInteraction reco

Definition at line 10 of file TruthMatch.cc.

10  {
11  TruthMatch match;
12  // Try to match against truth events
13  // Vertices
14  match.truth_vertex_distance = -1;
15  for (int truth_i = 0; truth_i < truth.size(); truth_i++) {
16  // find the closest vertex
17  if (match.truth_vertex_distance < 0 || (truth[truth_i].neutrino.position - reco.position).Mag() < match.truth_vertex_distance) {
18  match.truth_vertex_distance = (truth[truth_i].neutrino.position - reco.position).Mag();
19  match.mctruth_vertex_id = truth_i;
20  }
21  }
22 
23  // TODO: better way to do matching
24  // only keep match if distance is less than 5cm
25  if (match.truth_vertex_distance < 0. || match.truth_vertex_distance > 5.) {
26  std::cout << "Vertex not matched well.\n";
27  }
28  else {
29  std::cout << "Vertex matched to neutrino.\n";
30  }
31 
32  // Match to the primary track of the event
33  match.mctruth_track_id = -1;
34  if (reco.slice.primary_track_index >= 0) {
35  const numu::TrackTruthMatch &ptrack_match = reco_tracks.at(reco.slice.primary_track_index).match;
36  if (ptrack_match.has_match) {
37  for (int truth_i = 0; truth_i < truth.size(); truth_i++) {
38  if ((truth[truth_i].neutrino.position - ptrack_match.mctruth_vertex).Mag() < 1.) {
39  match.mctruth_track_id = truth_i;
40  }
41  }
42  }
43  std::cout << "Primary track index: " << reco.slice.primary_track_index << std::endl;
44 
45  // return the interaction mode
46  if (ptrack_match.has_match && ptrack_match.mctruth_origin == simb::Origin_t::kCosmicRay) {
47  match.tmode = numu::tmCosmic;
48  std::cout << "Track matched to cosmic.\n";
49  }
50  else if (ptrack_match.has_match && ptrack_match.mctruth_origin == simb::Origin_t::kBeamNeutrino) {
51  match.tmode = numu::tmNeutrino;
52  std::cout << "Track matched to neutrino.\n";
53  }
54  else if (ptrack_match.has_match) {
55  match.tmode = numu::tmOther;
56  std::cout << "Track matched to other: " << ptrack_match.mctruth_origin << std::endl;
57  std::cout << "Match PDG: " << ptrack_match.match_pdg << std::endl;
58  std::cout << "Match completion: " << ptrack_match.completion << std::endl;
59  std::cout << "Is primary: " << ptrack_match.is_primary << std::endl;
60  const numu::RecoTrack &track = reco_tracks.at(reco.slice.primary_track_index);
61  std::cout << "reco start: " << track.start.X() << " " << track.start.Y() << " " << track.start.Z() << std::endl;
62  std::cout << "reco end: " << track.end.X() << " " << track.end.Y() << " " << track.end.Z() << std::endl;
63  }
64  else {
65  match.tmode = numu::tmOther;
66  std::cout << "Track not matched.\n";
67  }
68  }
69  else {
70  match.tmode = numu::tmOther;
71  std::cout << "No track to match.\n";
72  }
73 
74  // Use the track to match the mode of the interaction
75  if (reco.slice.primary_track_index >= 0 &&
76  reco_tracks.at(reco.slice.primary_track_index).match.has_match) {
77  const numu::TrackTruthMatch &ptrack_match = reco_tracks.at(reco.slice.primary_track_index).match;
78  if (ptrack_match.mctruth_origin == simb::Origin_t::kBeamNeutrino && ptrack_match.mctruth_ccnc == simb::kCC) {
79  if (ptrack_match.is_primary) {
80  match.mode = numu::mCC;
81  }
82  else {
83  match.mode = numu::mCCNonPrimary;
84  }
85  }
86  else if (ptrack_match.mctruth_origin == simb::Origin_t::kBeamNeutrino && ptrack_match.mctruth_ccnc == simb::kNC) {
87  if (ptrack_match.is_primary) {
88  match.mode = numu::mNC;
89  }
90  else {
91  match.mode = numu::mNCNonPrimary;
92  }
93  }
94  else if (ptrack_match.mctruth_origin == simb::Origin_t::kCosmicRay) {
95  match.mode = numu::mCosmic;
96  }
97  else {
98  match.mode = numu::mOther;
99  }
100  }
101  // **shrug**
102  else {
103  match.mode = numu::mOther;
104  }
105  match.has_match = true;
106 
107  return match;
108 }
int mctruth_ccnc
CC (1) / NC (0) value of the MCTruth this object matches to.
Definition: Mode.h:11
int mctruth_origin
Value of Origin_t enum of the MCTruth object this track matches to.
int primary_track_index
Index of the primary track.
Definition: RecoEvent.h:27
process_name use argoneut_mc_hitfinder track
Charged-current interactions.
Definition: IPrediction.h:38
TVector3 position
location of the vertex
Definition: RecoEvent.h:40
Definition: Mode.h:9
TVector3 start
start position of track
Definition: RecoTrack.h:54
TVector3 mctruth_vertex
The interaction vertex of the MCTruth object this track matches to.
float completion
Fraction of energy deposits by true particle matched by this track.
const Cut kCosmicRay
bool has_match
Whether a track match exists.
bool is_primary
Whether this track was produced as the &quot;primary&quot; process.
RecoSlice slice
Particle content of the interaction.
Definition: RecoEvent.h:39
Neutral-current interactions.
Definition: IPrediction.h:39
TVector3 end
end position of track
Definition: RecoTrack.h:55
int match_pdg
PDG of the MCParticle this track matches to.
BEGIN_PROLOG could also be cout
float numu::MCSMomentum ( const numu::RecoTrack track)

Definition at line 20 of file TrackAlgo.cc.

20  {
21  return track.mcs_muon.fwd_mcs_momentum;
22 }
float fwd_mcs_momentum
MCS momentum calculation under hypothesis track is forward.
Definition: RecoTrack.h:14
MCSFitResult mcs_muon
MCS calculation result for Muon PID hypothesis.
Definition: RecoTrack.h:35
float numu::MeanTruncateddQdx ( const anab::Calorimetry calo)

Definition at line 25 of file TrackAlgo.cc.

25  {
26  // copy the dQdx list for partial sorting to find median
27  std::vector<float> dQdx = calo.dQdx();
28 
29  // if no or only 1 calo points, then no dQdx
30  if (dQdx.size() <= 1) return -9999.;
31 
32  // calculate the median
33  float median = -1;
34  if (dQdx.size() % 2 == 0 /* even */) {
35  const auto median_lo = dQdx.begin() + dQdx.size()/2 - 1;
36  const auto median_hi = dQdx.begin() + dQdx.size()/2;
37  std::nth_element(dQdx.begin(), median_lo, dQdx.end());
38  float median_lo_val = *median_lo;
39  std::nth_element(dQdx.begin(), median_hi, dQdx.end());
40  float median_hi_val = *median_hi;
41  median = (median_lo_val + median_hi_val) / 2.;
42  }
43  else /* odd */ {
44  const auto median_ptr = dQdx.begin() + dQdx.size()/2;
45  std::nth_element(dQdx.begin(), median_ptr, dQdx.end());
46  median = *median_ptr;
47  }
48 
49  // get the mean and variance
50  float mean = std::accumulate(dQdx.begin(), dQdx.end(), 0.) / dQdx.size();
51  float var = 0.;
52  for (float d: dQdx) {
53  var += (d - mean) * (d - mean);
54  }
55  var = var / dQdx.size();
56 
57  // truncated mean
58  float trunc_sum = 0.;
59  unsigned n_point = 0;
60  for (float d: dQdx) {
61  if (var < (d - mean) * (d - mean)) {
62  trunc_sum += d;
63  n_point ++;
64  }
65  }
66  // if no points available, return garbage
67  if (n_point == 0) return -9999.;
68 
69  return trunc_sum / n_point;
70 }
const std::vector< float > & dQdx() const
Definition: Calorimetry.h:102
double mean(const std::vector< short > &wf, size_t start, size_t nsample)
Definition: UtilFunc.cxx:13
std::vector< std::string > numu::MultiplyNames ( const std::vector< std::vector< std::string >> &  strings,
char  join = '_' 
)

Definition at line 25 of file DynamicSelector.cc.

25  {
26  std::vector<std::string> ret;
27  std::vector<unsigned> indices (strings.size(), 0);
28  do {
29  std::string this_str;
30  for (unsigned i = 0; i < indices.size(); i++) {
31  this_str += strings[i][indices[i]] + join;
32  }
33  ret.push_back(this_str);
34  Increment(indices, strings);
35  } while (NotZero(indices));
36  return ret;
37 }
void Increment(std::vector< unsigned > &indices, const std::vector< std::vector< T >> &lists)
S join(S const &sep, Coll const &s)
Returns a concatenation of strings in s separated by sep.
constexpr std::array< std::size_t, geo::vect::dimension< Vector >)> indices()
Returns a sequence of indices valid for a vector of the specified type.
bool NotZero(const std::vector< unsigned > &indices)
std::vector< numu::TrackSelector > numu::MultiplyTrackSelectors ( const std::vector< std::vector< std::string >> &  track_function_strings)

Definition at line 39 of file DynamicSelector.cc.

39  {
40  std::vector<numu::TrackSelector> ret;
41  std::vector<std::vector<numu::TrackFunction>> track_functions;
42  for (unsigned i = 0; i < track_function_strings.size(); i++) {
43  track_functions.emplace_back();
44  for (unsigned j = 0; j < track_function_strings[i].size(); j++) {
45  track_functions[i].push_back(uscript::compile<numu::RecoTrack, numu::TrueParticle, unsigned>("track", "particle", "mctype", track_function_strings[i][j].c_str()));
46  }
47  }
48  std::vector<unsigned> indices (track_functions.size(), 0);
49  do {
50  std::vector<numu::TrackFunction> functions;
51  for (unsigned i = 0; i < indices.size(); i++) {
52  functions.push_back(track_functions[i][indices[i]]);
53  }
54  ret.push_back(
55  [functions](const numu::RecoTrack &track, const numu::TrueParticle &part, unsigned mctype) {
56  for (const numu::TrackFunction &function: functions) {
57  uscript::Value ret = function(&track, &part, &mctype);
58  if (!ret) return false;
59  }
60  return true;
61  });
62  Increment(indices, track_functions);
63  } while (NotZero(indices));
64  return ret;
65 }
process_name use argoneut_mc_hitfinder track
void Increment(std::vector< unsigned > &indices, const std::vector< std::vector< T >> &lists)
std::function< uscript::Value(const numu::RecoTrack *, const numu::TrueParticle *, const unsigned *)> TrackFunction
constexpr std::array< std::size_t, geo::vect::dimension< Vector >)> indices()
Returns a sequence of indices valid for a vector of the specified type.
bool NotZero(const std::vector< unsigned > &indices)
float numu::RangeMomentum ( const numu::RecoTrack track)

Definition at line 15 of file TrackAlgo.cc.

15  {
16  return track.range_momentum_muon;
17 }
float range_momentum_muon
Range momentum calculation of the track using range under the assumption it is a muon [GeV]...
Definition: RecoTrack.h:33
int numu::SelectLongestIDdMuon ( const std::map< size_t, RecoTrack > &  tracks,
const RecoSlice &  slice 
)

Algorithm to select the primary track – get the longest one that has a better particle ID for muon than proton

Parameters
tracksAll reconstructed tracks in the event
sliceThe RecoSlice to select the primary track for
Returns
The ID of the primary track. Returns -1 if no such track exists.

Definition at line 22 of file PrimaryTrack.cc.

22  {
23  if (slice.primary_index < 0) return -1;
24 
25  const numu::RecoParticle &neutrino = slice.particles.at(slice.primary_index);
26  double max_len = -1.;
27  int ret = -1;
28  for (size_t pfp_index: neutrino.daughters) {
29  const numu::RecoParticle &daughter = slice.particles.at(pfp_index);
30  if (tracks.count(daughter.ID)) {
31  const numu::RecoTrack &track = tracks.at(daughter.ID);
32  // muon ID: exiting or chi2 hypothesis smaller than proton
33  if (!track.is_contained || track.chi2_muon < track.chi2_proton) {
34  if (ret < 0 || track.length > max_len) {
35  max_len = track.length;
36  ret = daughter.ID;
37  }
38  }
39  }
40  }
41 
42  return ret;
43 }
bool is_contained
is it contained in the &quot;containment volume&quot;?
Definition: RecoTrack.h:52
float chi2_muon
Chi2 of dE/dx to muon hypotheis. Combined agaisnt all planes.
Definition: RecoTrack.h:43
process_name use argoneut_mc_hitfinder track
size_t ID
ID of particle.
Definition: RecoParticle.h:18
std::vector< size_t > daughters
Daughters of the particle in the &quot;particle flow&quot;. Value represents index into pandora information...
Definition: RecoParticle.h:17
float length
Length of track.
Definition: RecoTrack.h:48
float chi2_proton
Chi2 of dE/dx to proton hypothesis. Combined against all planes.
Definition: RecoTrack.h:40
int numu::SelectLongestTrack ( const std::map< size_t, RecoTrack > &  tracks,
const RecoSlice &  slice 
)

Algorithm to select the primary track – get the longest one

Parameters
tracksAll reconstructed tracks in the event
sliceThe RecoSlice to select the primary track for
Returns
The ID of the primary track. Returns -1 if no such track exists.

Definition at line 3 of file PrimaryTrack.cc.

3  {
4  if (slice.primary_index < 0) return -1;
5 
6  const numu::RecoParticle &neutrino = slice.particles.at(slice.primary_index);
7  double max_len = -1.;
8  int ret = -1;
9  for (size_t pfp_index: neutrino.daughters) {
10  const numu::RecoParticle &daughter = slice.particles.at(pfp_index);
11  if (tracks.count(daughter.ID)) {
12  if (ret < 0 || tracks.at(daughter.ID).length > max_len) {
13  max_len = tracks.at(daughter.ID).length;
14  ret = daughter.ID;
15  }
16  }
17  }
18 
19  return ret;
20 }
size_t ID
ID of particle.
Definition: RecoParticle.h:18
std::vector< size_t > daughters
Daughters of the particle in the &quot;particle flow&quot;. Value represents index into pandora information...
Definition: RecoParticle.h:17
float numu::trackMatchCompletion ( unsigned  truth_index,
const numu::RecoEvent event 
)

Get the completon of a reconstructed track matching to a true track

Parameters
truth_indexThe index of the true particle to match
eventThe reconstructed event
Returns
The maximum completion of all reconstructed tracks to the true track. -1 if there is no matching reconstructed track

Definition at line 13 of file Derived.cc.

13  {
14  for (const numu::RecoInteraction &reco: event.reco) {
15  const numu::RecoTrack &primary_track = event.tracks.at(reco.slice.primary_track_index);
16  if (reco.match.mctruth_track_id == truth_index && primary_track.match.is_primary) {
17  return primary_track.match.completion;
18  }
19  }
20  return -1;
21 }
TruthMatch match
Info for mathing to truth.
Definition: RecoEvent.h:42
int primary_track_index
Index of the primary track.
Definition: RecoEvent.h:27
std::vector< RecoInteraction > reco
List of reconstructed vertices.
Definition: RecoEvent.h:50
process_name standard_reco_uboone reco
float completion
Fraction of energy deposits by true particle matched by this track.
bool is_primary
Whether this track was produced as the &quot;primary&quot; process.
TrackTruthMatch match
Truth matching information.
Definition: RecoTrack.h:57
RecoSlice slice
Particle content of the interaction.
Definition: RecoEvent.h:39
int mctruth_track_id
index of the primary track into the list of MCTruths. -1 if no match.
float numu::TrackMomentum ( const numu::RecoTrack track)

Definition at line 5 of file TrackAlgo.cc.

5  {
6  if (track.is_contained) {
7  return numu::RangeMomentum(track);
8  }
9  else {
10  return numu::MCSMomentum(track);
11  }
12 }
bool is_contained
is it contained in the &quot;containment volume&quot;?
Definition: RecoTrack.h:52
float RangeMomentum(const numu::RecoTrack &track)
Definition: TrackAlgo.cc:15
float MCSMomentum(const numu::RecoTrack &track)
Definition: TrackAlgo.cc:20
std::vector< numu::FlashTriggerPrimitive > numu::TriggerPrimitives ( const std::vector< raw::OpDetWaveform > &  waveforms,
double  tick_period,
std::pair< double, double > &  window,
int  thresh,
bool  is_sbnd 
)

Definition at line 7 of file PMTTrigger.cc.

7  {
8  std::vector<numu::FlashTriggerPrimitive> ret;
9 
10  opdet::sbndPDMapAlg mapalg;
11  for (const raw::OpDetWaveform &wvf: waveforms) {
12  // check if this is a PMT
13  bool is_pmt = (!is_sbnd) || mapalg.isPDType(wvf.ChannelNumber(), "pmt");
14 
15  if (!is_pmt) continue;
16 
17  // look for any overlap with the enable window
18  double waveform_start = wvf.TimeStamp();
19  int waveform_index_start = std::max((int)((window.first - waveform_start) / tick_period), 0);
20  int waveform_index_end = std::min((int)((window.second - waveform_start) / tick_period), (int) wvf.size());
21 
22  if (waveform_index_start < waveform_index_end) {
24  prim.channel = wvf.ChannelNumber();
25  for (int i = waveform_index_start; i < waveform_index_end; i++) {
26  if (wvf[i] <= thresh) { // PMT waveforms go down
27  FlashTriggerPrimitive::Trig this_trig {wvf[i], i + (int)((waveform_start - window.first) / tick_period)};
28  prim.triggers.push_back(this_trig);
29  }
30  }
31  ret.push_back(prim);
32  }
33  }
34  return ret;
35 }
bool isPDType(size_t ch, std::string pdname) const override
std::vector< Trig > triggers
Definition: DetInfo.h:26
int tick_period
Definition: constants.py:3
std::vector< int > numu::TriggerThresholds ( const std::vector< numu::FlashTriggerPrimitive > &  primitives,
unsigned  size 
)

Definition at line 37 of file PMTTrigger.cc.

37  {
38  std::map<int, std::vector<int>> trigs;
39  for (const numu::FlashTriggerPrimitive &primitive: primitives) {
40  for (const numu::FlashTriggerPrimitive::Trig &trig: primitive.triggers) {
41  trigs[trig.tdc].push_back(trig.adc);
42  }
43  }
44 
45  std::map<int, std::pair<unsigned, int>> time_sorted;
46  for (const auto &pair: trigs) {
47  time_sorted[pair.first].first = pair.second.size();
48  time_sorted[pair.first].second = *std::max_element(pair.second.begin(), pair.second.end());
49  }
50 
51  std::vector<int> ret;
52  for (unsigned i = 1; i <= size; i++) {
53  int min = 8000;
54  for (const auto &pair: time_sorted) {
55  if (pair.second.first >= i) {
56  if (pair.second.second < min) min = pair.second.second;
57  }
58  }
59  ret.push_back(min);
60  }
61 
62  return ret;
63 }
std::size_t size(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:561