Jpp  master_rocky-43-ge265d140c
the software that should make you happy
Functions
JDIS.cc File Reference
#include <string>
#include <iostream>
#include <iomanip>
#include <map>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TRandom3.h"
#include "JPhysics/JDIS.hh"
#include "JROOT/JManager.hh"
#include "JROOT/JRootToolkit.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"

Go to the source code of this file.

Functions

int main (int argc, char **argv)
 Example application to display DIS of muon. More...
 

Function Documentation

◆ main()

int main ( int  argc,
char **  argv 
)

Example application to display DIS of muon.

Author
mdejong

Definition at line 25 of file JDIS.cc.

26 {
27  using namespace std;
28  using namespace JPP;
29 
30  string outputFile;
31  int numberOfEvents;
32  UInt_t seed;
33  int debug;
34 
35  try {
36 
37  JParser<> zap("Example application to display DIS of muon.");
38 
39  zap['o'] = make_field(outputFile) = "dis.root";
40  zap['n'] = make_field(numberOfEvents) = 0;
41  zap['S'] = make_field(seed) = 0;
42  zap['d'] = make_field(debug) = 0;
43 
44  zap(argc, argv);
45  }
46  catch(const exception &error) {
47  FATAL(error.what() << endl);
48  }
49 
50  gRandom->SetSeed(seed);
51 
52 
53  const JDIS dis;
54 
55 
56  TH1D h0("sigma", NULL, 1000, -1.0, 9.0);
57 
58  JManager<double, TH1D> HA(new TH1D("v [% GeV]", NULL, 1000, -5.0, 0.0));
59  JManager<double, TH1D> H1(new TH1D("ran[% GeV]", NULL, 1000, -5.0, 0.0));
60  JManager<double, TH1D> HB(new TH1D("RMS[% GeV]", NULL, 1000, -5.0, 0.0));
61 
62  H1->Sumw2();
63 
64  const map<double, double> setups = { { 1.0e2, 0.7 },
65  { 1.0e3, 0.6 },
66  { 1.0e4, 0.5 } };
67 
68 
69  for (int i = 1; i <= h0.GetNbinsX(); ++i) {
70 
71  const double x = h0.GetBinCenter(i);
72  const double E = pow(10.0, x);
73  const double y = dis.getCrossSection(E);
74 
75  h0.SetBinContent(i, y * 1.0e+30); // ub
76  }
77 
78  for (const auto& setup : setups) {
79 
80  const double E = setup.first;
81 
82  TH1D* ha = HA[E];
83  TH1D* h1 = H1[E];
84  TH1D* hb = HB[E];
85 
86  for (int i = 1; i <= ha->GetNbinsX(); ++i) {
87 
88  const double x = ha->GetBinCenter(i);
89  const double v = pow(10.0, x);
90 
91  double y = dis.getP(E, v);
92 
93  if (numberOfEvents == 0) {
94  y *= setup.second / dis.getP(E, dis.getV(E));
95  }
96 
97  ha->SetBinContent(i, y);
98  }
99 
100  if (numberOfEvents != 0) {
101 
102  for (int i = 0; i != numberOfEvents; ++i) {
103 
104  const double Es = dis.getE(E);
105  const double v = Es/E;
106 
107  h1->Fill(log10(v));
108  }
109 
110  for (int i = 1; i <= h1->GetNbinsX(); ++i) {
111 
112  const double y = h1->GetBinContent(i);
113  const double z = h1->GetBinError (i);
114 
115  const double xmin = h1->GetXaxis()->GetBinLowEdge(i);
116  const double xmax = h1->GetXaxis()->GetBinUpEdge (i);
117 
118  const double w = (pow(10.0,xmax) - pow(10.0,xmin)) * numberOfEvents;
119 
120  h1->SetBinContent(i, y / w);
121  h1->SetBinError (i, z / w);
122  }
123  }
124 
125  for (int i = 1; i <= hb->GetNbinsX(); ++i) {
126 
127  const double x = hb->GetBinCenter(i);
128  const double v = pow(10.0, x);
129  const double y = dis.getThetaRMS(E, v);
130 
131  hb->SetBinContent(i, y);
132  }
133  }
134 
135 
136  TFile out(outputFile.c_str(), "recreate");
137 
138  out << h0 << HA << H1 << HB;
139 
140  out.Write();
141  out.Close();
142 }
string outputFile
#define FATAL(A)
Definition: JMessage.hh:67
int debug
debug level
Definition: JSirene.cc:69
#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
Deep-inelastic muon-nucleon scattering.
Definition: JDIS.hh:31
double getCrossSection(const double E) const
Get cross section.
Definition: JDIS.hh:46
double getP(const double E, const double v) const
Get probability of given energy fraction.
Definition: JDIS.hh:64
double getThetaRMS(const double E, const double v) const
Get RMS of scattering angle.
Definition: JDIS.hh:93
double getV(const double E) const
Get breakpoint.
Definition: JDIS.hh:110
double getE(const double E) const
Get shower energy.
Definition: JDIS.hh:78
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys.
Definition: JManager.hh:47
const double xmax
Definition: JQuadrature.cc:24
const double xmin
Definition: JQuadrature.cc:23
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).
data_type w[N+1][M+1]
Definition: JPolint.hh:867
data_type v[N+1][M+1]
Definition: JPolint.hh:866
Definition: JSTDTypes.hh:14