All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Private Attributes | List of all members
ana::SBNOsc::Flatten Class Reference

#include <Flatten.h>

Inheritance diagram for ana::SBNOsc::Flatten:
core::PostProcessorBase

Public Member Functions

void FileSetup (TFile *f, TTree *eventTree)
 
void Initialize (fhicl::ParameterSet *config)
 
void ProcessEvent (const event::Event *event)
 
void InitializeThread ()
 
void ProcessSubRun (const SubRun *subrun)
 
void ProcessFileMeta (const FileMeta *meta)
 
void Finalize ()
 
- Public Member Functions inherited from core::PostProcessorBase
 PostProcessorBase ()
 
virtual ~PostProcessorBase ()
 
void Run (std::vector< std::string > filelist)
 
void Initialize (char *config=NULL, const std::string &output_fname="", unsigned n_workers=1)
 

Private Attributes

Cuts fCuts
 
std::vector< TFile * > fOutputFiles
 
std::vector< TNtuple * > fNtuples
 
std::vector
< numu::flat::FlatInteraction
fInteractions
 
std::vector< numu::MCTypefMCTypes
 
std::vector< numu::RecoEvent * > fRecoEvents
 

Additional Inherited Members

- Protected Member Functions inherited from core::PostProcessorBase
virtual void FileCleanup (TTree *eventTree)
 
unsigned NWorkers ()
 
unsigned WorkerID ()
 
- Protected Attributes inherited from core::PostProcessorBase
ProviderManagerfProviderManager
 Interface for provider access. More...
 
int fConfigExperimentID
 

Detailed Description

Definition at line 30 of file Flatten.h.

Member Function Documentation

void ana::SBNOsc::Flatten::FileSetup ( TFile *  f,
TTree *  eventTree 
)
virtual

Setup anything needed per file

Parameters
fthe TFile being opened
eventTreethe TTree associated with the sbncode event. Use this TTree to set branch addresses for everything other than the sbncode event.

Files are guaranteed to be processed in the order they are specified on the command line for sbn-postprocess

Reimplemented from core::PostProcessorBase.

Definition at line 107 of file Flatten.cc.

107  {
108  unsigned ind = WorkerID();
109 
110  eventTree->SetBranchAddress("reco_event", &fRecoEvents[ind]);
111 
112  TParameter<int> *file_type = (TParameter<int> *)f->Get("MCType");
113  fMCTypes[ind] = (numu::MCType) file_type->GetVal();
114 
115  fInteractions[ind].meta.mc_type = fMCTypes[ind];
116  fInteractions[ind].meta.detector = fConfigExperimentID;
117 }
std::vector< numu::flat::FlatInteraction > fInteractions
Definition: Flatten.h:46
std::vector< numu::RecoEvent * > fRecoEvents
Definition: Flatten.h:48
std::vector< numu::MCType > fMCTypes
Definition: Flatten.h:47
MCType
Definition: MCType.h:5
void ana::SBNOsc::Flatten::Finalize ( )
virtual

Perform user-level finalization. Called after all events have been processed.

Reimplemented from core::PostProcessorBase.

Definition at line 198 of file Flatten.cc.

198  {
199  TList *merge = new TList;
200  fOutputFiles[0]->cd();
201  for (unsigned i = 1; i < NWorkers(); i++) {
202  merge->Add(fNtuples[i]);
203  }
204  fNtuples[0]->Merge(merge);
205  fNtuples[0]->Write();
206  delete merge;
207 
208  for (unsigned i = 1; i < NWorkers(); i++) {
209  const char *name = fOutputFiles[i]->GetName();
210  fOutputFiles[i]->Close();
211  gSystem->Unlink(name);
212  }
213 
214 }
std::vector< TFile * > fOutputFiles
Definition: Flatten.h:44
then echo fcl name
std::vector< TNtuple * > fNtuples
Definition: Flatten.h:45
void ana::SBNOsc::Flatten::Initialize ( fhicl::ParameterSet *  config)
virtual

Perform user-level initialization.

Parameters
configA configuration, as a JSON object.

Implements core::PostProcessorBase.

Definition at line 82 of file Flatten.cc.

