-
Notifications
You must be signed in to change notification settings - Fork 1
/
TSSelection.cxx
708 lines (598 loc) · 22 KB
/
TSSelection.cxx
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
#include <cassert>
#include <iostream>
#include <map>
#include <sstream>
#include <string>
#include <vector>
// includes for random draws from gaussian
#include <random>
#include "canvas/Utilities/InputTag.h"
#include "canvas/Persistency/Common/FindMany.h"
#include "nusimdata/SimulationBase/MCTruth.h"
#include "nusimdata/SimulationBase/GTruth.h"
#include "lardataobj/MCBase/MCShower.h"
#include "lardataobj/MCBase/MCTrack.h"
#include "uboone/EventWeight/MCEventWeight.h"
#include "gallery/Event.h"
#include <TDatabasePDG.h>
#include <TFile.h>
#include <TKey.h>
#include <TLorentzVector.h>
#include <TH2F.h>
#include <TNtuple.h>
#include <TTree.h>
#include "TSUtil.h"
#include "TSSelection.h"
namespace galleryfmwk {
void TSSelection::setShowerEnergyResolution(float res, bool by_percent) {
_shower_energy_resolution = res;
_shower_energy_by_percent = by_percent;
_shower_energy_distribution = std::normal_distribution<float>(0., res);
}
void TSSelection::setTrackEnergyResolution(float res, bool by_percent) {
_track_energy_resolution = res;
_track_energy_by_percent = by_percent;
_track_energy_distribution = std::normal_distribution<float>(0., res);
}
void TSSelection::setShowerAngleResolution(float res, bool by_percent) {
_shower_angle_resolution = res;
_shower_angle_by_percent = by_percent;
_shower_angle_distribution = std::normal_distribution<float>(0., res);
}
void TSSelection::setTrackAngleResolution(float res, bool by_percent) {
_track_angle_resolution = res;
_track_angle_by_percent = by_percent;
_track_angle_distribution = std::normal_distribution<float>(0., res);
}
void TSSelection::setTrackShowerConfusion(float percent) {
_track_shower_confusion_distribution = std::bernoulli_distribution(percent);
}
void TSSelection::setAcceptP(bool b, int n_protons) {
if (n_protons == 1)
_accept_1p = b;
else if (n_protons > 1)
_accept_np = b;
}
float TSSelection::nextTrackEnergyDistortion(float this_energy=0.) {
if (_track_energy_resolution < 1e-4)
return 0.;
if (_track_energy_by_percent) {
return _track_energy_distribution( _gen ) * this_energy;
}
else {
return _track_energy_distribution( _gen );
}
}
float TSSelection::nextShowerEnergyDistortion(float this_energy=0.) {
if (_shower_energy_resolution < 1e-4)
return 0.;
if (_shower_energy_by_percent) {
return _shower_energy_distribution( _gen) * this_energy;
}
else {
return _shower_energy_distribution( _gen );
}
}
float TSSelection::nextTrackAngleDistortion(float this_angle=0.) {
if (_track_angle_resolution < 1e-4)
return 0.;
if (_track_angle_by_percent) {
return _track_angle_distribution(_gen) * this_angle;
}
else {
return _track_angle_distribution(_gen);
}
}
bool TSSelection::nextTrackShowerConfusion() {
if (_track_shower_confusion < 1e-4)
return false;
else
return _track_shower_confusion_distribution(_gen);
}
float TSSelection::nextShowerAngleDistortion(float this_angle=0.) {
if (_shower_angle_resolution < 1e-4)
return 0.;
if (_shower_angle_by_percent) {
return _shower_angle_distribution(_gen) * this_angle;
}
else {
return _shower_angle_distribution(_gen);
}
}
bool TSSelection::pass_selection(std::vector<PIDParticle>& p, int lpdg) {
// Count protons and the chosen lepton type
size_t np = 0;
size_t nl = 0;
size_t n_trk = 0;
for (size_t i=0; i<p.size(); i++) {
if (p[i].pdg == 2212) {
np++;
n_trk ++;
}
// test for charged pions (the only other possible track?)
if (abs(p[i].pdg) == 211) {
n_trk++;
}
if (p[i].pdg == lpdg) {
nl++;
}
}
bool pass_1l1p = nl == 1 && np == 1 && nl + np == p.size();
bool pass_1lnp = nl == 1 && nl + np == p.size();
bool pass_1lntrk = nl == 1 && nl + n_trk == p.size();
return (pass_1l1p && _accept_1p)
|| (pass_1lnp && _accept_np)
|| (pass_1lntrk && _accept_ntrk);
}
bool TSSelection::initialize(std::vector<std::string> input_files) {
// Load track dE/dx distributions from file
_pdf_file = TFile::Open("./dedx_pdfs.root");
assert(_pdf_file->IsOpen());
TIter next(_pdf_file->GetListOfKeys());
TKey* key;
while ((key = (TKey*)next())) {
const char* name = key->GetName();
TObjArray* tokens = TString(name).Tokenize("_");
TString htype = ((TObjString*)tokens->At(0))->GetString();
int pdg = atoi(((TObjString*)tokens->At(1))->GetString());
if (htype.Contains("htrackdedx")) {
TH2F* h = (TH2F*) _pdf_file->Get(name);
if (h->Integral() == 0 || pdg < 0 || pdg > 10000) { continue; }
// Ignore a few low bins
for (int i=0; i<h->GetNbinsX(); i++) {
for (int j=0; j<h->GetNbinsY(); j++) {
if (i < 2 || j < 2) {
h->SetBinContent(i, j, 0);
}
}
}
_trackdedxs[pdg] = h;
}
}
// initialize input files
_input_files = input_files;
// Initialize dataset identifier
_dataset_id = -1;
// Initialize event counters
true_1e1p = 0;
good_1e1p = 0;
miss_1e1p = 0;
true_1m1p = 0;
good_1m1p = 0;
miss_1m1p = 0;
// set producer to null
_ew_producer = "";
// initialize resolutions to 0 (perfect resolution)
_shower_energy_resolution = 0.;
_track_energy_resolution = 0.;
_shower_energy_by_percent = false;
_track_energy_by_percent = false;
_shower_energy_distribution = std::normal_distribution<float>(0.0, 0.0);
_track_energy_distribution = std::normal_distribution<float>(0.0, 0.0);
_shower_angle_resolution = 0.;
_track_angle_resolution = 0.;
_shower_angle_by_percent = false;
_track_angle_by_percent = false;
_shower_angle_distribution = std::normal_distribution<float>(0.0, 0.0);
_track_angle_distribution = std::normal_distribution<float>(0.0, 0.0);
_track_shower_confusion = 0.;
_track_shower_confusion_distribution = std::bernoulli_distribution();
_accept_1p = true;
_accept_ntrk = false;
_accept_np = false;
// setting up random # stuff
std::random_device rd;
_gen = std::mt19937( rd() );;
// Set up the output trees
assert(_fout);
assert(_fout->IsOpen());
_fout->cd();
_data = new OutputData;
_tree = new TTree("data", "");
_tree->Branch("nupdg", &_data->nupdg);
_tree->Branch("enu", &_data->enu);
_tree->Branch("q2", &_data->q2);
_tree->Branch("w", &_data->w);
_tree->Branch("q0", &_data->q0);
_tree->Branch("q3", &_data->q3);
_tree->Branch("int", &_data->int_type);
_tree->Branch("mode", &_data->int_mode);
_tree->Branch("ccnc", &_data->ccnc);
_tree->Branch("eccqe", &_data->eccqe);
_tree->Branch("eps", &_data->eps);
_tree->Branch("ppdgs", &_data->ppdgs);
_tree->Branch("elep", &_data->elep);
_tree->Branch("thetalep", &_data->thetalep);
_tree->Branch("philep", &_data->philep);
_tree->Branch("lpdg", &_data->lpdg);
_tree->Branch("lpid", &_data->lpid);
_tree->Branch("llen", &_data->llen);
_tree->Branch("lexit", &_data->lexit);
_tree->Branch("bnbweight", &_data->bnbweight);
_tree->Branch("dataset", &_data->dataset);
_tree->Branch("weights", &_data->weights);
_truthtree = new TNtuple("truth", "", "nupdg:enu:ccnc:int:mode:w:q2:lpdg:elep:tlep:npip:npim:npi0:np:nn:fw:ttrk:rtrk:texit:tshr:rshr:sexit");
_mectree = new TNtuple("mec", "", "nupdg:enu:ccnc:mode:w:q2:lpdg:tlep:ep0:ep1:ep2:ep3:ep4");
_record_truth = true;
_record_mec = true;
return true;
}
bool TSSelection::run() {
bool ret = true;
for (gallery::Event ev(_input_files) ; !ev.atEnd(); ev.next()) {
ret = ret && analyze(&ev);
}
return ret;
}
bool TSSelection::analyze(gallery::Event* ev) {
// Get handles for event data
art::InputTag gtruth_tag(_mct_producer);
auto const& gtruth_list = \
(*ev->getValidHandle<std::vector<simb::GTruth> >(gtruth_tag));
std::vector<evwgh::MCEventWeight> eventweights_list;
if (_ew_producer.size() > 0) {
art::InputTag eventweight_tag(_ew_producer);
eventweights_list = \
(*ev->getValidHandle<std::vector<evwgh::MCEventWeight> >(eventweight_tag));
}
art::InputTag mctruth_tag(_mct_producer);
auto const& mctruth_list = \
(*ev->getValidHandle<std::vector<simb::MCTruth> >(mctruth_tag));
art::InputTag mcshower_tag(_mcshw_producer);
auto const& mcshower_list = \
(*ev->getValidHandle<std::vector<sim::MCShower> >(mcshower_tag));
art::InputTag mctrack_tag(_mctrk_producer);
auto const& mctrack_list = \
(*ev->getValidHandle<std::vector<sim::MCTrack> >(mctrack_tag));
// Sanity checks
assert(_fout);
assert(mctruth_list.size() == gtruth_list.size());
// BNB flux weight
double wbnb = 1.0;
if (!eventweights_list.empty() &&
eventweights_list[0].fWeight.find("bnbcorrection_FluxHist") != eventweights_list[0].fWeight.end()) {
wbnb = eventweights_list[0].fWeight.at("bnbcorrection_FluxHist")[0];
}
// Loop through MC truth interactions
for(size_t i=0; i<mctruth_list.size(); i++) {
const simb::MCTruth& mctruth = mctruth_list.at(i);
const simb::GTruth& gtruth = gtruth_list.at(i);
size_t ntracks = 0, nshowers = 0;
// Keep track of event particle content (currently a little redundant)
std::vector<PIDParticle> particles_found;
std::vector<PIDParticle> particles_true;
// Get vertex-associated contained tracks
for (size_t j=0; j<mctrack_list.size(); j++) {
const sim::MCTrack& mct = mctrack_list.at(j);
float this_angle = mct.Start().Momentum().Theta();
float this_energy = mct.Start().E() - tsutil::get_pdg_mass(mct.PdgCode());
float energy_distortion = nextTrackEnergyDistortion( this_energy );
float angle_distortion = nextTrackAngleDistortion(this_angle);
// Apply track cuts
if (!goodTrack(mct, mctruth, energy_distortion, angle_distortion)) {
continue;
}
// Track length
double s = 0;
TLorentzVector pos = mct.End().Position();
for (long k=mct.size()-2; k>=0; k--) {
s += (pos.Vect() - mct[k].Position().Vect()).Mag();
pos = mct[k].Position();
}
// don't apply energy distortion to the "true" particle data
particles_true.push_back({
mct.PdgCode(),
mct.PdgCode(),
mct.Start().Momentum(),
mct.Start().E() - tsutil::get_pdg_mass(mct.PdgCode()),
tsutil::eccqe(mct.Start().Momentum()),
s,
!tsutil::inFV(mct)
});
ntracks++;
// KS test on dE/dx vs. range to pick the best match for PID
TH2F* htemp = tsutil::HistDEdx(mct);
double ks_best = 0;
int pdg_best = -999;
for (auto const& it : _trackdedxs) {
if (htemp->Integral() == 0 || it.second->Integral() == 0) { continue; }
double ks = it.second->KolmogorovTest(htemp);
if (ks > ks_best) {
ks_best = ks;
pdg_best = it.first;
}
}
delete htemp;
// Un-PID "protons" that are too long or short
if (pdg_best == 2212 && (s > 80 || s < 12)) {
pdg_best = -888;
}
// Call all remaining unmatched tracks protons
if (pdg_best == -999) {
pdg_best = 2212;
}
auto new_momentum = TLorentzVector(mct.Start().Momentum());
new_momentum.SetTheta(this_angle + angle_distortion);
particles_found.push_back({
pdg_best,
mct.PdgCode(),
new_momentum,
mct.Start().E() - tsutil::get_pdg_mass(mct.PdgCode()) + energy_distortion,
tsutil::eccqe(mct.Start().Momentum(), energy_distortion, angle_distortion),
s,
!tsutil::inFV(mct)
});
}
// Get vertex-associated contained showers
for (size_t j=0; j<mcshower_list.size(); j++) {
const sim::MCShower& mcs = mcshower_list.at(j);
float this_energy = mcs.Start().E() - tsutil::get_pdg_mass(mcs.PdgCode());
float energy_distortion = nextShowerEnergyDistortion( this_energy );
// Apply shower cuts
if (!goodShower(mcs, mctruth, energy_distortion)) {
continue;
}
// don't apply energy distortion to the "true" particle data
particles_true.push_back({
mcs.PdgCode(),
mcs.PdgCode(),
mcs.Start().Momentum(),
mcs.Start().E() - tsutil::get_pdg_mass(mcs.PdgCode()),
tsutil::eccqe(mcs.Start().Momentum()),
-1,
!tsutil::inFV(mcs)
});
nshowers++;
// Guess the PDG based on shower dE/dx (cut at 3.5)
int pdg_best = (mcs.dEdx() < 3.5 ? 11 : 22);
particles_found.push_back({
pdg_best,
mcs.PdgCode(),
mcs.Start().Momentum(),
// @ANDY: IS THIS A BUG???
// previous lines have this energy as:
// mcs.Start().E() - tsutil::get_pdg_mass(mcs.PdgCode())
// note: fixed
mcs.Start().E() + energy_distortion - tsutil::get_pdg_mass(mcs.PdgCode()),
tsutil::eccqe(mcs.Start().Momentum(), energy_distortion),
-1,
!tsutil::inFV(mcs)
});
}
// Classify the event (found/true 1lip/1m1p)
// "True" good_event here means there are one true l and one true p that pass the
// track/shower cuts (i.e. are in within this specific signal definition).
bool f_1e1p = pass_selection(particles_found, 11);
bool t_1e1p = pass_selection(particles_true, 11);
bool f_1m1p = pass_selection(particles_found, 13);
bool t_1m1p = pass_selection(particles_true, 13);
// Where have all the muons gone?
//if (t_1m1p && !f_1m1p) {
// std::cout << "1m1p missed" << std::endl;
// std::cout << "true: n=" << particles_true.size() << ": ";
// for (size_t i=0; i<particles_true.size(); i++) {
// std::cout << particles_true[i].pdg << "/" << particles_true[i].pdgtrue << "(" << particles_true[i].evis << ") ";
// }
// std::cout << std::endl;
// std::cout << "found: n=" << particles_found.size() << ": ";
// for (size_t i=0; i<particles_found.size(); i++) {
// std::cout << particles_found[i].pdg << "/" << particles_found[i].pdgtrue << "(" << particles_found[i].evis << ") ";
// }
// std::cout << std::endl;
//}
if (t_1e1p) true_1e1p++;
if (f_1e1p && t_1e1p) good_1e1p++;
if (f_1e1p && !t_1e1p) miss_1e1p++;
if (t_1m1p) true_1m1p++;
if (f_1m1p && t_1m1p) good_1m1p++;
if (f_1m1p && !t_1m1p) miss_1m1p++;
// Print out PID information mis-IDs
if ((f_1e1p && !t_1e1p) || (f_1m1p && !t_1m1p)) {
std::cout << "true: " << mctruth.GetNeutrino().Nu().E() * 1000
<< "[" << mctruth.GetNeutrino().InteractionType() << "] ";
for (size_t k=0; k<particles_true.size(); k++) {
std::cout << particles_true[k] << " ";
}
std::cout << std::endl;
std::cout << "est: " << ntracks << " tracks, "
<< nshowers << " showers; ";
for (size_t k=0; k<particles_found.size(); k++) {
std::cout << particles_found[k] << " ";
}
std::cout << std::endl;
}
// Write event to output tree for found 1l1p events
if (f_1e1p || f_1m1p) {
double eccqe=-1, elep=-1, thetalep=-1, philep=-1, lpdg=-1, lpid=-1, llen=-1, lexit=-1;
std::vector<double> eps;
std::vector<int> ppdgs;
for (size_t k=0; k<particles_found.size(); k++) {
if (particles_found[k].pdg == 2212) {
eps.push_back( particles_found[k].evis );
ppdgs.push_back( particles_found[k].pdgtrue );
}
else {
eccqe = particles_found[k].eccqe;
elep = particles_found[k].evis;
thetalep = particles_found[k].p.Theta();
philep = particles_found[k].p.Phi();
lpid = particles_found[k].pdg;
lpdg = particles_found[k].pdgtrue;
llen = particles_found[k].len;
lexit = particles_found[k].exiting;
}
}
const simb::MCNeutrino& nu = mctruth.GetNeutrino();
const simb::MCParticle& pnu = nu.Nu();
const simb::MCParticle& plep = nu.Lepton();
TLorentzVector xp = (pnu.Momentum() - plep.Momentum());
std::map<std::string, std::vector<double> > wgh;
if (!eventweights_list.empty()) {
wgh = eventweights_list[0].fWeight;
}
_data->nupdg = nu.Nu().PdgCode();
_data->enu = nu.Nu().E();
_data->q2 = nu.QSqr();
_data->w = nu.W();
_data->q0 = xp.E();
_data->q3 = xp.Vect().Mag();
_data->int_type = nu.InteractionType();
_data->int_mode = nu.Mode();
_data->ccnc = nu.CCNC();
_data->eccqe = eccqe;
_data->eps = eps;
_data->ppdgs = ppdgs;
_data->elep = elep;
_data->thetalep = thetalep;
_data->philep = philep;
_data->lpdg = lpdg;
_data->lpid = lpid;
_data->llen = llen;
_data->lexit = lexit;
_data->bnbweight = wbnb;
_data->dataset = _dataset_id;
_data->weights = &wgh;
_tree->Fill();
}
// Fill the event truth tree
if (_record_truth) {
float vtt[22] = {
(float) mctruth.GetNeutrino().Nu().PdgCode(),
(float) mctruth.GetNeutrino().Nu().E(),
(float) mctruth.GetNeutrino().CCNC(),
(float) mctruth.GetNeutrino().InteractionType(),
(float) mctruth.GetNeutrino().Mode(),
(float) mctruth.GetNeutrino().W(),
(float) mctruth.GetNeutrino().QSqr(),
(float) mctruth.GetNeutrino().Lepton().PdgCode(),
(float) mctruth.GetNeutrino().Lepton().E(),
(float) gtruth.fGint,
(float) gtruth.fNumPiPlus,
(float) gtruth.fNumPiMinus,
(float) gtruth.fNumPi0,
(float) gtruth.fNumProton,
(float) gtruth.fNumNeutron,
(float) wbnb,
(float) mctrack_list.size(),
(float) ntracks,
(float) 0,
(float) mcshower_list.size(),
(float) nshowers,
(float) 0
};
_truthtree->Fill(vtt);
}
if (_record_mec) {
// Fill tree for MEC events with sorted proton energies
if (mctruth.GetNeutrino().Mode() == simb::kMEC) {
std::vector<float> epmec(5, -1);
for (size_t k=0; k<mctrack_list.size(); k++) {
const sim::MCTrack& t = mctrack_list[k];
if (t.PdgCode() == 2212 && t.Process() == "primary") {
double ke = t.Start().E() - tsutil::get_pdg_mass(t.PdgCode());
epmec.push_back(ke);
}
}
std::sort(epmec.begin(), epmec.end(), std::greater<>());
float v2[13] = {
(float) mctruth.GetNeutrino().Nu().PdgCode(),
(float) mctruth.GetNeutrino().Nu().E(),
(float) mctruth.GetNeutrino().CCNC(),
(float) mctruth.GetNeutrino().Mode(),
(float) mctruth.GetNeutrino().W(),
(float) mctruth.GetNeutrino().QSqr(),
(float) mctruth.GetNeutrino().Lepton().PdgCode(),
(float) mctruth.GetNeutrino().Lepton().E(),
(float) epmec[0],
(float) epmec[1],
(float) epmec[2],
(float) epmec[3],
(float) epmec[4]
};
_mectree->Fill(v2);
}
}
}
return true;
}
bool TSSelection::finalize() {
// Print out statistics
std::cout << "1e1p true: " << true_1e1p
<< ", good: " << good_1e1p
<< ", miss: " << miss_1e1p
<< std::endl;
std::cout << "1e1p eff: "
<< 1.0 * (good_1e1p + miss_1e1p) / true_1e1p
<< std::endl;
std::cout << "1e1p pur: "
<< 1.0 * good_1e1p / (good_1e1p + miss_1e1p)
<< std::endl;
std::cout << "1m1p true: " << true_1m1p
<< ", good: " << good_1m1p
<< ", miss: " << miss_1m1p << std::endl;
std::cout << "1m1p eff: "
<< 1.0 * (good_1m1p + miss_1m1p) / true_1m1p
<< std::endl;
std::cout << "1m1p pur: "
<< 1.0 * good_1m1p / (good_1m1p + miss_1m1p)
<< std::endl;
std::cout << "SHOWER, TRACK ENERGY RESOLUTION: " << _shower_energy_resolution << " "<< _track_energy_resolution << std::endl;
// record header data
_fout->cd();
HeaderData header;
HeaderData *to_header = &header;
TTree *header_tree = new TTree("header", "");
header_tree->Branch("track_producer", &to_header->track_producer);
header_tree->Branch("fw_producer", &to_header->fw_producer);
header_tree->Branch("ew_producer", &to_header->ew_producer);
header_tree->Branch("mct_producer", &to_header->mct_producer);
header_tree->Branch("mcf_producer", &to_header->mcf_producer);
header_tree->Branch("mctrk_producer", &to_header->mctrk_producer);
header_tree->Branch("mcshw_producer", &to_header->mcshw_producer);
header_tree->Branch("shower_energy_resolution", &to_header->shower_energy_resolution);
header_tree->Branch("shower_energy_by_percent", &to_header->shower_energy_by_percent);
header_tree->Branch("track_energy_resolution", &to_header->track_energy_resolution);
header_tree->Branch("track_energy_by_percent", &to_header->track_energy_by_percent);
header_tree->Branch("shower_angle_resolution", &to_header->shower_angle_resolution);
header_tree->Branch("shower_angle_by_percent", &to_header->shower_angle_by_percent);
header_tree->Branch("track_angle_resolution", &to_header->track_angle_resolution);
header_tree->Branch("track_angle_by_percent", &to_header->track_angle_by_percent);
header_tree->Branch("accept_1p", &to_header->accept_1p);
header_tree->Branch("accept_np", &to_header->accept_np);
header_tree->Branch("accept_ntrk", &to_header->accept_ntrk);
header_tree->Branch("input_files", &to_header->input_files);
header.track_producer = _track_producer;
header.fw_producer = _fw_producer;
header.ew_producer = _ew_producer;
header.mct_producer = _mct_producer;
header.mcf_producer = _mcf_producer;
header.mctrk_producer = _mctrk_producer;
header.mcshw_producer = _mcshw_producer;
header.shower_energy_resolution = _shower_energy_resolution;
header.shower_energy_by_percent = _shower_energy_by_percent;
header.track_energy_resolution = _track_energy_resolution;
header.track_energy_by_percent = _track_energy_by_percent;
header.shower_angle_resolution = _shower_angle_resolution;
header.shower_angle_by_percent = _shower_angle_by_percent;
header.track_angle_resolution = _track_angle_resolution;
header.track_angle_by_percent = _track_angle_by_percent;
header.accept_1p = _accept_1p;
header.accept_np = _accept_np;
header.accept_ntrk = _accept_ntrk;
header.input_files = _input_files;
header_tree->Fill();
// Write output to ROOT file
if (_fout) {
_fout->cd();
_tree->Write();
_truthtree->Write();
_mectree->Write();
header_tree->Write();
}
return true;
}
std::ostream& operator<<(std::ostream& os, const TSSelection::PIDParticle& dt) {
os << dt.pdg << "(" << dt.evis << ")";
return os;
}
} // namespace gallery_fmwk