Loading...
Searching...
No Matches
Public Member Functions | Private Member Functions | Private Attributes | List of all members
RadioBio::LET Class Reference

#include <Doxymodules_medical.h>

Inheritance diagram for RadioBio::LET:
RadioBio::VRadiobiologicalQuantity

Public Member Functions

 LET ()
 
 ~LET () override
 
void AddFromAccumulable (G4VAccumulable *) override
 
void Initialize () override
 
void Compute () override
 
void Reset () override
 
void Store () override
 
void PrintParameters () override
 
- Public Member Functions inherited from RadioBio::VRadiobiologicalQuantity
 VRadiobiologicalQuantity ()
 
virtual ~VRadiobiologicalQuantity ()
 
G4bool IsCalculationEnabled () const
 
G4bool HasBeenCalculated () const
 
void SetCalculationEnabled (G4bool enabled)
 
void SetVerboseLevel (G4int level)
 
G4int GetVerboseLevel () const
 
void SetPath (G4String fN)
 
void PrintVirtualParameters ()
 

Private Member Functions

void SetNTotalLETT (const array_type NTotalLETT)
 
void SetNTotalLETD (const array_type NTotalLETD)
 
void SetDTotalLETT (const array_type DTotalLETT)
 
void SetDTotalLETD (const array_type DTotalLETD)
 
void AddNTotalLETT (const array_type NTotalLETT)
 
void AddNTotalLETD (const array_type NTotalLETD)
 
void AddDTotalLETT (const array_type DTotalLETT)
 
void AddDTotalLETD (const array_type DTotalLETD)
 
void AddIon (const IonLet ion)
 

Private Attributes

array_type fNTotalLETT = {}
 
array_type fDTotalLETT = {}
 
array_type fNTotalLETD = {}
 
array_type fDTotalLETD = {}
 
array_type fTotalLETT = {}
 
array_type fTotalLETD = {}
 
std::vector< IonLetfIonLetStore
 

Additional Inherited Members

- Public Types inherited from RadioBio::VRadiobiologicalQuantity
using array_type = std::valarray< G4double >
 
- Protected Attributes inherited from RadioBio::VRadiobiologicalQuantity
G4UImessengerfMessenger = nullptr
 
G4int fVerboseLevel = 1
 
G4bool fInitialized = false
 
G4bool fSaved = false
 
G4bool fCalculated = false
 
G4String fPath
 
G4bool fCalculationEnabled = false
 

Detailed Description

Definition at line 134 of file Doxymodules_medical.h.

Constructor & Destructor Documentation

◆ LET()

RadioBio::LET::LET ( )

Definition at line 48 of file LET.cc.

49{
50 fPath = "RadioBio_LET.out"; // Default output filename
51
52 fMessenger = new LETMessenger(this);
53
54 Initialize();
55}
void Initialize() override
Definition LET.cc:66

◆ ~LET()

RadioBio::LET::~LET ( )
override

Definition at line 59 of file LET.cc.

60{
61 delete fMessenger;
62}

Member Function Documentation

◆ AddFromAccumulable()

void RadioBio::LET::AddFromAccumulable ( G4VAccumulable GenAcc)
overridevirtual

Implements RadioBio::VRadiobiologicalQuantity.

Definition at line 253 of file LET.cc.

254{
255 LETAccumulable* acc = (LETAccumulable*)GenAcc;
256
257 AddNTotalLETT(acc->GetTotalLETT());
258 AddDTotalLETT(acc->GetDTotalLETT());
259 AddNTotalLETD(acc->GetTotalLETD());
260 AddDTotalLETD(acc->GetDTotalLETD());
261
262 // Merges ion counters
263 // Loop over rhs ions
264 for (unsigned int l = 0; l < acc->GetIonLetStore().size(); l++) {
265 G4int PDGencoding = acc->GetIonLetStore()[l].GetPDGencoding();
266 size_t q;
267 // Loop over lhs ions to find the one
268 for (q = 0; q < fIonLetStore.size(); q++) {
269 // If he found it, sums values
270 if (fIonLetStore[q].GetPDGencoding() == PDGencoding) {
271 if (acc->GetIonLetStore()[l].IsPrimary() == fIonLetStore[q].IsPrimary()) break;
272 }
273 }
274 // If ion is missing, copy it
275 if (q == fIonLetStore.size())
276 fIonLetStore.push_back(acc->GetIonLetStore()[l]);
277 else // Merge rhs data with lhs ones
278 fIonLetStore[q].Merge(&(acc->GetIonLetStore()[l]));
279 }
280 fCalculated = false;
281}
std::vector< IonLet > fIonLetStore
Definition LET.hh:80
void AddNTotalLETT(const array_type NTotalLETT)
Definition LET.hh:102
void AddNTotalLETD(const array_type NTotalLETD)
Definition LET.hh:106
void AddDTotalLETD(const array_type DTotalLETD)
Definition LET.hh:114
void AddDTotalLETT(const array_type DTotalLETT)
Definition LET.hh:110