82  {
83  fCuts.Initialize(config->get<fhicl::ParameterSet>("Cuts"), fProviderManager->GetGeometryProvider());
84  std::string output_name = config->get<std::string>("OutputFile", "output.root");
85 
86  fOutputFiles = std::vector<TFile *>(NWorkers());
87  fNtuples = std::vector<TNtuple *>(NWorkers());
88  fInteractions = std::vector<numu::flat::FlatInteraction>(NWorkers());
89  fRecoEvents = std::vector<numu::RecoEvent *>(NWorkers(), 0);
90  fMCTypes = std::vector<numu::MCType>(NWorkers());
91  std::string ntuple = ClasstoNTupleName("numu::flat::FlatInteraction", "interaction");
92 
93  for (unsigned i = 0; i < NWorkers(); i++) {
94  std::string name = output_name;
95  if (i > 0) name += std::to_string(i);
96  fOutputFiles[i] = new TFile(name.c_str(), "RECREATE");
97  fOutputFiles[i]->cd();
98  fNtuples[i] = new TNtuple("interaction", "interaction", ntuple.c_str());
99  }
100 }
std::vector< numu::flat::FlatInteraction > fInteractions
Definition: Flatten.h:46
const geo::GeometryCore * GetGeometryProvider() const
std::vector< TFile * > fOutputFiles
Definition: Flatten.h:44
std::vector< numu::RecoEvent * > fRecoEvents
Definition: Flatten.h:48
BEGIN_PROLOG ntuple
Definition: genie.fcl:8
std::vector< numu::MCType > fMCTypes
Definition: Flatten.h:47
std::string to_string(WindowPattern const &pattern)
std::string ClasstoNTupleName(const std::string &classname, const std::string &prefix)
Definition: Flatten.cc:71
then echo fcl name
std::vector< TNtuple * > fNtuples
Definition: Flatten.h:45
ProviderManager * fProviderManager
Interface for provider access.
void Initialize(const fhicl::ParameterSet &cfg, const geo::GeometryCore *geometry)
Definition: Cuts.cc:10
void ana::SBNOsc::Flatten::InitializeThread ( )
virtual

Perform user-level initialization per-thread. Only called when multiple workers are specified on the command line.

Reimplemented from core::PostProcessorBase.

Definition at line 102 of file Flatten.cc.

102  {
103  unsigned i = WorkerID();
104  fOutputFiles[i]->cd();
105 }
std::vector< TFile * > fOutputFiles
Definition: Flatten.h:44
void ana::SBNOsc::Flatten::ProcessEvent ( const event::Event event)
virtual

Process one event.

Parameters
eventThe sbncode event for the current event

Implements core::PostProcessorBase.

Definition at line 119 of file Flatten.cc.

