source: trunk/examples/advanced/medical_linac/src/MedLinacPhantomROGeometry.cc @ 1292

Last change on this file since 1292 was 807, checked in by garnier, 16 years ago

update

File size: 6.9 KB
Line 
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// $Id: MedLinacPhantomROGeometry.cc,v 1.5 2006/06/29 16:04:35 gunter Exp $
27//
28//
29// Code developed by: M. Piergentili
30//
31#include "MedLinacPhantomROGeometry.hh"
32#include "MedLinacDummySD.hh"
33
34#include "G4LogicalVolume.hh"
35#include "G4VPhysicalVolume.hh"
36#include "G4PVPlacement.hh"
37#include "G4PVReplica.hh"
38#include "G4SDManager.hh"
39#include "G4Box.hh"
40#include "G4Tubs.hh"
41#include "G4SubtractionSolid.hh"
42#include "G4ThreeVector.hh"
43#include "G4Material.hh"
44
45MedLinacPhantomROGeometry::MedLinacPhantomROGeometry(G4String aString,
46                                                 G4double phantomDim_x,
47                                                 G4double phantomDim_y,
48                                                 G4double phantomDim_z,
49                                                 G4int numberOfVoxelsX,
50                                                 G4int numberOfVoxelsY,
51                                                 G4int numberOfVoxelsZ):
52  G4VReadOutGeometry(aString),
53  PhantomDimensionX(phantomDim_x),
54  PhantomDimensionY(phantomDim_y),
55  PhantomDimensionZ(phantomDim_z),
56  NumberOfVoxelsAlongX(numberOfVoxelsX),
57  NumberOfVoxelsAlongY(numberOfVoxelsY),
58  NumberOfVoxelsAlongZ(numberOfVoxelsZ)
59{}
60MedLinacPhantomROGeometry::~MedLinacPhantomROGeometry()
61{
62}
63
64G4VPhysicalVolume* MedLinacPhantomROGeometry::Build()
65{
66 
67  G4cout << " MedLinacPhantomROGeometry::Build()-----------------------" << G4endl;
68 
69  //G4Material* dummyMat = new G4Material(name="dummyMat", 1., 1.*g/mole, 1.*g/cm3);
70  G4double a;  // atomic mass
71  G4double z;  // atomic number
72  G4int natoms;
73  G4String name;
74  G4String symbol;
75  a = 16.00*g/mole;
76  G4Element* elO = new G4Element(name="Oxygen",symbol="O", z=8., a);
77
78  a = 1.00794*g/mole;
79  G4Element* elH = new G4Element(name="Hydrogen",symbol="H", z=1., a);
80
81  G4double density;
82  G4int ncomponents;
83  density = 1.000*g/cm3;
84  G4Material* H2O = new G4Material(name="Water",density, ncomponents=2);
85  H2O->AddElement(elH, natoms=2);
86  H2O->AddElement(elO, natoms=1);
87
88  G4double expHall_x = 3.0*m;
89  G4double expHall_y = 3.0*m;
90  G4double expHall_z = 3.0*m;
91
92  // variables for  division ...
93  G4double halfXVoxelDimensionX = PhantomDimensionX/NumberOfVoxelsAlongX;
94  G4double halfXVoxelDimensionZ = PhantomDimensionZ;
95  G4double voxelXThickness = 2*halfXVoxelDimensionX;
96
97
98  // world volume of ROGeometry ...
99  G4Box *ROWorld = new G4Box("ROWorld",
100                             expHall_x,
101                             expHall_y,
102                             expHall_z);
103
104  G4LogicalVolume *ROWorldLog = new G4LogicalVolume(ROWorld,
105                                                    H2O,
106                                                    "ROWorldLog",
107                                                    0,0,0);
108
109  G4VPhysicalVolume *ROWorldPhys = new G4PVPlacement(0,
110                                                     G4ThreeVector(),
111                                                     "ROWorldPhys",
112                                                     ROWorldLog,
113                                                     0,false,0);
114 
115  // phantom ROGeometry ...
116  G4Box *ROPhantom = new G4Box("ROPhantom", 
117                               PhantomDimensionX, 
118                               PhantomDimensionY, 
119                               PhantomDimensionZ);
120
121  G4LogicalVolume *ROPhantomLog = new G4LogicalVolume(ROPhantom,
122                                                      H2O,
123                                                      "ROPhantomLog",
124                                                      0,0,0);
125
126  G4VPhysicalVolume *ROPhantomPhys = new G4PVPlacement(0,
127                                                       G4ThreeVector(),
128                                                       "PhantomPhys",
129                                                       ROPhantomLog,
130                                                       ROWorldPhys,
131                                                       false,0);
132 
133  // ROGeomtry: Voxel division
134 
135  // X division first...
136  G4Box *ROPhantomXDivision = new G4Box("ROPhantomXDivision",
137                                        halfXVoxelDimensionX,
138                                        halfXVoxelDimensionZ,
139                                        halfXVoxelDimensionZ);
140
141  G4LogicalVolume *ROPhantomXDivisionLog = new G4LogicalVolume(ROPhantomXDivision,
142                                                               H2O,
143                                                               "ROPhantomXDivisionLog",
144                                                               0,0,0);
145
146  G4VPhysicalVolume *ROPhantomXDivisionPhys = new G4PVReplica("ROPhantomXDivisionPhys",
147                                                              ROPhantomXDivisionLog,
148                                                              ROPhantomPhys,
149                                                              kXAxis,
150                                                              NumberOfVoxelsAlongX,
151                                                              voxelXThickness);
152  // ...then Z division
153 
154  G4Box *ROPhantomZDivision = new G4Box("ROPhantomZDivision",
155                                        halfXVoxelDimensionX,
156                                        halfXVoxelDimensionZ, 
157                                        halfXVoxelDimensionX);
158
159  G4LogicalVolume *ROPhantomZDivisionLog = new G4LogicalVolume(ROPhantomZDivision,
160                                                               H2O,
161                                                               "ROPhantomZDivisionLog",
162                                                               0,0,0);
163
164  G4VPhysicalVolume *ROPhantomZDivisionPhys = new G4PVReplica("ROPhantomZDivisionPhys",
165                                                              ROPhantomZDivisionLog,
166                                                              ROPhantomXDivisionPhys,
167                                                              kZAxis,
168                                                              NumberOfVoxelsAlongZ,
169                                                              voxelXThickness);
170  // ...then Y  division
171
172  G4Box *ROPhantomYDivision = new G4Box("ROPhantomYDivision",
173                                        halfXVoxelDimensionX, 
174                                        halfXVoxelDimensionX,
175                                        halfXVoxelDimensionX);
176
177  G4LogicalVolume *ROPhantomYDivisionLog = new G4LogicalVolume(ROPhantomYDivision,
178                                                               H2O,
179                                                               "ROPhantomYDivisionLog",
180                                                               0,0,0);
181 
182  ROPhantomYDivisionPhys = new G4PVReplica("ROPhantomYDivisionPhys",
183                                                              ROPhantomYDivisionLog,
184                                                              ROPhantomZDivisionPhys,
185                                                              kYAxis,
186                                                              NumberOfVoxelsAlongZ,
187                                                              voxelXThickness);
188
189  MedLinacDummySD *dummySD = new  MedLinacDummySD;
190  ROPhantomYDivisionLog->SetSensitiveDetector(dummySD);
191
192  return ROWorldPhys;
193}
194
195
196
Note: See TracBrowser for help on using the repository browser.