◆ Initialize()

void RadioBio::LET::Initialize ( )
overridevirtual

Implements RadioBio::VRadiobiologicalQuantity.

Definition at line 66 of file LET.cc.

67{
69
70 // Instances for Total LET
71 fNTotalLETT.resize(VoxelNumber);
72 fNTotalLETD.resize(VoxelNumber);
73 fDTotalLETT.resize(VoxelNumber);
74 fDTotalLETD.resize(VoxelNumber);
75
76 fTotalLETD.resize(VoxelNumber);
77 fTotalLETT.resize(VoxelNumber);
78
79 fInitialized = true;
80}
array_type fDTotalLETT
Definition LET.hh:70
array_type fTotalLETT
Definition LET.hh:77
array_type fNTotalLETT
Definition LET.hh:69
array_type fNTotalLETD
Definition LET.hh:73
array_type fTotalLETD
Definition LET.hh:78
array_type fDTotalLETD
Definition LET.hh:74
static VoxelizedSensitiveDetector * GetInstance()
Static method to retrieve a pointer to the only object existing in the simulation.
G4int GetTotalVoxelNumber() const
Method to get the total voxel number.

◆ Compute()

void RadioBio::LET::Compute ( )
overridevirtual

Implements RadioBio::VRadiobiologicalQuantity.

Definition at line 84 of file LET.cc.

85{
86 // Skip LET computation if calculation not enabled.
88 if (fVerboseLevel > 0) {
89 G4cout << "LET::Compute() called but skipped as calculation not enabled" << G4endl;
90 }
91 return;
92 }
93 if (fVerboseLevel > 0) G4cout << "LET::Compute()" << G4endl;
94
96 G4Exception("LET::ComputeLET", "PointerNotAvailable", FatalException,
97 "Computing LET without voxelized geometry pointer!");
98
100
101 for (G4int v = 0; v < VoxelNumber; v++) {
102 if (fVerboseLevel > 1) G4cout << "COMPUTING LET of voxel " << v << G4endl;
103
104 // Compute total LET
105 if (fDTotalLETD[v] > 0.) fTotalLETD[v] = fNTotalLETD[v] / fDTotalLETD[v];
106 // G4cout << "ERASE ME. DEBUG: voxel = " << v << " fNTotalLETD[v] = " << fNTotalLETD[v]
107 // << " fDTotalLETD[v] = " << fDTotalLETD[v] << G4endl;
108 if (fDTotalLETT[v] > 0.) fTotalLETT[v] = fNTotalLETT[v] / fDTotalLETT[v];
109 }
110
111 // Sort ions by A and then by Z ...
112 std::sort(fIonLetStore.begin(), fIonLetStore.end());
113
114 // Compute LET Track and LET Dose for ions
115 for (size_t ion = 0; ion < fIonLetStore.size(); ion++) {
116 fIonLetStore[ion].Calculate();
117 }
118
119 fCalculated = true;
120}

◆ Reset()

void RadioBio::LET::Reset ( )
overridevirtual

Implements RadioBio::VRadiobiologicalQuantity.

Definition at line 233 of file LET.cc.

234{
235 if (fVerboseLevel > 1) {
236 G4cout << "LET::Reset(): ";
237 }
238 fNTotalLETT = 0.0;
239 fNTotalLETD = 0.0;
240 fDTotalLETT = 0.0;
241 fDTotalLETD = 0.0;
242
243 fTotalLETD = 0.0;
244 fTotalLETT = 0.0;
245
246 fIonLetStore.clear();
247 fCalculated = false;
248}

◆ Store()

void RadioBio::LET::Store ( )
overridevirtual

Implements RadioBio::VRadiobiologicalQuantity.

Definition at line 125 of file LET.cc.

