Jpp  18.3.0
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JGinneken.cc
Go to the documentation of this file.
1 #include <string>
2 #include <iostream>
3 #include <iomanip>
4 
5 #include "TROOT.h"
6 #include "TFile.h"
7 #include "TH1D.h"
8 
9 #include "JPhysics/JRadiation.hh"
10 
11 #include "JROOT/JManager.hh"
12 #include "JROOT/JRootToolkit.hh"
13 
14 #include "Jeep/JParser.hh"
15 #include "Jeep/JMessage.hh"
16 
17 /**
18  * Example application to display theta RMS of muon energy loss.
19  */
20 int main(int argc, char **argv)
21 {
22  using namespace std;
23  using namespace JPP;
24 
25  string outputFile;
26  int debug;
27 
28  try {
29 
30  JParser<> zap;
31 
32  zap['o'] = make_field(outputFile) = "ginneken.root";
33  zap['d'] = make_field(debug) = 0;
34 
35  zap(argc, argv);
36  }
37  catch(const exception &error) {
38  FATAL(error.what() << endl);
39  }
40 
41 
42  const JRadiation radiation(4.0, 8.0, 40, 0.01, 0.1, 0.1); // Be
43 
44  JManager<double, TH1D> HA(new TH1D("Brems1 [% GeV]", NULL, 10000, -7.5, log10(0.5)));
45  JManager<double, TH1D> HB(new TH1D("Brems2 [% GeV]", NULL, 10000, -3.2, log10(0.5)));
46  JManager<double, TH1D> HC(new TH1D("EErad [% GeV]", NULL, 10000, -7.8, 0.0));
47 
48  for (const double E : { 1.0e2, 1.0e3 }) {
49 
50  TH1D* ha = HA[E];
51  TH1D* hb = HB[E];
52  TH1D* hc = HC[E];
53 
54  for (int i = 1; i <= ha->GetNbinsX(); ++i) {
55 
56  const double x = ha->GetBinCenter(i);
57  const double v = pow(10.0, x);
58  const double y = radiation.ThetaRMSfromBrems(E, v);
59 
60  ha->SetBinContent(i, y);
61  }
62 
63  for (int i = 1; i <= hb->GetNbinsX(); ++i) {
64 
65  const double x = hb->GetBinCenter(i);
66  const double v = 1.0 - pow(10.0, x);
67  const double y = radiation.ThetaRMSfromBrems(E, v);
68 
69  hb->SetBinContent(i, y);
70  }
71 
72  for (int i = 1; i <= hc->GetNbinsX(); ++i) {
73 
74  const double x = hc->GetBinCenter(i);
75  const double v = pow(10.0, x);
76  const double y = radiation.ThetaRMSfromEErad(E, v);
77 
78  hc->SetBinContent(i, y);
79  }
80  }
81 
82  TFile out(outputFile.c_str(), "recreate");
83 
84  out << HA << HB << HC;
85 
86  out.Write();
87  out.Close();
88 }
Utility class to parse command line options.
Definition: JParser.hh:1514
int main(int argc, char *argv[])
Definition: Main.cc:15
then usage $script< input file >[option[primary[working directory]]] nWhere option can be E
Definition: JMuonPostfit.sh:40
Dynamic ROOT object management.
Muon radiative cross sections.
string outputFile
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1989
set_variable E_E log10(E_{fit}/E_{#mu})"
T pow(const T &x, const double y)
Power .
Definition: JMath.hh:97
General purpose messaging.
#define FATAL(A)
Definition: JMessage.hh:67
Utility class to parse command line options.
data_type v[N+1][M+1]
Definition: JPolint.hh:866
int debug
debug level