Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAMultipleIonisationManager Class Reference

#include <G4DNAMultipleIonisationManager.hh>

Public Member Functions

 G4DNAMultipleIonisationManager ()=default
 
 ~G4DNAMultipleIonisationManager ()=default
 
void CreateMultipleIonisedWaterMolecule (MultipleIonisedModification mod, G4int *shell_level, const G4Track *incoming_track)
 
G4bool CheckShellEnergy (MultipleIonisedModification mod, G4double *shell_energy)
 
void LoadAlphaParam (const G4String &filepath, G4double Z, G4double A)
 
G4double GetAlphaParam (G4double energy)
 

Detailed Description

Definition at line 46 of file G4DNAMultipleIonisationManager.hh.

Constructor & Destructor Documentation

◆ G4DNAMultipleIonisationManager()

G4DNAMultipleIonisationManager::G4DNAMultipleIonisationManager ( )
default

◆ ~G4DNAMultipleIonisationManager()

G4DNAMultipleIonisationManager::~G4DNAMultipleIonisationManager ( )
default

Member Function Documentation

◆ CheckShellEnergy()

G4bool G4DNAMultipleIonisationManager::CheckShellEnergy ( MultipleIonisedModification mod,
G4double * shell_energy )

Definition at line 78 of file G4DNAMultipleIonisationManager.cc.

80{
81 G4int num_shells{0};
82 switch (mod) {
84 num_shells = 2;
85 break;
87 num_shells = 3;
88 break;
90 num_shells = 4;
91 break;
92 default: // never happen
93 break;
94 }
95
96 G4bool stop_process{false};
97 for (G4int i = 0; i < num_shells; i++) {
98 if (shell_energy[i] < 0.0) {
99 stop_process = true;
100 break;
101 }
102 }
103
104 return stop_process;
105}
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85

◆ CreateMultipleIonisedWaterMolecule()

void G4DNAMultipleIonisationManager::CreateMultipleIonisedWaterMolecule ( MultipleIonisedModification mod,
G4int * shell_level,
const G4Track * incoming_track )

Definition at line 42 of file G4DNAMultipleIonisationManager.cc.

45{
46 if (!G4DNAChemistryManager::IsActivated()) { return; }
47
48 G4int num_shells{0};
49 switch (mod) {
51 num_shells = 2;
52 break;
54 num_shells = 3;
55 break;
57 num_shells = 4;
58 break;
59 default: // never happen
60 return;
61 }
62
63 auto* H2O = new G4Molecule(G4H2O::Definition());
64 for (G4int i = 0; i < num_shells; i++) {
65 H2O->IonizeMolecule(4 - shell_level[i]);
66 }
67
68 constexpr G4double kT0 = 1.0 * picosecond;
69 auto* H2O_track = H2O->BuildTrack(kT0, incoming_track->GetPosition());
70 H2O_track->SetParentID(incoming_track->GetTrackID());
71 H2O_track->SetTrackStatus(fStopButAlive);
72 H2O_track->SetKineticEnergy(0.0);
73
74 G4VITTrackHolder::Instance()->Push(H2O_track);
75}
@ fStopButAlive
double G4double
Definition G4Types.hh:83
static G4H2O * Definition()
Definition G4H2O.cc:42
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
static G4VITTrackHolder * Instance()
virtual void Push(G4Track *)

◆ GetAlphaParam()

G4double G4DNAMultipleIonisationManager::GetAlphaParam ( G4double energy)

Definition at line 137 of file G4DNAMultipleIonisationManager.cc.

138{
139 auto find_lower_bound = [this](G4double e) {
140 auto low = 0;
141 auto upp = num_node_ - 1;
142 if (e < Eion_[0]) { return low; }
143 while (low <= upp) {
144 const auto mid = static_cast<int>((low + upp) * 0.5);
145 if (e < Eion_[mid]) { upp = mid - 1; }
146 else { low = mid + 1; }
147 if (upp < 0) { upp = 0; }
148 }
149 return upp;
150 };
151
152 auto interp_log_log = [this](G4int bin1, G4double e) {
153 if (e < Eion_[0]) { return alpha_[0]; }
154 const auto num_bin = num_node_ - 1;
155 const auto bin2 = bin1 + 1;
156 G4double value{0.0};
157 if (bin2 <= num_bin) {
158 auto log10_e = std::log10(e);
159 auto log10_e1 = Eion_[bin1];
160 auto log10_e2 = Eion_[bin2];
161 auto log10_a1 = alpha_[bin1];
162 auto log10_a2 = alpha_[bin2];
163 if (log10_a1 != 0.0 && log10_a2 != 0.0) {
164 log10_e1 = std::log10(log10_e1);
165 log10_e2 = std::log10(log10_e2);
166 log10_a1 = std::log10(log10_a1);
167 log10_a2 = std::log10(log10_a2);
168 value = log10_a1 + (log10_a2 - log10_a1)
169 * (log10_e - log10_e1) / (log10_e2 - log10_e1);
170 value = std::pow(10.0, value);
171 }
172 } else {
173 value = alpha_[num_bin];
174 }
175 return value;
176 };
177
178 const auto bin1 = find_lower_bound(energy);
179 return interp_log_log(bin1, energy);
180}

◆ LoadAlphaParam()

void G4DNAMultipleIonisationManager::LoadAlphaParam ( const G4String & filepath,
G4double Z,
G4double A )

Definition at line 108 of file G4DNAMultipleIonisationManager.cc.

110{
111 const char* path = G4FindDataDir("G4LEDATA");
112 if (path == nullptr) {
113 G4Exception("G4DNAMultipleIonisationManager::LoadAlphaParam","em0006",
114 FatalException,"G4LEDATA environment variable not set.");
115 }
116
117 std::stringstream fullpath;
118 fullpath << path << "/" << filepath;
119
120 std::fstream fin(fullpath.str());
121
122 G4double e, a;
123 std::string line = "";
124 while (getline(fin, line)) {
125 std::stringstream ss;
126 ss << line;
127 ss >> e >> a;
128 Eion_.push_back(e * Z * A * MeV);
129 alpha_.push_back(a);
130 }
131
132 num_node_ = (G4int)Eion_.size();
133 fin.close();
134}
const char * G4FindDataDir(const char *)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
const G4double A[17]

The documentation for this class was generated from the following files: