Jpp
Functions
JGeanz.cc File Reference
#include <string>
#include <iostream>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TGraph.h"
#include "TRandom3.h"
#include "JPhysics/JGeanz.hh"
#include "JGizmo/JGizmoToolkit.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 longitudinal shower profile using JPHYSICS::JGeanz.

Author
mdejong

Definition in file JGeanz.cc.

Function Documentation

◆ main()

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

Definition at line 23 of file JGeanz.cc.

24 {
25  using namespace std;
26 
27  string outputFile;
28  double E_GeV;
29  int numberOfEvents;
30  bool debug;
31 
32  try {
33 
34  JParser<> zap("Example program to histogram longitudinal shower profile.");
35 
36  zap['o'] = make_field(outputFile) = "geanz.root";
37  zap['E'] = make_field(E_GeV);
38  zap['n'] = make_field(numberOfEvents) = 0;
39  zap['d'] = make_field(debug);
40 
41  zap(argc, argv);
42  }
43  catch(const exception &error) {
44  FATAL(error.what() << endl);
45  }
46 
47  using namespace JPP;
48 
49 
50  TFile out(outputFile.c_str(), "recreate");
51 
52  TH1D h1("Probability", NULL, 1000, 0.0, 25.0);
53  TH1D h2("Integral", NULL, 1000, 0.0, 25.0);
54  TH1D h3("MC", NULL, 100, 0.0, 25.0);
55 
56  for(int i = 1; i <= h1.GetNbinsX(); ++i) {
57 
58  const double x = h1.GetBinCenter(i);
59 
60  h1.Fill(x, geanz.getProbability(E_GeV, x));
61  h2.Fill(x, geanz.getIntegral (E_GeV, x));
62  }
63 
64  if (numberOfEvents != 0) {
65 
66  for (int i = 0; i != numberOfEvents; ++i) {
67 
68  const double x = gRandom->Rndm();
69 
70  h3.Fill(geanz.getLength(E_GeV, x));
71  }
72 
73  convertToPDF(h3, "NW");
74  }
75 
76  const Double_t x = geanz.getMaximum(E_GeV);
77  const Double_t y = geanz.getProbability(E_GeV, x);
78 
79  TGraph g1(1, &x, &y);
80 
81  g1.SetName("Maximum");
82  g1.Write();
83 
84  out.Write();
85  out.Close();
86 }
g1
Double_t g1(const Double_t x)
Function.
Definition: JQuantiles.cc:25
JGIZMO::convertToPDF
void convertToPDF(TH1 &h1, const std::string &option="NW", const double factor=1.0)
Convert 1D histogram to PDF.
Definition: JGizmoToolkit.hh:366
JPARSER::JParser
Utility class to parse command line options.
Definition: JParser.hh:1493
JPHYSICS::JGeanz::getIntegral
double getIntegral(const double E, const double z) const
Integral of PDF (starting from 0).
Definition: JGeanz.hh:95
JPHYSICS::JGeanz::getProbability
double getProbability(const double E, const double z) const
Probability Density Function.
Definition: JGeanz.hh:56
JPP
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Definition: JAAnetToolkit.hh:37
JPHYSICS::JGeanz::getMaximum
double getMaximum(const double E) const
Get depth of shower maximum.
Definition: JGeanz.hh:162
debug
int debug
debug level
Definition: JSirene.cc:59
make_field
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1954
JPHYSICS::JGeanz::getLength
double getLength(const double E, const double P, const double eps=1.0e-3) const
Get shower length for a given integrated probability.
Definition: JGeanz.hh:122
std
Definition: jaanetDictionary.h:36
FATAL
#define FATAL(A)
Definition: JMessage.hh:67
outputFile
string outputFile
Definition: JDAQTimesliceSelector.cc:37
JPHYSICS::geanz
static const JGeanz geanz(1.85, 0.62, 0.54)
Function object for longitudinal EM-shower profile.