Jpp 19.3.0-rc.3
the software that should make you happy
Loading...
Searching...
No Matches
JGeane.cc File Reference

Example program to histogram muon energy loss. More...

#include <string>
#include <iostream>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2D.h"
#include "JPhysics/JGeane.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

Example program to histogram muon energy loss.

Author
mdejong

Definition in file JGeane.cc.

Function Documentation

◆ main()

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

Definition at line 21 of file JGeane.cc.

22{
23 using namespace std;
24 using namespace JPP;
25
26 string outputFile;
27 int debug;
28
29 try {
30
31 JParser<> zap("Example program to histogram muon energy loss.");
32
33 zap['o'] = make_field(outputFile) = "geane.root";
34 zap['d'] = make_field(debug) = 2;
35
36 zap(argc, argv);
37 }
38 catch(const exception &error) {
39 FATAL(error.what() << endl);
40 }
41
42
43 TFile out(outputFile.c_str(), "recreate");
44
45 TH1D h1("h1", NULL, 900, -1.0, +8.0);
46 TH2D h2("h2", NULL, 90, -1.0, +8.0, 100, 0.0, 5.0);
47 TH2D h3("h3", NULL, 130, -1.0, +12.0, 90, -1.0, +8.0);
48
49 for(int ix = 1; ix <= h1.GetNbinsX(); ++ix) {
50
51 const double x = h1.GetBinCenter(ix);
52 const double E = pow(10.0, x);
53
54 h1.SetBinContent(ix, gWater(E));
55 }
56
57 for(int ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
58 for(int iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
59
60 const double x = h2.GetXaxis()->GetBinCenter(ix);
61 const double y = h2.GetYaxis()->GetBinCenter(iy);
62
63 const double E = pow(10.0, x);
64 const double dx = pow(10.0, y);
65
66 h2.SetBinContent(ix, iy, gWater(E, dx));
67 }
68 }
69
70 const double precision = 1.0e-3;
71
72 for(int ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
73 for(int iy = 1; iy <= ix; ++iy) {
74
75 const double x = h3.GetXaxis()->GetBinCenter(ix);
76 const double y = h3.GetYaxis()->GetBinCenter(iy);
77
78 const double E1 = pow(10.0, x);
79 const double E2 = pow(10.0, y);
80
81 double Rmin = 0.0;
82 double Rmax = gWater(E1);
83
84 if (E2 < E1 && Rmax > Rmin) {
85
86 double R = 0.0;
87
88 for ( ; ; ) {
89
90 R = 0.5 * (Rmax + Rmin);
91
92 const double E = gWater(E1, R);
93
94 if (fabs(E - E2) <= precision) {
95 break;
96 }
97
98 if (E > E2)
99 Rmin = R;
100 else
101 Rmax = R;
102 }
103
104 h3.SetBinContent(ix, iy, R);
105 }
106 }
107 }
108
109 out.Write();
110 out.Close();
111}
string outputFile
#define FATAL(A)
Definition JMessage.hh:67
int debug
debug level
Definition JSirene.cc:72
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
Utility class to parse command line options.
Definition JParser.hh:1698
T pow(const T &x, const double y)
Power .
Definition JMath.hh:97
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).