Jpp  17.3.1
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JEarth.cc File Reference

Mickey-mouse program to propagate neutrinos through the Earth. More...

#include <string>
#include <iostream>
#include <cmath>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2D.h"
#include "TRandom.h"
#include "JPhysics/JNeutrino.hh"
#include "JPhysics/JConstants.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Detailed Description

Mickey-mouse program to propagate neutrinos through the Earth.

Author
mdejong

Definition in file JEarth.cc.

Function Documentation

int main ( int  argc,
char *  argv[] 
)

Definition at line 33 of file JEarth.cc.

34 {
35  using namespace std;
36  using namespace JPP;
37 
38  typedef long long int counter_type;
39 
40  string outputFile;
41  double E_GeV;
42  counter_type numberOfEvents;
43  UInt_t seed;
44  char option;
45  int debug;
46 
47  try {
48 
49  JParser<> zap("Mickey-mouse program to propagate neutrinos through the Earth.");
50 
51  zap['o'] = make_field(outputFile, "output file") = "earth.root";
52  zap['E'] = make_field(E_GeV, "neutrino energy at source");
53  zap['n'] = make_field(numberOfEvents, "number of generated events");
54  zap['S'] = make_field(seed) = 0;
55  zap['O'] = make_field(option, endl
56  << nc_enabled << " -> " << "Neutral-current interactions as simulated." << endl
57  << nc_no_energy_loss << " -> " << "Neutral-current interactions with Bjorken-Y set to zero." << endl
58  << nc_no_scattering << " -> " << "Neutral-current interactions with transverse energy set to zero." << endl
59  << nc_disabled << " -> " << "Neutral-current interactions which have no effect." << endl
60  << nc_absorption << " -> " << "Neutral-current interactions lead to absorption." << endl) =
61  nc_enabled,
62  nc_no_energy_loss,
63  nc_no_scattering,
64  nc_disabled,
65  nc_absorption;
66 
67  zap['d'] = make_field(debug) = 0;
68 
69  zap(argc, argv);
70  }
71  catch(const exception &error) {
72  FATAL(error.what() << endl);
73  }
74 
75 
76  gRandom->SetSeed(seed);
77 
78 
79  const double R_SOURCE_M = 40000.0; // [m]
80  const double R_TARGET_M = 500.0; // [m]
81  const double EMIN_GEV = 10.0; // [GeV]
82 
83 
84  const double Zmin = 0.0; // [m]
85  const double Zmax = R_EARTH_KM * 2.0e3; // [m]
86 
87  enum {
88  CC = 0, // charged-current interaction
89  NC, // neutral-current interaction
90  N // number of interactions
91  };
92 
93  const JCrossSection* sigma[N] = { NULL }; // cross-sections
94 
95  counter_type number_of_events[2][2] = { 0 }; // [source inside/outside][target inside/outside]
96 
97  TFile out(outputFile.c_str(), "recreate");
98 
99  TH2D hs("hs", NULL,
100  101, -R_SOURCE_M, +R_SOURCE_M,
101  101, -R_SOURCE_M, +R_SOURCE_M);
102 
103  TH2D ht("ht", NULL,
104  101, -R_SOURCE_M, +R_SOURCE_M,
105  101, -R_SOURCE_M, +R_SOURCE_M);
106 
107  TH1D hn("hn", NULL, 11, -0.5, 10.5);
108 
109  TH1D ha("ha", NULL, 100, 0.0, 1.0);
110  TH1D hb("hb", NULL, 100, 0.0, 1.0);
111 
112  TH2D h2("h2", NULL,
113  101, -0.01, +0.01,
114  101, -0.01, +0.01);
115 
116 
117  for (counter_type count = 0; count != numberOfEvents; ++count) {
118 
119  if (count%10000 == 0) {
120  STATUS("event " << setw(9) << count << '\r'); DEBUG(endl);
121  }
122 
123  const bool neutrino = (count%2 == 0); // [anti-]neutrino ratio 1:1
124 
125  if (neutrino) {
126  sigma[CC] = &cc_nu;
127  sigma[NC] = &nc_nu;
128  } else {
129  sigma[CC] = &cc_nubar;
130  sigma[NC] = &nc_nubar;
131  }
132 
133  const double r2 = gRandom->Uniform(0.0, R_SOURCE_M * R_SOURCE_M);
134  const double r1 = sqrt(r2);
135  const double phi = gRandom->Uniform(-PI, +PI);
136 
137  // start values
138 
139  double x = r1 * cos(phi);
140  double y = r1 * sin(phi);
141  double z = Zmin;
142  double Tx = 0.0;
143  double Ty = 0.0;
144  double E = E_GeV;
145 
146  hs.Fill(x,y);
147 
148  const int i1 = (sqrt(x*x + y*y) <= R_TARGET_M ? 0 : 1); // inside -> 0; outside -> 1
149 
150  int ni[N] = { 0 }; // number of interactions
151 
152  while (ni[CC] == 0 && E > EMIN_GEV) {
153 
154  double li[N]; // inverse interaction length [m^-1]
155  double ls = 0.0;
156 
157  for (int i = 0; i != N; ++i) {
158  ls += li[i] = 1.0e+2 * (*sigma[i])(E) * DENSITY_EARTH * AVOGADRO;
159  }
160 
161  double dz = gRandom->Exp(1.0) / ls;
162 
163  if (dz > Zmax - z) {
164  dz = Zmax - z;
165  }
166 
167  // step
168 
169  if (option != nc_disabled && option != nc_no_scattering) {
170  x += dz * Tx;
171  y += dz * Ty;
172  }
173  z += dz;
174 
175  if (z >= Zmax) {
176 
177  break; // ready
178 
179  } else {
180 
181  double y = gRandom->Uniform(ls);
182 
183  if ((y-= li[CC]) < 0.0) { // charged-current interaction
184 
185  ++ni[CC];
186 
187  break; // stop
188  }
189 
190  if ((y-= li[NC]) < 0.0) { // neutral-current interaction
191 
192  ++ni[NC];
193 
194  if (option == nc_absorption) {
195  break;
196  }
197 
198  // simplified neutrino scattering
199 
200  double s = E * (MASS_PROTON + MASS_NEUTRON);
201  double ET = 0.5 * sqrt(s) * gRandom->Uniform(0.0, 1.0);
202  double phi = gRandom->Uniform(-PI, +PI);
203 
204  double By = gRandom->Uniform(0.0, 1.0);
205 
206  if (!neutrino) {
207  while (gRandom->Rndm() > (1.0 - By) * (1.0 - By)) {
208  By = gRandom->Uniform(0.0, 1.0);
209  }
210  }
211 
212  if (option == nc_disabled || option == nc_no_scattering) {
213  ET = 0.0;
214  phi = 0.0;
215  }
216 
217  if (option == nc_disabled || option == nc_no_energy_loss) {
218  By = 0.0;
219  }
220 
221  if (neutrino)
222  ha.Fill(By);
223  else
224  hb.Fill(By);
225 
226  Tx = ET * cos(phi) / E;
227  Ty = ET * sin(phi) / E;
228  E = (1.0 - By) * E;
229  }
230  }
231  }
232 
233  if (z >= Zmax) { // reach target
234  hn.Fill((double) ni[NC]);
235  ht.Fill(x,y);
236  h2.Fill(Tx,Ty);
237  }
238 
239  if (z >= Zmax && sqrt(x*x + y*y) <= R_TARGET_M) { // hit target
240 
241  number_of_events[i1][0] += 1; // count events
242 
243  if (ni[NC] == 0) {
244  number_of_events[i1][1] += 1; // count events
245  }
246  }
247  }
248  STATUS(endl);
249 
250  // summary
251 
252  const double n00 = number_of_events[0][0]; // number of events in target from source area covered by target
253  const double n10 = number_of_events[1][0]; // number of events in target from source area not covered by target
254  const double n01 = number_of_events[0][1]; // number of events in target from source area covered by target, without neutral-current interaction
255 
256  if (option == nc_enabled) {
257  cout << " " << SCIENTIFIC(7,1) << E_GeV << flush;
258  }
259  {
260  cout << " & " << FIXED(7,0) << n00 << flush;
261  }
262  if (option != nc_disabled && option != nc_no_scattering) {
263  cout << " & " << FIXED(7,0) << n10 << flush;
264  }
265  if (option == nc_enabled) {
266  cout << " & " << FIXED(7,0) << n01 << " & " << FIXED(5,2) << (n00 + n10) / n01 << flush;
267  }
268  if (option == nc_disabled) {
269  cout << " \\\\" << endl;
270  }
271 
272  if (debug >= status_t) {
273 
274  cout << ' ' << setw(8) << right << "inside" << ' ' << setw(8) << right << "outside" << endl;
275 
276  for (int row = 0; row != N; ++row) {
277  for (int col = 0; col != N; ++col) {
278  cout << ' ' << setw(8) << number_of_events[row][col];
279  }
280  cout << endl;
281  }
282  }
283 
284  out.Write();
285  out.Close();
286 }
Utility class to parse command line options.
Definition: JParser.hh:1517
static const double R_EARTH_KM
Geophysics constants.
static const JNCnu nc_nu
Function object for neutral current neutrino cross section [cm^2] as a function of neutrino energy [G...
Definition: JNeutrino.hh:620
then usage E
Definition: JMuonPostfit.sh:35
#define STATUS(A)
Definition: JMessage.hh:63
then JShowerPostfit f $INPUT_FILE o $OUTPUT_FILE N
static const double AVOGADRO
Avogadro&#39;s number [gr^-1].
Long64_t counter_type
Type definition for counter.
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:446
string outputFile
static const JCCnubar cc_nubar
Function object for charged current anti-neutrino cross section [cm^2] as a function of neutrino ener...
Definition: JNeutrino.hh:621
const double sigma[]
Definition: JQuadrature.cc:74
static const double MASS_NEUTRON
neutron mass [GeV]
then usage set_variable ACOUSTICS_WORKDIR $WORKDIR set_variable FORMULA sin([0]+2 *$PI *([1]+[2]*x)*x)" set_variable DY 1.0e-8 mkdir $WORKDIR for DETECTOR in $DETECTORS[*]
static const JNCnubar nc_nubar
Function object for neutral current anti-neutrino cross section [cm^2] as a function of neutrino ener...
Definition: JNeutrino.hh:622
status
Definition: JMessage.hh:30
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1993
static const double PI
Mathematical constants.
static const JCCnu cc_nu
Function object for charged current neutrino cross section [cm^2] as a function of neutrino energy [G...
Definition: JNeutrino.hh:618
#define FATAL(A)
Definition: JMessage.hh:67
static const double MASS_PROTON
proton mass [GeV]
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:484
int debug
debug level
static const double DENSITY_EARTH
Average density of the Earth [gr/cm³].
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62