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