Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Analysis.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
27 
28 #include "Analysis.hh"
29 #include "G4SystemOfUnits.hh"
30 
31 // ROOT headers
32 #include "TROOT.h"
33 #include "TFile.h"
34 #include "TH1.h"
35 #include "TH2.h"
36 
37 Analysis* Analysis::myanalysis = NULL;
38 
39 // --------------------------------------------------------------------------
40 Analysis::Analysis()
41 {
42  // ROOT style
43  gROOT-> Reset();
44 
45  // define histograms
46  incident_map = new TH2D("incident map", "Incident Distributuon",
47  50, -5., 5.,
48  50, -5., 5.);
49  incident_map-> GetXaxis()-> SetTitle("X (cm)");
50  incident_map-> GetYaxis()-> SetTitle("Y (cm)");
51  incident_map-> SetStats(0);
52 
53  incident_x_hist = new TH1D("incident x", "Incident X", 100, -5., 5.);
54  incident_x_hist-> GetXaxis()-> SetTitle("X (cm)");
55  incident_x_hist-> SetFillColor(kRed);
56 
57  dose_map = new TH2D("dose map", "Dose Distribution",
58  500, 0., 50.,
59  200, -10., 10.);
60  dose_map-> GetXaxis()-> SetTitle("Z (cm)");
61  dose_map-> GetYaxis()-> SetTitle("X (cm)");
62  dose_map-> SetStats(0);
63 
64  dose_hist = new TH1D("dose", "Dose Distribution", 500, 0., 50.);
65  dose_hist-> GetXaxis()-> SetTitle("Z (cm)");
66  dose_hist-> GetYaxis()-> SetTitle("Dose (GeV)");
67  dose_hist-> SetFillColor(kBlue);
68  dose_hist-> SetStats(0);
69 
70 }
71 
72 // --------------------------------------------------------------------------
74 {
75  delete incident_map;
76  delete incident_x_hist;
77  delete dose_map;
78  delete dose_hist;
79 
80  myanalysis = NULL;
81 }
82 
83 // --------------------------------------------------------------------------
85 {
86  if ( myanalysis == NULL ) {
87  myanalysis = new Analysis();
88  }
89 
90  return myanalysis;
91 }
92 
93 // --------------------------------------------------------------------------
95 {
96  return;
97 }
98 
99 // --------------------------------------------------------------------------
101 {
102  incident_map-> Reset();
103  incident_x_hist-> Reset();
104  dose_map-> Reset();
105  dose_hist-> Reset();
106 
107  return;
108 }
109 
110 // --------------------------------------------------------------------------
112 {
113  TFile* file = new TFile(fname.c_str(),
114  "RECREATE", "Geant4 ROOT analysis");
115  incident_map-> Write();
116  incident_x_hist-> Write();
117  dose_map-> Write();
118  dose_hist-> Write();
119 
120  file-> Close();
121 
122  delete file;
123 
124  return;
125 }
126 
127 // --------------------------------------------------------------------------
129 {
130  if ( ! incidentFlag ) {
131  incident_map-> Fill(p.x()/cm, p.y()/cm);
132  incident_x_hist-> Fill(p.x()/cm);
133 
134  incidentFlag = true;
135  }
136 }
137 
138 // --------------------------------------------------------------------------
140 {
141  const G4double Z0 = 25.*cm;
142  const G4double dxy = 10.*mm;
143 
144  if(std::abs(p.y()) < dxy ) {
145  dose_map-> Fill((p.z()+Z0)/cm, p.x()/cm, dedx/GeV);
146 
147  if(std::abs(p.x()) < dxy) {
148  dose_hist-> Fill((p.z()+Z0)/cm, dedx/GeV);
149  }
150  }
151 }
152