Geant4_10
RunAction.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 // ********************************************************************
25 //
26 // Please cite the following paper if you use this software
27 // Nucl.Instrum.Meth.B260:20-27, 2007
28 
29 // #define MATRIX_BOUND_CHECK
30 
31 #include "RunAction.hh"
32 #include "Analysis.hh"
33 
34 using namespace std;
35 
36 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
37 
39 :fDetector(det),fPrimary(pri)
40 {;}
41 
42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
43 
45 {}
46 
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48 
49 void RunAction::BeginOfRunAction(const G4Run* /*aRun*/)
50 {
51  // Vector initialization
52 
53  fRow=0;
54 
55  fXVector = CLHEP::HepVector(32);
56  fYVector = CLHEP::HepVector(32);
57  fThetaVector = CLHEP::HepVector(32);
58  fPhiVector = CLHEP::HepVector(32);
59 
60  // Histograms
61 
62  // Get/create analysis manager
63  G4cout << "##### Create analysis manager " << " " << this << G4endl;
64 
65  G4AnalysisManager* man = G4AnalysisManager::Instance();
66 
67  G4cout << "Using " << man->GetType() << " analysis manager" << G4endl;
68 
69 
70  // Open an output file
71  man->OpenFile("nanobeam");
72  man->SetFirstHistoId(1);
73  man->SetFirstNtupleId(1);
74 
75  // Create 1st ntuple (id = 1)
76  //
77  man->CreateNtuple("ntuple0", "BeamProfile");
78  man->CreateNtupleDColumn("xIn");
79  man->CreateNtupleDColumn("yIn");
80  man->CreateNtupleDColumn("zIn");
81  man->FinishNtuple();
82  G4cout << "Ntuple-1 created" << G4endl;
83 
84  // Create 2nd htuple (id = 2)
85  //
86  man->CreateNtuple("ntuple1","Grid");
87  man->CreateNtupleDColumn("xIn");
88  man->CreateNtupleDColumn("yIn");
89  man->CreateNtupleDColumn("e");
90  man->FinishNtuple();
91  G4cout << "Ntuple-2 created" << G4endl;
92 
93  // Create 3rd ntuple (id = 3)
94  man->CreateNtuple("ntuple2","Coef");
95  man->CreateNtupleDColumn("xIn");
96  man->CreateNtupleDColumn("yIn");
97  man->CreateNtupleDColumn("thetaIn");
98  man->CreateNtupleDColumn("phiIn");
99  man->FinishNtuple();
100  G4cout << "Ntuple-3 created" << G4endl;
101 
102  return;
103 
104 }
105 
106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107 
108 void RunAction::EndOfRunAction(const G4Run* /*aRun*/)
109 {
110 
111 if (fDetector->GetCoef()==1)
112 {
113  CLHEP::HepMatrix m;
114 
115  // VECTOR READING
116  // VECTOR READING
117 
118  m = CLHEP::HepMatrix(32,32);
119  m = fPrimary->GetMatrix();
120 
121  G4cout << G4endl;
122  G4cout << "===> NANOBEAM LINE INTRINSIC ABERRATION COEFFICIENTS (units of micrometer and mrad) :" << G4endl;
123  G4cout << G4endl;
124 
125  int inv;
126 
127  m.invert(inv);
128  CLHEP::HepVector tmp(32,0);
129  tmp=m*fXVector;
130  CLHEP::HepVector b;
131  b=tmp.sub(2,2); G4cout << "<x|theta>=" << b << G4endl;
132  b=tmp.sub(8,8); G4cout << "<x|theta*delta>=" << b << G4endl;
133  b=tmp.sub(10,10); G4cout << "<x|theta^3>=" << b << G4endl;
134  b=tmp.sub(12,12); G4cout << "<x|theta*phi^2>=" << b << G4endl;
135  m.invert(inv);
136 
137  m.invert(inv);
138  tmp = m*fThetaVector;
139  m.invert(inv);
140  b=tmp.sub(2,2); G4cout << "<x|x>=" << b << G4endl;
141 
142  m.invert(inv);
143  tmp=m*fYVector;
144  b=tmp.sub(3,3); G4cout << "<y|phi>=" << b << G4endl;
145  b=tmp.sub(9,9); G4cout << "<y|phi*delta>=" << b << G4endl;
146  b=tmp.sub(11,11); G4cout << "<y|theta^2*phi>=" << b << G4endl;
147  b=tmp.sub(13,13); G4cout << "<y|phi^3>=" << b << G4endl;
148  m.invert(inv);
149 
150  m.invert(inv);
151  tmp = m*fPhiVector;
152  m.invert(inv);
153  b=tmp.sub(3,3); G4cout << "<y|y>=" << b << G4endl;
154 
155 }
156 
157  // Save histograms
158 
159  G4AnalysisManager* man = G4AnalysisManager::Instance();
160  man->Write();
161  man->CloseFile();
162 
163  // Complete clean-up
164 
165  delete G4AnalysisManager::Instance();
166 
167 }
Float_t tmp
Definition: plot.C:37
void BeginOfRunAction(const G4Run *)
Definition: RunAction.cc:57
tuple b
Definition: test.py:12
G4GLOB_DLL std::ostream G4cout
void EndOfRunAction(const G4Run *)
Definition: RunAction.cc:260
Definition: G4Run.hh:46
ExG4HbookAnalysisManager G4AnalysisManager
Definition: g4hbook_defs.hh:46
~RunAction()
Definition: RunAction.cc:52
#define G4endl
Definition: G4ios.hh:61