126{
127 // Skip LET storing if calculation not enabled.
128 if (!fCalculationEnabled) {
129 if (fVerboseLevel > 0) {
130 G4cout << "LET::Store() called but skipped as calculation not enabled" << G4endl;
131 }
132 return;
133 }
134
135 if (fSaved == true)
136 G4Exception("LET::StoreLET", "NtuplesAlreadySaved", FatalException,
137 "Overwriting LET file. For multiple runs, change filename.");
138
139 Compute();
140#define width 15L
141 if (fVerboseLevel > 0) G4cout << "LET::StoreLET" << G4endl;
142
143 // If there is at least one ion
144 if (fIonLetStore.size()) {
145 std::ofstream ofs(fPath);
146 // ofs.open(fPath, std::ios::out);
147 if (ofs.is_open()) {
148 // Write the voxels index and total LETs and the list of particles/ions
149 ofs << std::setprecision(6) << std::left << "i\tj\tk\t";
150 ofs << std::setw(width) << "LDT";
151 ofs << std::setw(width) << "LTT";
152
153 for (size_t l = 0; l < fIonLetStore.size(); l++) // Write ions name
154 {
155 G4String a = (fIonLetStore[l].IsPrimary()) ? "_1_D" : "_D";
156 ofs << std::setw(width) << fIonLetStore[l].GetName() + a;
157 G4String b = (fIonLetStore[l].IsPrimary()) ? "_1_T" : "_T";
158 ofs << std::setw(width) << fIonLetStore[l].GetName() + b;
159 }
160 ofs << G4endl;
161
162 // Write data
163 G4AnalysisManager* LETFragmentTuple = G4AnalysisManager::Instance();
164
165 LETFragmentTuple->SetVerboseLevel(1);
166 LETFragmentTuple->SetFirstHistoId(1);
167 LETFragmentTuple->SetFirstNtupleId(1);
168
169 LETFragmentTuple->SetDefaultFileType("xml");
170 LETFragmentTuple->OpenFile("LET");
171
172 LETFragmentTuple->CreateNtuple("coordinate", "LET");
173
174 LETFragmentTuple->CreateNtupleIColumn("i"); // 0
175 LETFragmentTuple->CreateNtupleIColumn("j"); // 1
176 LETFragmentTuple->CreateNtupleIColumn("k"); // 2
177 LETFragmentTuple->CreateNtupleDColumn("TotalLETD"); // 3
178 LETFragmentTuple->CreateNtupleDColumn("TotalLETT"); // 4
179 LETFragmentTuple->CreateNtupleIColumn("A"); // 5
180 LETFragmentTuple->CreateNtupleIColumn("Z"); // 6
181 LETFragmentTuple->CreateNtupleDColumn("IonLetD"); // 7
182 LETFragmentTuple->CreateNtupleDColumn("IonLetT"); // 8
183 LETFragmentTuple->FinishNtuple();
184
185 auto voxSensDet = VoxelizedSensitiveDetector::GetInstance();
186
187 for (G4int i = 0; i < voxSensDet->GetVoxelNumberAlongX(); i++)
188 for (G4int j = 0; j < voxSensDet->GetVoxelNumberAlongY(); j++)
189 for (G4int k = 0; k < voxSensDet->GetVoxelNumberAlongZ(); k++) {
190 LETFragmentTuple->FillNtupleIColumn(1, 0, i);
191 LETFragmentTuple->FillNtupleIColumn(1, 1, j);
192 LETFragmentTuple->FillNtupleIColumn(1, 2, k);
193
194 G4int v = voxSensDet->GetThisVoxelNumber(i, j, k);
195
196 // Write total LETs and voxels index
197 ofs << G4endl;
198 ofs << i << '\t' << j << '\t' << k << '\t';
199 ofs << std::setw(width) << fTotalLETD[v] / (keV / um);
200 ofs << std::setw(width) << fTotalLETT[v] / (keV / um);
201
202 // Write ions LETs
203 for (size_t l = 0; l < fIonLetStore.size(); l++) {
204 // Write ions LETs
205 ofs << std::setw(width) << fIonLetStore[l].GetLETD()[v] / (keV / um);
206 ofs << std::setw(width) << fIonLetStore[l].GetLETT()[v] / (keV / um);
207 }
208
209 LETFragmentTuple->FillNtupleDColumn(1, 3, fTotalLETD[v] / (keV / um));
210 LETFragmentTuple->FillNtupleDColumn(1, 4, fTotalLETT[v] / (keV / um));
211
212 for (size_t ll = 0; ll < fIonLetStore.size(); ll++) {
213 LETFragmentTuple->FillNtupleIColumn(1, 5, fIonLetStore[ll].GetA());
214 LETFragmentTuple->FillNtupleIColumn(1, 6, fIonLetStore[ll].GetZ());
215
216 LETFragmentTuple->FillNtupleDColumn(1, 7, fIonLetStore[ll].GetLETD()[v] / (keV / um));
217 LETFragmentTuple->FillNtupleDColumn(1, 8, fIonLetStore[ll].GetLETT()[v] / (keV / um));
218 LETFragmentTuple->AddNtupleRow(1);
219 }
220 }
221 ofs.close();
222
223 // LETFragmentTuple->Write();
224 LETFragmentTuple->CloseFile();
225 }
226 }
227
228 fSaved = true;
229}
#define width
Definition Dose.cc:44
std::vector< ExP01TrackerHit * > a
void Compute() override
Definition LET.cc:84

