Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4MIRDRightLegBone.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 //
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 
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 
59 
61  const G4String& colourName, G4bool wireFrame,G4bool)
62 {
63 
65 
66  G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
67 
68  G4Material* skeleton = material -> GetMaterial("skeleton");
69 
70  delete material;
71 
72  G4double dz = 79.8 * cm;
73  G4double rmin1 = 0.0 * cm;
74  G4double rmin2 = 0.0 * cm;
75  G4double rmax1 = 1. * cm;
76  G4double rmax2 = 3.5 * cm;
77  G4double startphi = 0. * degree;
78  G4double deltaphi = 360. * degree;
79 
80  G4Cons* leg_bone = new G4Cons("OneRightLegBone",
81  rmin1, rmax1,
82  rmin2, rmax2, dz/2.,
83  startphi, deltaphi);
84 
85  //G4RotationMatrix* rm_relative = new G4RotationMatrix();
86  //rm_relative -> rotateY(-12.5 * degree);
87 
88  // G4UnionSolid* legs_bones = new G4UnionSolid("RightLegBone",
89  // leg_bone, leg_bone,
90  // 0,
91  // G4ThreeVector(20.* cm, 0.0,0. * cm));
92 
93 
94  G4LogicalVolume* logicRightLegBone = new G4LogicalVolume(leg_bone, skeleton,"logical" + volumeName,
95  0, 0, 0);
96 
97 
98  // Define rotation and position here!
99  G4VPhysicalVolume* physRightLegBone = new G4PVPlacement(0,
100  G4ThreeVector(0.0 * cm, 0.0, 0.1*cm),
101  "physicalRightLegBone",
102  logicRightLegBone,
103  mother,
104  false,
105  0, true);
106 
107 
108 
109  // Visualization Attributes
110  G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour();
111  G4Colour colour = colourPointer -> GetColour(colourName);
112  G4VisAttributes* RightLegBoneVisAtt = new G4VisAttributes(colour);
113 
114  RightLegBoneVisAtt->SetForceSolid(wireFrame);
115  logicRightLegBone->SetVisAttributes(RightLegBoneVisAtt);
116 
117  G4cout << "RightLegBone created !!!!!!" << G4endl;
118 
119  // Testing RightLegBone Volume
120  G4double RightLegBoneVol = logicRightLegBone->GetSolid()->GetCubicVolume();
121  G4cout << "Volume of RightLegBone = " << RightLegBoneVol/cm3 << " cm^3" << G4endl;
122 
123  // Testing RightLegBone Material
124  G4String RightLegBoneMat = logicRightLegBone->GetMaterial()->GetName();
125  G4cout << "Material of RightLegBone = " << RightLegBoneMat << G4endl;
126 
127  // Testing Density
128  G4double RightLegBoneDensity = logicRightLegBone->GetMaterial()->GetDensity();
129  G4cout << "Density of Material = " << RightLegBoneDensity*cm3/g << " g/cm^3" << G4endl;
130 
131  // Testing Mass
132  G4double RightLegBoneMass = (RightLegBoneVol)*RightLegBoneDensity;
133  G4cout << "Mass of RightLegBone = " << RightLegBoneMass/gram << " g" << G4endl;
134 
135 
136  return physRightLegBone;
137 }
CLHEP::Hep3Vector G4ThreeVector
#define G4endl
Definition: G4ios.hh:61
G4VPhysicalVolume * Construct(const G4String &, G4VPhysicalVolume *, const G4String &, G4bool, G4bool)
G4Material * GetMaterial() const
const G4String & GetName() const
Definition: G4Material.hh:179
static constexpr double g
Definition: G4SIunits.hh:183
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
Definition: G4Cons.hh:83
void SetVisAttributes(const G4VisAttributes *pVA)
void SetForceSolid(G4bool=true)
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:188
static constexpr double degree
Definition: G4SIunits.hh:144
static constexpr double cm
Definition: G4SIunits.hh:119
G4GLOB_DLL std::ostream G4cout
G4VSolid * GetSolid() const
static constexpr double cm3
Definition: G4SIunits.hh:121
const G4String & GetName() const
static constexpr double gram
Definition: G4SIunits.hh:178
G4double GetDensity() const
Definition: G4Material.hh:181