Geant4  10.00.p02
HadrontherapyLet.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 
28 #include "HadrontherapyLet.hh"
32 #include "HadrontherapyMatrix.hh"
33 #include "G4RunManager.hh"
34 #include "G4SystemOfUnits.hh"
35 #include <cmath>
36 
39 
41 {
42  if (instance) delete instance;
43  instance = new HadrontherapyLet(pDet);
44  return instance;
45 }
46 
48 {
49  return instance;
50 }
51 
53  :filename("Let.out") // Default output filename
54 {
55 
57 
58  if (!matrix) G4Exception("HadrontherapyLet::HadrontherapyLet","Hadrontherapy0005", FatalException, "HadrontherapyMatrix not found. Firstly create an instance of it.");
59 
60  nVoxels = matrix -> GetNvoxel();
61 
62  numberOfVoxelAlongX = matrix -> GetNumberOfVoxelAlongX();
63  numberOfVoxelAlongY = matrix -> GetNumberOfVoxelAlongY();
64  numberOfVoxelAlongZ = matrix -> GetNumberOfVoxelAlongZ();
65 
67  pPGA = (HadrontherapyPrimaryGeneratorAction*)runManager -> GetUserPrimaryGeneratorAction();
68  // Pointer to the detector material
69  detectorMat = pDet -> GetDetectorLogicalVolume() -> GetMaterial();
70  density = detectorMat -> GetDensity();
71 // Instances for Total LET
72  totalLetD = new G4double[nVoxels];
74 
75 }
77 {
78  Clear();
79  delete [] totalLetD;
80  delete [] DtotalLetD;
81 }
82 
83 // Fill energy spectrum for every voxel (local energy spectrum)
85 {
86  for(G4int v=0; v < nVoxels; v++) totalLetD[v] = DtotalLetD[v] = 0.;
87  Clear();
88 }
93 {
94  for (size_t i=0; i < ionLetStore.size(); i++)
95  {
96  delete [] ionLetStore[i].letDN;
97  delete [] ionLetStore[i].letDD;
98  }
99  ionLetStore.clear();
100 }
102  G4ParticleDefinition* particleDef,
103  /*G4double kinEnergy,*/
104  G4double DE,
105  G4double DX,
106  G4int i, G4int j, G4int k)
107 {
108  if (DE <= 0. || DX <=0.) return;
109  if (!doCalculation) return;
110  G4int Z = particleDef -> GetAtomicNumber();
111 
112 
113  G4int PDGencoding = particleDef -> GetPDGEncoding();
114  PDGencoding -= PDGencoding%10;
115 
116  G4int voxel = matrix -> Index(i,j,k);
117  // Total LET calculation...
118  totalLetD[voxel] += DE*(DE/DX);
119  DtotalLetD[voxel] += DE;
120  // Single ion LET
121  if (Z>=1)
122  {
123  // Search for already allocated data...
124  size_t l;
125  for (l=0; l < ionLetStore.size(); l++)
126  {
127  if (ionLetStore[l].PDGencoding == PDGencoding)
128  if ( ((trackID ==1) && (ionLetStore[l].isPrimary)) || ((trackID !=1) && (!ionLetStore[l].isPrimary)))
129  break;
130  }
131 
132  if (l == ionLetStore.size()) // Just another type of ion/particle for our store...
133  {
134 
135  G4int A = particleDef -> GetAtomicMass();
136 
137  G4String fullName = particleDef -> GetParticleName();
138  G4String name = fullName.substr (0, fullName.find("[") ); // cut excitation energy [x.y]
139 
140  ionLet ion =
141  {
142  (trackID == 1) ? true:false, // is it the primary particle?
143  PDGencoding,
144  fullName,
145  name,
146  Z,
147  A,
148  new G4double[nVoxels], // Let Dose Numerator
149  new G4double[nVoxels] // Let Dose Denominator
150  };
151 
152  // Initialize let
153  for(G4int v=0; v < nVoxels; v++) ion.letDN[v] = ion.letDD[v] = 0.;
154  ionLetStore.push_back(ion);
155  //G4cout << "Allocated LET data for " << ion.name << G4endl;
156 
157  }
158  ionLetStore[l].letDN[voxel] += DE*(DE/DX);
159  ionLetStore[l].letDD[voxel] += DE;
160  }
161 }
162 
163 
164 
165 
167 {
168  for(G4int v=0; v < nVoxels; v++) if (DtotalLetD[v]>0.) totalLetD[v] = totalLetD[v]/DtotalLetD[v];
169  // Sort ions by A and then by Z ...
170  std::sort(ionLetStore.begin(), ionLetStore.end());
171  // Compute Let Track and Let Dose for any single ion
172 
173  for(G4int v=0; v < nVoxels; v++)
174  {
175  for (size_t ion=0; ion < ionLetStore.size(); ion++)
176  {
177  if (ionLetStore[ion].letDD[v] >0.) ionLetStore[ion].letDN[v] = ionLetStore[ion].letDN[v] / ionLetStore[ion].letDD[v];
178 
179  }// end loop over ions
180  }
181 
182 }// end loop over voxels
183 
185 {
186 #define width 15L
187  if(ionLetStore.size())
188  {
189  ofs.open(filename, std::ios::out);
190  if (ofs.is_open())
191  {
192 
193  // Write the voxels index and the list of particles/ions
194  ofs << std::setprecision(6) << std::left <<
195  "i\tj\tk\t";
196  ofs << std::setw(width) << "LDT";
197  for (size_t l=0; l < ionLetStore.size(); l++)
198  {
199  G4String a = (ionLetStore[l].isPrimary) ? "_1":"";
200  ofs << std::setw(width) << ionLetStore[l].name + a ;
201  }
202  ofs << G4endl;
203 
204  // Write data
205  for(G4int i = 0; i < numberOfVoxelAlongX; i++)
206  for(G4int j = 0; j < numberOfVoxelAlongY; j++)
207  for(G4int k = 0; k < numberOfVoxelAlongZ; k++)
208  {
209  G4int v = matrix -> Index(i, j, k);
210  // row write
211  for (size_t l=0; l < ionLetStore.size(); l++)
212  {
213  // Write only not identically null data lines
214  if(ionLetStore[l].letDN)
215  {
216  ofs << G4endl;
217  ofs << i << '\t' << j << '\t' << k << '\t';
218 
219  ofs << std::setw(width) << totalLetD[v]/(keV/um);
220  for (size_t ll=0; ll < ionLetStore.size(); ll++)
221  {
222  ofs << std::setw(width) << ionLetStore[ll].letDN[v]/(keV/um) ;
223  }
224  break;
225  }
226  }
227  }
228  ofs.close();
229  G4cout << "Let is being written to " << filename << G4endl;
230  }
231 
232  }
233 }
234 
236 {
237 #ifdef G4ANALYSIS_USE_ROOT
238 
240 
241  for(G4int i = 0; i < numberOfVoxelAlongX; i++)
242  for(G4int j = 0; j < numberOfVoxelAlongY; j++)
243  for(G4int k = 0; k < numberOfVoxelAlongZ; k++)
244  {
245  G4int v = matrix -> Index(i, j, k);
246  for (size_t ion=0; ion < ionLetStore.size(); ion++)
247  {
248 
249  analysis -> FillLetFragmentTuple( i, j, k, ionLetStore[ion].A, ionLetStore[ion].Z, ionLetStore[ion].letDN[v]);
250 
251 
252  }
253  }
254 
255 #endif
256 }
257 
static G4bool doCalculation
void FillEnergySpectrum(G4int trackID, G4ParticleDefinition *particleDef, G4double DE, G4double DX, G4int i, G4int j, G4int k)
G4double * letDD
static HadrontherapyAnalysisManager * GetInstance()
Get the pointer to the analysis manager.
HadrontherapyMatrix * matrix
static HadrontherapyLet * GetInstance()
HadrontherapyPrimaryGeneratorAction * pPGA
G4String name
Definition: TRTMaterials.hh:40
#define width
G4double a
Definition: TRTMaterials.hh:39
int G4int
Definition: G4Types.hh:78
static HadrontherapyMatrix * GetInstance()
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
HadrontherapyLet(HadrontherapyDetectorConstruction *)
G4Material * detectorMat
std::vector< ionLet > ionLetStore
void Clear()
Clear all stored data.
static const G4double A[nN]
std::ofstream ofs
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:74
A class for connecting the simulation to an analysis package.
#define G4endl
Definition: G4ios.hh:61
static HadrontherapyLet * instance
G4double * letDN
static const double keV
Definition: G4SIunits.hh:195
double G4double
Definition: G4Types.hh:76