◆ PrintParameters()

void RadioBio::LET::PrintParameters ( )
overridevirtual

Reimplemented from RadioBio::VRadiobiologicalQuantity.

Definition at line 285 of file LET.cc.

286{
287 G4cout << "*******************************************" << G4endl
288 << "******* Parameters of the class LET *******" << G4endl
289 << "*******************************************" << G4endl;
291}

◆ SetNTotalLETT()

void RadioBio::LET::SetNTotalLETT ( const array_type  NTotalLETT)
inlineprivate

Definition at line 85 of file LET.hh.

86 {
87 fNTotalLETT = NTotalLETT;
88 }

◆ SetNTotalLETD()

void RadioBio::LET::SetNTotalLETD ( const array_type  NTotalLETD)
inlineprivate

Definition at line 89 of file LET.hh.

90 {
91 fNTotalLETD = NTotalLETD;
92 }

◆ SetDTotalLETT()

void RadioBio::LET::SetDTotalLETT ( const array_type  DTotalLETT)
inlineprivate

Definition at line 93 of file LET.hh.

94 {
95 fDTotalLETT = DTotalLETT;
96 }

◆ SetDTotalLETD()

void RadioBio::LET::SetDTotalLETD ( const array_type  DTotalLETD)
inlineprivate

Definition at line 97 of file LET.hh.

98 {
99 fDTotalLETD = DTotalLETD;
100 }

◆ AddNTotalLETT()

void RadioBio::LET::AddNTotalLETT ( const array_type  NTotalLETT)
inlineprivate

Definition at line 102 of file LET.hh.

103 {
104 fNTotalLETT += NTotalLETT;
105 }

◆ AddNTotalLETD()

void RadioBio::LET::AddNTotalLETD ( const array_type  NTotalLETD)
inlineprivate

Definition at line 106 of file LET.hh.

107 {
108 fNTotalLETD += NTotalLETD;
109 }

◆ AddDTotalLETT()

void RadioBio::LET::AddDTotalLETT ( const array_type  DTotalLETT)
inlineprivate

Definition at line 110 of file LET.hh.

111 {
112 fDTotalLETT += DTotalLETT;
113 }

◆ AddDTotalLETD()

void RadioBio::LET::AddDTotalLETD ( const array_type  DTotalLETD)
inlineprivate

Definition at line 114 of file LET.hh.

115 {
116 fDTotalLETD += DTotalLETD;
117 }

◆ AddIon()

void RadioBio::LET::AddIon ( const IonLet  ion)
private

Member Data Documentation

◆ fNTotalLETT

array_type RadioBio::LET::fNTotalLETT = {}
private

Definition at line 69 of file LET.hh.

69{};

◆ fDTotalLETT

array_type RadioBio::LET::fDTotalLETT = {}
private

Definition at line 70 of file LET.hh.

70{};

◆ fNTotalLETD

array_type RadioBio::LET::fNTotalLETD = {}
private

Definition at line 73 of file LET.hh.

73{};

◆ fDTotalLETD

array_type RadioBio::LET::fDTotalLETD = {}
private

Definition at line 74 of file LET.hh.

74{};

◆ fTotalLETT

array_type RadioBio::LET::fTotalLETT = {}
private

Definition at line 77 of file LET.hh.

77{};

◆ fTotalLETD

array_type RadioBio::LET::fTotalLETD = {}
private

Definition at line 78 of file LET.hh.

78{};

◆ fIonLetStore

std::vector<IonLet> RadioBio::LET::fIonLetStore
private

Definition at line 80 of file LET.hh.


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

Applications | User Support | Publications | Collaboration