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

#include <G4MIRDLeftBreast.hh>

Inheritance diagram for G4MIRDLeftBreast:
G4VOrgan

Public Member Functions

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

Detailed Description

Definition at line 43 of file G4MIRDLeftBreast.hh.

Constructor & Destructor Documentation

G4MIRDLeftBreast::G4MIRDLeftBreast ( )

Definition at line 53 of file G4MIRDLeftBreast.cc.

54 {
55 }
G4MIRDLeftBreast::~G4MIRDLeftBreast ( )

Definition at line 57 of file G4MIRDLeftBreast.cc.

58 {
59 }

Member Function Documentation

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

Implements G4VOrgan.

Definition at line 62 of file G4MIRDLeftBreast.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(), CLHEP::HepRotation::rotateZ(), 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 ax= 4.95* cm;
72  G4double by= 4.35* cm;
73  G4double cz= 4.15*cm;
74 
75  G4Ellipsoid* oneLeftBreast = new G4Ellipsoid("OneLeftBreast",
76  ax, by, cz);
77 
78  G4double dx= 20.* cm;
79  G4double dy= 10.* cm;
80  G4double dz= 35.*cm;
81 
82  G4EllipticalTube* Trunk = new G4EllipticalTube("Trunk",dx, dy, dz );
83 
84 
85 
86  G4RotationMatrix* rm_relative = new G4RotationMatrix();
87  rm_relative -> rotateX(90. * degree);
88 
89  G4SubtractionSolid* breast = new G4SubtractionSolid("LeftBreast",
90  oneLeftBreast,
91  Trunk,
92  rm_relative,
93  G4ThreeVector(-10.*cm,
94  0.0*cm,
95  -8.66*cm));
96 
97 
98  G4LogicalVolume* logicLeftBreast = new G4LogicalVolume(breast, soft,
99  "logical" + volumeName, 0, 0,0);
100 
101 
102  // Define rotation and position here!
104  rm->rotateX(90.*degree);
105  rm->rotateY(0.*degree);
106  rm->rotateZ(16.*degree);
107  G4VPhysicalVolume* physLeftBreast = new G4PVPlacement(rm,G4ThreeVector(10*cm, 9.1 *cm, 52.* cm),
108  "physicalLeftBreast",
109  logicLeftBreast,
110  mother,
111  false,
112  0, true);
113 
114  // Visualization Attributes
115  //G4VisAttributes* LeftBreastVisAtt = new G4VisAttributes(G4Colour(1.0,0.41,0.71));
116  G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour();
117  G4Colour colour = colourPointer -> GetColour(colourName);
118  G4VisAttributes* LeftBreastVisAtt = new G4VisAttributes(colour);
119  LeftBreastVisAtt->SetForceSolid(wireFrame);
120 
121  logicLeftBreast->SetVisAttributes(LeftBreastVisAtt);
122 
123  G4cout << "LeftBreast created !!!!!!" << G4endl;
124 
125  // Testing LeftBreast Volume
126  G4double LeftBreastVol = logicLeftBreast->GetSolid()->GetCubicVolume();
127  G4cout << "Volume of LeftBreast = " << LeftBreastVol/cm3 << " cm^3" << G4endl;
128 
129  // Testing LeftBreast Material
130  G4String LeftBreastMat = logicLeftBreast->GetMaterial()->GetName();
131  G4cout << "Material of LeftBreast = " << LeftBreastMat << G4endl;
132 
133  // Testing Density
134  G4double LeftBreastDensity = logicLeftBreast->GetMaterial()->GetDensity();
135  G4cout << "Density of Material = " << LeftBreastDensity*cm3/g << " g/cm^3" << G4endl;
136 
137  // Testing Mass
138  G4double LeftBreastMass = (LeftBreastVol)*LeftBreastDensity;
139  G4cout << "Mass of LeftBreast = " << LeftBreastMass/gram << " g" << G4endl;
140 
141 
142  return physLeftBreast;
143 }
CLHEP::Hep3Vector G4ThreeVector
HepRotation & rotateX(double delta)
Definition: Rotation.cc:66
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
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
HepRotation & rotateZ(double delta)
Definition: Rotation.cc:92
#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: