Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PhysChemIO.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/*
27 * G4PhysChemIO.cc
28 *
29 * Created on: 3 févr. 2017
30 * Author: matkara
31 */
32
33#include "G4PhysChemIO.hh"
34#include "G4SystemOfUnits.hh"
35#include "G4Track.hh"
36#include "G4VAnalysisManager.hh"
37
38using namespace std;
39
40//------------------------------------------------------------------------------
41
42namespace G4PhysChemIO{
43
45 fRunID = -1;
46 fEventID = -1;
47 fFileInitialized = false;
48}
49
50//------------------------------------------------------------------------------
51
53 CloseFile();
54}
55
56//------------------------------------------------------------------------------
57
59{
60 if(fFileInitialized) return;
61
62 fOfstream << std::setprecision(6) << std::scientific;
63 fOfstream << setw(11) << left << "#Parent ID" << setw(10) << "Molecule"
64 << setw(14) << "Elec Modif" << setw(13) << "Energy (eV)"
65 << setw(22) << "X pos of parent [nm]" << setw(22)
66 << "Y pos of parent [nm]" << setw(22) << "Z pos of parent [nm]"
67 << setw(14) << "X pos [nm]" << setw(14) << "Y pos [nm]"
68 << setw(14) << "Z pos [nm]" << G4endl<< setw(21) << "#"
69 << setw(13) << "1)io/ex=0/1"
70 << G4endl
71 << setw(21) << "#"
72 << setw(13) << "2)level=0...5"
73 << G4endl;
74
75 fFileInitialized = true;
76}
77
78//------------------------------------------------------------------------------
79
81 ios_base::openmode mode)
82{
83 fOfstream.open(output.data(), mode);
84 fFileInitialized = false;
85}
86
87//------------------------------------------------------------------------------
88
90{
92}
93
94//------------------------------------------------------------------------------
95
97{
98 if (fFileInitialized == false) return;
99
100 if (fOfstream.is_open())
101 {
102 fOfstream.close();
103 }
104}
105
106//------------------------------------------------------------------------------
107
109 G4int electronicLevel,
110 G4double energy,
111 const G4Track* theIncomingTrack)
112{
114
115 fOfstream << setw(11) << left << theIncomingTrack->GetTrackID()
116 << setw(10) << "H2O" << left << modification << internal
117 << ":" << right << electronicLevel << left << setw(11) << ""
118 << std::setprecision(2) << std::fixed << setw(13)
119 << energy / eV << std::setprecision(6) << std::scientific
120 << setw(22)
121 << (theIncomingTrack->GetPosition().x()) / nanometer
122 << setw(22)
123 << (theIncomingTrack->GetPosition().y()) / nanometer
124 << setw(22)
125 << (theIncomingTrack->GetPosition().z()) / nanometer
126 << G4endl;
127}
128
129//------------------------------------------------------------------------------
130
132 G4ThreeVector* finalPosition)
133{
135
136 fOfstream << setw(11) << theIncomingTrack->GetTrackID() << setw(10)
137 << "e_aq" << setw(14) << -1 << std::setprecision(2)
138 << std::fixed << setw(13)
139 << theIncomingTrack->GetKineticEnergy() / eV
140 << std::setprecision(6) << std::scientific << setw(22)
141 << (theIncomingTrack->GetPosition().x()) / nanometer
142 << setw(22)
143 << (theIncomingTrack->GetPosition().y()) / nanometer
144 << setw(22)
145 << (theIncomingTrack->GetPosition().z()) / nanometer;
146
147 if (finalPosition != 0)
148 {
149 fOfstream << setw(14) << (finalPosition->x()) / nanometer << setw(14)
150 << (finalPosition->y()) / nanometer << setw(14)
151 << (finalPosition->z()) / nanometer;
152 }
153
154 fOfstream << G4endl;
155}
156
157//------------------------------------------------------------------------------
158//
159// Using G4analysis
160//
161
163fpAnalysisManager(analysisManager)
164{
165 fFileInitialized = false;
166 fNtupleID = -1;
167}
168
169//------------------------------------------------------------------------------
170
172{
174}
175
176//------------------------------------------------------------------------------
177
179{
180 if (fFileInitialized) return;
181
182 fNtupleID = fpAnalysisManager->CreateNtuple("PhysChem","PhysChem");
185
186 //----------------------------------------------------------------------------
187 // valid for H2O only
189 // ionization = 0 / excitation = 1 / diss att = 2
191 // valid for ion and exc only
193 // valid for ion and exc only
194
195 //----------------------------------------------------------------------------
203
204 fFileInitialized = true;
205}
206
207//------------------------------------------------------------------------------
208
210 ios_base::openmode)
211{
213 fFileInitialized = false;
214}
215
216//------------------------------------------------------------------------------
217
219{
220// fpAnalysisManager->Write();
221// fpAnalysisManager->CloseFile();
222}
223
224//------------------------------------------------------------------------------
225
227 G4int electronicLevel,
228 G4double energy,
229 const G4Track* theIncomingTrack)
230{
232
233 // parent ID
235 theIncomingTrack->GetTrackID());
236
237 // molecule type
239
240 //----------------------------------------------------------------------------
241 // valid for H2O only
242
243 // electronic modif
245 // ionization = 0 / excitation = 1 / diss att = 2
246 fpAnalysisManager->FillNtupleIColumn(fNtupleID, 3, electronicLevel);
248
249 //----------------------------------------------------------------------------
250 const G4ThreeVector& parentPos = theIncomingTrack->GetPosition();
251
252 fpAnalysisManager->FillNtupleDColumn(fNtupleID,5,(parentPos.x())/nanometer);
253 fpAnalysisManager->FillNtupleDColumn(fNtupleID,6,(parentPos.y())/nanometer);
254 fpAnalysisManager->FillNtupleDColumn(fNtupleID,7,(parentPos.z())/nanometer);
255
256 fpAnalysisManager->FillNtupleDColumn(fNtupleID,8,(parentPos.x())/nanometer);
257 fpAnalysisManager->FillNtupleDColumn(fNtupleID,9,(parentPos.y())/nanometer);
258 fpAnalysisManager->FillNtupleDColumn(fNtupleID,10,(parentPos.z())/nanometer);
260}
261
262//------------------------------------------------------------------------------
263
265 G4ThreeVector* finalPosition)
266{
268
269 // parent ID
271 electronTrack->GetTrackID());
272
273 // molecule type
275
276 //----------------------------------------------------------------------------
277 // valid for H2O only
278
279 // electronic modif
280 fpAnalysisManager->FillNtupleIColumn(fNtupleID, 2, -1); // electronic modif
281 fpAnalysisManager->FillNtupleIColumn(fNtupleID, 3, -1); // electronic level
283 electronTrack->GetKineticEnergy() / eV);
284
285 //----------------------------------------------------------------------------
286 const G4ThreeVector& parentPos = electronTrack->GetPosition();
287 const double i_nm = 1./nanometer;
288
289 fpAnalysisManager->FillNtupleDColumn(fNtupleID,5, parentPos.x() *i_nm);
290 fpAnalysisManager->FillNtupleDColumn(fNtupleID,6, parentPos.y() *i_nm);
291 fpAnalysisManager->FillNtupleDColumn(fNtupleID,7, parentPos.z() *i_nm);
292
293 if (finalPosition != 0)
294 {
295 fpAnalysisManager->FillNtupleDColumn(fNtupleID,8, finalPosition->x()*i_nm);
296 fpAnalysisManager->FillNtupleDColumn(fNtupleID,9, finalPosition->y()*i_nm);
297 fpAnalysisManager->FillNtupleDColumn(fNtupleID,10, finalPosition->z()*i_nm);
298 }
299 else
300 {
301 fpAnalysisManager->FillNtupleDColumn(fNtupleID,8, parentPos.x() *i_nm);
302 fpAnalysisManager->FillNtupleDColumn(fNtupleID,9, parentPos.y() *i_nm);
303 fpAnalysisManager->FillNtupleDColumn(fNtupleID,10, parentPos.z() *i_nm);
304 }
305
307}
308}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
double z() const
double x() const
double y() const
virtual void CreateWaterMolecule(G4int electronicModif, G4int, G4double energy, const G4Track *)
virtual void CreateSolvatedElectron(const G4Track *, G4ThreeVector *finalPosition=0)
virtual void InitializeFile()
Definition: G4PhysChemIO.cc:58
virtual void AddEmptyLineInOutputFile()
Definition: G4PhysChemIO.cc:89
virtual void WriteInto(const G4String &, std::ios_base::openmode mode=std::ios_base::out)
Definition: G4PhysChemIO.cc:80
virtual void CreateSolvatedElectron(const G4Track *, G4ThreeVector *finalPosition=0)
virtual void CreateWaterMolecule(G4int electronicModif, G4int, G4double energy, const G4Track *)
G4Analysis(G4VAnalysisManager *)
virtual void InitializeFile()
virtual void CloseFile()
virtual void WriteInto(const G4String &, std::ios_base::openmode mode=std::ios_base::out)
G4VAnalysisManager * fpAnalysisManager
const char * data() const
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
G4double GetKineticEnergy() const
G4bool FillNtupleSColumn(G4int id, const G4String &value)
G4int CreateNtupleIColumn(const G4String &name)
G4int CreateNtupleDColumn(const G4String &name)
G4bool FillNtupleDColumn(G4int id, G4double value)
G4bool FillNtupleIColumn(G4int id, G4int value)
G4bool OpenFile(const G4String &fileName="")
G4int CreateNtupleSColumn(const G4String &name)
G4int CreateNtuple(const G4String &name, const G4String &title)