Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4MIRDRightAdrenal.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 "G4MIRDRightAdrenal.hh"
36 
37 #include "globals.hh"
38 #include "G4SystemOfUnits.hh"
39 #include "G4SDManager.hh"
40 #include "G4VisAttributes.hh"
42 #include "G4SDManager.hh"
43 #include "G4PVPlacement.hh"
44 #include "G4SubtractionSolid.hh"
45 #include "G4Ellipsoid.hh"
46 #include "G4ThreeVector.hh"
47 #include "G4VPhysicalVolume.hh"
48 #include "G4RotationMatrix.hh"
49 #include "G4Material.hh"
50 #include "G4EllipticalTube.hh"
51 #include "G4Box.hh"
52 #include "G4UnionSolid.hh"
53 #include "G4HumanPhantomColour.hh"
54 
56 {
57 }
58 
60 {
61 }
62 
63 
65  const G4String& colourName, G4bool wireFrame, G4bool)
66 {
67  G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
68 
69 
71  G4Material* soft = material -> GetMaterial("soft_tissue");
72  delete material;
73 
74  G4double ax= 1.5 *cm; //a
75  G4double by= 0.5 *cm; //b
76  G4double cz= 5.0 *cm; //c
77 
78  G4VSolid* rightAdrenal = new G4Ellipsoid("OneRightAdrenal",ax, by, cz, 0. *cm, cz);
79 
80 
81  G4LogicalVolume* logicRightAdrenal = new G4LogicalVolume(rightAdrenal,
82  soft,
83  "logical" + volumeName,
84  0, 0, 0);
85 
86  G4VPhysicalVolume* physRightAdrenal = new G4PVPlacement(0 ,G4ThreeVector(-4.5*cm, // xo
87  6.5 *cm, //yo
88  3. *cm),//zo
89  "physicalRightAdrenal", logicRightAdrenal,
90  mother,
91  false,
92  0, true);
93 
94  // Visualization Attributes
95  // G4VisAttributes* RightAdrenalVisAtt = new G4VisAttributes(G4Colour(0.72,0.52,0.04));
96  G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour();
97  G4Colour colour = colourPointer -> GetColour(colourName);
98  G4VisAttributes* RightAdrenalVisAtt = new G4VisAttributes(colour);
99  RightAdrenalVisAtt->SetForceSolid(wireFrame);
100  logicRightAdrenal->SetVisAttributes(RightAdrenalVisAtt);
101 
102  G4cout << "Right RightAdrenal created !!!!!!" << G4endl;
103 
104  // Testing RightAdrenal Volume
105  G4double RightAdrenalVol = logicRightAdrenal->GetSolid()->GetCubicVolume();
106  G4cout << "Volume of RightAdrenal = " << RightAdrenalVol/cm3 << " cm^3" << G4endl;
107 
108  // Testing RightAdrenal Material
109  G4String RightAdrenalMat = logicRightAdrenal->GetMaterial()->GetName();
110  G4cout << "Material of RightAdrenal = " << RightAdrenalMat << G4endl;
111 
112  // Testing Density
113  G4double RightAdrenalDensity = logicRightAdrenal->GetMaterial()->GetDensity();
114  G4cout << "Density of Material = " << RightAdrenalDensity*cm3/g << " g/cm^3" << G4endl;
115 
116  // Testing Mass
117  G4double RightAdrenalMass = (RightAdrenalVol)*RightAdrenalDensity;
118  G4cout << "Mass of RightAdrenal = " << RightAdrenalMass/gram << " g" << G4endl;
119 
120 
121  return physRightAdrenal;
122 }
CLHEP::Hep3Vector G4ThreeVector
#define G4endl
Definition: G4ios.hh:61
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
G4VPhysicalVolume * Construct(const G4String &, G4VPhysicalVolume *, const G4String &, G4bool, G4bool)
void SetVisAttributes(const G4VisAttributes *pVA)
void SetForceSolid(G4bool=true)
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:188
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