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

#include <G4MIRDLeftLeg.hh>

Inheritance diagram for G4MIRDLeftLeg:
G4VOrgan

Public Member Functions

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

Constructor & Destructor Documentation

G4MIRDLeftLeg::G4MIRDLeftLeg ( )

Definition at line 51 of file G4MIRDLeftLeg.cc.

52 {
53 }
G4MIRDLeftLeg::~G4MIRDLeftLeg ( )

Definition at line 55 of file G4MIRDLeftLeg.cc.

56 {
57 }

Member Function Documentation

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

Implements G4VOrgan.

Definition at line 60 of file G4MIRDLeftLeg.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().

62 {
63 
64  G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
65 
67  G4Material* soft = material -> GetMaterial("soft_tissue");
68 
69  G4double rmin1 = 0.* cm;
70  G4double rmin2 = 0.* cm;
71  G4double dz= 80.0 * cm;
72  G4double rmax1= 2.0 * cm;
73  G4double rmax2= 10. * cm;
74  G4double startphi= 0.* degree;
75  G4double deltaphi= 360. * degree;
76 
77  G4Cons* leg1 = new G4Cons("Leg1",
78  rmin1, rmax1,
79  rmin2, rmax2, dz/2.,
80  startphi, deltaphi);
81 
82  G4LogicalVolume* logicLeftLeg = new G4LogicalVolume(leg1,
83  soft,
84  "logical" + volumeName,
85  0, 0, 0);
86 
88  rm->rotateX(180.*degree);
89  rm->rotateY(180.*degree);
90  G4VPhysicalVolume* physLeftLeg = new G4PVPlacement(rm,
91  //G4ThreeVector(10. * cm, 0. * cm, -47. *cm), //FA
92  G4ThreeVector(10. * cm, 0. * cm, -40. *cm),
93  "physicalLeftLeg",
94  logicLeftLeg,
95  mother,
96  false,
97  0, true);
98 
99 
100  // Visualization Attributes
101  //G4VisAttributes* LeftLegVisAtt = new G4VisAttributes(G4Colour(0.94,0.5,0.5));
102  G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour();
103  G4Colour colour = colourPointer -> GetColour(colourName);
104  G4VisAttributes* LeftLegVisAtt = new G4VisAttributes(colour);
105  LeftLegVisAtt->SetForceSolid(wireFrame);
106  logicLeftLeg->SetVisAttributes(LeftLegVisAtt);
107 
108  G4cout << "LeftLeg created !!!!!!" << G4endl;
109 
110  // Testing LeftLeg Volume
111  G4double LeftLegVol = logicLeftLeg->GetSolid()->GetCubicVolume();
112  G4cout << "Volume of LeftLeg = " << LeftLegVol/cm3 << " cm^3" << G4endl;
113 
114  // Testing LeftLeg Material
115  G4String LeftLegMat = logicLeftLeg->GetMaterial()->GetName();
116  G4cout << "Material of LeftLeg = " << LeftLegMat << G4endl;
117 
118  // Testing Density
119  G4double LeftLegDensity = logicLeftLeg->GetMaterial()->GetDensity();
120  G4cout << "Density of Material = " << LeftLegDensity*cm3/g << " g/cm^3" << G4endl;
121 
122  // Testing Mass
123  G4double LeftLegMass = (LeftLegVol)*LeftLegDensity;
124  G4cout << "Mass of LeftLeg = " << LeftLegMass/gram << " g" << G4endl;
125 
126 
127  return physLeftLeg;
128 }
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
Definition: G4Cons.hh:82
#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: