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

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

update

File size: 6.4 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 "G4MIRDLowerLargeIntestine.hh"
35#include "globals.hh"
36#include "G4SDManager.hh"
37#include "G4VisAttributes.hh"
38#include "G4EllipticalTube.hh"
39#include "G4UnionSolid.hh"
40#include "G4RotationMatrix.hh"
41#include "G4ThreeVector.hh"
42#include "G4VPhysicalVolume.hh"
43#include "G4PVPlacement.hh"
44#include "G4LogicalVolume.hh"
45#include "G4Torus.hh"
46#include "G4HumanPhantomMaterial.hh"
47#include "G4HumanPhantomColour.hh"
48G4MIRDLowerLargeIntestine::G4MIRDLowerLargeIntestine()
49{
50}
51
52G4MIRDLowerLargeIntestine::~G4MIRDLowerLargeIntestine()
53{
54
55}
56
57G4VPhysicalVolume* G4MIRDLowerLargeIntestine::Construct(const G4String& volumeName,
58 G4VPhysicalVolume* mother,
59 const G4String& colourName, G4bool wireFrame,G4bool sensitivity)
60{
61 G4cout << "Construct "<< volumeName <<G4endl;
62
63 G4HumanPhantomMaterial* material = new G4HumanPhantomMaterial();
64 G4Material* soft = material -> GetMaterial("soft_tissue");
65 delete material;
66
67 G4double dx = 1.88 * cm; //a
68 G4double dy = 2.13 *cm; //b
69 G4double dz = 7.64 *cm; //(z1-z2)/2
70
71 G4EllipticalTube* DescendingColonLowerLargeIntestine = new G4EllipticalTube("DiscendingColon",dx, dy, dz);
72
73
74 G4double rmin= 0.0 *cm;
75 G4double rmax = 1.88 * cm;//a
76 G4double rtor= 5.72*cm; //R1
77 G4double startphi= 0. * degree;
78 G4double deltaphi= 90. * degree;
79
80 G4Torus* SigmoidColonUpLowerLargeIntestine = new G4Torus("SigmoidColonUpLowerLargeIntestine",
81 rmin, rmax,rtor,
82 startphi, deltaphi);
83
84 rtor = 3. * cm;//R2
85 G4VSolid* SigmoidColonDownLowerLargeIntestine = new G4Torus("SigmoidColonDownLowerLargeIntestine",
86 rmin, rmax,
87 rtor,startphi,deltaphi);
88
89 G4RotationMatrix* relative_rm = new G4RotationMatrix();
90 relative_rm -> rotateY(180. * degree);
91 relative_rm -> rotateZ(90. * degree);
92
93 G4UnionSolid* SigmoidColonLowerLargeIntestine = new G4UnionSolid( "SigmoidColonLowerLargeIntestine",
94 SigmoidColonUpLowerLargeIntestine,
95 SigmoidColonDownLowerLargeIntestine,
96 relative_rm,
97 G4ThreeVector(0.0,8.72*cm,0.0));
98 // R1 + R2
99
100 G4RotationMatrix* relative_rm_2 = new G4RotationMatrix();
101 relative_rm_2 -> rotateX(90. * degree);
102
103 G4UnionSolid* LowerLargeIntestine = new G4UnionSolid( "LowerLargeIntestine",
104 DescendingColonLowerLargeIntestine,
105 SigmoidColonLowerLargeIntestine,
106 relative_rm_2,
107 G4ThreeVector(-5.72*cm,0.0*cm, -7.64*cm)
108 ); // -rtor,0, -dz
109
110
111 G4LogicalVolume* logicLowerLargeIntestine = new G4LogicalVolume( LowerLargeIntestine, soft,
112 "logical" + volumeName,
113 0, 0, 0);
114
115 G4VPhysicalVolume* physLowerLargeIntestine = new G4PVPlacement(0, // R1+ R2, -2.36 (y0), z0
116 G4ThreeVector(8.72*cm, -2.36*cm,-18.64 *cm),
117 "physicalLowerLargeIntestine",
118 logicLowerLargeIntestine,
119 mother,
120 false,
121 0, true);
122 // Sensitive Body Part
123 if (sensitivity==true)
124 {
125 G4SDManager* SDman = G4SDManager::GetSDMpointer();
126 logicLowerLargeIntestine->SetSensitiveDetector( SDman->FindSensitiveDetector("BodyPartSD") );
127 }
128
129 // Visualization Attributes
130 //G4VisAttributes* LowerLargeIntestineVisAtt = new G4VisAttributes(G4Colour(1.0,1.0,0.0));
131 G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour();
132 G4Colour colour = colourPointer -> GetColour(colourName);
133 G4VisAttributes* LowerLargeIntestineVisAtt = new G4VisAttributes(colour);
134 LowerLargeIntestineVisAtt->SetForceSolid(wireFrame);
135 logicLowerLargeIntestine->SetVisAttributes(LowerLargeIntestineVisAtt);
136
137 G4cout << "LowerLargeIntestine created !!!!!!" << G4endl;
138
139 // Testing LowerLargeIntestine Volume
140 G4double LowerLargeIntestineVol = logicLowerLargeIntestine->GetSolid()->GetCubicVolume();
141 G4cout << "Volume of LowerLargeIntestine = " << LowerLargeIntestineVol/cm3 << " cm^3" << G4endl;
142
143 // Testing LowerLargeIntestine Material
144 G4String LowerLargeIntestineMat = logicLowerLargeIntestine->GetMaterial()->GetName();
145 G4cout << "Material of LowerLargeIntestine = " << LowerLargeIntestineMat << G4endl;
146
147 // Testing Density
148 G4double LowerLargeIntestineDensity = logicLowerLargeIntestine->GetMaterial()->GetDensity();
149 G4cout << "Density of Material = " << LowerLargeIntestineDensity*cm3/g << " g/cm^3" << G4endl;
150
151 // Testing Mass
152 G4double LowerLargeIntestineMass = (LowerLargeIntestineVol)*LowerLargeIntestineDensity;
153 G4cout << "Mass of LowerLargeIntestine = " << LowerLargeIntestineMass/gram << " g" << G4endl;
154
155
156 return physLowerLargeIntestine;
157}
Note: See TracBrowser for help on using the repository browser.