Geant4_10
PhantomConfiguration.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 // This example is provided by the Geant4-DNA collaboration
27 // Any report or published results obtained using the Geant4-DNA software
28 // shall cite the following Geant4-DNA collaboration publication:
29 // Med. Phys. 37 (2010) 4692-4708
30 // The Geant4-DNA web site is available at http://geant4-dna.org
31 //
32 // If you use this example, please cite the following publication:
33 // Rad. Prot. Dos. 133 (2009) 2-11
34 
35 #include "PhantomConfiguration.hh"
36 #include "G4SystemOfUnits.hh"
37 
38 G4int PhantomConfiguration::fPhantomTotalPixels = 0;
39 G4int PhantomConfiguration::fNucleusTotalPixels = 0;
40 G4int PhantomConfiguration::fCytoplasmTotalPixels = 0;
41 G4float PhantomConfiguration::fDx = 0;
42 G4float PhantomConfiguration::fDy = 0;
43 G4float PhantomConfiguration::fDz = 0;
44 G4float PhantomConfiguration::fNucleusMass = 0;
45 G4float PhantomConfiguration::fCytoplasmMass = 0;
46 
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48 
50 {
51  Initialize();
52 }
53 
55 {
56  delete[] fVoxelThreeVector;
57 }
58 
59 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60 
62 {
63  G4float vx, vy, vz, tmp, density;
64  G4int den, mat;
65  G4float denCyto1, denCyto2, denCyto3, denNucl1, denNucl2, denNucl3;
66 
67  FILE* fMap;
68 
69  density=0;
70 
71  fPhantomTotalPixels=0;
72  fNucleusTotalPixels=0;
73  fCytoplasmTotalPixels=0;
74  fDx=0;
75  fDy=0;
76  fDz=0;
77  fNucleusMass=0;
78  fCytoplasmMass=0;
79 
80  // READ PHANTOM PARAMETERS
81  fMap = fopen("phantom.dat","r");
82 
83  fscanf(fMap,"%i %i %i",&fPhantomTotalPixels, &fNucleusTotalPixels, &fCytoplasmTotalPixels);
84  fscanf(fMap,"%f %f %f",&fDx, &fDy, &fDz);
85  fscanf(fMap,"%f %f %f",&tmp, &tmp, &tmp);
86  fscanf(fMap,"%f %f %f",&denCyto1, &denCyto2, &denCyto3);
87  fscanf(fMap,"%f %f %f",&denNucl1, &denNucl2, &denNucl3);
88  fDx = fDx * micrometer;
89  fDy = fDy * micrometer;
90  fDz = fDz * micrometer;
91 
92  fVoxelThreeVector = new G4ThreeVector [fPhantomTotalPixels];
93 
94  for (G4int i=0; i<fPhantomTotalPixels; i++)
95  {
96  fscanf(fMap,"%f %f %f %i %i %f",&vx, &vy, &vz, &mat, &den, &tmp);
97 
98  if (std::abs(mat-2)<1.e-30) // NUCLEUS
99  {
100  if (std::abs(den-1)<1.e-30) density = denNucl1*(g/cm3);
101  if (std::abs(den-2)<1.e-30) density = denNucl2*(g/cm3);
102  if (std::abs(den-3)<1.e-30) density = denNucl3*(g/cm3);
103  fNucleusMass = fNucleusMass + density * fDx * fDy * fDz ;
104  }
105 
106  if (std::abs(mat-1)<1.e-30) // CYTOPLASM
107  {
108  if (std::abs(den-1)<1e-30) density = denCyto1*(g/cm3);
109  if (std::abs(den-2)<1e-30) density = denCyto2*(g/cm3);
110  if (std::abs(den-3)<1e-30) density = denCyto3*(g/cm3);
111  fCytoplasmMass = fCytoplasmMass + density * fDx * fDy * fDz ;
112  }
113 
114  G4ThreeVector v(vx,vy,vz);
115  fVoxelThreeVector[i] = v;
116  }
117 
118  fclose(fMap);
119 
120  return 0;
121 }
122 
Float_t den
Definition: plot.C:37
Float_t tmp
Definition: plot.C:37
float G4float
Definition: G4Types.hh:77
int G4int
Definition: G4Types.hh:78
Float_t mat
Definition: plot.C:40
G4double density
Definition: TRTMaterials.hh:39
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
tuple v
Definition: test.py:18
fclose(fp)
int micrometer
Definition: hepunit.py:34