Geant4  10.00.p01
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 , M. G. Pia, INFN Genova and F. Ambroglini INFN Perugia, 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"
50 #include "G4HumanPhantomColour.hh"
51 
52 
54 {
55 }
56 
58 {
59 }
60 
61 
63  G4VPhysicalVolume* mother,
64  const G4String& colourName,
65  G4bool wireFrame,G4bool)
66 
67 {
68  G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
70  G4Material* adipose = material -> GetMaterial("adipose");
71  G4Material* adipose_glandular = material -> GetMaterial("adipose_glandular");
72  delete material;
73 
74  G4double rmin = 0.* cm;
75  G4double rmax = 6. * cm;
76  G4double zz = 6. * cm;
77  G4double startPhi = 0.* degree;
78  G4double spanningPhi = 180. * degree;
79 
80  G4Tubs* breast = new G4Tubs("out_breast", rmin, rmax,
81  zz/2., startPhi, spanningPhi);
82 
83  G4LogicalVolume* breast_log = new G4LogicalVolume(breast,
84  adipose,
85  "logicalOut"+ volumeName,
86  0, 0, 0);
87  rmax = 5.5 *cm;
88  zz = 5. *cm;
89  G4Tubs* innerBreast = new G4Tubs("inner_breast",
90  rmin, rmax,
91  zz/2., startPhi, spanningPhi);
92 
93  G4LogicalVolume* innerBreast_log = new G4LogicalVolume(innerBreast,
94  adipose_glandular,
95  "logical"+ volumeName,
96  0, 0, 0);
97 
98  G4RotationMatrix* matrix = new G4RotationMatrix();
99  matrix -> rotateX(0.* degree);
100  matrix -> rotateY(0.* degree);
101  matrix -> rotateZ(-18. *degree);
102 
103 
104  G4VPhysicalVolume* physBreast = new G4PVPlacement(matrix,G4ThreeVector(-10.*cm, 8.7 *cm, 52.* cm),
105  "physicalVoxelRightBreast",
106  breast_log,
107  mother,
108  false,
109  0, true);
110 
111  G4VPhysicalVolume* physInnerBreast = new G4PVPlacement(0,G4ThreeVector(),
112  "RightBreast",
113  innerBreast_log,
114  physBreast,
115  false,
116  0, true);
117 
118  G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour();
119  G4Colour colour = colourPointer -> GetColour(colourName);
120  delete colourPointer;
121 
122  G4VisAttributes* BreastVisAtt = new G4VisAttributes(colour);
123  BreastVisAtt -> SetForceSolid(wireFrame);
124  breast_log -> SetVisAttributes(BreastVisAtt);
125 
126  G4VisAttributes* innerBreastVisAtt = new G4VisAttributes(colour);
127  innerBreastVisAtt -> SetForceSolid(false);
128  innerBreast_log -> SetVisAttributes(innerBreastVisAtt);
129 
130 
131  G4cout << "Voxel Right Breast created !!!!!! This model must be refined!" << G4endl;
132 
133 
134  return physInnerBreast;
135 }
static const double cm
Definition: G4SIunits.hh:106
CLHEP::Hep3Vector G4ThreeVector
CLHEP::HepRotation G4RotationMatrix
Definition: G4Tubs.hh:84
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
bool G4bool
Definition: G4Types.hh:79
G4VPhysicalVolume * Construct(const G4String &, G4VPhysicalVolume *, const G4String &, G4bool, G4bool)
static const double degree
Definition: G4SIunits.hh:125
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76