Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4MIRDRightLegBone.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 
35 #include "G4MIRDRightLegBone.hh"
36 
37 #include "globals.hh"
38 #include "G4SystemOfUnits.hh"
39 #include "G4SDManager.hh"
40 #include "G4VisAttributes.hh"
42 #include "G4EllipticalTube.hh"
43 #include "G4RotationMatrix.hh"
44 #include "G4ThreeVector.hh"
45 #include "G4VPhysicalVolume.hh"
46 #include "G4PVPlacement.hh"
47 #include "G4Cons.hh"
48 #include "G4UnionSolid.hh"
49 #include "G4HumanPhantomColour.hh"
50 
52 {
53 }
54 
56 {
57 }
58 
60  const G4String& colourName, G4bool wireFrame,G4bool sensitivity)
61 {
62 
64 
65  G4cout << "Construct " << volumeName << G4endl;
66 
67  G4Material* skeleton = material -> GetMaterial("skeleton");
68 
69  delete material;
70 
71  G4double dz = 79.8 * cm;
72  G4double rmin1 = 0.0 * cm;
73  G4double rmin2 = 0.0 * cm;
74  G4double rmax1 = 1. * cm;
75  G4double rmax2 = 3.5 * cm;
76  G4double startphi = 0. * degree;
77  G4double deltaphi = 360. * degree;
78 
79  G4Cons* leg_bone = new G4Cons("OneRightLegBone",
80  rmin1, rmax1,
81  rmin2, rmax2, dz/2.,
82  startphi, deltaphi);
83 
84  //G4RotationMatrix* rm_relative = new G4RotationMatrix();
85  //rm_relative -> rotateY(-12.5 * degree);
86 
87  // G4UnionSolid* legs_bones = new G4UnionSolid("RightLegBone",
88  // leg_bone, leg_bone,
89  // 0,
90  // G4ThreeVector(20.* cm, 0.0,0. * cm));
91 
92 
93  G4LogicalVolume* logicRightLegBone = new G4LogicalVolume(leg_bone, skeleton,"logical" + volumeName,
94  0, 0, 0);
95 
96 
97  // Define rotation and position here!
98  G4VPhysicalVolume* physRightLegBone = new G4PVPlacement(0,
99  G4ThreeVector(0.0 * cm, 0.0, 0.1*cm),
100  "physicalRightLegBone",
101  logicRightLegBone,
102  mother,
103  false,
104  0, true);
105 
106  // Sensitive Body Part
107  if (sensitivity==true)
108  {
110  logicRightLegBone->SetSensitiveDetector( SDman->FindSensitiveDetector("BodyPartSD") );
111  }
112 
113  // Visualization Attributes
114  G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour();
115  G4Colour colour = colourPointer -> GetColour(colourName);
116  G4VisAttributes* RightLegBoneVisAtt = new G4VisAttributes(colour);
117 
118  RightLegBoneVisAtt->SetForceSolid(wireFrame);
119  logicRightLegBone->SetVisAttributes(RightLegBoneVisAtt);
120 
121  G4cout << "RightLegBone created !!!!!!" << G4endl;
122 
123  // Testing RightLegBone Volume
124  G4double RightLegBoneVol = logicRightLegBone->GetSolid()->GetCubicVolume();
125  G4cout << "Volume of RightLegBone = " << RightLegBoneVol/cm3 << " cm^3" << G4endl;
126 
127  // Testing RightLegBone Material
128  G4String RightLegBoneMat = logicRightLegBone->GetMaterial()->GetName();
129  G4cout << "Material of RightLegBone = " << RightLegBoneMat << G4endl;
130 
131  // Testing Density
132  G4double RightLegBoneDensity = logicRightLegBone->GetMaterial()->GetDensity();
133  G4cout << "Density of Material = " << RightLegBoneDensity*cm3/g << " g/cm^3" << G4endl;
134 
135  // Testing Mass
136  G4double RightLegBoneMass = (RightLegBoneVol)*RightLegBoneDensity;
137  G4cout << "Mass of RightLegBone = " << RightLegBoneMass/gram << " g" << G4endl;
138 
139 
140  return physRightLegBone;
141 }