Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4MIRDRightKidney Class Reference

#include <G4MIRDRightKidney.hh>

Inheritance diagram for G4MIRDRightKidney:
Collaboration diagram for G4MIRDRightKidney:

Public Member Functions

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

Constructor & Destructor Documentation

G4MIRDRightKidney::G4MIRDRightKidney ( )

Definition at line 55 of file G4MIRDRightKidney.cc.

56 {
57 }
G4MIRDRightKidney::~G4MIRDRightKidney ( )

Definition at line 59 of file G4MIRDRightKidney.cc.

60 {
61 }

Member Function Documentation

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

Implements G4VOrgan.

Definition at line 64 of file G4MIRDRightKidney.cc.

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 ax= 4.5 *cm; //a
76  G4double by= 1.5 *cm; //b
77  G4double cz= 5.5 *cm; //c
78 
79  G4VSolid* oneRightKidney = new G4Ellipsoid("OneRightKidney",ax, by, cz);
80 
81  G4double xx = 6. * cm;
82  G4double yy = 12.00*cm;
83  G4double zz = 12.00*cm;
84  G4VSolid* subtrRightKidney = new G4Box("SubtrRightKidney",xx/2., yy/2., zz/2.);
85 
86  G4SubtractionSolid* kidney = new G4SubtractionSolid("RightKidney",
87  oneRightKidney,
88  subtrRightKidney,
89  0,
90  G4ThreeVector(6. *cm, // x0
91  0.0 *cm,
92  0.0 * cm));
93 
94  G4LogicalVolume* logicRightKidney = new G4LogicalVolume(kidney,
95  soft,
96  "logical" + volumeName,
97  0, 0, 0);
98 
99  G4VPhysicalVolume* physRightKidney = new G4PVPlacement(0 ,G4ThreeVector(-6.*cm, // xo
100  6. *cm, //yo
101  -2.50 *cm),//zo
102  "physicalRightKidney", logicRightKidney,
103  mother,
104  false,
105  0, true);
106 
107 
108 
109  // Visualization Attributes
110  //G4VisAttributes* RightKidneyVisAtt = new G4VisAttributes(G4Colour(0.72,0.52,0.04));
111  G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour();
112  G4Colour colour = colourPointer -> GetColour(colourName);
113  G4VisAttributes* RightKidneyVisAtt = new G4VisAttributes(colour);
114  RightKidneyVisAtt->SetForceSolid(wireFrame);
115  logicRightKidney->SetVisAttributes(RightKidneyVisAtt);
116 
117  G4cout << "RightKidney created !!!!!!" << G4endl;
118 
119  // Testing RightKidney Volume
120  G4double RightKidneyVol = logicRightKidney->GetSolid()->GetCubicVolume();
121  G4cout << "Volume of RightKidney = " << RightKidneyVol/cm3 << " cm^3" << G4endl;
122 
123  // Testing RightKidney Material
124  G4String RightKidneyMat = logicRightKidney->GetMaterial()->GetName();
125  G4cout << "Material of RightKidney = " << RightKidneyMat << G4endl;
126 
127  // Testing Density
128  G4double RightKidneyDensity = logicRightKidney->GetMaterial()->GetDensity();
129  G4cout << "Density of Material = " << RightKidneyDensity*cm3/g << " g/cm^3" << G4endl;
130 
131  // Testing Mass
132  G4double RightKidneyMass = (RightKidneyVol)*RightKidneyDensity;
133  G4cout << "Mass of RightKidney = " << RightKidneyMass/gram << " g" << G4endl;
134 
135 
136  return physRightKidney;
137 }
G4Material * GetMaterial() const
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Box.hh:64
const G4String & GetName() const
Definition: G4Material.hh:178
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:189
G4double GetDensity() const
Definition: G4Material.hh:180
G4VSolid * GetSolid() const
void SetForceSolid(G4bool=true)
static constexpr double gram
Definition: G4SIunits.hh:178
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
static constexpr double cm
Definition: G4SIunits.hh:119
static constexpr double cm3
Definition: G4SIunits.hh:121
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void SetVisAttributes(const G4VisAttributes *pVA)

Here is the call graph for this function:


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