Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4VoxelRightBreast.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 // Authors: S. Guatelli and M. G. Pia, INFN Genova, Italy
27 //
28 // Based on code developed by the undergraduate student G. Guerrieri
29 // Note: this is a preliminary beta-version of the code; an improved
30 // version will be distributed in the next Geant4 public release, compliant
31 // with the design in a forthcoming publication, and subject to a
32 // design and code review.
33 //
34 #include "G4VoxelRightBreast.hh"
35 
36 #include "globals.hh"
37 #include "G4SystemOfUnits.hh"
38 #include "G4SDManager.hh"
39 #include "G4VisAttributes.hh"
40 #include "G4Tubs.hh"
41 #include "G4ThreeVector.hh"
42 #include "G4VPhysicalVolume.hh"
43 #include "G4RotationMatrix.hh"
44 #include "G4Material.hh"
45 #include "G4LogicalVolume.hh"
47 #include "G4VPhysicalVolume.hh"
48 #include "G4PVPlacement.hh"
49 #include "G4PVReplica.hh"
51 #include "G4VoxelRightBreastSD.hh"
52 #include "G4HumanPhantomColour.hh"
53 
55 {
56 }
57 
59 {
60 }
61 
63  G4VPhysicalVolume* mother,
64  const G4String& colourName,
65  G4bool wireFrame,G4bool sensitivity)
66 
67 {
69  G4Material* adipose = material -> GetMaterial("adipose");
70  G4Material* adipose_glandular = material -> GetMaterial("adipose_glandular");
71  delete material;
72 
73  G4double rmin = 0.* cm;
74  G4double rmax = 6. * cm;
75  G4double zz = 6. * cm;
76  G4double startPhi = 0.* degree;
77  G4double spanningPhi = 180. * degree;
78 
79  G4Tubs* breast = new G4Tubs("out_breast", rmin, rmax,
80  zz/2., startPhi, spanningPhi);
81 
82  G4LogicalVolume* breast_log = new G4LogicalVolume(breast,
83  adipose,
84  "logicalOut"+ volumeName,
85  0, 0, 0);
86  rmax = 5.5 *cm;
87  zz = 5. *cm;
88  G4Tubs* innerBreast = new G4Tubs("inner_breast",
89  rmin, rmax,
90  zz/2., startPhi, spanningPhi);
91 
92  G4LogicalVolume* innerBreast_log = new G4LogicalVolume(innerBreast,
93  adipose_glandular,
94  "logical"+ volumeName,
95  0, 0, 0);
96 
97  G4RotationMatrix* matrix = new G4RotationMatrix();
98  matrix -> rotateX(-90.* degree);
99  matrix -> rotateY(180.* degree);
100  matrix -> rotateZ(18. * degree);
101 
102 
103  G4VPhysicalVolume* physBreast = new G4PVPlacement(matrix,G4ThreeVector(-10.*cm, 52.* cm, 8.7 *cm),
104  "physicalVoxelRightBreast",
105  breast_log,
106  mother,
107  false,
108  0, true);
109 
110  G4VPhysicalVolume* physInnerBreast = new G4PVPlacement(0,G4ThreeVector(),
111  "RightBreast",
112  innerBreast_log,
113  physBreast,
114  false,
115  0, true);
116 
117 
118  // Parameterisation with Replicas
119  /*
120  G4double rmin_voxelz = 0.* cm;
121  G4double rmax_voxelz = 5.5 * cm;
122  G4double zz_voxels = 5. * cm;
123  G4double startPhi_voxelz = 0.* degree;
124  G4double spanningPhi_voxelz = 180. * degree;
125  G4int nslice = 10;
126 
127  G4Tubs* voxelz = new G4Tubs("voxel_z", rmin_voxelz, rmax_voxelz, zz_voxels/(2*nslice),startPhi_voxelz, spanningPhi_voxelz);
128  G4LogicalVolume* voxelz_log = new G4LogicalVolume(voxelz, adipose_glandular, "voxelz_log",0 , 0, 0);
129  G4VPhysicalVolume* voxelz_phys = new G4PVReplica("LinearArray", voxelz_log, physInnerBreast,
130  kZAxis, nslice, zz_voxels/nslice);
131 
132  G4double voxel_height = (zz_voxels/nslice);
133 
134  G4int n_voxels = 10;
135 
136  G4double voxel_phi = spanningPhi_voxelz/n_voxels;
137 
138  G4cout <<"voxel_phi: " <<voxel_phi/degree << " deg"<< G4endl;
139 
140  G4Tubs* voxel = new G4Tubs("voxel", rmin_voxelz, rmax_voxelz, voxel_height/2.,
141  startPhi_voxelz,voxel_phi);
142 
143  G4LogicalVolume* voxel_log = new G4LogicalVolume(voxel, adipose_glandular, "voxel_log", 0,0,0);
144 
145  G4VPhysicalVolume* voxel_phys = new G4PVReplica("RZPhiSlices",
146  voxel_log,
147  voxelz_phys,
148  kPhi, n_voxels, voxel_phi,- (voxel_phi/2.));
149 
150 
151  */
152 
153  G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour();
154  G4Colour colour = colourPointer -> GetColour(colourName);
155  delete colourPointer;
156 
157  G4VisAttributes* BreastVisAtt = new G4VisAttributes(colour);
158  BreastVisAtt -> SetForceSolid(wireFrame);
159  breast_log -> SetVisAttributes(BreastVisAtt);
160 
161  G4VisAttributes* innerBreastVisAtt = new G4VisAttributes(colour);
162  innerBreastVisAtt -> SetForceSolid(true);
163  innerBreast_log -> SetVisAttributes(innerBreastVisAtt);
164 
165  // voxelz_log -> SetVisAttributes(innerBreastVisAtt);
166  // voxel_log -> SetVisAttributes(innerBreastVisAtt);
167 
168  G4cout << "Voxel Right Breast created !!!!!!" << G4endl;
169 
170 if (sensitivity==true)
171  {
173 
174  G4String sensitiveDetectorName = "VoxelRightBreast";
175 
176  G4VoxelRightBreastSD* breastSD = new G4VoxelRightBreastSD(sensitiveDetectorName);
177 
178  G4String ROGeometryName = "VoxelRightROGeometry";
179 
180  G4VoxelRightBreastROGeometry* phantomROGeometry = new G4VoxelRightBreastROGeometry(ROGeometryName);
181  phantomROGeometry -> BuildROGeometry();
182  breastSD -> SetROgeometry(phantomROGeometry);
183  SDman -> AddNewDetector(breastSD);
184  innerBreast_log -> SetSensitiveDetector(breastSD);
185  }
186  return physInnerBreast;
187 
188 }