ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4MIRDLeftLegBone.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4MIRDLeftLegBone.cc
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 // Previous authors: G. Guerrieri, S. Guatelli and M. G. Pia, INFN Genova, Italy
26 // Authors (since 2007): S. Guatelli, University of Wollongong, Australia
27 //
28 
29 #include "G4MIRDLeftLegBone.hh"
30 
31 #include "globals.hh"
32 #include "G4SystemOfUnits.hh"
33 #include "G4SDManager.hh"
34 #include "G4VisAttributes.hh"
36 #include "G4EllipticalTube.hh"
37 #include "G4RotationMatrix.hh"
38 #include "G4ThreeVector.hh"
39 #include "G4VPhysicalVolume.hh"
40 #include "G4PVPlacement.hh"
41 #include "G4Cons.hh"
42 #include "G4UnionSolid.hh"
43 #include "G4HumanPhantomColour.hh"
44 
46 {
47 }
48 
50 {
51 }
52 
53 
55  const G4String& colourName, G4bool wireFrame,G4bool)
56 {
57 
59 
60  G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
61 
62  G4Material* skeleton = material -> GetMaterial("skeleton");
63 
64  delete material;
65 
66  G4double dz = 79.8 * cm;
67  G4double rmin1 = 0.0 * cm;
68  G4double rmin2 = 0.0 * cm;
69  G4double rmax1 = 1. * cm;
70  G4double rmax2 = 3.5 * cm;
71  G4double startphi = 0. * degree;
72  G4double deltaphi = 360. * degree;
73 
74  G4Cons* leg_bone = new G4Cons("OneLeftLegBone",
75  rmin1, rmax1,
76  rmin2, rmax2, dz/2.,
77  startphi, deltaphi);
78 
79  G4LogicalVolume* logicLeftLegBone = new G4LogicalVolume(leg_bone, skeleton,"logical" + volumeName,
80  0, 0, 0);
81 
82 
83  // Define rotation and position here!
84  G4VPhysicalVolume* physLeftLegBone = new G4PVPlacement(0,
85  G4ThreeVector(0.0 * cm, 0.0, 0.1*cm),
86  "physicalLeftLegBone",
87  logicLeftLegBone,
88  mother,
89  false,
90  0, true);
91 
92 
93  // Visualization Attributes
94  //G4VisAttributes* LeftLegBoneVisAtt = new G4VisAttributes(G4Colour(0.46,0.53,0.6));
95  G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour();
96  G4Colour colour = colourPointer -> GetColour(colourName);
97  G4VisAttributes* LeftLegBoneVisAtt = new G4VisAttributes(colour);
98 
99  LeftLegBoneVisAtt->SetForceSolid(wireFrame);
100  logicLeftLegBone->SetVisAttributes(LeftLegBoneVisAtt);
101 
102  G4cout << "LeftLegBone created !!!!!!" << G4endl;
103 
104  // Testing LeftLegBone Volume
105  G4double LeftLegBoneVol = logicLeftLegBone->GetSolid()->GetCubicVolume();
106  G4cout << "Volume of LeftLegBone = " << LeftLegBoneVol/cm3 << " cm^3" << G4endl;
107 
108  // Testing LeftLegBone Material
109  G4String LeftLegBoneMat = logicLeftLegBone->GetMaterial()->GetName();
110  G4cout << "Material of LeftLegBone = " << LeftLegBoneMat << G4endl;
111 
112  // Testing Density
113  G4double LeftLegBoneDensity = logicLeftLegBone->GetMaterial()->GetDensity();
114  G4cout << "Density of Material = " << LeftLegBoneDensity*cm3/g << " g/cm^3" << G4endl;
115 
116  // Testing Mass
117  G4double LeftLegBoneMass = (LeftLegBoneVol)*LeftLegBoneDensity;
118  G4cout << "Mass of LeftLegBone = " << LeftLegBoneMass/gram << " g" << G4endl;
119 
120 
121  return physLeftLegBone;
122 }