INukeHadroData2015.cxx
Go to the documentation of this file.
1 //____________________________________________________________________________
2 /*
3  Copyright (C) 2003-2015, GENIE Neutrino MC Generator Collaboration
4  For the full text of the license visit http://copyright.genie-mc.org
5  or see $GENIE/LICENSE
6 
7  Author: Costas Andreopoulos <costas.andreopoulos \at stfc.ac.uk>, Rutherford Lab.
8  Steve Dytman <dytman+@pitt.edu>, Pittsburgh Univ.
9  Aaron Meyer <asm58@pitt.edu>, Pittsburgh Univ.
10  Alex Bell, Pittsburgh Univ.
11  February 01, 2007
12 
13  For the class documentation see the corresponding header file.
14 
15  Important revisions after version 2.0.0 :
16  @ Dec 06, 2008 - CA
17  Tweak dtor so as not to clutter the output if GENIE exits in err so as to
18  spot the fatal mesg immediately.
19  @ Jul 15, 2010 - AM
20  MeanFreePath now distinguishes between protons and neutrons. To account for
21  this, additional splines were added for each nucleon. Absorption fates
22  condensed into a single fate, splines updated accordingly. Added gamma and kaon
23  splines.Outdated splines were removed. Function IntBounce implemented to calculate
24  a CM scattering angle for given probe, target, product, and fate. AngleAndProduct
25  similar to IntBounce, but also determines the target nucleon.
26  @ May 01, 2012 - CA
27  Pick data from $GENIE/data/evgen/intranuke/
28  @ Jan 9, 2015 - SD, Nick Geary, Tomek Golan
29  Added 2014 version of INTRANUKE codes for independent development.
30  @ Oct, 2015 - TG
31  Added 2015 version of INTRANUKE codes for independent development.
32  Include Oset data files.
33  @ Apr, 2016 - Flor Blasczyk
34  Added K+ cex data files
35 */
36 //____________________________________________________________________________
37 
38 #include <cassert>
39 #include <string>
40 
41 #include <TSystem.h>
42 #include <TNtupleD.h>
43 #include <TGraph2D.h>
44 #include <TTree.h>
45 #include <TMath.h>
46 #include <TFile.h>
47 #include <iostream>
48 
49 #include "Conventions/Constants.h"
50 #include "GHEP/GHepParticle.h"
52 #include "Messenger/Messenger.h"
53 #include "Numerical/RandomGen.h"
54 #include "Numerical/Spline.h"
55 #include "PDG/PDGCodes.h"
56 
57 using std::ostringstream;
58 using std::ios;
59 
60 using namespace genie;
61 using namespace genie::constants;
62 
63 //____________________________________________________________________________
65 //____________________________________________________________________________
66 double INukeHadroData2015::fMinKinEnergy = 1.0; // MeV
67 double INukeHadroData2015::fMaxKinEnergyHA = 999.0; // MeV
68 double INukeHadroData2015::fMaxKinEnergyHN = 1799.0; // MeV
69 //____________________________________________________________________________
71 {
72  this->LoadCrossSections();
73  fInstance = 0;
74 }
75 //____________________________________________________________________________
77 {
78  // pi+n/p hA x-section splines
79  delete fXSecPipn_Tot;
80  delete fXSecPipn_CEx;
81  delete fXSecPipn_Elas;
82  delete fXSecPipn_Reac;
83  delete fXSecPipp_Tot;
84  delete fXSecPipp_CEx;
85  delete fXSecPipp_Elas;
86  delete fXSecPipp_Reac;
87  delete fXSecPipd_Abs;
88 
89  delete fXSecPp_Cmp; //added to fix memory leak; no noticeable effect, but good convention.
90  delete fXSecPn_Cmp;
91  delete fXSecNn_Cmp;
92  delete fXSecPp_Tot;
93  delete fXSecPp_Elas;
94  delete fXSecPp_Reac;
95  delete fXSecPn_Tot;
96  delete fXSecPn_Elas;
97  delete fXSecPn_Reac;
98  delete fXSecNn_Tot;
99  delete fXSecNn_Elas;
100  delete fXSecNn_Reac;
101 
102  // pi0n/p hA x-section splines
103  delete fXSecPi0n_Tot;
104  delete fXSecPi0n_CEx;
105  delete fXSecPi0n_Elas;
106  delete fXSecPi0n_Reac;
107  delete fXSecPi0p_Tot;
108  delete fXSecPi0p_CEx;
109  delete fXSecPi0p_Elas;
110  delete fXSecPi0p_Reac;
111  delete fXSecPi0d_Abs;
112 
113  // K+N x-section splines
114  delete fXSecKpn_Elas;
115  delete fXSecKpp_Elas;
116  delete fXSecKpN_Abs;
117  delete fXSecKpN_Tot;
118 
119  // gamma x-section splines (inelastic only)
120  delete fXSecGamp_fs;
121  delete fXSecGamn_fs;
122  delete fXSecGamN_Tot;
123 
124  // N+A x-section splines
125  delete fFracPA_Tot;
126  delete fFracPA_Elas;
127  delete fFracPA_Inel;
128  delete fFracPA_CEx;
129  delete fFracPA_Abs;
130  delete fFracPA_PiPro;
131  delete fFracNA_Tot;
132  delete fFracNA_Elas;
133  delete fFracNA_Inel;
134  delete fFracNA_CEx;
135  delete fFracNA_Abs;
136  delete fFracNA_PiPro;
137 
138  delete fFracPA_Cmp; // cmp - add support later
139  delete fFracNA_Cmp;
140 
141  // hN data
142  delete fhN2dXSecPP_Elas;
143  delete fhN2dXSecNP_Elas;
144  delete fhN2dXSecPipN_Elas;
145  delete fhN2dXSecPi0N_Elas;
146  delete fhN2dXSecPimN_Elas;
147  delete fhN2dXSecKpN_Elas;
148  delete fhN2dXSecKpP_Elas;
149  delete fhN2dXSecPiN_CEx;
150  delete fhN2dXSecPiN_Abs;
151  delete fhN2dXSecGamPi0P_Inelas;
152  delete fhN2dXSecGamPi0N_Inelas;
153  delete fhN2dXSecGamPipN_Inelas;
154  delete fhN2dXSecGamPimP_Inelas;
155 }
156 //____________________________________________________________________________
158 {
159  if(fInstance == 0) {
160  LOG("INukeData", pINFO) << "INukeHadroData2015 late initialization";
161  static INukeHadroData2015::Cleaner cleaner;
163  fInstance = new INukeHadroData2015;
164  }
165  return fInstance;
166 }
167 //____________________________________________________________________________
169 {
170 // Loads hadronic x-section data
171 
172  //-- Get the top-level directory with input hadron cross-section data
173  // (search for $GINUKEHADRONDATA or use default location)
174  string data_dir = (gSystem->Getenv("GINUKEHADRONDATA")) ?
175  string(gSystem->Getenv("GINUKEHADRONDATA")) :
176  string(gSystem->Getenv("GENIE")) + string("/data/evgen/intranuke");
177 
178  LOG("INukeData", pINFO)
179  << "Loading INTRANUKE hadron data from: " << data_dir;
180 
181  //-- Build filenames
182 
183  string datafile_NN = data_dir + "/tot_xsec/intranuke-xsections-NN2014.dat";
184  string datafile_pipN = data_dir + "/tot_xsec/intranuke-xsections-pi+N.dat";
185  string datafile_pi0N = data_dir + "/tot_xsec/intranuke-xsections-pi0N.dat";
186  string datafile_NA = data_dir + "/tot_xsec/intranuke-fractions-NA.dat";
187  string datafile_KA = data_dir + "/tot_xsec/intranuke-fractions-KA.dat";
188  string datafile_gamN = data_dir + "/tot_xsec/intranuke-xsections-gamN.dat";
189  string datafile_kN = data_dir + "/tot_xsec/intranuke-xsections-kaonN2015.dat";
190 
191  //-- Make sure that all data files are available
192 
193  assert( ! gSystem->AccessPathName(datafile_NN. c_str()) );
194  assert( ! gSystem->AccessPathName(datafile_pipN.c_str()) );
195  assert( ! gSystem->AccessPathName(datafile_pi0N.c_str()) );
196  assert( ! gSystem->AccessPathName(datafile_KA. c_str()) );
197  assert( ! gSystem->AccessPathName(datafile_gamN.c_str()) );
198  assert( ! gSystem->AccessPathName(datafile_kN. c_str()) );
199 
200  LOG("INukeData", pINFO) << "Found all necessary data files...";
201 
202  //-- Load data files
203 
204  TTree data_NN;
205  TTree data_pipN;
206  TTree data_pi0N;
207  TTree data_NA;
208  TTree data_KA;
209  TTree data_gamN;
210  TTree data_kN;
211 
212  data_NN.ReadFile(datafile_NN.c_str(),"ke/D:pp_tot/D:pp_elas/D:pp_reac/D:pn_tot/D:pn_elas/D:pn_reac/D:nn_tot/D:nn_elas/D:nn_reac/D:pp_cmp/D:pn_cmp/D:nn_cmp/D");
213  data_pipN.ReadFile(datafile_pipN.c_str(),
214  "ke/D:pipn_tot/D:pipn_cex/D:pipn_elas/D:pipn_reac/D:pipp_tot/D:pipp_cex/D:pipp_elas/D:pipp_reac/D:pipd_abs");
215  data_pi0N.ReadFile(datafile_pi0N.c_str(),
216  "ke/D:pi0n_tot/D:pi0n_cex/D:pi0n_elas/D:pi0n_reac/D:pi0p_tot/D:pi0p_cex/D:pi0p_elas/D:pi0p_reac/D:pi0d_abs");
217  //data_NA.ReadFile(datafile_NA.c_str(),
218  //"ke/D:pA_tot/D:pA_elas/D:pA_inel/D:pA_cex/D:pA_abs/D:pA_pipro/D"); // add support for cmp here (?)
219  data_NA.ReadFile(datafile_NA.c_str(),
220  "ke/D:pA_tot/D:pA_elas/D:pA_inel/D:pA_cex/D:pA_abs/D:pA_pipro/D:pA_cmp/D"); // add support for cmp here (?)
221  data_gamN.ReadFile(datafile_gamN.c_str(),
222  "ke/D:pi0p_tot/D:pipn_tot/D:pimp_tot/D:pi0n_tot/D:gamp_fs/D:gamn_fs/D:gamN_tot/D");
223  data_kN.ReadFile(datafile_kN.c_str(),
224  "ke/D:kpp_elas/D:kpn_elas/D:kpn_cex/D:kp_abs/D:kpN_tot/D");
225  data_KA.ReadFile(datafile_KA.c_str(),
226  "ke/D:KA_tot/D:KA_elas/D:KA_inel/D:KA_abs/D");
227 
228  LOG("INukeData", pDEBUG) << "Number of data rows in NN : " << data_NN.GetEntries();
229  LOG("INukeData", pDEBUG) << "Number of data rows in pipN : " << data_pipN.GetEntries();
230  LOG("INukeData", pDEBUG) << "Number of data rows in pi0N : " << data_pi0N.GetEntries();
231  LOG("INukeData", pDEBUG) << "Number of data rows in NA : " << data_NA.GetEntries();
232  LOG("INukeData", pDEBUG) << "Number of data rows in KA : " << data_KA.GetEntries();
233  LOG("INukeData", pDEBUG) << "Number of data rows in gamN : " << data_gamN.GetEntries();
234  LOG("INukeData", pDEBUG) << "Number of data rows in kN : " << data_kN.GetEntries();
235 
236  LOG("INukeData", pINFO) << "Done loading all x-section files...";
237 
238  //-- Build x-section splines
239 
240  // p/n+p/n hA x-section splines
241  fXSecPp_Tot = new Spline(&data_NN, "ke:pp_tot");
242  fXSecPp_Elas = new Spline(&data_NN, "ke:pp_elas");
243  fXSecPp_Reac = new Spline(&data_NN, "ke:pp_reac");
244  fXSecPn_Tot = new Spline(&data_NN, "ke:pn_tot");
245  fXSecPn_Elas = new Spline(&data_NN, "ke:pn_elas");
246  fXSecPn_Reac = new Spline(&data_NN, "ke:pn_reac");
247  fXSecNn_Tot = new Spline(&data_NN, "ke:nn_tot");
248  fXSecNn_Elas = new Spline(&data_NN, "ke:nn_elas");
249  fXSecNn_Reac = new Spline(&data_NN, "ke:nn_reac");
250  fXSecPp_Cmp = new Spline(&data_NN, "ke:pp_cmp"); //for compound nucleus fate
251  fXSecPn_Cmp = new Spline(&data_NN, "ke:pn_cmp");
252  fXSecNn_Cmp = new Spline(&data_NN, "ke:nn_cmp");
253 
254  // pi+n/p hA x-section splines
255  fXSecPipn_Tot = new Spline(&data_pipN, "ke:pipn_tot");
256  fXSecPipn_CEx = new Spline(&data_pipN, "ke:pipn_cex");
257  fXSecPipn_Elas = new Spline(&data_pipN, "ke:pipn_elas");
258  fXSecPipn_Reac = new Spline(&data_pipN, "ke:pipn_reac");
259  fXSecPipp_Tot = new Spline(&data_pipN, "ke:pipp_tot");
260  fXSecPipp_CEx = new Spline(&data_pipN, "ke:pipp_cex");
261  fXSecPipp_Elas = new Spline(&data_pipN, "ke:pipp_elas");
262  fXSecPipp_Reac = new Spline(&data_pipN, "ke:pipp_reac");
263  fXSecPipd_Abs = new Spline(&data_pipN, "ke:pipd_abs");
264 
265  // pi0n/p hA x-section splines
266  fXSecPi0n_Tot = new Spline(&data_pi0N, "ke:pi0n_tot");
267  fXSecPi0n_CEx = new Spline(&data_pi0N, "ke:pi0n_cex");
268  fXSecPi0n_Elas = new Spline(&data_pi0N, "ke:pi0n_elas");
269  fXSecPi0n_Reac = new Spline(&data_pi0N, "ke:pi0n_reac");
270  fXSecPi0p_Tot = new Spline(&data_pi0N, "ke:pi0p_tot");
271  fXSecPi0p_CEx = new Spline(&data_pi0N, "ke:pi0p_cex");
272  fXSecPi0p_Elas = new Spline(&data_pi0N, "ke:pi0p_elas");
273  fXSecPi0p_Reac = new Spline(&data_pi0N, "ke:pi0p_reac");
274  fXSecPi0d_Abs = new Spline(&data_pi0N, "ke:pi0d_abs");
275 
276  // K+N x-section splines
277  fXSecKpn_Elas = new Spline(&data_kN, "ke:kpn_elas");
278  fXSecKpp_Elas = new Spline(&data_kN, "ke:kpp_elas");
279  fXSecKpn_CEx = new Spline(&data_kN, "ke:kpn_cex");
280  fXSecKpN_Abs = 0; // new Spline(&data_kN, "ke:kp_abs"); why not used?
281  fXSecKpN_Tot = new Spline(&data_kN, "ke:kpN_tot");
282 
283  // gamma x-section splines
284  fXSecGamp_fs = new Spline(&data_gamN, "ke:gamp_fs");
285  fXSecGamn_fs = new Spline(&data_gamN, "ke:gamn_fs");
286  fXSecGamN_Tot = new Spline(&data_gamN, "ke:gamN_tot");
287 
288  // N+A x-section fraction splines
289  fFracPA_Tot = new Spline(&data_NA, "ke:pA_tot");
290  fFracPA_Elas = new Spline(&data_NA, "ke:pA_elas");
291  fFracPA_Inel = new Spline(&data_NA, "ke:pA_inel");
292  fFracPA_CEx = new Spline(&data_NA, "ke:pA_cex");
293  fFracPA_Abs = new Spline(&data_NA, "ke:pA_abs");
294  fFracPA_PiPro = new Spline(&data_NA, "ke:pA_pipro");
295  fFracNA_Tot = new Spline(&data_NA, "ke:pA_tot"); // assuming nA same as pA
296  fFracNA_Elas = new Spline(&data_NA, "ke:pA_elas");
297  fFracNA_Inel = new Spline(&data_NA, "ke:pA_inel");
298  fFracNA_CEx = new Spline(&data_NA, "ke:pA_cex");
299  fFracNA_Abs = new Spline(&data_NA, "ke:pA_abs");
300  fFracNA_PiPro = new Spline(&data_NA, "ke:pA_pipro");
301 
302  fFracPA_Cmp = new Spline(&data_NA, "ke:pA_cmp"); //cmp - add support later
303  fFracNA_Cmp = new Spline(&data_NA, "ke:pA_cmp");
304 
305  // K+A x-section fraction splines
306  fFracKA_Tot = new Spline(&data_KA, "ke:KA_tot");
307  fFracKA_Elas = new Spline(&data_KA, "ke:KA_elas");
308  fFracKA_CEx = 0; // new Spline(&data_KA, "ke:KA_cex"); //Added, but needs to be computed
309  fFracKA_Inel = new Spline(&data_KA, "ke:KA_inel");
310  fFracKA_Abs = new Spline(&data_KA, "ke:KA_abs");
311  //
312  // hN stuff
313  //
314 
315 
316  // kIHNFtElas
317  // pp, nn --> read from pp/pp%.txt
318  // pn, np --> read from pp/pn%.txt
319  // pi+ N --> read from pip/pip%.txt
320  // pi0 N --> read from pip/pip%.txt
321  // pi- N --> read from pim/pim%.txt
322  // K+ N --> read from kpn/kpn%.txt
323  // K+ P --> read from kpp/kpp%.txt
324  // kIHNFtCEx
325  // pi+, pi0, pi- --> read from pie/pie%.txt (using pip+n->pi0+p data)
326  // kIHNFtAbs
327  // pi+, pi0, pi- --> read from pid2p/pid2p%.txt (using pip+D->2p data)
328  // kIHNFtInelas
329  // gamma p -> p pi0 --> read from gampi0p/%-pi0p.txt
330  // gamma p -> n pi+ --> read from gampi+n/%-pi+n.txt
331  // gamma n -> n pi0 --> read from gampi0n/%-pi0n.txt
332  // gamma n -> p pi- --> read from gampi-p/%-pi-p.txt
333 
334 
335  // kIHNFtElas, pp&nn :
336  {
337  const int hN_ppelas_nfiles = 20;
338  const int hN_ppelas_points_per_file = 21;
339  const int hN_ppelas_npoints = hN_ppelas_points_per_file * hN_ppelas_nfiles;
340 
341  double hN_ppelas_energies[hN_ppelas_nfiles] = {
342  50, 100, 150, 200, 250, 300, 350, 400, 450, 500,
343  550, 600, 650, 700, 750, 800, 850, 900, 950, 1000
344  };
345 
346  double hN_ppelas_costh [hN_ppelas_points_per_file];
347  double hN_ppelas_xsec [hN_ppelas_npoints];
348 
349  int ipoint=0;
350 
351  for(int ifile = 0; ifile < hN_ppelas_nfiles; ifile++) {
352  // build filename
353  ostringstream hN_datafile;
354  double ke = hN_ppelas_energies[ifile];
355  hN_datafile << data_dir << "/diff_ang/pp/pp" << ke << ".txt";
356  // read data
357  ReadhNFile(
358  hN_datafile.str(), ke, hN_ppelas_points_per_file,
359  ipoint, hN_ppelas_costh, hN_ppelas_xsec,2);
360  }//loop over files
361 
362  fhN2dXSecPP_Elas = new BLI2DNonUnifGrid(hN_ppelas_nfiles,hN_ppelas_points_per_file,
363  hN_ppelas_energies,hN_ppelas_costh,hN_ppelas_xsec);
364  }
365 
366  // kIHNFtElas, pn&np :
367  {
368  const int hN_npelas_nfiles = 20;
369  const int hN_npelas_points_per_file = 21;
370  const int hN_npelas_npoints = hN_npelas_points_per_file * hN_npelas_nfiles;
371 
372  double hN_npelas_energies[hN_npelas_nfiles] = {
373  50, 100, 150, 200, 250, 300, 350, 400, 450, 500,
374  550, 600, 650, 700, 750, 800, 850, 900, 950, 1000
375  };
376 
377  double hN_npelas_costh [hN_npelas_points_per_file];
378  double hN_npelas_xsec [hN_npelas_npoints];
379 
380  int ipoint=0;
381 
382  for(int ifile = 0; ifile < hN_npelas_nfiles; ifile++) {
383  // build filename
384  ostringstream hN_datafile;
385  double ke = hN_npelas_energies[ifile];
386  hN_datafile << data_dir << "/diff_ang/pn/pn" << ke << ".txt";
387  // read data
388  ReadhNFile(
389  hN_datafile.str(), ke, hN_npelas_points_per_file,
390  ipoint, hN_npelas_costh, hN_npelas_xsec,2);
391  }//loop over files
392 
393  fhN2dXSecNP_Elas = new BLI2DNonUnifGrid(hN_npelas_nfiles,hN_npelas_points_per_file,
394  hN_npelas_energies,hN_npelas_costh,hN_npelas_xsec);
395  }
396 
397  // kIHNFtElas, pipN :
398  {
399  const int hN_pipNelas_nfiles = 60;
400  const int hN_pipNelas_points_per_file = 21;
401  const int hN_pipNelas_npoints = hN_pipNelas_points_per_file * hN_pipNelas_nfiles;
402 
403  double hN_pipNelas_energies[hN_pipNelas_nfiles] = {
404  10, 20, 30, 40, 50, 60, 70, 80, 90,
405  100, 110, 120, 130, 140, 150, 160, 170, 180, 190,
406  200, 210, 220, 230, 240, 250, 260, 270, 280, 290,
407  300, 340, 380, 420, 460, 500, 540, 580, 620, 660,
408  700, 740, 780, 820, 860, 900, 940, 980,
409  1020, 1060, 1100, 1140, 1180, 1220, 1260,
410  1300, 1340, 1380, 1420, 1460, 1500
411  };
412 
413  double hN_pipNelas_costh [hN_pipNelas_points_per_file];
414  double hN_pipNelas_xsec [hN_pipNelas_npoints];
415 
416  int ipoint=0;
417 
418  for(int ifile = 0; ifile < hN_pipNelas_nfiles; ifile++) {
419  // build filename
420  ostringstream hN_datafile;
421  double ke = hN_pipNelas_energies[ifile];
422  hN_datafile << data_dir << "/diff_ang/pip/pip" << ke << ".txt";
423  // read data
424  ReadhNFile(
425  hN_datafile.str(), ke, hN_pipNelas_points_per_file,
426  ipoint, hN_pipNelas_costh, hN_pipNelas_xsec,2);
427  }//loop over files
428 
429  fhN2dXSecPipN_Elas = new BLI2DNonUnifGrid(hN_pipNelas_nfiles,hN_pipNelas_points_per_file,
430  hN_pipNelas_energies,hN_pipNelas_costh,hN_pipNelas_xsec);
431  }
432 
433  // kIHNFtElas, pi0N :
434  {
435  const int hN_pi0Nelas_nfiles = 60;
436  const int hN_pi0Nelas_points_per_file = 21;
437  const int hN_pi0Nelas_npoints = hN_pi0Nelas_points_per_file * hN_pi0Nelas_nfiles;
438 
439  double hN_pi0Nelas_energies[hN_pi0Nelas_nfiles] = {
440  10, 20, 30, 40, 50, 60, 70, 80, 90,
441  100, 110, 120, 130, 140, 150, 160, 170, 180, 190,
442  200, 210, 220, 230, 240, 250, 260, 270, 280, 290,
443  300, 340, 380, 420, 460, 500, 540, 580, 620, 660,
444  700, 740, 780, 820, 860, 900, 940, 980,
445  1020, 1060, 1100, 1140, 1180, 1220, 1260,
446  1300, 1340, 1380, 1420, 1460, 1500
447  };
448 
449  double hN_pi0Nelas_costh [hN_pi0Nelas_points_per_file];
450  double hN_pi0Nelas_xsec [hN_pi0Nelas_npoints];
451 
452  int ipoint=0;
453 
454  for(int ifile = 0; ifile < hN_pi0Nelas_nfiles; ifile++) {
455  // build filename
456  ostringstream hN_datafile;
457  double ke = hN_pi0Nelas_energies[ifile];
458  hN_datafile << data_dir << "/diff_ang/pip/pip" << ke << ".txt";
459  // read data
460  ReadhNFile(
461  hN_datafile.str(), ke, hN_pi0Nelas_points_per_file,
462  ipoint, hN_pi0Nelas_costh, hN_pi0Nelas_xsec,2);
463  }//loop over files
464 
465  fhN2dXSecPi0N_Elas = new BLI2DNonUnifGrid(hN_pi0Nelas_nfiles,hN_pi0Nelas_points_per_file,
466  hN_pi0Nelas_energies,hN_pi0Nelas_costh,hN_pi0Nelas_xsec);
467  }
468 
469  // kIHNFtElas, pimN :
470  {
471  const int hN_pimNelas_nfiles = 60;
472  const int hN_pimNelas_points_per_file = 21;
473  const int hN_pimNelas_npoints = hN_pimNelas_points_per_file * hN_pimNelas_nfiles;
474 
475  double hN_pimNelas_energies[hN_pimNelas_nfiles] = {
476  10, 20, 30, 40, 50, 60, 70, 80, 90,
477  100, 110, 120, 130, 140, 150, 160, 170, 180, 190,
478  200, 210, 220, 230, 240, 250, 260, 270, 280, 290,
479  300, 340, 380, 420, 460, 500, 540, 580, 620, 660,
480  700, 740, 780, 820, 860, 900, 940, 980,
481  1020, 1060, 1100, 1140, 1180, 1220, 1260,
482  1300, 1340, 1380, 1420, 1460, 1500
483  };
484 
485  double hN_pimNelas_costh [hN_pimNelas_points_per_file];
486  double hN_pimNelas_xsec [hN_pimNelas_npoints];
487 
488  int ipoint=0;
489 
490  for(int ifile = 0; ifile < hN_pimNelas_nfiles; ifile++) {
491  // build filename
492  ostringstream hN_datafile;
493  double ke = hN_pimNelas_energies[ifile];
494  hN_datafile << data_dir << "/diff_ang/pim/pim" << ke << ".txt";
495  // read data
496  ReadhNFile(
497  hN_datafile.str(), ke, hN_pimNelas_points_per_file,
498  ipoint, hN_pimNelas_costh, hN_pimNelas_xsec,2);
499  }//loop over files
500 
501  fhN2dXSecPimN_Elas = new BLI2DNonUnifGrid(hN_pimNelas_nfiles,hN_pimNelas_points_per_file,
502  hN_pimNelas_energies,hN_pimNelas_costh,hN_pimNelas_xsec);
503  }
504 
505  // kIHNFtElas, kpn :
506  {
507  const int hN_kpNelas_nfiles = 18;
508  const int hN_kpNelas_points_per_file = 37;
509  const int hN_kpNelas_npoints = hN_kpNelas_points_per_file * hN_kpNelas_nfiles;
510 
511  double hN_kpNelas_energies[hN_kpNelas_nfiles] = {
512  100, 200, 300, 400, 500, 600, 700, 800, 900, 1000,
513  1100, 1200, 1300, 1400, 1500, 1600, 1700, 1800
514  };
515 
516  double hN_kpNelas_costh [hN_kpNelas_points_per_file];
517  double hN_kpNelas_xsec [hN_kpNelas_npoints];
518 
519  int ipoint=0;
520 
521  for(int ifile = 0; ifile < hN_kpNelas_nfiles; ifile++) {
522  // build filename
523  ostringstream hN_datafile;
524  double ke = hN_kpNelas_energies[ifile];
525  hN_datafile << data_dir << "/diff_ang/kpn/kpn" << ke << ".txt";
526  // read data
527  ReadhNFile(
528  hN_datafile.str(), ke, hN_kpNelas_points_per_file,
529  ipoint, hN_kpNelas_costh, hN_kpNelas_xsec,2);
530  }//loop over files
531 
532  fhN2dXSecKpN_Elas = new BLI2DNonUnifGrid(hN_kpNelas_nfiles,hN_kpNelas_points_per_file,
533  hN_kpNelas_energies,hN_kpNelas_costh,hN_kpNelas_xsec);
534  }
535  // kIHNFtCEx, kpn :
536  {
537  const int hN_kpNcex_nfiles = 18;
538  const int hN_kpNcex_points_per_file = 37;
539  const int hN_kpNcex_npoints = hN_kpNcex_points_per_file * hN_kpNcex_nfiles;
540 
541  double hN_kpNcex_energies[hN_kpNcex_nfiles] = {
542  100, 200, 300, 400, 500, 600, 700, 800, 900, 1000,
543  1100, 1200, 1300, 1400, 1500, 1600, 1700, 1800
544  };
545 
546  double hN_kpNcex_costh [hN_kpNcex_points_per_file];
547  double hN_kpNcex_xsec [hN_kpNcex_npoints];
548 
549  int ipoint=0;
550 
551  for(int ifile = 0; ifile < hN_kpNcex_nfiles; ifile++) {
552  // build filename
553  ostringstream hN_datafile;
554  double ke = hN_kpNcex_energies[ifile];
555  hN_datafile << data_dir << "/diff_ang/kpncex/kpcex" << ke << ".txt";
556  // read data
557  ReadhNFile(
558  hN_datafile.str(), ke, hN_kpNcex_points_per_file,
559  ipoint, hN_kpNcex_costh, hN_kpNcex_xsec,2);
560  }//loop over files
561 
562  /*double hN_kpNcex_costh_cond [hN_kpNcex_points_per_file];
563  for (int ient = 0; ient < hN_kpNcex_points_per_file; ient++) {
564  hN_kpNcex_costh_cond[ient] = hN_kpNcex_costh[ient];
565  }*/
566 
567  fhN2dXSecKpN_CEx = new BLI2DNonUnifGrid(hN_kpNcex_nfiles,hN_kpNcex_points_per_file,
568  hN_kpNcex_energies,hN_kpNcex_costh,hN_kpNcex_xsec);
569  }
570 //----------------------------------------------------------------------------------------
571 
572 
573  // kIHNFtElas, kpp :
574  {
575  const int hN_kpPelas_nfiles = 18;
576  const int hN_kpPelas_points_per_file = 37;
577  const int hN_kpPelas_npoints = hN_kpPelas_points_per_file * hN_kpPelas_nfiles;
578 
579  double hN_kpPelas_energies[hN_kpPelas_nfiles] = {
580  100, 200, 300, 400, 500, 600, 700, 800, 900, 1000,
581  1100, 1200, 1300, 1400, 1500, 1600, 1700, 1800
582  };
583 
584  double hN_kpPelas_costh [hN_kpPelas_points_per_file];
585  double hN_kpPelas_xsec [hN_kpPelas_npoints];
586 
587  int ipoint=0;
588 
589  for(int ifile = 0; ifile < hN_kpPelas_nfiles; ifile++) {
590  // build filename
591  ostringstream hN_datafile;
592  double ke = hN_kpPelas_energies[ifile];
593  hN_datafile << data_dir << "/diff_ang/kpp/kpp" << ke << ".txt";
594  // read data
595  ReadhNFile(
596  hN_datafile.str(), ke, hN_kpPelas_points_per_file,
597  ipoint, hN_kpPelas_costh, hN_kpPelas_xsec,2);
598  }//loop over files
599 
600  fhN2dXSecKpP_Elas = new BLI2DNonUnifGrid(hN_kpPelas_nfiles,hN_kpPelas_points_per_file,
601  hN_kpPelas_energies,hN_kpPelas_costh,hN_kpPelas_xsec);
602  }
603 
604  // kIHNFtCEx, (pi+, pi0, pi-) N
605  {
606  const int hN_piNcex_nfiles = 60;
607  const int hN_piNcex_points_per_file = 21;
608  const int hN_piNcex_npoints = hN_piNcex_points_per_file * hN_piNcex_nfiles;
609 
610  double hN_piNcex_energies[hN_piNcex_nfiles] = {
611  10, 20, 30, 40, 50, 60, 70, 80, 90,
612  100, 110, 120, 130, 140, 150, 160, 170, 180, 190,
613  200, 210, 220, 230, 240, 250, 260, 270, 280, 290,
614  300, 340, 380, 420, 460, 500, 540, 580, 620, 660,
615  700, 740, 780, 820, 860, 900, 940, 980,
616  1020, 1060, 1100, 1140, 1180, 1220, 1260,
617  1300, 1340, 1380, 1420, 1460, 1500
618  };
619 
620  double hN_piNcex_costh [hN_piNcex_points_per_file];
621  double hN_piNcex_xsec [hN_piNcex_npoints];
622 
623  int ipoint=0;
624 
625  for(int ifile = 0; ifile < hN_piNcex_nfiles; ifile++) {
626  // build filename
627  ostringstream hN_datafile;
628  double ke = hN_piNcex_energies[ifile];
629  hN_datafile << data_dir << "/diff_ang/pie/pie" << ke << ".txt";
630  // read data
631  ReadhNFile(
632  hN_datafile.str(), ke, hN_piNcex_points_per_file,
633  ipoint, hN_piNcex_costh, hN_piNcex_xsec,2);
634  }//loop over files
635 
636  fhN2dXSecPiN_CEx = new BLI2DNonUnifGrid(hN_piNcex_nfiles,hN_piNcex_points_per_file,
637  hN_piNcex_energies,hN_piNcex_costh,hN_piNcex_xsec);
638  }
639 
640  // kIHNFtAbs, (pi+, pi0, pi-) N
641  {
642  const int hN_piNabs_nfiles = 19;
643  const int hN_piNabs_points_per_file = 21;
644  const int hN_piNabs_npoints = hN_piNabs_points_per_file * hN_piNabs_nfiles;
645 
646  double hN_piNabs_energies[hN_piNabs_nfiles] = {
647  50, 75, 100, 125, 150, 175, 200, 225, 250, 275,
648  300, 325, 350, 375, 400, 425, 450, 475, 500
649  };
650 
651  double hN_piNabs_costh [hN_piNabs_points_per_file];
652  double hN_piNabs_xsec [hN_piNabs_npoints];
653 
654  int ipoint=0;
655 
656  for(int ifile = 0; ifile < hN_piNabs_nfiles; ifile++) {
657  // build filename
658  ostringstream hN_datafile;
659  double ke = hN_piNabs_energies[ifile];
660  hN_datafile << data_dir << "/diff_ang/pid2p/pid2p" << ke << ".txt";
661  // read data
662  ReadhNFile(
663  hN_datafile.str(), ke, hN_piNabs_points_per_file,
664  ipoint, hN_piNabs_costh, hN_piNabs_xsec,2);
665  }//loop over files
666 
667  fhN2dXSecPiN_Abs = new BLI2DNonUnifGrid(hN_piNabs_nfiles,hN_piNabs_points_per_file,
668  hN_piNabs_energies,hN_piNabs_costh,hN_piNabs_xsec);
669  }
670 
671  // kIHNFtInelas, gamma p -> p pi0
672  {
673  const int hN_gampi0pInelas_nfiles = 29;
674  const int hN_gampi0pInelas_points_per_file = 37;
675  const int hN_gampi0pInelas_npoints = hN_gampi0pInelas_points_per_file * hN_gampi0pInelas_nfiles;
676 
677  double hN_gampi0pInelas_energies[hN_gampi0pInelas_nfiles] = {
678  160, 180, 200, 220, 240, 260, 280, 300, 320, 340,
679  360, 380, 400, 450, 500, 550, 600, 650, 700, 750,
680  800, 850, 900, 950, 1000, 1050, 1100, 1150, 1200
681  };
682 
683  double hN_gampi0pInelas_costh [hN_gampi0pInelas_points_per_file];
684  double hN_gampi0pInelas_xsec [hN_gampi0pInelas_npoints];
685 
686  int ipoint=0;
687 
688  for(int ifile = 0; ifile < hN_gampi0pInelas_nfiles; ifile++) {
689  // build filename
690  ostringstream hN_datafile;
691  double ke = hN_gampi0pInelas_energies[ifile];
692  hN_datafile << data_dir << "/diff_ang/gampi0p/" << ke << "-pi0p.txt";
693  // read data
694  ReadhNFile(
695  hN_datafile.str(), ke, hN_gampi0pInelas_points_per_file,
696  ipoint, hN_gampi0pInelas_costh, hN_gampi0pInelas_xsec,3);
697  }//loop over files
698 
699  fhN2dXSecGamPi0P_Inelas = new BLI2DNonUnifGrid(hN_gampi0pInelas_nfiles,hN_gampi0pInelas_points_per_file,
700  hN_gampi0pInelas_energies,hN_gampi0pInelas_costh,hN_gampi0pInelas_xsec);
701  }
702 
703  // kIHNFtInelas, gamma n -> n pi0
704  {
705  const int hN_gampi0nInelas_nfiles = 29;
706  const int hN_gampi0nInelas_points_per_file = 37;
707  const int hN_gampi0nInelas_npoints = hN_gampi0nInelas_points_per_file * hN_gampi0nInelas_nfiles;
708 
709  double hN_gampi0nInelas_energies[hN_gampi0nInelas_nfiles] = {
710  160, 180, 200, 220, 240, 260, 280, 300, 320, 340,
711  360, 380, 400, 450, 500, 550, 600, 650, 700, 750,
712  800, 850, 900, 950, 1000, 1050, 1100, 1150, 1200
713  };
714 
715  double hN_gampi0nInelas_costh [hN_gampi0nInelas_points_per_file];
716  double hN_gampi0nInelas_xsec [hN_gampi0nInelas_npoints];
717  int ipoint=0;
718 
719  for(int ifile = 0; ifile < hN_gampi0nInelas_nfiles; ifile++) {
720  // build filename
721  ostringstream hN_datafile;
722  double ke = hN_gampi0nInelas_energies[ifile];
723  hN_datafile << data_dir << "/diff_ang/gampi0n/" << ke << "-pi0n.txt";
724  // read data
725  ReadhNFile(
726  hN_datafile.str(), ke, hN_gampi0nInelas_points_per_file,
727  ipoint, hN_gampi0nInelas_costh, hN_gampi0nInelas_xsec,3);
728  }//loop over files
729 
730  fhN2dXSecGamPi0N_Inelas = new BLI2DNonUnifGrid(hN_gampi0nInelas_nfiles,hN_gampi0nInelas_points_per_file,
731  hN_gampi0nInelas_energies,hN_gampi0nInelas_costh,hN_gampi0nInelas_xsec);
732  }
733 
734  // kIHNFtInelas, gamma p -> n pi+
735  {
736  const int hN_gampipnInelas_nfiles = 29;
737  const int hN_gampipnInelas_points_per_file = 37;
738  const int hN_gampipnInelas_npoints = hN_gampipnInelas_points_per_file * hN_gampipnInelas_nfiles;
739 
740  double hN_gampipnInelas_energies[hN_gampipnInelas_nfiles] = {
741  160, 180, 200, 220, 240, 260, 280, 300, 320, 340,
742  360, 380, 400, 450, 500, 550, 600, 650, 700, 750,
743  800, 850, 900, 950, 1000, 1050, 1100, 1150, 1200
744  };
745 
746  double hN_gampipnInelas_costh [hN_gampipnInelas_points_per_file];
747  double hN_gampipnInelas_xsec [hN_gampipnInelas_npoints];
748 
749  int ipoint=0;
750 
751  for(int ifile = 0; ifile < hN_gampipnInelas_nfiles; ifile++) {
752  // build filename
753  ostringstream hN_datafile;
754  double ke = hN_gampipnInelas_energies[ifile];
755  hN_datafile << data_dir << "/diff_ang/gampi+n/" << ke << "-pi+n.txt";
756  // read data
757  ReadhNFile(
758  hN_datafile.str(), ke, hN_gampipnInelas_points_per_file,
759  ipoint, hN_gampipnInelas_costh, hN_gampipnInelas_xsec,3);
760  }//loop over files
761 
762  fhN2dXSecGamPipN_Inelas = new BLI2DNonUnifGrid(hN_gampipnInelas_nfiles,hN_gampipnInelas_points_per_file,
763  hN_gampipnInelas_energies,hN_gampipnInelas_costh,hN_gampipnInelas_xsec);
764  }
765 
766  // kIHNFtInelas, gamma n -> p pi-
767  {
768  const int hN_gampimpInelas_nfiles = 29;
769  const int hN_gampimpInelas_points_per_file = 37;
770  const int hN_gampimpInelas_npoints = hN_gampimpInelas_points_per_file * hN_gampimpInelas_nfiles;
771 
772  double hN_gampimpInelas_energies[hN_gampimpInelas_nfiles] = {
773  160, 180, 200, 220, 240, 260, 280, 300, 320, 340,
774  360, 380, 400, 450, 500, 550, 600, 650, 700, 750,
775  800, 850, 900, 950, 1000, 1050, 1100, 1150, 1200
776  };
777 
778  double hN_gampimpInelas_costh [hN_gampimpInelas_points_per_file];
779  double hN_gampimpInelas_xsec [hN_gampimpInelas_npoints];
780 
781  int ipoint=0;
782 
783  for(int ifile = 0; ifile < hN_gampimpInelas_nfiles; ifile++) {
784  // build filename
785  ostringstream hN_datafile;
786  double ke = hN_gampimpInelas_energies[ifile];
787  hN_datafile << data_dir << "/diff_ang/gampi-p/" << ke << "-pi-p.txt";
788  // read data
789  ReadhNFile(
790  hN_datafile.str(), ke, hN_gampimpInelas_points_per_file,
791  ipoint, hN_gampimpInelas_costh, hN_gampimpInelas_xsec,3);
792  }//loop over files
793 
794  fhN2dXSecGamPimP_Inelas = new BLI2DNonUnifGrid(hN_gampimpInelas_nfiles,hN_gampimpInelas_points_per_file,
795  hN_gampimpInelas_energies,hN_gampimpInelas_costh,hN_gampimpInelas_xsec);
796  }
797 
798 
799 
800 
801  TFile TGraphs_file;
802  bool saveTGraphsToFile = true;
803 
804  if (saveTGraphsToFile) {
805  string filename = "TGraphs.root";
806  LOG("INukeHadroData2015", pNOTICE) << "Saving INTRANUKE hadron x-section data to ROOT file: " << filename;
807  TGraphs_file.Open(filename.c_str(), "RECREATE");
808  }
809 
810  // kIHNFtTot, pip + A PipA_Tot
811  {
812  const int pipATot_nfiles = 22;
813  const int pipATot_nuclei[pipATot_nfiles] = {1, 2, 3, 4, 6, 7, 9, 12, 16, 27, 28, 32, 40, 48, 56, 58, 63, 93, 120, 165, 181, 209};
814  const int pipATot_npoints = 203;
815 
816  TPipA_Tot = new TGraph2D(pipATot_npoints);
817 
818  int ipoint=0;
819  double x, y;
820 
821  for(int ifile=0; ifile < pipATot_nfiles; ifile++) {
822  ostringstream ADep_datafile;
823  int nucleus = pipATot_nuclei[ifile];
824  ADep_datafile << data_dir << "/tot_xsec/pipA_tot/pip" << nucleus << "_tot.txt";
825  TGraph * buff = new TGraph(ADep_datafile.str().c_str());
826  for(int i=0; i < buff->GetN(); i++) {
827  buff -> GetPoint(i,x,y);
828  TPipA_Tot -> SetPoint(ipoint,(double)nucleus,x,y);
829  ipoint++;
830  }
831  delete buff;
832  }
833 
834  if (saveTGraphsToFile) {
835  TPipA_Tot -> Write("TPipA_Tot");
836  }
837 
838  }
839 
840 
841  // kIHNFtAbs, pip + A PipA_Abs_frac
842  {
843  const int pipAAbs_f_nfiles = 18;
844  const int pipAAbs_f_nuclei[pipAAbs_f_nfiles] = {1, 2, 3, 4, 7, 9, 12, 16, 27, 48, 56, 58, 63, 93, 120, 165, 181, 209};
845  const int pipAAbs_f_npoints = 111;
846 
847  TfracPipA_Abs = new TGraph2D(pipAAbs_f_npoints);
848 
849  int ipoint=0;
850  double x, y;
851  for(int ifile=0; ifile < pipAAbs_f_nfiles; ifile++) {
852  ostringstream ADep_datafile;
853  int nucleus = pipAAbs_f_nuclei[ifile];
854  ADep_datafile << data_dir << "/tot_xsec/pipA_abs_frac/pip" << nucleus << "_abs_frac.txt";
855  TGraph * buff = new TGraph(ADep_datafile.str().c_str());
856  for(int i=0; i < buff->GetN(); i++) {
857  buff -> GetPoint(i,x,y);
858  TfracPipA_Abs -> SetPoint(ipoint,(double)nucleus,x,y);
859  ipoint++;
860  }
861  delete buff;
862  }
863  if (saveTGraphsToFile) {
864  TfracPipA_Abs -> Write("TfracPipA_Abs");
865  }
866 
867  }
868 
869 
870  // kIHNFtCEx, pip + A PipA_CEx_frac
871  {
872  const int pipACEx_f_nfiles = 18;
873  const int pipACEx_f_nuclei[pipACEx_f_nfiles] = {1, 2, 3, 4, 7, 9, 12, 16, 27, 48, 56, 58, 63, 93, 120, 165, 181, 209};
874  const int pipACEx_f_npoints = 129;
875 
876  TfracPipA_CEx = new TGraph2D(pipACEx_f_npoints);
877 
878  int ipoint=0;
879  double x, y;
880 
881  for(int ifile=0; ifile < pipACEx_f_nfiles; ifile++) {
882  ostringstream ADep_datafile;
883  int nucleus = pipACEx_f_nuclei[ifile];
884  ADep_datafile << data_dir << "/tot_xsec/pipA_cex_frac/pip" << nucleus << "_cex_frac.txt";
885  TGraph * buff = new TGraph(ADep_datafile.str().c_str());
886  for(int i=0; i < buff->GetN(); i++) {
887  buff -> GetPoint(i,x,y);
888  TfracPipA_CEx -> SetPoint(ipoint,(double)nucleus,x,y);
889  ipoint++;
890  }
891  delete buff;
892  }
893 
894  if (saveTGraphsToFile) {
895  TfracPipA_CEx -> Write("TfracPipA_CEx");
896  }
897 
898  }
899 
900 
901 
902  // kIHNFtCEx, pip + A PipA_CEx (just for developmental purposes)
903  {
904  TGraph2D * TPipA_CEx;
905 
906  const int pipACEx_nfiles = 18;
907  const int pipACEx_nuclei[pipACEx_nfiles] = {1, 2, 3, 4, 7, 9, 12, 16, 27, 48, 56, 58, 63, 93, 120, 165, 181, 209};
908  const int pipACEx_npoints = 129;
909 
910  TPipA_CEx = new TGraph2D(pipACEx_npoints);
911 
912  int ipoint=0;
913  double x, y;
914 
915  for(int ifile=0; ifile < pipACEx_nfiles; ifile++) {
916  ostringstream ADep_datafile;
917  int nucleus = pipACEx_nuclei[ifile];
918  ADep_datafile << data_dir << "/tot_xsec/pipA_cex/pip" << nucleus << "_cex.txt";
919  TGraph * buff = new TGraph(ADep_datafile.str().c_str());
920  for(int i=0; i < buff->GetN(); i++) {
921  buff -> GetPoint(i,x,y);
922  TPipA_CEx -> SetPoint(ipoint,(double)nucleus,x,y);
923  ipoint++;
924  }
925  delete buff;
926  }
927 
928  if (saveTGraphsToFile) {
929  TPipA_CEx -> Write("TPipA_CEx");
930  }
931 
932  }
933 
934  // kIHNFtAbs, pip + A PipA_Abs (just for developmental purposes)
935  {
936  TGraph2D * TPipA_Abs;
937 
938  const int pipAAbs_nfiles = 18;
939  const int pipAAbs_nuclei[pipAAbs_nfiles] = {1, 2, 3, 4, 7, 9, 12, 16, 27, 48, 56, 58, 63, 93, 120, 165, 181, 209};
940  const int pipAAbs_npoints = 111;
941 
942  TPipA_Abs = new TGraph2D(pipAAbs_npoints);
943 
944  int ipoint=0;
945  double x, y;
946 
947  for(int ifile=0; ifile < pipAAbs_nfiles; ifile++) {
948  ostringstream ADep_datafile;
949  int nucleus = pipAAbs_nuclei[ifile];
950  ADep_datafile << data_dir << "/tot_xsec/pipA_abs/pip" << nucleus << "_abs.txt";
951  TGraph * buff = new TGraph(ADep_datafile.str().c_str());
952  for(int i=0; i < buff->GetN(); i++) {
953  buff -> GetPoint(i,x,y);
954  TPipA_Abs -> SetPoint(ipoint,(double)nucleus,x,y);
955  ipoint++;
956  }
957  delete buff;
958  }
959 
960  if (saveTGraphsToFile) {
961  TPipA_Abs -> Write("TPipA_Abs");
962  }
963 
964  }
965 
966  // kIHNFtElas, pip + A PipA_Elas (just for developmental purposes)
967  {
968  TGraph2D * TPipA_Elas;
969 
970  const int pipAElas_nfiles = 18;
971  const int pipAElas_nuclei[pipAElas_nfiles] = {1, 2, 3, 4, 7, 9, 12, 16, 27, 48, 56, 58, 63, 93, 120, 165, 181, 209};
972  const int pipAElas_npoints = 125;
973 
974  TPipA_Elas = new TGraph2D(pipAElas_npoints);
975 
976  int ipoint=0;
977  double x, y;
978 
979  for(int ifile=0; ifile < pipAElas_nfiles; ifile++) {
980  ostringstream ADep_datafile;
981  int nucleus = pipAElas_nuclei[ifile];
982  ADep_datafile << data_dir << "/tot_xsec/pipA_elas/pip" << nucleus << "_elas.txt";
983  TGraph * buff = new TGraph(ADep_datafile.str().c_str());
984  for(int i=0; i < buff->GetN(); i++) {
985  buff -> GetPoint(i,x,y);
986  TPipA_Elas -> SetPoint(ipoint,(double)nucleus,x,y);
987  ipoint++;
988  }
989  delete buff;
990  }
991 
992  if (saveTGraphsToFile) {
993  TPipA_Elas -> Write("TPipA_Elas");
994  }
995 
996  }
997 
998  // kIHNFtInelas, pip + A PipA_Inelas (just for developmental purposes)
999  {
1000  TGraph2D * TPipA_Inelas;
1001 
1002  const int pipAInelas_nfiles = 20;
1003  const int pipAInelas_nuclei[pipAInelas_nfiles] = {1, 2, 3, 4, 7, 9, 12, 16, 27, 40, 48, 56, 58, 63, 93, 120, 165, 181, 208, 209};
1004  const int pipAInelas_npoints = 118;
1005 
1006  TPipA_Inelas = new TGraph2D(pipAInelas_npoints);
1007 
1008  int ipoint=0;
1009  double x, y;
1010 
1011  for(int ifile=0; ifile < pipAInelas_nfiles; ifile++) {
1012  ostringstream ADep_datafile;
1013  int nucleus = pipAInelas_nuclei[ifile];
1014  ADep_datafile << data_dir << "/tot_xsec/pipA_inelas/pip" << nucleus << "_inelas.txt";
1015  TGraph * buff = new TGraph(ADep_datafile.str().c_str());
1016  for(int i=0; i < buff->GetN(); i++) {
1017  buff -> GetPoint(i,x,y);
1018  TPipA_Inelas -> SetPoint(ipoint,(double)nucleus,x,y);
1019  ipoint++;
1020  }
1021  delete buff;
1022  }
1023 
1024  if (saveTGraphsToFile) {
1025  TPipA_Inelas -> Write("TPipA_Inelas");
1026  }
1027 
1028  }
1029 
1030 
1031 
1032  // kIHNFtElas, pip + A PipA_Elas_frac
1033  {
1034  const int pipAElas_f_nfiles = 18;
1035  const int pipAElas_f_nuclei[pipAElas_f_nfiles] = {1, 2, 3, 4, 7, 9, 12, 16, 27, 48, 56, 58, 63, 93, 120, 165, 181, 209};
1036  const int pipAElas_f_npoints = 125;
1037 
1038  TfracPipA_Elas = new TGraph2D(pipAElas_f_npoints);
1039 
1040  int ipoint=0;
1041  double x, y;
1042 
1043  for(int ifile=0; ifile < pipAElas_f_nfiles; ifile++) {
1044  ostringstream ADep_datafile;
1045  int nucleus = pipAElas_f_nuclei[ifile];
1046  ADep_datafile << data_dir << "/tot_xsec/pipA_elas_frac/pip" << nucleus << "_elas_frac.txt";
1047  TGraph * buff = new TGraph(ADep_datafile.str().c_str());
1048  for(int i=0; i < buff->GetN(); i++) {
1049  buff -> GetPoint(i,x,y);
1050  TfracPipA_Elas -> SetPoint(ipoint,(double)nucleus,x,y);
1051  ipoint++;
1052  }
1053  delete buff;
1054  }
1055 
1056  if (saveTGraphsToFile) {
1057  TfracPipA_Elas -> Write("TfracPipA_Elas");
1058  }
1059 
1060  }
1061 
1062 
1063  // kIHNFtInelas, pip + A PipA_Inelas_frac
1064  {
1065  const int pipAInelas_f_nfiles = 20;
1066  const int pipAInelas_f_nuclei[pipAInelas_f_nfiles] = {1, 2, 3, 4, 7, 9, 12, 16, 27, 40, 48, 56, 58, 63, 93, 120, 165, 181, 208, 209};
1067  const int pipAInelas_f_npoints = 118;
1068 
1069  TfracPipA_Inelas = new TGraph2D(pipAInelas_f_npoints);
1070 
1071  int ipoint=0;
1072  double x, y;
1073 
1074  for(int ifile=0; ifile < pipAInelas_f_nfiles; ifile++) {
1075  ostringstream ADep_datafile;
1076  int nucleus = pipAInelas_f_nuclei[ifile];
1077  ADep_datafile << data_dir << "/tot_xsec/pipA_inelas_frac/pip" << nucleus << "_inelas_frac.txt";
1078  TGraph * buff = new TGraph(ADep_datafile.str().c_str());
1079  for(int i=0; i < buff->GetN(); i++) {
1080  buff -> GetPoint(i,x,y);
1081  TfracPipA_Inelas -> SetPoint(ipoint,(double)nucleus,x,y);
1082  ipoint++;
1083  }
1084  delete buff;
1085  }
1086 
1087  if (saveTGraphsToFile) {
1088  TfracPipA_Inelas -> Write("TfracPipA_Inelas");
1089  }
1090 
1091  }
1092 
1093 
1094  // kIHNFtPiPro, pip + A PipA_PiPro_frac
1095  {
1096  const int pipAPiPro_f_nfiles = 17;
1097  const int pipAPiPro_f_nuclei[pipAPiPro_f_nfiles] = {1, 2, 3, 4, 7, 9, 12, 16, 48, 56, 58, 63, 93, 120, 165, 181, 209};
1098  const int pipAPiPro_f_npoints = 76;
1099 
1100  TfracPipA_PiPro = new TGraph2D(pipAPiPro_f_npoints);
1101 
1102  int ipoint=0;
1103  double x, y;
1104 
1105  for(int ifile=0; ifile < pipAPiPro_f_nfiles; ifile++) {
1106  ostringstream ADep_datafile;
1107  int nucleus = pipAPiPro_f_nuclei[ifile];
1108  ADep_datafile << data_dir << "/tot_xsec/pipA_pipro_frac/pip" << nucleus << "_pipro_frac.txt";
1109  TGraph * buff = new TGraph(ADep_datafile.str().c_str());
1110  for(int i=0; i < buff->GetN(); i++) {
1111  buff -> GetPoint(i,x,y);
1112  TfracPipA_PiPro -> SetPoint(ipoint,(double)nucleus,x,y);
1113  ipoint++;
1114  }
1115  delete buff;
1116  }
1117 
1118  if (saveTGraphsToFile) {
1119  TfracPipA_PiPro -> Write("TfracPipA_PiPro");
1120  }
1121  }
1122 
1123  TGraphs_file.Close();
1124 
1125  LOG("INukeData", pINFO) << "Done building x-section splines...";
1126 
1127 }
1128 //____________________________________________________________________________
1130  string filename, double ke, int npoints, int & curr_point,
1131  double * costh_array, double * xsec_array, int cols)
1132 {
1133  // open
1134  std::ifstream hN_stream(filename.c_str(), ios::in);
1135  if(!hN_stream.good()) {
1136  LOG("INukeData", pERROR)
1137  << "Error reading INTRANUKE/hN data from: " << filename;
1138  return;
1139  }
1140 
1141  if(cols<2) {
1142  LOG("INukeData", pERROR)
1143  << "Error reading INTRANUKE/hN data from: " << filename;
1144  LOG("INukeData", pERROR)
1145  << "Too few columns: " << cols;
1146  return;
1147  }
1148 
1149  LOG("INukeData", pINFO)
1150  << "Reading INTRANUKE/hN data from: " << filename;
1151 
1152  // skip initial comments
1153  char cbuf[501];
1154  hN_stream.getline(cbuf,400);
1155  hN_stream.getline(cbuf,400);
1156  hN_stream.getline(cbuf,400);
1157 
1158  // read
1159  double angle = 0;
1160  double xsec = 0;
1161  double trash = 0;
1162 
1163  for(int ip = 0; ip < npoints; ip++) {
1164  hN_stream >> angle >> xsec;
1165 
1166  for(int ic = 0; ic < (cols-2); ic++) {
1167  hN_stream >> trash;
1168  }
1169 
1170  LOG("INukeData", pDEBUG)
1171  << "Adding data point: (KE = " << ke << " MeV, angle = "
1172  << angle << ", sigma = " << xsec << " mbarn)";
1173  costh_array[ip] = TMath::Cos(angle*kPi/180.);
1174  xsec_array [curr_point] = xsec;
1175  curr_point++;
1176  }
1177 }
1178 //____________________________________________________________________________
1180  int hpdgc, int tgtpdgc, int nppdgc, INukeFateHN_t fate, double ke, double costh) const
1181 {
1182 // inputs
1183 // fate : h+N fate code
1184 // hpdgc : h PDG code
1185 // tgtpdgc : N PDG code
1186 // nppdgc : product N PDG code
1187 // ke : kinetic energy (MeV)
1188 // costh : cos(scattering angle)
1189 // returns
1190 // xsec : mbarn
1191 
1192  double ke_eval = ke;
1193  double costh_eval = costh;
1194 
1195  costh_eval = TMath::Min(costh, 1.);
1196  costh_eval = TMath::Max(costh_eval, -1.);
1197 
1198  if(fate==kIHNFtElas) {
1199 
1200  if( (hpdgc==kPdgProton && tgtpdgc==kPdgProton) ||
1201  (hpdgc==kPdgNeutron && tgtpdgc==kPdgNeutron) )
1202  {
1203  ke_eval = TMath::Min(ke_eval, 999.);
1204  ke_eval = TMath::Max(ke_eval, 50.);
1205  return fhN2dXSecPP_Elas->Evaluate(ke_eval, costh_eval);
1206  }
1207  else
1208  if( (hpdgc==kPdgProton && tgtpdgc==kPdgNeutron) ||
1209  (hpdgc==kPdgNeutron && tgtpdgc==kPdgProton) )
1210  {
1211  ke_eval = TMath::Min(ke_eval, 999.);
1212  ke_eval = TMath::Max(ke_eval, 50.);
1213  return fhN2dXSecNP_Elas->Evaluate(ke_eval, costh_eval);
1214  }
1215  else
1216  if(hpdgc==kPdgPiP)
1217  {
1218  ke_eval = TMath::Min(ke_eval, 1499.);
1219  ke_eval = TMath::Max(ke_eval, 10.);
1220  return fhN2dXSecPipN_Elas->Evaluate(ke_eval, costh_eval);
1221  }
1222  else
1223  if(hpdgc==kPdgPi0)
1224  {
1225  ke_eval = TMath::Min(ke_eval, 1499.);
1226  ke_eval = TMath::Max(ke_eval, 10.);
1227  return fhN2dXSecPi0N_Elas->Evaluate(ke_eval, costh_eval);
1228  }
1229  else
1230  if(hpdgc==kPdgPiM)
1231  {
1232  ke_eval = TMath::Min(ke_eval, 1499.);
1233  ke_eval = TMath::Max(ke_eval, 10.);
1234  return fhN2dXSecPimN_Elas->Evaluate(ke_eval, costh_eval);
1235  }
1236  else
1237  if(hpdgc==kPdgKP && tgtpdgc==kPdgNeutron)
1238  {
1239  ke_eval = TMath::Min(ke_eval, 1799.);
1240  ke_eval = TMath::Max(ke_eval, 100.);
1241  return fhN2dXSecKpN_Elas->Evaluate(ke_eval, costh_eval);
1242  }
1243  else
1244  if(hpdgc==kPdgKP && tgtpdgc==kPdgProton)
1245  {
1246  ke_eval = TMath::Min(ke_eval, 1799.);
1247  ke_eval = TMath::Max(ke_eval, 100.);
1248  return fhN2dXSecKpP_Elas->Evaluate(ke_eval, costh_eval);
1249  }
1250  }
1251 
1252  else if(fate == kIHNFtCEx) {
1253  if( (hpdgc==kPdgPiP || hpdgc==kPdgPi0 || hpdgc==kPdgPiM) &&
1254  (tgtpdgc==kPdgProton || tgtpdgc==kPdgNeutron) )
1255  {
1256  ke_eval = TMath::Min(ke_eval, 1499.);
1257  ke_eval = TMath::Max(ke_eval, 10.);
1258  return fhN2dXSecPiN_CEx->Evaluate(ke_eval, costh_eval);
1259  }
1260  else if( (hpdgc == kPdgProton && tgtpdgc == kPdgProton) ||
1261  (hpdgc == kPdgNeutron && tgtpdgc == kPdgNeutron) )
1262  {
1263  LOG("INukeData", pWARN) << "Inelastic pp does not exist!";
1264  ke_eval = TMath::Min(ke_eval, 999.);
1265  ke_eval = TMath::Max(ke_eval, 50.);
1266  return fhN2dXSecPP_Elas->Evaluate(ke_eval, costh_eval);
1267  }
1268  else if( (hpdgc == kPdgProton && tgtpdgc == kPdgNeutron) ||
1269  (hpdgc == kPdgNeutron && tgtpdgc == kPdgProton) )
1270  {
1271  ke_eval = TMath::Min(ke_eval, 999.);
1272  ke_eval = TMath::Max(ke_eval, 50.);
1273  return fhN2dXSecNP_Elas->Evaluate(ke_eval, costh_eval);
1274  }
1275  else if(hpdgc == kPdgKP && tgtpdgc == kPdgNeutron) {
1276  ke_eval = TMath::Min(ke_eval, 1799.);
1277  ke_eval = TMath::Max(ke_eval, 100.);
1278  return fhN2dXSecKpN_CEx->Evaluate(ke_eval, costh_eval);
1279  }
1280  }
1281 
1282  else if(fate == kIHNFtAbs) {
1283  if( (hpdgc==kPdgPiP || hpdgc==kPdgPi0 || hpdgc==kPdgPiM) &&
1284  (tgtpdgc==kPdgProton || tgtpdgc==kPdgNeutron) )
1285  {
1286  ke_eval = TMath::Min(ke_eval, 499.);
1287  ke_eval = TMath::Max(ke_eval, 50.);
1288  return fhN2dXSecPiN_Abs->Evaluate(ke_eval, costh_eval);
1289  }
1290  if(hpdgc==kPdgKP) return 1.; //isotropic since no data ???
1291  }
1292 
1293  else if(fate == kIHNFtInelas) {
1294  if( hpdgc==kPdgGamma && tgtpdgc==kPdgProton &&nppdgc==kPdgProton )
1295  {
1296  ke_eval = TMath::Min(ke_eval, 1199.);
1297  ke_eval = TMath::Max(ke_eval, 160.);
1298  return fhN2dXSecGamPi0P_Inelas->Evaluate(ke_eval, costh_eval);
1299  }
1300  else
1301  if( hpdgc==kPdgGamma && tgtpdgc==kPdgProton && nppdgc==kPdgNeutron )
1302  {
1303  ke_eval = TMath::Min(ke_eval, 1199.);
1304  ke_eval = TMath::Max(ke_eval, 160.);
1305  return fhN2dXSecGamPipN_Inelas->Evaluate(ke_eval, costh_eval);
1306  }
1307  else
1308  if( hpdgc==kPdgGamma && tgtpdgc==kPdgNeutron && nppdgc==kPdgProton )
1309  {
1310  ke_eval = TMath::Min(ke_eval, 1199.);
1311  ke_eval = TMath::Max(ke_eval, 160.);
1312  return fhN2dXSecGamPimP_Inelas->Evaluate(ke_eval, costh_eval);
1313  }
1314  else
1315  if( hpdgc==kPdgGamma && tgtpdgc==kPdgNeutron && nppdgc==kPdgNeutron )
1316  {
1317  ke_eval = TMath::Min(ke_eval, 1199.);
1318  ke_eval = TMath::Max(ke_eval, 160.);
1319  return fhN2dXSecGamPi0N_Inelas->Evaluate(ke_eval, costh_eval);
1320  }
1321  }
1322 
1323  return 0;
1324 }
1325 //____________________________________________________________________________
1326 double INukeHadroData2015::FracADep(int hpdgc, INukeFateHA_t fate, double ke, int targA) const
1327 {
1328 // return the x-section fraction for the input fate for the particle with the input pdg
1329 // code and the target with the input mass number at the input kinetic energy
1330 
1331  ke = TMath::Max(fMinKinEnergy, ke); // ke >= 1 MeV
1332  ke = TMath::Min(fMaxKinEnergyHA, ke); // ke <= 999 MeV
1333 
1334  targA = TMath::Min(208, targA); // A <= 208
1335 
1336  LOG("INukeData", pDEBUG) << "Querying hA cross section at ke = " << ke << " and target " << targA;
1337 
1338  if (hpdgc == kPdgPiP) {
1339  // handle pi+
1340  if (fate == kIHAFtCEx ) return TMath::Max(0., TfracPipA_CEx -> Interpolate (targA,ke));
1341  else if (fate == kIHAFtElas ) return TMath::Max(0., TfracPipA_Elas -> Interpolate (targA,ke));
1342  else if (fate == kIHAFtInelas ) return TMath::Max(0., TfracPipA_Inelas -> Interpolate (targA,ke));
1343  else if (fate == kIHAFtAbs ) return TMath::Max(0., TfracPipA_Abs -> Interpolate (targA,ke));
1344  else if (fate == kIHAFtPiProd ) return TMath::Max(0., TfracPipA_PiPro -> Interpolate (targA,ke));
1345  else {
1346  LOG("INukeData", pWARN)
1347  << "Pi+'s don't have this fate: " << INukeHadroFates::AsString(fate);
1348  return 0;
1349  }
1350 
1351  } else if (hpdgc == kPdgPiM) {
1352  // handle pi-
1353  if (fate == kIHAFtCEx ) return TMath::Max(0., TfracPipA_CEx -> Interpolate (targA,ke));
1354  else if (fate == kIHAFtElas ) return TMath::Max(0., TfracPipA_Elas -> Interpolate (targA,ke));
1355  else if (fate == kIHAFtInelas ) return TMath::Max(0., TfracPipA_Inelas -> Interpolate (targA,ke));
1356  else if (fate == kIHAFtAbs ) return TMath::Max(0., TfracPipA_Abs -> Interpolate (targA,ke));
1357  else if (fate == kIHAFtPiProd ) return TMath::Max(0., TfracPipA_PiPro -> Interpolate (targA,ke));
1358  else {
1359  LOG("INukeData", pWARN)
1360  << "Pi-'s don't have this fate: " << INukeHadroFates::AsString(fate);
1361  return 0;
1362  }
1363 
1364  } else if (hpdgc == kPdgPi0) {
1365  // handle pi0
1366  if (fate == kIHAFtCEx ) return TMath::Max(0., TfracPipA_CEx -> Interpolate (targA,ke));
1367  else if (fate == kIHAFtElas ) return TMath::Max(0., TfracPipA_Elas -> Interpolate (targA,ke));
1368  else if (fate == kIHAFtInelas ) return TMath::Max(0., TfracPipA_Inelas -> Interpolate (targA,ke));
1369  else if (fate == kIHAFtAbs ) return TMath::Max(0., TfracPipA_Abs -> Interpolate (targA,ke));
1370  else if (fate == kIHAFtPiProd ) return TMath::Max(0., TfracPipA_PiPro -> Interpolate (targA,ke));
1371  else {
1372  LOG("INukeData", pWARN)
1373  << "Pi0's don't have this fate: " << INukeHadroFates::AsString(fate);
1374  return 0;
1375  }
1376  }
1377 
1378  LOG("INukeData", pWARN)
1379  << "Can't handle particles with pdg code = " << hpdgc;
1380  return 0;
1381 }
1382 //____________________________________________________________________________
1383 double INukeHadroData2015::FracAIndep(int hpdgc, INukeFateHA_t fate, double ke) const
1384 {
1385 // return the x-section fraction for the input fate for the particle with the input pdg
1386 // code at the input kinetic energy
1387 
1388  ke = TMath::Max(fMinKinEnergy, ke);
1389  ke = TMath::Min(fMaxKinEnergyHA, ke);
1390 
1391  LOG("INukeData", pDEBUG) << "Querying hA cross section at ke = " << ke;
1392 
1393  if(hpdgc == kPdgProton) {
1394  /* handle protons */
1395  if (fate == kIHAFtCEx ) return TMath::Max(0., fFracPA_CEx -> Evaluate (ke));
1396  else if (fate == kIHAFtElas ) return TMath::Max(0., fFracPA_Elas -> Evaluate (ke));
1397  else if (fate == kIHAFtInelas ) return TMath::Max(0., fFracPA_Inel -> Evaluate (ke));
1398  else if (fate == kIHAFtAbs ) return TMath::Max(0., fFracPA_Abs -> Evaluate (ke));
1399  else if (fate == kIHAFtPiProd ) return TMath::Max(0., fFracPA_PiPro -> Evaluate (ke));
1400  else if (fate == kIHAFtCmp ) return TMath::Max(0., fFracPA_Cmp -> Evaluate (ke)); // cmp - add support for this later
1401  else {
1402  LOG("INukeData", pWARN)
1403  << "Protons don't have this fate: " << INukeHadroFates::AsString(fate);
1404  return 0;
1405  }
1406 
1407  } else if (hpdgc == kPdgNeutron) {
1408  /* handle neutrons */
1409  if (fate == kIHAFtCEx ) return TMath::Max(0., fFracNA_CEx -> Evaluate (ke));
1410  else if (fate == kIHAFtElas ) return TMath::Max(0., fFracNA_Elas -> Evaluate (ke));
1411  else if (fate == kIHAFtInelas ) return TMath::Max(0., fFracNA_Inel -> Evaluate (ke));
1412  else if (fate == kIHAFtAbs ) return TMath::Max(0., fFracNA_Abs -> Evaluate (ke));
1413  else if (fate == kIHAFtPiProd ) return TMath::Max(0., fFracNA_PiPro -> Evaluate (ke));
1414  else if (fate == kIHAFtCmp ) return TMath::Max(0., fFracNA_Cmp -> Evaluate (ke)); // cmp - add support for this later
1415  else {
1416  LOG("INukeData", pWARN)
1417  << "Neutrons don't have this fate: " << INukeHadroFates::AsString(fate);
1418  return 0;
1419  }
1420 
1421  } else if (hpdgc == kPdgKP) {
1422  /* handle K+ */
1423  if (fate == kIHAFtInelas ) return TMath::Max(0., fFracKA_Inel -> Evaluate (ke));
1424  else if (fate == kIHAFtAbs ) return TMath::Max(0., fFracKA_Abs -> Evaluate (ke));
1425  // else if (fate == kIHAFtElas ) return TMath::Max(0., fFracKA_Elas -> Evaluate (ke));
1426  //else if (fate == kIHAFtCEx ) return TMath::Max(0., fFracKA_CEx -> Evaluate(ke)); //why is this not used?
1427  else {
1428  LOG("INukeData", pWARN)
1429  << "K+'s don't have this fate: " << INukeHadroFates::AsString(fate);
1430  return 0;
1431  }
1432  }
1433  LOG("INukeData", pWARN)
1434  << "Can't handle particles with pdg code = " << hpdgc;
1435  return 0;
1436 }
1437 //____________________________________________________________________________
1438 double INukeHadroData2015::XSec(int hpdgc, INukeFateHN_t fate, double ke, int targA, int targZ) const
1439 {
1440 // return the x-section for the input fate for the particle with the input pdg
1441 // code at the input kinetic energy
1442 //
1443  ke = TMath::Max(fMinKinEnergy, ke);
1444  ke = TMath::Min(fMaxKinEnergyHN, ke);
1445 
1446  LOG("INukeData", pDEBUG) << "Querying hN cross section at ke = " << ke;
1447 
1448  double xsec=0;
1449 
1450  if (hpdgc == kPdgPiP) {
1451  /* handle pi+ */
1452  if (fate == kIHNFtCEx ) {xsec = TMath::Max(0., fXSecPipp_CEx -> Evaluate(ke)) * targZ;
1453  xsec+= TMath::Max(0., fXSecPipn_CEx -> Evaluate(ke)) * (targA-targZ);
1454  return xsec;}
1455  else if (fate == kIHNFtElas ) {xsec = TMath::Max(0., fXSecPipp_Elas -> Evaluate(ke)) * targZ;
1456  xsec+= TMath::Max(0., fXSecPipn_Elas -> Evaluate(ke)) * (targA-targZ);
1457  return xsec;}
1458  else if (fate == kIHNFtInelas) {xsec = TMath::Max(0., fXSecPipp_Reac -> Evaluate(ke)) * targZ;
1459  xsec+= TMath::Max(0., fXSecPipn_Reac -> Evaluate(ke)) * (targA-targZ);
1460  return xsec;}
1461  else if (fate == kIHNFtAbs ) {xsec = TMath::Max(0., fXSecPipd_Abs -> Evaluate(ke)) * targA;
1462  return xsec;}
1463  else {
1464  LOG("INukeData", pWARN)
1465  << "Pi+'s don't have this fate: " << INukeHadroFates::AsString(fate);
1466  return 0;
1467  }
1468 
1469  } else if (hpdgc == kPdgPiM) {
1470  /* handle pi- */
1471  if (fate == kIHNFtCEx ) {xsec = TMath::Max(0., fXSecPipn_CEx -> Evaluate(ke)) * targZ;
1472  xsec+= TMath::Max(0., fXSecPipp_CEx -> Evaluate(ke)) * (targA-targZ);
1473  return xsec;}
1474  else if (fate == kIHNFtElas ) {xsec = TMath::Max(0., fXSecPipn_Elas -> Evaluate(ke)) * targZ;
1475  xsec+= TMath::Max(0., fXSecPipp_Elas -> Evaluate(ke)) * (targA-targZ);
1476  return xsec;}
1477  else if (fate == kIHNFtInelas) {xsec = TMath::Max(0., fXSecPipn_Reac -> Evaluate(ke)) * targZ;
1478  xsec+= TMath::Max(0., fXSecPipp_Reac -> Evaluate(ke)) * (targA-targZ);
1479  return xsec;}
1480  else if (fate == kIHNFtAbs ) {xsec = TMath::Max(0., fXSecPipd_Abs -> Evaluate(ke)) * targA;
1481  return xsec;}
1482  else {
1483  LOG("INukeData", pWARN)
1484  << "Pi-'s don't have this fate: " << INukeHadroFates::AsString(fate);
1485  return 0;
1486  }
1487 
1488  } else if (hpdgc == kPdgPi0) {
1489  /* handle pi0 */
1490  if (fate == kIHNFtCEx ) {xsec = TMath::Max(0., fXSecPi0p_CEx -> Evaluate(ke)) * targZ;
1491  xsec+= TMath::Max(0., fXSecPi0n_CEx -> Evaluate(ke)) * (targA-targZ);
1492  return xsec;}
1493  else if (fate == kIHNFtElas ) {xsec = TMath::Max(0., fXSecPi0p_Elas -> Evaluate(ke)) * targZ;
1494  xsec+= TMath::Max(0., fXSecPi0n_Elas -> Evaluate(ke)) * (targA-targZ);
1495  return xsec;}
1496  else if (fate == kIHNFtInelas) {xsec = TMath::Max(0., fXSecPi0p_Reac -> Evaluate(ke)) * targZ;
1497  xsec+= TMath::Max(0., fXSecPi0n_Reac -> Evaluate(ke)) * (targA-targZ);
1498  return xsec;}
1499  else if (fate == kIHNFtAbs ) {xsec = TMath::Max(0., fXSecPi0d_Abs -> Evaluate(ke)) * targA;
1500  return xsec;}
1501  else {
1502  LOG("INukeData", pWARN)
1503  << "Pi0's don't have this fate: " << INukeHadroFates::AsString(fate);
1504  return 0;
1505  }
1506 
1507  } else if (hpdgc == kPdgProton) {
1508  /* handle protons */
1509  if (fate == kIHNFtElas ) {xsec = TMath::Max(0., fXSecPp_Elas -> Evaluate(ke)) * targZ;
1510  xsec+= TMath::Max(0., fXSecPn_Elas -> Evaluate(ke)) * (targA-targZ);
1511  return xsec;}
1512  else if (fate == kIHNFtInelas) {xsec = TMath::Max(0., fXSecPp_Reac -> Evaluate(ke)) * targZ;
1513  xsec+= TMath::Max(0., fXSecPn_Reac -> Evaluate(ke)) * (targA-targZ);
1514  return xsec;}
1515  else if (fate == kIHNFtCmp) {xsec = TMath::Max(0., fXSecPp_Cmp -> Evaluate(ke)) * targZ;
1516  xsec+= TMath::Max(0., fXSecPn_Cmp -> Evaluate(ke)) * (targA-targZ);
1517  return xsec;}
1518  else {
1519  LOG("INukeData", pWARN)
1520  << "Protons don't have this fate: " << INukeHadroFates::AsString(fate);
1521  return 0;
1522  }
1523 
1524  } else if (hpdgc == kPdgNeutron) {
1525  /* handle protons */
1526  if (fate == kIHNFtElas ) {xsec = TMath::Max(0., fXSecPn_Elas -> Evaluate(ke)) * targZ;
1527  xsec+= TMath::Max(0., fXSecNn_Elas -> Evaluate(ke)) * (targA-targZ);
1528  return xsec;}
1529  else if (fate == kIHNFtInelas) {xsec = TMath::Max(0., fXSecPn_Reac -> Evaluate(ke)) * targZ;
1530  xsec+= TMath::Max(0., fXSecNn_Reac -> Evaluate(ke)) * (targA-targZ);
1531  return xsec;}
1532  else if (fate == kIHNFtCmp) {xsec = TMath::Max(0., fXSecPp_Cmp -> Evaluate(ke)) * targZ;
1533  xsec+= TMath::Max(0., fXSecPn_Cmp -> Evaluate(ke)) * (targA-targZ);
1534  return xsec;}
1535  else {
1536  LOG("INukeData", pWARN)
1537  << "Neutrons don't have this fate: " << INukeHadroFates::AsString(fate);
1538  return 0;
1539  }
1540  //Adding here kaons, why elastic only on protons? hA or hN? No _Reac for kaons...
1541  } else if (hpdgc == kPdgKP) {
1542  /* handle K+ */
1543  if (fate == kIHNFtCEx ) {xsec = TMath::Max(0., fXSecKpn_CEx -> Evaluate(ke)) * targZ;
1544  xsec+= TMath::Max(0., fXSecKpn_CEx -> Evaluate(ke)) * (targA-targZ);
1545  return xsec;}
1546  else if (fate == kIHNFtElas ) {xsec = TMath::Max(0., fXSecKpn_Elas -> Evaluate(ke)) * targZ;
1547  xsec+= TMath::Max(0., fXSecKpn_Elas -> Evaluate(ke)) * (targA-targZ);
1548  return xsec;}
1549  /*else if (fate == kIHNFtAbs ) {xsec = TMath::Max(0., fXSecKpd_Abs -> Evaluate(ke)) * targA;
1550  return xsec;}*/
1551  else {
1552  LOG("INukeData", pWARN)
1553  << "K+'s don't have this fate: " << INukeHadroFates::AsString(fate);
1554  return 0;
1555  }
1556  //------------------------------------------------
1557  /* } else if (hpdgc == kPdgGamma) {
1558  / * handle gamma * /
1559  if (fate == kIHNFtInelas) {xsec = TMath::Max(0., fXSecGamp_fs -> Evaluate(ke)) * targZ;
1560  xsec+= TMath::Max(0., fXSecGamn_fs -> Evaluate(ke)) * (targA-targZ);
1561  return xsec;}
1562  else {
1563  LOG("INukeData", pWARN)
1564  << "Gamma's don't have this fate: " << INukeHadroFates::AsString(fate);
1565  return 0;
1566  }*/
1567  }
1568  LOG("INukeData", pWARN)
1569  << "Can't handle particles with pdg code = " << hpdgc;
1570 
1571  return 0;
1572 }
1573 
1574 double INukeHadroData2015::Frac(int hpdgc, INukeFateHN_t fate, double ke, int targA, int targZ) const
1575 {
1576 // return the x-section fraction for the input fate for the particle with the
1577 // input pdg code at the input kinetic energy
1578 
1579  ke = TMath::Max(fMinKinEnergy, ke);
1580  ke = TMath::Min(fMaxKinEnergyHN, ke);
1581 
1582  // get x-section
1583  double xsec = this->XSec(hpdgc,fate,ke,targA,targZ);
1584 
1585  // get max x-section
1586  double xsec_tot = 0;
1587  if (hpdgc == kPdgPiP ){xsec_tot = TMath::Max(0., fXSecPipp_Tot -> Evaluate(ke)) * targZ;
1588  xsec_tot+= TMath::Max(0., fXSecPipn_Tot -> Evaluate(ke)) * (targA-targZ);}
1589  else if (hpdgc == kPdgPiM ){xsec_tot = TMath::Max(0., fXSecPipn_Tot -> Evaluate(ke)) * targZ;
1590  xsec_tot+= TMath::Max(0., fXSecPipp_Tot -> Evaluate(ke)) * (targA-targZ);}
1591  else if (hpdgc == kPdgPi0 ){xsec_tot = TMath::Max(0., fXSecPi0p_Tot -> Evaluate(ke)) * targZ;
1592  xsec_tot+= TMath::Max(0., fXSecPi0n_Tot -> Evaluate(ke)) * (targA-targZ);}
1593  else if (hpdgc == kPdgProton ){xsec_tot = TMath::Max(0., fXSecPp_Tot -> Evaluate(ke)) * targZ;
1594  xsec_tot+= TMath::Max(0., fXSecPn_Tot -> Evaluate(ke)) * (targA-targZ);}
1595  else if (hpdgc == kPdgNeutron){xsec_tot = TMath::Max(0., fXSecPn_Tot -> Evaluate(ke)) * targZ;
1596  xsec_tot+= TMath::Max(0., fXSecNn_Tot -> Evaluate(ke)) * (targA-targZ);}
1597  else if (hpdgc == kPdgGamma ) xsec_tot = TMath::Max(0., fXSecGamN_Tot -> Evaluate(ke));
1598  else if (hpdgc == kPdgKP ) xsec_tot = TMath::Max(0., fXSecKpN_Tot -> Evaluate(ke));
1599 
1600  // compute fraction
1601  double frac = (xsec_tot>0) ? xsec/xsec_tot : 0.;
1602  return frac;
1603 }
1604 //____________________________________________________________________________
1606 {
1607  // This method returns a random cos(ang) according to a distribution
1608  // based upon the particle and fate. The sampling uses the
1609  // Accept/Reject method, whereby a distribution is bounded above by
1610  // an envelope, or in this case, a number of envelopes, which can be
1611  // easily sampled (here, we use uniform distributions).
1612  // To get a random value, first the envelope is sampled to
1613  // obtain an x-coordinate (cos(ang)), and then another random value
1614  // is obtained uniformally in the range [0,h(j,0)], where h(j,0)
1615  // is the height of the j-th envelope. If the point is beneath the
1616  // distribution, the x-coordinate is accepted, otherwise, we try
1617  // again.
1618 
1619  RandomGen * rnd = RandomGen::Instance();
1620 
1621  // numEnv is the number of envelopes in the total envelope,
1622  // that is, the number of seperate simple uniform distributions
1623  // that will be fit against the distribution in question in the
1624  // Accept/Reject process of sampling
1625  int numEnv = 4;
1626  int numPoints = 1000; // The number of points to be evaluated
1627  // for the purpose of finding the max
1628  // value of the distribution
1629  assert((numPoints%numEnv)==0); // numPoints/numEnv has to be an integer
1630  double sr = 2.0 / numEnv; // Subrange, i.e., range of an envelope
1631  double cstep = 2.0 / (numPoints); // Magnitude of the step between eval. points
1632 
1633  double ke = (p->E() - p->Mass()) * 1000.0; // ke in MeV
1634  if (TMath::Abs((int)ke-ke)<.01) ke+=.3; // make sure ke isn't an integer,
1635  // otherwise sometimes gives weird results
1636  // due to ROOT's Interpolate() function
1637  double avg = 0.0; // average value in envelop
1638 
1639  // Matrices to hold data; buff holds the distribution
1640  // data per envelope from which the max value is
1641  // obtained. That value is then recorded in dist, where
1642  // the integral of the envelope to that point is
1643  // also recorded
1644 
1645  double * buff = new double[numPoints/numEnv + 1];
1646  double ** dist = new double*[numEnv];
1647  for(int ih=0;ih<numEnv;ih++)
1648  {
1649  dist[ih] = new double[3];
1650  }
1651 
1652  // Acc-Rej Sampling Method
1653  // -- Starting at the beginning of each envelope,
1654  // this loop evaluates (numPoints) amount of points
1655  // on the distribution and holds them in buff;
1656  // then takes the max and places it in the first row
1657  // of h. The second row of h contains the interval of
1658  // the total envelope, up to the current envelope.
1659  // Thus, when properly normalized, the last value
1660  // in the second row of h should be 1.
1661  double totxsec = 0.0;
1662  for(int i=0;i<numEnv;i++)
1663  {
1664  double lbound = -1 + i*sr;
1665 
1666  for(int j=0;j<=numPoints / numEnv; j++)
1667  {
1668  buff[j] = this->XSec(p->Pdg(),target,scode,fate,ke,lbound+j*cstep);
1669  avg += buff[j];
1670  }
1671 
1672  totxsec+=avg;
1673  avg/= (double(numPoints)/double(numEnv));
1674  dist[i][0] = TMath::MaxElement(numPoints/numEnv+1,buff);
1675  dist[i][1] = avg;
1676  dist[i][2] = dist[i][1] + ((i==0)?0.0:dist[i-1][2]);
1677  avg=0.0;
1678  }
1679 
1680 
1681  delete [] buff;
1682 
1683  int iter=1; // keep track of iterations
1684  int env=0; // envelope index
1685  double rval = 0.0; // random value
1686  double val = 0.0; // angle value
1687 
1688  // Get a random point, see if its in the distribution, and if not
1689  // then try again.
1690 
1691  rval = rnd->RndFsi().Rndm()*dist[numEnv-1][2];
1692 
1693  env=0;
1694  // Check which envelope it's in, to
1695  // get proper height
1696  while(env<numEnv)
1697  {
1698  if(rval<=dist[env][2]) break;
1699  else env++;
1700  }
1701  if(env==numEnv) env=numEnv - 1;
1702 
1703 while(iter)
1704  {
1705 
1706  // Obtain the correct x-coordinate from the random sample
1707  val = rnd->RndFsi().Rndm()*sr;
1708  val += sr*env-1;
1709  rval = rnd->RndFsi().Rndm()*dist[env][0];
1710 
1711  // Test to see if point is in distribution, if it is, stop and return
1712  if(rval < this->XSec(p->Pdg(),target,scode,fate,ke,val)) break;
1713 
1714  // Possibly an extremely long loop, don't want to
1715  // hold up the program
1716  if(iter==1000)
1717  {
1718  int NUM_POINTS=2000;
1719  int pvalues=0;
1720  double points[200]={0};
1721  for(int k=0;k<NUM_POINTS;k++)
1722  {
1723  points[int(k/10)]=this->XSec(p->Pdg(),target,scode,fate,ke,-1+(2.0/NUM_POINTS)*k);
1724  if(points[int(k/10)]>0) pvalues++;
1725  }
1726  if(pvalues<(.05*NUM_POINTS))
1727  {
1728  // if it reaches here, one more test...if momenta of particle is
1729  // extremely low, just give it an angle from a uniform distribution
1730  if(p->P4()->P()<.005) // 5 MeV
1731  {
1732  val = 2*rnd->RndFsi().Rndm()-1;
1733  break;
1734  }
1735  else
1736  {
1737  LOG("Intranuke", pWARN) << "Hung-up in IntBounce method - Exiting";
1738  LOG("Intranuke", pWARN) << (*p);
1739  LOG("Intranuke", pWARN) << "Target: " << target << ", Scode: " << scode << ", fate: " << INukeHadroFates::AsString(fate);
1740  for(int ie=0;ie<200;ie+=10) {
1741  LOG("Intranuke", pWARN) << points[ie+0] << ", " << points[ie+1] << ", " << points[ie+2] << ", "
1742  << points[ie+3] << ", " << points[ie+4] << ", " << points[ie+5] << ", " << points[ie+6] << ", "
1743  << points[ie+7] << ", " << points[ie+8] << ", " << points[ie+9];
1744  }
1745  for(int ih=0;ih<numEnv;ih++)
1746  {
1747  delete [] dist[ih];
1748  }
1749  delete [] dist;
1750 
1751  return -2.;
1752  }
1753  }
1754  }
1755  iter++;
1756  }
1757 
1758  for(int ih=0;ih<numEnv;ih++)
1759  {
1760  delete [] dist[ih];
1761  }
1762  delete [] dist;
1763 
1764  return val;
1765 }
1766 //___________________________________________________________________________
enum genie::EINukeFateHA_t INukeFateHA_t
Basic constants.
TRandom3 & RndFsi(void) const
rnd number generator used by intranuclear cascade monte carlos
Definition: RandomGen.h:60
static INukeHadroData2015 * fInstance
#include "Numerical/GSFunc.h"
Definition: AlgCmp.h:26
#define pERROR
Definition: Messenger.h:50
static INukeHadroData2015 * Instance(void)
double E(void) const
Get energy.
Definition: GHepParticle.h:90
double IntBounce(const GHepParticle *p, int target, int s1, INukeFateHN_t fate)
static RandomGen * Instance()
Access instance.
Definition: RandomGen.cxx:79
A numeric analysis tool class for interpolating 1-D functions.
Definition: Spline.h:47
double FracAIndep(int hpdgc, INukeFateHA_t fate, double ke) const
double Mass(void) const
Mass that corresponds to the PDG code.
double Frac(int hpdgc, INukeFateHN_t fate, double ke, int targA=0, int targZ=0) const
A singleton holding random number generator classes. All random number generation in GENIE should tak...
Definition: RandomGen.h:30
int Pdg(void) const
Definition: GHepParticle.h:64
double y
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:87
double XSec(int hpdgc, int tgt, int nprod, INukeFateHN_t rxnType, double ke, double costh) const
const int kPdgGamma
Definition: PDGCodes.h:163
const int kPdgKP
Definition: PDGCodes.h:146
const int kPdgPiP
Definition: PDGCodes.h:132
const int kPdgPi0
Definition: PDGCodes.h:134
def Interpolate(x1, y1, x2, y2, yvalue)
Definition: HandyFuncs.py:200
#define pINFO
Definition: Messenger.h:53
enum genie::EINukeFateHN_t INukeFateHN_t
#define pWARN
Definition: Messenger.h:51
p
Definition: test.py:228
static string AsString(INukeFateHN_t fate)
const int kPdgPiM
Definition: PDGCodes.h:133
void ReadhNFile(string filename, double ke, int npoints, int &curr_point, double *costh_array, double *xsec_array, int cols)
const int kPdgProton
Definition: PDGCodes.h:62
#define pNOTICE
Definition: Messenger.h:52
TLorentzVector * P4(void) const
Definition: GHepParticle.h:79
static const double sr
Definition: Units.h:167
const int kPdgNeutron
Definition: PDGCodes.h:64
static const double kPi
Definition: Constants.h:38
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils...
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:40
double FracADep(int hpdgc, INukeFateHA_t fate, double ke, int targA) const
#define pDEBUG
Definition: Messenger.h:54