source: trunk/examples/advanced/human_phantom/src/G4MIRDRibCage.cc@ 1230

Last change on this file since 1230 was 807, checked in by garnier, 17 years ago

update

File size: 8.5 KB
RevLine 
[807]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 "G4MIRDRibCage.hh"
35
36#include "globals.hh"
37#include "G4SDManager.hh"
38#include "G4VisAttributes.hh"
39#include "G4HumanPhantomMaterial.hh"
40#include "G4SubtractionSolid.hh"
41#include "G4EllipticalTube.hh"
42#include "G4PVReplica.hh"
43#include "G4Box.hh"
44#include "G4PVPlacement.hh"
45#include "G4HumanPhantomColour.hh"
46
47G4MIRDRibCage::G4MIRDRibCage()
48{
49}
50
51G4MIRDRibCage::~G4MIRDRibCage()
52{
53 }
54
55G4VPhysicalVolume* G4MIRDRibCage::Construct(const G4String& volumeName, G4VPhysicalVolume* mother,
56 const G4String& colourName, G4bool wireFrame, G4bool sensitivity)
57{
58 G4HumanPhantomMaterial* material = new G4HumanPhantomMaterial();
59
60 G4cout << "Construct "<< volumeName <<G4endl;
61
62 G4Material* skeleton = material -> GetMaterial("skeleton");
63 G4Material* soft = material -> GetMaterial("soft_tissue");
64
65 delete material;
66
67 G4double dx= 17. *cm; // a2
68 G4double dy= 9.8 * cm; //b2
69 G4double thickness= 32.4 * cm; // z2/2 of cage
70
71 G4EllipticalTube* outCage = new G4EllipticalTube("outCage",dx, dy, thickness/2.);
72
73 dx = 16.4 * cm; // a1
74 dy = 9.2 * cm; // b1
75 G4double dz = 34. *cm; // z2/2
76
77 G4EllipticalTube* inCage = new G4EllipticalTube("inCage",dx, dy, dz/2.);
78
79 G4SubtractionSolid* cage = new G4SubtractionSolid("Cage",
80 outCage,
81 inCage, 0, G4ThreeVector(0.*cm, 0.*cm, 0. * cm));
82
83
84 G4LogicalVolume* logicRibCage = new G4LogicalVolume(cage, soft, "LogicalCage", 0, 0, 0);
85
86 G4VPhysicalVolume* physRibCage = new G4PVPlacement(0,G4ThreeVector(0.0, 0.0, thickness/2. + 0.1 * cm),
87 // with respect to the trunk
88 "physicalRibCage",
89 logicRibCage,
90 mother,
91 false,
92 0, true);
93
94
95 G4double xx = 17.*cm;
96 G4double yy = 9.8*cm;
97 G4double ribThickness = 1.4*cm;
98 G4EllipticalTube* rib_out = new G4EllipticalTube("rib_out",xx, yy, ribThickness/2.);
99
100 xx = 16.5 *cm;
101 yy = 9.3 * cm;
102 G4double zz = 1.5 * cm;
103 G4EllipticalTube* rib_in = new G4EllipticalTube("rib_in",xx, yy, zz/2.);
104 G4SubtractionSolid* rib = new G4SubtractionSolid("rib",rib_out, rib_in);
105
106 G4LogicalVolume* logicRib= new G4LogicalVolume(rib, skeleton, "logical" + volumeName, 0, 0, 0);
107
108 physRib1 = new G4PVPlacement(0,G4ThreeVector(0.0, 0.0, (- 32.2*cm/2. + 0.8 *cm)),
109 // with respect to the trunk
110 "physicalRib",
111 logicRib,
112 physRibCage,
113 false,
114 0, true);
115
116 physRib2 = new G4PVPlacement(0,G4ThreeVector(0.0, 0.0, ( - 32.2*cm/2. + 0.8 *cm + 2.8 *cm)),
117 // with respect to the trunk
118 "physicalRib",
119 logicRib,
120 physRibCage,
121 false,
122 0, true);
123
124 physRib3 = new G4PVPlacement(0,G4ThreeVector(0.0, 0.0, (-thickness/2. + 0.8 * cm + 5.6 *cm)),
125 // with respect to the trunk
126 "physicalRib",
127 logicRib,
128 physRibCage,
129 false,
130 0, true);
131
132 physRib4 = new G4PVPlacement(0,G4ThreeVector(0.0, 0.0, (-thickness/2. + 0.8 * cm + 8.4 *cm)),
133 // with respect to the trunk
134 "physicalRib",
135 logicRib,
136 physRibCage,
137 false,
138 0, true);
139
140 physRib5 = new G4PVPlacement(0,G4ThreeVector(0.0, 0.0, (-thickness/2. + 0.8 * cm + 11.2 *cm)),
141 // with respect to the trunk
142 "physicalRib",
143 logicRib,
144 physRibCage,
145 false,
146 0, true);
147
148 physRib6 = new G4PVPlacement(0,G4ThreeVector(0.0, 0.0, (-thickness/2. + 0.8 * cm + 14. *cm)),
149 // with respect to the trunk
150 "physicalRib",
151 logicRib,
152 physRibCage,
153 false,
154 0, true);
155
156 physRib7 = new G4PVPlacement(0,G4ThreeVector(0.0, 0.0, (-thickness/2. + 0.8 *cm + 16.8 *cm)),
157 // with respect to the trunk
158 "physicalRib",
159 logicRib,
160 physRibCage,
161 false,
162 0, true);
163
164 physRib8 = new G4PVPlacement(0,G4ThreeVector(0.0, 0.0, (-thickness/2. + 0.8 *cm + 19.6 *cm)),
165 // with respect to the trunk
166 "physicalRib",
167 logicRib,
168 physRibCage,
169 false,
170 0, true);
171
172 physRib9 = new G4PVPlacement(0,G4ThreeVector(0.0, 0.0, (-thickness/2. + 0.8*cm + 22.4 *cm)),
173 // with respect to the trunk
174 "physicalRib",
175 logicRib,
176 physRibCage,
177 false,
178 0, true);
179
180 physRib10 = new G4PVPlacement(0,G4ThreeVector(0.0, 0.0, (-thickness/2. + 0.8*cm + 25.2 *cm)),
181 // with respect to the trunk
182 "physicalRib",
183 logicRib,
184 physRibCage,
185 false,
186 0, true);
187
188 physRib11 = new G4PVPlacement(0,G4ThreeVector(0.0, 0.0, (-thickness/2. + 0.8*cm + 28. *cm)),
189 // with respect to the trunk
190 "physicalRib",
191 logicRib,
192 physRibCage,
193 false,
194 0, true);
195
196 physRib12 = new G4PVPlacement(0,G4ThreeVector(0.0, 0.0, (-thickness/2. + 0.8*cm + 30.8 *cm)),
197 // with respect to the trunk
198 "physicalRib",
199 logicRib,
200 physRibCage,
201 false,
202 0, true);
203
204 // Sensitive Body Part
205 if (sensitivity==true)
206 {
207 G4SDManager* SDman = G4SDManager::GetSDMpointer();
208 // logicRibCage->SetSensitiveDetector( SDman->FindSensitiveDetector("BodyPartSD") );
209 logicRib->SetSensitiveDetector( SDman->FindSensitiveDetector("BodyPartSD") );
210 }
211
212 // Visualization Attributes
213 logicRibCage -> SetVisAttributes(G4VisAttributes::Invisible);
214
215 //G4VisAttributes* RibCageVisAtt = new G4VisAttributes(G4Colour(0.46,0.53,0.6));
216
217 G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour();
218 G4Colour colour = colourPointer -> GetColour(colourName);
219 G4VisAttributes* RibCageVisAtt = new G4VisAttributes(colour);
220 RibCageVisAtt->SetForceSolid(wireFrame);
221 logicRib->SetVisAttributes(RibCageVisAtt);
222
223 G4cout << "RibCage created !!!!!!" << G4endl;
224 // Testing Pelvis Volume
225 G4double RibCageVol = logicRib->GetSolid()->GetCubicVolume();
226 G4cout << "Volume of RibCage = " << ((RibCageVol)*12.)/cm3 << " cm^3" << G4endl;
227
228 // Testing RibCage Material
229 G4String RibCageMat = logicRib->GetMaterial()->GetName();
230 G4cout << "Material of RibCage = " << RibCageMat << G4endl;
231
232 // Testing Density
233 G4double RibCageDensity = logicRib->GetMaterial()->GetDensity();
234 G4cout << "Density of Material = " << RibCageDensity*cm3/g << " g/cm^3" << G4endl;
235
236 // Testing Mass
237 G4double RibCageMass = (RibCageVol)* RibCageDensity * 12;// 12 is the total number of ribs;
238 G4cout << "Mass of RibCage = " << (RibCageMass)/gram << " g" << G4endl;
239
240 return physRibCage;
241}
Note: See TracBrowser for help on using the repository browser.