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

#include <G4MIRDSmallIntestine.hh>

Inheritance diagram for G4MIRDSmallIntestine:
G4VOrgan

Public Member Functions

 G4MIRDSmallIntestine ()
 
 ~G4MIRDSmallIntestine ()
 
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 G4MIRDSmallIntestine.hh.

Constructor & Destructor Documentation

G4MIRDSmallIntestine::G4MIRDSmallIntestine ( )

Definition at line 55 of file G4MIRDSmallIntestine.cc.

56 {
57 }
G4MIRDSmallIntestine::~G4MIRDSmallIntestine ( )

Definition at line 59 of file G4MIRDSmallIntestine.cc.

60 {
61 
62 }

Member Function Documentation

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

Implements G4VOrgan.

Definition at line 65 of file G4MIRDSmallIntestine.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, CLHEP::HepRotation::rotateX(), CLHEP::HepRotation::rotateY(), G4VisAttributes::SetForceSolid(), and G4LogicalVolume::SetVisAttributes().

68 {
69  G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
70 
72  G4Material* soft = material -> GetMaterial("soft_tissue");
73  delete material;
74 
75  G4double boxX = 11.*cm;
76  G4double boxY = 3.53*cm;
77  G4double boxZ = 5*cm;
78 
79 
80  G4Box* smallIntestineBox = new G4Box("smallIntestineBox",boxX,boxY,boxZ);
81 
82  G4double tubsRmin = 0*cm;
83  G4double tubsRmax = 11.*cm;
84  G4double tubsZ = 5*cm;
85  G4double tubsSphi = 0*degree;
86  G4double tubsDphi = 360*degree;
87 
88 
89  G4Tubs* smallIntestineTubs = new G4Tubs("smallIntestineTubs",tubsRmin,tubsRmax,tubsZ,tubsSphi,tubsDphi);
90 
91  //G4IntersectionSolid* SmallIntestine = new G4IntersectionSolid("SmallIntestine",smallIntestineTubs,smallIntestineBox,
92  G4IntersectionSolid* filledSmallIntestine1 = new G4IntersectionSolid("filledSmallIntestine1",smallIntestineTubs,smallIntestineBox,
93  0,G4ThreeVector(0*cm,-1.33*cm, 0*cm));
94 
95  G4IntersectionSolid* filledSmallIntestine = new G4IntersectionSolid("filledSmallIntestine",filledSmallIntestine1,smallIntestineTubs,
96  0,G4ThreeVector(0*cm,0.8*cm, 0*cm));
97 
98  G4double dx = 2.50*cm; // aU
99  G4double dy = 2.50*cm; //bU
100  G4double dz = 4.775*cm; //dzU
101 
102  G4VSolid* AscendingColonUpperLargeIntestine = new G4EllipticalTube("AscendingColon",dx, dy, dz);
103 
104  dx = 2.50 * cm;//bt
105  dy = 1.50 *cm;//ct
106  dz = 10.50* cm;//x1t
107 
108  G4VSolid* TraverseColonUpperLargeIntestine = new G4EllipticalTube("TraverseColon",dx, dy, dz);
109 
110  G4RotationMatrix* relative_rm = new G4RotationMatrix();
111  relative_rm -> rotateX(90. * degree);
112  //relative_rm -> rotateZ(180. * degree);
113  relative_rm -> rotateY(90. * degree);
114  G4UnionSolid* upperLargeIntestine = new G4UnionSolid("UpperLargeIntestine",
115  AscendingColonUpperLargeIntestine,
116  TraverseColonUpperLargeIntestine,
117  relative_rm,
118  G4ThreeVector(-8.0 *cm, 0.0*cm,6.275* cm)); //,0,dzU + ct transverse
119 
120  dx = 1.88 * cm; //a
121  dy = 2.13 *cm; //b
122  dz = 7.64 *cm; //(z1-z2)/2
123 
124  G4EllipticalTube* DescendingColonLowerLargeIntestine = new G4EllipticalTube("DiscendingColon",dx, dy, dz);
125 
126  G4UnionSolid* upperlowerLargeIntestine = new G4UnionSolid("UpperLowerLargeIntestine",
127  upperLargeIntestine,
128  DescendingColonLowerLargeIntestine,
129  0,
130  G4ThreeVector(-16.72*cm, 0.0*cm,-2.865* cm)); //,0,dzU + ct t
131 
132 
133  G4SubtractionSolid* SmallIntestine = new G4SubtractionSolid("SmallIntestine",
134  filledSmallIntestine,
135  upperlowerLargeIntestine,
136  0,
137  G4ThreeVector(8.0*cm,-0.3*cm,-2.775*cm));
138 
139 
140  G4LogicalVolume* logicSmallIntestine = new G4LogicalVolume( SmallIntestine,
141  soft,
142  "logical"+volumeName,
143  0, 0, 0);
145  rm->rotateX(180.*degree);
146  rm->rotateY(180.*degree);
147  G4VPhysicalVolume* physSmallIntestine = new G4PVPlacement(rm,
148  G4ThreeVector(0*cm, -2.66*cm, -13*cm), // Xcheck the spina position the correct placement shuod be this one
149  //G4ThreeVector(0*cm, -5.13*cm, -13*cm), // Xcheck the spina position the correct placement shuod be this one
150  //G4ThreeVector(0*cm, -6*cm, -13*cm),
151  "physical"+volumeName,
152  logicSmallIntestine,
153  mother,
154  false,
155  0, true);
156 
157 
158  // Visualization Attributes
159  //G4VisAttributes* SmallIntestineVisAtt = new G4VisAttributes(G4Colour(1.0,1.0,0.0));
160  G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour();
161  G4Colour colour = colourPointer -> GetColour(colourName);
162  G4VisAttributes* SmallIntestineVisAtt = new G4VisAttributes(colour);
163  SmallIntestineVisAtt->SetForceSolid(wireFrame);
164  logicSmallIntestine->SetVisAttributes(SmallIntestineVisAtt);
165 
166  G4cout << "SmallIntestine created !!!!!!" << G4endl;
167 
168  // Testing SmallIntestine Volume
169  G4double SmallIntestineVol = logicSmallIntestine->GetSolid()->GetCubicVolume();
170  G4cout << "Volume of SmallIntestine = " << SmallIntestineVol/cm3 << " cm^3" << G4endl;
171 
172  // Testing SmallIntestine Material
173  G4String SmallIntestineMat = logicSmallIntestine->GetMaterial()->GetName();
174  G4cout << "Material of SmallIntestine = " << SmallIntestineMat << G4endl;
175 
176  // Testing Density
177  G4double SmallIntestineDensity = logicSmallIntestine->GetMaterial()->GetDensity();
178  G4cout << "Density of Material = " << SmallIntestineDensity*cm3/g << " g/cm^3" << G4endl;
179 
180  // Testing Mass
181  G4double SmallIntestineMass = (SmallIntestineVol)*SmallIntestineDensity;
182  G4cout << "Mass of SmallIntestine = " << SmallIntestineMass/gram << " g" << G4endl;
183 
184 
185  return physSmallIntestine;
186 }
CLHEP::Hep3Vector G4ThreeVector
HepRotation & rotateX(double delta)
Definition: Rotation.cc:66
CLHEP::HepRotation G4RotationMatrix
G4Material * GetMaterial() const
Definition: G4Box.hh:63
const G4String & GetName() const
Definition: G4Material.hh:176
Definition: G4Tubs.hh:84
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:188
G4double GetDensity() const
Definition: G4Material.hh:178
HepRotation & rotateY(double delta)
Definition: Rotation.cc:79
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: