-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathwarp_gsimple.C
283 lines (241 loc) · 9.16 KB
/
warp_gsimple.C
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
//========================================================================
//
// This is a script for creating a modified GSimpleNtpFlux file from
// an original file, but modifying entries. Assumes it will be passed
// the name of a single file.
//
// User must supply two functions "warp_entry()" and "warp_meta()"
//
//========================================================================
#include "Numerical/RandomGen.h"
#include "FluxDrivers/GSimpleNtpFlux.h"
#include "Utils/UnitUtils.h"
#include "TSystem.h"
#include "TStopwatch.h"
#include "TLorentzVector.h"
#include "TNtuple.h"
#include "TFile.h"
#include "TChain.h"
#include <iostream>
#include <iomanip>
#include <string>
#include <sstream>
#include <set>
#include <algorithm>
using namespace std;
using namespace genie;
using namespace genie::flux;
// function prototype
double pick_energy(double, double);
bool gDoWarp = true;
//========================================================================
// this is what the USER must supply
void warp_entry(GSimpleNtpEntry* entry,
GSimpleNtpAux* aux,
GSimpleNtpNuMI* numi)
{
if ( ! gDoWarp ) return;
// this is a dumb weighting scheme
//if ( entry->E > 1.5 ) {
//// don't set absolute weight in case we started with weighted entries
//// scale whatever is there from GSimpleNtpFlux driver
//entry->wgt *= 0.5;
//}
//// change some flavors
//if ( gRandom->Rndm() < 0.25 ) {
//entry->pdg = ( entry->pdg > 0 ) ? 16 : -16;
//}
// record the 4 momentum before warp
TLorentzVector fmb4 = TLorentzVector(entry->px, entry->py, entry->pz,
entry->E);
// sample an energy value from the 1/E distribution
// min energy 0.01 GeV, where GENIE spline starts
// max energy 100 GeV, by eye
entry->E = pick_energy(0.01, 20);
// scale momentum so that the new momentum is in the same direction as the
// old one but obeys the energy-momentum relation
double k = sqrt(1+(entry->E*entry->E-fmb4.E()*fmb4.E())/fmb4.Vect().Mag2());
entry->px *= k;
entry->py *= k;
entry->pz *= k;
// use all branches so there are no complaints from compiler
static double trash = 0;
if ( numi ) trash = numi->entryno;
if ( aux ) trash = aux->auxint.size();
if ( trash == -1 ) cout << "trash was -1" << endl;
}
void warp_meta(GSimpleNtpMeta* meta, string fnamein)
{
meta->infiles.push_back("THIS IS MDC WARPED FLUX FROM:");
meta->infiles.push_back(fnamein);
if ( ! gDoWarp ) {
meta->infiles.push_back("NO ACTUAL WARPING APPLIED");
} else {
string msg ="USING WARP FUNCTION CODENAMED WHISKEY-TANGO-FOXTROT";
meta->infiles.push_back(msg);
}
}
double pick_energy(double emin, double emax)
{
// pick a distribution of shape 1/E
// between the values [emin, emax]
// must have low energy cut to avoid infrared divergence
double c = log(emin);
double A = 1. / (log(emax) - c);
double r = gRandom->Rndm(); // flat variate between (0,1]
return exp(r/A + c);
}
void update_flavors(GSimpleNtpMeta* meta, std::set<int>& gFlavors);
//========================================================================
// main routine
// assumes a single input file (meta-data handling)
void warp_gsimple(string fnameout="gsimple_output.root",
string fnamein="gsimple_input.root",
bool dowarp=true)
{
gDoWarp = dowarp;
string fnameinlong(gSystem->ExpandPathName(fnamein.c_str()));
TFile* fin = TFile::Open(fnameinlong.c_str(),"READONLY");
TTree* etreein = (TTree*)fin->Get("flux");
TTree* mtreein = (TTree*)fin->Get("meta");
genie::flux::GSimpleNtpEntry* entry_in = new genie::flux::GSimpleNtpEntry;
genie::flux::GSimpleNtpNuMI* numi_in = new genie::flux::GSimpleNtpNuMI;
genie::flux::GSimpleNtpAux* aux_in = new genie::flux::GSimpleNtpAux;
genie::flux::GSimpleNtpMeta* meta_in = new genie::flux::GSimpleNtpMeta;
long int nentries = etreein->GetEntries();
int sba_status[4] = { -999, -999, -999, -999 };
sba_status[0] = etreein->SetBranchAddress("entry",&entry_in);
sba_status[1] = etreein->SetBranchAddress("numi",&numi_in);
sba_status[2] = etreein->SetBranchAddress("aux",&aux_in);
sba_status[3] = mtreein->SetBranchAddress("meta",&meta_in);
cout << "SetBranchAddress results "
<< sba_status[0] << ","
<< sba_status[1] << ","
<< sba_status[2] << ","
<< sba_status[3]
<< endl;
bool donumi = ( sba_status[1] == 0 );
bool doaux = ( sba_status[2] == 0 );
int nindices = mtreein->BuildIndex("metakey"); // tie metadata to entry
cout << "saw " << nindices << " metakey indices" << endl;
cout << "Creating: " << fnameout << endl;
cout << "Input file: " << fnamein << endl;
cout << "Branches: entry "
<< ( (doaux) ? "aux ":"" )
<< ( (donumi) ? "numi ":"" )
<< endl;
if ( ! gDoWarp ) cout << "++++++++ NO ACTUAL WARP APPLIED" << endl;
TFile* fout = TFile::Open(fnameout.c_str(),"RECREATE");
TTree* fluxntp = new TTree("flux","a simple flux n-tuple");
TTree* metantp = new TTree("meta","metadata for flux n-tuple");
genie::flux::GSimpleNtpEntry* fentry = new genie::flux::GSimpleNtpEntry;
genie::flux::GSimpleNtpAux* faux = new genie::flux::GSimpleNtpAux;
genie::flux::GSimpleNtpNuMI* fnumi = new genie::flux::GSimpleNtpNuMI;
genie::flux::GSimpleNtpMeta* fmeta = new genie::flux::GSimpleNtpMeta;
fluxntp->Branch("entry",&fentry);
if ( doaux ) fluxntp->Branch("aux",&faux);
if ( donumi ) fluxntp->Branch("numi",&fnumi);
metantp->Branch("meta",&fmeta);
cout << "=========================== Start " << endl;
TStopwatch sw;
sw.Start();
const double large = 1.0e300;
double minwgt = large, maxwgt = -large, maxenergy = -large;
std::set<int> gFlavors;
for ( long int ientry = 0; ientry < nentries; ++ientry ) {
// reset what's been read in anticipation of a new entry
entry_in->Reset();
aux_in->Reset();
numi_in->Reset();
// read the next entry, get metadata if it's different
//int nbytes =
etreein->GetEntry(ientry);
UInt_t metakey = entry_in->metakey;
if ( fmeta->metakey != metakey ) {
// UInt_t oldkey = meta_in->metakey;
int nbmeta = mtreein->GetEntryWithIndex(metakey);
cout << "on entry " << ientry << " fetched metadata w/ key "
<< metakey << " read " << nbmeta << " bytes" << endl;
}
// reset what we're writing
fentry->Reset();
fnumi->Reset();
faux->Reset();
// copy read in objects to output objects
*fentry = *entry_in;
*fnumi = *numi_in;
*faux = *aux_in;
// mess with the values to your hearts content
warp_entry(fentry,faux,fnumi);
// keep track of any weights for the meta data
if ( fentry->wgt < minwgt ) minwgt = fentry->wgt;
if ( fentry->wgt > maxwgt ) maxwgt = fentry->wgt;
// user might have changed an energy larger than any in original file
if ( fentry->E > maxenergy ) maxenergy = fentry->E;
// user might have changed a flavor keep track of all we see
gFlavors.insert(fentry->pdg);
// process currently held metadata after transition to metadata key
if ( fmeta->metakey != meta_in->metakey ) {
// new meta data found
cout << "new meta found " << fmeta->metakey << "/" << meta_in->metakey << endl;
if ( fmeta->metakey != 0 ) {
// have metadata that needs adjustment and writing
// before processing new metadata
warp_meta(fmeta,fnamein);
fmeta->minWgt = minwgt;
fmeta->maxWgt = maxwgt;
fmeta->maxEnergy = maxenergy;
update_flavors(fmeta,gFlavors);
gFlavors.clear();
cout << "metantp->Fill() " << *fmeta << endl;
metantp->Fill();
// next meta-data restarts weight range (and energy max)
minwgt = large;
maxwgt = -large;
maxenergy = -large;
}
// get a copy of metadata for accumulating adjustments
*fmeta = *meta_in;
//cout << " copy meta_in " << *meta_in << endl << *fmeta << endl;
}
fluxntp->Fill();
}
// write last set of meta-data after all is done
if ( fmeta->metakey != 0 ) {
cout << "final meta flush " << fmeta->metakey << endl;
// have metadata that needs adjustment and writing
warp_meta(fmeta,fnamein);
fmeta->minWgt = minwgt;
fmeta->maxWgt = maxwgt;
fmeta->maxEnergy = maxenergy;
update_flavors(fmeta,gFlavors);
gFlavors.clear();
cout << "metantp->Fill() " << *fmeta << endl;
metantp->Fill();
}
cout << "=========================== Complete " << endl;
sw.Stop();
cout << "Generated " << nentries
<< endl
<< "Time to generate: " << endl;
sw.Print();
fout->cd();
// write ntuples out
fluxntp->Write();
metantp->Write();
fout->Close();
cout << endl << endl;
}
void update_flavors(GSimpleNtpMeta* meta, std::set<int>& gFlavors)
{
// add any new flavors here, if necessary
std::set<int>::const_iterator flitr = gFlavors.begin();
std::vector<Int_t>& flvlist = meta->pdglist;
for ( ; flitr != gFlavors.end(); ++flitr ) {
int seen_pdg = *flitr;
std::vector<Int_t>::iterator knwitr =
std::find(flvlist.begin(),flvlist.end(),seen_pdg);
bool known_pdg = ( knwitr != flvlist.end() );
if ( ! known_pdg ) meta->pdglist.push_back(seen_pdg);
}
}