Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4MIRDLeftArmBone.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 "G4MIRDLeftArmBone.hh"
36 #include "globals.hh"
37 #include "G4SystemOfUnits.hh"
38 #include "G4SDManager.hh"
39 #include "G4VisAttributes.hh"
41 #include "G4EllipticalTube.hh"
42 #include "G4RotationMatrix.hh"
43 #include "G4ThreeVector.hh"
44 #include "G4VPhysicalVolume.hh"
45 #include "G4PVPlacement.hh"
46 #include "G4UnionSolid.hh"
47 #include "G4EllipticalCone.hh"
48 #include "G4HumanPhantomColour.hh"
49 
51 {
52 }
53 
55 {
56 }
57 
59  const G4String& colourName, G4bool wireFrame, G4bool sensitivity)
60 {
61  // Remind! the elliptical cone gives problems! Intersections of volumes,
62  // wrong calculation of the volume!
63 
65 
66  G4cout << "Construct " << volumeName <<G4endl;
67 
68  G4Material* skeleton = material -> GetMaterial("skeleton");
69 
70  delete material;
71 
72  G4double dx = 1.4 * cm;//a
73  G4double dy = 2.7 * cm;//b
74  // G4double dz= 46. * cm;//z0
75 
76  //G4EllipticalCone* arm = new G4EllipticalCone("OneLeftArmBone",dx/2.,dy/2.,dz, 34.5 *cm);
77  G4EllipticalTube* leftArm = new G4EllipticalTube("OneLeftArmBone",dx,dy,34.5 *cm);
78 
79  G4LogicalVolume* logicLeftArmBone = new G4LogicalVolume(leftArm,
80  skeleton,
81  "logical" + volumeName,
82  0, 0,0);
83 
84  G4RotationMatrix* matrix = new G4RotationMatrix();
85  matrix -> rotateX(180. * degree);
86 
87  G4VPhysicalVolume* physLeftArmBone = new G4PVPlacement(matrix,
88  G4ThreeVector(18.4 * cm, 0.0, -0.5*cm),
89  //-x0
90  "physicalLeftArmBone",
91  logicLeftArmBone,
92  mother,
93  false,0,true);
94 
95 
96  // Sensitive Body Part
97  if (sensitivity==true)
98  {
100  logicLeftArmBone->SetSensitiveDetector( SDman->FindSensitiveDetector("BodyPartSD") );
101  }
102 
103 
104  // Visualization Attributes
105 
106  G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour();
107  G4Colour colour = colourPointer -> GetColour(colourName);
108  G4VisAttributes* LeftArmBoneVisAtt = new G4VisAttributes(colour);
109  LeftArmBoneVisAtt->SetForceSolid(wireFrame);
110  logicLeftArmBone->SetVisAttributes(LeftArmBoneVisAtt);
111 
112  G4cout << "LeftArmBone created !!!!!!" << G4endl;
113 
114  // Testing LeftArmBone Volume
115  G4double LeftArmBoneVol = logicLeftArmBone->GetSolid()->GetCubicVolume();
116  G4cout << "Volume of LeftArmBone = " << LeftArmBoneVol/cm3 << " cm^3" << G4endl;
117 
118  // Testing LeftArmBone Material
119  G4String LeftArmBoneMat = logicLeftArmBone->GetMaterial()->GetName();
120  G4cout << "Material of LeftArmBone = " << LeftArmBoneMat << G4endl;
121 
122  // Testing Density
123  G4double LeftArmBoneDensity = logicLeftArmBone->GetMaterial()->GetDensity();
124  G4cout << "Density of Material = " << LeftArmBoneDensity*cm3/g << " g/cm^3" << G4endl;
125 
126  // Testing Mass
127  G4double LeftArmBoneMass = (LeftArmBoneVol)*LeftArmBoneDensity;
128  G4cout << "Mass of LeftArmBone = " << LeftArmBoneMass/gram << " g" << G4endl;
129 
130  return physLeftArmBone;
131 }