Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4MIRDLowerLargeIntestine Class Reference

#include <G4MIRDLowerLargeIntestine.hh>

Inheritance diagram for G4MIRDLowerLargeIntestine:
G4VOrgan

Public Member Functions

 G4MIRDLowerLargeIntestine ()
 
 ~G4MIRDLowerLargeIntestine ()
 
G4VPhysicalVolumeConstruct (const G4String &, G4VPhysicalVolume *, const G4String &, G4bool, G4bool)
 
- Public Member Functions inherited from G4VOrgan
 G4VOrgan ()
 
virtual ~G4VOrgan ()
 

Detailed Description

Definition at line 42 of file G4MIRDLowerLargeIntestine.hh.

Constructor & Destructor Documentation

G4MIRDLowerLargeIntestine::G4MIRDLowerLargeIntestine ( )

Definition at line 51 of file G4MIRDLowerLargeIntestine.cc.

52 {
53 }
G4MIRDLowerLargeIntestine::~G4MIRDLowerLargeIntestine ( )

Definition at line 55 of file G4MIRDLowerLargeIntestine.cc.

56 {
57 
58 }

Member Function Documentation

G4VPhysicalVolume * G4MIRDLowerLargeIntestine::Construct ( const G4String volumeName,
G4VPhysicalVolume mother,
const G4String colourName,
G4bool  wireFrame,
G4bool   
)
virtual

Implements G4VOrgan.

Definition at line 61 of file G4MIRDLowerLargeIntestine.cc.

References python.hepunit::cm, python.hepunit::cm3, python.hepunit::degree, g(), G4cout, G4endl, G4VSolid::GetCubicVolume(), G4Material::GetDensity(), G4LogicalVolume::GetMaterial(), G4VPhysicalVolume::GetName(), G4Material::GetName(), G4LogicalVolume::GetSolid(), python.hepunit::gram, eplot::material, G4VisAttributes::SetForceSolid(), and G4LogicalVolume::SetVisAttributes().

64 {
65  G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
66 
68  G4Material* soft = material -> GetMaterial("soft_tissue");
69  delete material;
70 
71  G4double dx = 1.88 * cm; //a
72  G4double dy = 2.13 *cm; //b
73  G4double dz = 7.64 *cm; //(z1-z2)/2
74 
75  G4EllipticalTube* DescendingColonLowerLargeIntestine = new G4EllipticalTube("DiscendingColon",dx, dy, dz);
76 
77 
78  G4double rmin= 0.0 *cm;
79  G4double rmax = 1.88 * cm;//a
80  G4double rtor= 5.72*cm; //R1
81  G4double startphi= 0. * degree;
82  G4double deltaphi= 90. * degree;
83 
84  G4Torus* SigmoidColonUpLowerLargeIntestine = new G4Torus("SigmoidColonUpLowerLargeIntestine",
85  rmin, rmax,rtor,
86  startphi, deltaphi);
87 
88  rtor = 3. * cm;//R2
89  G4VSolid* SigmoidColonDownLowerLargeIntestine = new G4Torus("SigmoidColonDownLowerLargeIntestine",
90  rmin, rmax,
91  rtor,startphi,deltaphi);
92 
93  G4RotationMatrix* relative_rm = new G4RotationMatrix();
94  relative_rm -> rotateY(180. * degree);
95  relative_rm -> rotateZ(90. * degree);
96 
97  G4UnionSolid* SigmoidColonLowerLargeIntestine = new G4UnionSolid( "SigmoidColonLowerLargeIntestine",
98  SigmoidColonUpLowerLargeIntestine,
99  SigmoidColonDownLowerLargeIntestine,
100  relative_rm,
101  G4ThreeVector(0.0,8.72*cm,0.0));
102  // R1 + R2
103 
104  G4RotationMatrix* relative_rm_2 = new G4RotationMatrix();
105  relative_rm_2 -> rotateX(90. * degree);
106 
107  G4UnionSolid* LowerLargeIntestine = new G4UnionSolid( "LowerLargeIntestine",
108  DescendingColonLowerLargeIntestine,
109  SigmoidColonLowerLargeIntestine,
110  relative_rm_2,
111  G4ThreeVector(-5.72*cm,0.0*cm, -7.64*cm)
112  ); // -rtor,0, -dz
113 
114 
115  G4LogicalVolume* logicLowerLargeIntestine = new G4LogicalVolume( LowerLargeIntestine, soft,
116  "logical" + volumeName,
117  0, 0, 0);
118 
119  G4VPhysicalVolume* physLowerLargeIntestine = new G4PVPlacement(0, // R1+ R2, -2.36 (y0), z0
120  G4ThreeVector(8.72*cm, -2.36*cm,-18.64 *cm),
121  "physicalLowerLargeIntestine",
122  logicLowerLargeIntestine,
123  mother,
124  false,
125  0, true);
126 
127 
128  // Visualization Attributes
129  //G4VisAttributes* LowerLargeIntestineVisAtt = new G4VisAttributes(G4Colour(1.0,1.0,0.0));
130  G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour();
131  G4Colour colour = colourPointer -> GetColour(colourName);
132  G4VisAttributes* LowerLargeIntestineVisAtt = new G4VisAttributes(colour);
133  LowerLargeIntestineVisAtt->SetForceSolid(wireFrame);
134  logicLowerLargeIntestine->SetVisAttributes(LowerLargeIntestineVisAtt);
135 
136  G4cout << "LowerLargeIntestine created !!!!!!" << G4endl;
137 
138  // Testing LowerLargeIntestine Volume
139  G4double LowerLargeIntestineVol = logicLowerLargeIntestine->GetSolid()->GetCubicVolume();
140  G4cout << "Volume of LowerLargeIntestine = " << LowerLargeIntestineVol/cm3 << " cm^3" << G4endl;
141 
142  // Testing LowerLargeIntestine Material
143  G4String LowerLargeIntestineMat = logicLowerLargeIntestine->GetMaterial()->GetName();
144  G4cout << "Material of LowerLargeIntestine = " << LowerLargeIntestineMat << G4endl;
145 
146  // Testing Density
147  G4double LowerLargeIntestineDensity = logicLowerLargeIntestine->GetMaterial()->GetDensity();
148  G4cout << "Density of Material = " << LowerLargeIntestineDensity*cm3/g << " g/cm^3" << G4endl;
149 
150  // Testing Mass
151  G4double LowerLargeIntestineMass = (LowerLargeIntestineVol)*LowerLargeIntestineDensity;
152  G4cout << "Mass of LowerLargeIntestine = " << LowerLargeIntestineMass/gram << " g" << G4endl;
153 
154 
155  return physLowerLargeIntestine;
156 }
CLHEP::Hep3Vector G4ThreeVector
CLHEP::HepRotation G4RotationMatrix
G4Material * GetMaterial() const
const G4String & GetName() const
Definition: G4Material.hh:176
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:188
G4double GetDensity() const
Definition: G4Material.hh:178
void SetForceSolid(G4bool)
string material
Definition: eplot.py:19
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
tuple degree
Definition: hepunit.py:69
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void SetVisAttributes(const G4VisAttributes *pVA)
G4VSolid * GetSolid() const

The documentation for this class was generated from the following files: