Jpp  test_elongated_shower_pde
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JEarth.cc
Go to the documentation of this file.
1 #include <string>
2 #include <iostream>
3 #include <cmath>
4 
5 #include "TROOT.h"
6 #include "TFile.h"
7 #include "TH1D.h"
8 #include "TH2D.h"
9 #include "TRandom.h"
10 
11 #include "JPhysics/JNeutrino.hh"
12 #include "JPhysics/JConstants.hh"
13 
14 #include "Jeep/JParser.hh"
15 #include "Jeep/JMessage.hh"
16 
17 namespace {
18 
19  const char nc_enabled = 'A'; //!< Neutral-current interactions as simulated.
20  const char nc_no_energy_loss = 'C'; //!< Neutral-current interactions with Bjorken-Y set to zero.
21  const char nc_no_scattering = 'D'; //!< Neutral-current interactions with transverse energy set to zero.
22  const char nc_disabled = 'E'; //!< Neutral-current interactions which have no effect.
23  const char nc_absorption = 'F'; //!< Neutral-current interactions lead to absorption.
24 }
25 
26 
27 /**
28  * \file
29  *
30  * Mickey-mouse program to propagate neutrinos through the Earth.
31  * \author mdejong
32  */
33 int main(int argc, char* argv[])
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, MAKE_CSTRING(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 DIAMETER_EARTH_KM = 12724.0; // [km]
80  const double DENSITY_EARTH = 5.5; // [gr/cm³]
81  const double R_SOURCE_M = 40000.0; // [m]
82  const double R_TARGET_M = 500.0; // [m]
83  const double EMIN_GEV = 10.0; // [GeV]
84 
85 
86  const double Zmin = 0.0; // [m]
87  const double Zmax = DIAMETER_EARTH_KM * 1.0e3; // [m]
88 
89  enum {
90  CC = 0, // charged-current interaction
91  NC, // neutral-current interaction
92  N // number of interactions
93  };
94 
95  const JCrossSection* sigma[N] = { NULL }; // cross-sections
96 
97  counter_type number_of_events[2][2] = { 0 }; // [source inside/outside][target inside/outside]
98 
99  TFile out(outputFile.c_str(), "recreate");
100 
101  TH2D hs("hs", NULL,
102  101, -R_SOURCE_M, +R_SOURCE_M,
103  101, -R_SOURCE_M, +R_SOURCE_M);
104 
105  TH2D ht("ht", NULL,
106  101, -R_SOURCE_M, +R_SOURCE_M,
107  101, -R_SOURCE_M, +R_SOURCE_M);
108 
109  TH1D hn("hn", NULL, 11, -0.5, 10.5);
110 
111  TH1D ha("ha", NULL, 100, 0.0, 1.0);
112  TH1D hb("hb", NULL, 100, 0.0, 1.0);
113 
114  TH2D h2("h2", NULL,
115  101, -0.01, +0.01,
116  101, -0.01, +0.01);
117 
118 
119  for (counter_type count = 0; count != numberOfEvents; ++count) {
120 
121  if (count%10000 == 0) {
122  STATUS("event " << setw(9) << count << '\r'); DEBUG(endl);
123  }
124 
125  const bool neutrino = (count%2 == 0); // [anti-]neutrino ratio 1:1
126 
127  if (neutrino) {
128  sigma[CC] = &cc_nu;
129  sigma[NC] = &nc_nu;
130  } else {
131  sigma[CC] = &cc_nubar;
132  sigma[NC] = &nc_nubar;
133  }
134 
135  const double r2 = gRandom->Uniform(0.0, R_SOURCE_M * R_SOURCE_M);
136  const double r1 = sqrt(r2);
137  const double phi = gRandom->Uniform(-PI, +PI);
138 
139  // start values
140 
141  double x = r1 * cos(phi);
142  double y = r1 * sin(phi);
143  double z = Zmin;
144  double Tx = 0.0;
145  double Ty = 0.0;
146  double E = E_GeV;
147 
148  hs.Fill(x,y);
149 
150  const int i1 = (sqrt(x*x + y*y) <= R_TARGET_M ? 0 : 1); // inside -> 0; outside -> 1
151 
152  int ni[N] = { 0 }; // number of interactions
153 
154  while (ni[CC] == 0 && E > EMIN_GEV) {
155 
156  double li[N]; // inverse interaction length [m^-1]
157  double ls = 0.0;
158 
159  for (int i = 0; i != N; ++i) {
160  ls += li[i] = 1.0e+2 * (*sigma[i])(E) * DENSITY_EARTH * AVOGADRO;
161  }
162 
163  double dz = gRandom->Exp(1.0) / ls;
164 
165  if (dz > Zmax - z) {
166  dz = Zmax - z;
167  }
168 
169  // step
170 
171  if (option != nc_disabled && option != nc_no_scattering) {
172  x += dz * Tx;
173  y += dz * Ty;
174  }
175  z += dz;
176 
177  if (z >= Zmax) {
178 
179  break; // ready
180 
181  } else {
182 
183  double y = gRandom->Uniform(ls);
184 
185  if ((y-= li[CC]) < 0.0) { // charged-current interaction
186 
187  ++ni[CC];
188 
189  break; // stop
190  }
191 
192  if ((y-= li[NC]) < 0.0) { // neutral-current interaction
193 
194  ++ni[NC];
195 
196  if (option == nc_absorption) {
197  break;
198  }
199 
200  // simplified neutrino scattering
201 
202  double s = E * (MASS_PROTON + MASS_NEUTRON);
203  double ET = 0.5 * sqrt(s) * gRandom->Uniform(0.0, 1.0);
204  double phi = gRandom->Uniform(-PI, +PI);
205 
206  double By = gRandom->Uniform(0.0, 1.0);
207 
208  if (!neutrino) {
209  while (gRandom->Rndm() > (1.0 - By) * (1.0 - By)) {
210  By = gRandom->Uniform(0.0, 1.0);
211  }
212  }
213 
214  if (option == nc_disabled || option == nc_no_scattering) {
215  ET = 0.0;
216  phi = 0.0;
217  }
218 
219  if (option == nc_disabled || option == nc_no_energy_loss) {
220  By = 0.0;
221  }
222 
223  if (neutrino)
224  ha.Fill(By);
225  else
226  hb.Fill(By);
227 
228  Tx = ET * cos(phi) / E;
229  Ty = ET * sin(phi) / E;
230  E = (1.0 - By) * E;
231  }
232  }
233  }
234 
235  if (z >= Zmax) { // reach target
236  hn.Fill((double) ni[NC]);
237  ht.Fill(x,y);
238  h2.Fill(Tx,Ty);
239  }
240 
241  if (z >= Zmax && sqrt(x*x + y*y) <= R_TARGET_M) { // hit target
242 
243  number_of_events[i1][0] += 1; // count events
244 
245  if (ni[NC] == 0) {
246  number_of_events[i1][1] += 1; // count events
247  }
248  }
249  }
250  STATUS(endl);
251 
252  // summary
253 
254  const double n00 = number_of_events[0][0]; // number of events in target from source area covered by target
255  const double n10 = number_of_events[1][0]; // number of events in target from source area not covered by target
256  const double n01 = number_of_events[0][1]; // number of events in target from source area covered by target, without neutral-current interaction
257 
258  if (option == nc_enabled) {
259  cout << " " << SCIENTIFIC(7,1) << E_GeV << flush;
260  }
261  {
262  cout << " & " << FIXED(7,0) << n00 << flush;
263  }
264  if (option != nc_disabled && option != nc_no_scattering) {
265  cout << " & " << FIXED(7,0) << n10 << flush;
266  }
267  if (option == nc_enabled) {
268  cout << " & " << FIXED(7,0) << n01 << " & " << FIXED(5,2) << (n00 + n10) / n01 << flush;
269  }
270  if (option == nc_disabled) {
271  cout << " \\\\" << endl;
272  }
273 
274  if (debug >= status_t) {
275 
276  cout << ' ' << setw(8) << right << "inside" << ' ' << setw(8) << right << "outside" << endl;
277 
278  for (int row = 0; row != N; ++row) {
279  for (int col = 0; col != N; ++col) {
280  cout << ' ' << setw(8) << number_of_events[row][col];
281  }
282  cout << endl;
283  }
284  }
285 
286  out.Write();
287  out.Close();
288 }
Utility class to parse command line options.
Definition: JParser.hh:1500
int main(int argc, char *argv[])
Definition: Main.cc:15
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
#define MAKE_CSTRING(A)
Make C-string.
Definition: JPrint.hh:151
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
static const double MASS_NEUTRON
neutron mass [GeV]
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:1961
Physics constants.
static const double PI
Mathematical constants.
int debug
debug level
Definition: JSirene.cc:68
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
General purpose messaging.
#define FATAL(A)
Definition: JMessage.hh:67
std::vector< int > count
Definition: JAlgorithm.hh:180
Utility class to parse command line options.
static const double MASS_PROTON
proton mass [GeV]
then usage $script< input_file >< detector_file > fi set_variable OUTPUT_DIR set_variable SELECTOR JDAQTimesliceL1 set_variable DEBUG case set_variable DEBUG
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:484