119  {
120  unsigned index = WorkerID();
121  ana::SBNOsc::SetEvent(*fRecoEvents[index], *ev, fCuts, fMCTypes[index], true);
122  for (const numu::RecoInteraction &interaction: fRecoEvents[index]->reco) {
123  const numu::RecoTrack &primary_track = fRecoEvents[index]->tracks.at(interaction.slice.primary_track_index);
124  // set stuff
125  //
126  // Primary Track
128  track.length = primary_track.length;
129  track.costh = primary_track.costh;
130  track.range_momentum = numu::RangeMomentum(primary_track);
131  track.mcs_momentum = numu::MCSMomentum(primary_track);
132  track.crt_hit_distance = (primary_track.crt_match.hit_match.present) ? primary_track.crt_match.hit_match.distance : -1;
133  track.crt_hit_time = primary_track.crt_match.hit_match.time;
134  track.crt_track_angle = (primary_track.crt_match.track.present) ? primary_track.crt_match.track.angle: -1;
135  track.crt_track_time = primary_track.crt_match.track.time;
136  primary_track.start.GetXYZ(track.start);
137  primary_track.end.GetXYZ(track.end);
138  if (primary_track.match.has_match) {
139  const numu::TrueParticle &truth = fRecoEvents[index]->particles.at(primary_track.match.mcparticle_id);
140  truth.start_momentum.GetXYZ(track.truth.momentum);
141  track.truth.wall_enter = truth.wall_enter;
142  track.truth.wall_exit = truth.wall_exit;
143  track.truth.time = truth.start_time;
144  track.truth.is_cosmic = truth.is_cosmic;
145  track.truth.pdgid = truth.pdgid;
146  track.truth.is_contained = truth.is_contained;
147  }
148  else track.truth.pdgid = -1;
149 
150  // Slice
151  numu::flat::Slice slice;
152  slice.flash_score = interaction.slice.flash_match.score;
153  slice.flash_time = interaction.slice.flash_match.time;
154 
155  // True Neutrino
156  numu::flat::TrueNeutrino neutrino;
157  if (interaction.match.has_match) {
158  neutrino.E = ev->truth[interaction.match.mctruth_track_id].neutrino.energy;
159  neutrino.Q2 = ev->truth[interaction.match.mctruth_track_id].neutrino.Q2;
160  ev->truth[interaction.match.mctruth_track_id].neutrino.position.GetXYZ(neutrino.vertex);
161  }
162  else neutrino.E = -1;
163 
164  // Event Info
165  numu::flat::EventInfo event_info;
166  for (unsigned i = 0; i < 10; i++) event_info.crt_hit_pes[i] = -1;
167  // Event Info
168  unsigned i = 0;
169  for (const numu::CRTHit &hit: fRecoEvents[index]->in_time_crt_hits) {
170  assert(i < 10);
171  event_info.crt_hit_pes[i] = hit.pes;
172  event_info.crt_hit_times[i] = hit.time;
173  i += 1;
174  }
175  event_info.pass_trig = fCuts.PassFlashTrigger(*fRecoEvents[index]);
176 
177  // Event Meta data -- set previously
178 
179  fNtuples[index]->Fill((float *)&fInteractions[index]);
180 
181  // overwrite normalization stuff so it isn't duplicated
182  fInteractions[index].meta.pot = 0;
183  fInteractions[index].meta.n_gen_events = 0;
184  }
185 
186 }
Track track
CRT Track match.
Definition: CRTMatch.h:30
TVector3 start_momentum
Particle directional momentum for first trajectory point inside TPC AV [GeV].
Definition: TrueParticle.h:24
float distance
//!&lt; Distance from projected track to CRT Hit. Nonsense if present is false.
Definition: CRTMatch.h:26
bool is_cosmic
Whether this particle is of cosmic origin.
Definition: TrueParticle.h:45
CRTMatch crt_match
CRTMatch.
Definition: RecoTrack.h:59
bool present
Whether this CRTMatch has a matching track.
Definition: CRTMatch.h:16
int mcparticle_id
MCParticle ID of the particle this track matches to (same as the ID of the RecoTrack of that particle...
std::vector< numu::flat::FlatInteraction > fInteractions
Definition: Flatten.h:46
process_name use argoneut_mc_hitfinder track
process_name hit
Definition: cheaterreco.fcl:51
bool is_contained
is it contained in a single cryostat active volume
Definition: TrueParticle.h:40
std::vector< numu::RecoEvent * > fRecoEvents
Definition: Flatten.h:48
process_name standard_reco_uboone reco
float RangeMomentum(const numu::RecoTrack &track)
Definition: TrackAlgo.cc:15
TVector3 start
start position of track
Definition: RecoTrack.h:54
for($it=0;$it< $RaceTrack_number;$it++)
float length
Length of track.
Definition: RecoTrack.h:48
float MCSMomentum(const numu::RecoTrack &track)
Definition: TrackAlgo.cc:20
std::vector< numu::MCType > fMCTypes
Definition: Flatten.h:47
bool present
Whether this CRTMatch has a matching hit.
Definition: CRTMatch.h:25
bool PassFlashTrigger(const numu::RecoEvent &event) const
Definition: Cuts.cc:121
bool has_match
Whether a track match exists.
TrackTruthMatch match
Truth matching information.
Definition: RecoTrack.h:57
int pdgid
Particle ID code.
Definition: TrueParticle.h:43
float angle
Angle between TPC track and CRT track.
Definition: CRTMatch.h:18
TVector3 end
end position of track
Definition: RecoTrack.h:55
std::vector< TNtuple * > fNtuples
Definition: Flatten.h:45
void SetEvent(numu::RecoEvent &event, const event::Event &core, const ana::SBNOsc::Cuts &cuts, numu::MCType file_type, bool use_calorimetry=true)
Definition: SetEvent.cc:5
Wall wall_enter
the face of the TPC that the particle crosses on enter
Definition: TrueParticle.h:35
Wall wall_exit
the face of the TPC that the particle crosses on exit
Definition: TrueParticle.h:36
HitMatch hit_match
CRT Hit match.
Definition: CRTMatch.h:31
float costh
cosine of angle to z axis
Definition: RecoTrack.h:49
float start_time
start time of track
Definition: TrueParticle.h:32
float time
Matching time [us] of track. T==0 is set to beam spill start time.
Definition: CRTMatch.h:17
void ana::SBNOsc::Flatten::ProcessFileMeta ( const FileMeta meta)
virtual

Reimplemented from core::PostProcessorBase.

Definition at line 193 of file Flatten.cc.

193  {
194  unsigned index = WorkerID();
195  fInteractions[index].meta.n_gen_events += meta->n_gen_events;
196 }
std::vector< numu::flat::FlatInteraction > fInteractions
Definition: Flatten.h:46
unsigned n_gen_events
Definition: FileMeta.hh:25
void ana::SBNOsc::Flatten::ProcessSubRun ( const SubRun subrun)
virtual

Reimplemented from core::PostProcessorBase.

Definition at line 188 of file Flatten.cc.

188  {
189  unsigned index = WorkerID();
190  fInteractions[index].meta.pot += subrun->totgoodpot;
191 }
std::vector< numu::flat::FlatInteraction > fInteractions
Definition: Flatten.h:46
float totgoodpot
Definition: SubRun.hh:40

Member Data Documentation

Cuts ana::SBNOsc::Flatten::fCuts
private

Definition at line 43 of file Flatten.h.

std::vector<numu::flat::FlatInteraction> ana::SBNOsc::Flatten::fInteractions
private

Definition at line 46 of file Flatten.h.

std::vector<numu::MCType> ana::SBNOsc::Flatten::fMCTypes
private

Definition at line 47 of file Flatten.h.

std::vector<TNtuple *> ana::SBNOsc::Flatten::fNtuples
private

Definition at line 45 of file Flatten.h.

std::vector<TFile *> ana::SBNOsc::Flatten::fOutputFiles
private

Definition at line 44 of file Flatten.h.

std::vector<numu::RecoEvent *> ana::SBNOsc::Flatten::fRecoEvents
private

Definition at line 48 of file Flatten.h.


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