84 lines
3.3 KiB
C++
84 lines
3.3 KiB
C++
#include "G4HumanPhantomColour.h"
|
|
#include "G4HumanPhantomMaterial.h"
|
|
#include "G4MIRDRightKidney.h"
|
|
|
|
#include "G4Box.hh"
|
|
#include "G4Ellipsoid.hh"
|
|
#include "G4EllipticalTube.hh"
|
|
#include "G4Material.hh"
|
|
#include "G4PVPlacement.hh"
|
|
#include "G4RotationMatrix.hh"
|
|
#include "G4SDManager.hh"
|
|
#include "G4SubtractionSolid.hh"
|
|
#include "G4SystemOfUnits.hh"
|
|
#include "G4ThreeVector.hh"
|
|
#include "G4UnionSolid.hh"
|
|
#include "G4VPhysicalVolume.hh"
|
|
#include "G4VisAttributes.hh"
|
|
#include "globals.hh"
|
|
|
|
G4MIRDRightKidney::G4MIRDRightKidney() {}
|
|
|
|
G4MIRDRightKidney::~G4MIRDRightKidney() {}
|
|
|
|
G4VPhysicalVolume* G4MIRDRightKidney::Construct(const G4String& volumeName, G4VPhysicalVolume* mother,
|
|
const G4String& colourName, G4bool wireFrame, G4bool) {
|
|
G4cout << "Construct " << volumeName << " with mother volume " << mother->GetName() << G4endl;
|
|
|
|
G4HumanPhantomMaterial* material = new G4HumanPhantomMaterial();
|
|
G4Material* soft = material->GetMaterial("soft_tissue");
|
|
delete material;
|
|
|
|
G4double ax = 4.5 * cm; // a
|
|
G4double by = 1.5 * cm; // b
|
|
G4double cz = 5.5 * cm; // c
|
|
|
|
G4VSolid* oneRightKidney = new G4Ellipsoid("OneRightKidney", ax, by, cz);
|
|
|
|
G4double xx = 6. * cm;
|
|
G4double yy = 12.00 * cm;
|
|
G4double zz = 12.00 * cm;
|
|
G4VSolid* subtrRightKidney = new G4Box("SubtrRightKidney", xx / 2., yy / 2., zz / 2.);
|
|
|
|
G4SubtractionSolid* kidney = new G4SubtractionSolid("RightKidney", oneRightKidney, subtrRightKidney, 0,
|
|
G4ThreeVector(6. * cm, // x0
|
|
0.0 * cm, 0.0 * cm));
|
|
|
|
G4LogicalVolume* logicRightKidney = new G4LogicalVolume(kidney, soft, "logical" + volumeName, 0, 0, 0);
|
|
|
|
G4VPhysicalVolume* physRightKidney =
|
|
new G4PVPlacement(0,
|
|
G4ThreeVector(-6. * cm, // xo
|
|
6. * cm, // yo
|
|
-2.50 * cm), // zo
|
|
"physicalRightKidney", logicRightKidney, mother, false, 0, true);
|
|
|
|
// Visualization Attributes
|
|
// G4VisAttributes* RightKidneyVisAtt = new G4VisAttributes(G4Colour(0.72,0.52,0.04));
|
|
G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour();
|
|
G4Colour colour = colourPointer->GetColour(colourName);
|
|
G4VisAttributes* RightKidneyVisAtt = new G4VisAttributes(colour);
|
|
RightKidneyVisAtt->SetForceSolid(wireFrame);
|
|
logicRightKidney->SetVisAttributes(RightKidneyVisAtt);
|
|
|
|
G4cout << "RightKidney created !!!!!!" << G4endl;
|
|
|
|
// Testing RightKidney Volume
|
|
G4double RightKidneyVol = logicRightKidney->GetSolid()->GetCubicVolume();
|
|
G4cout << "Volume of RightKidney = " << RightKidneyVol / cm3 << " cm^3" << G4endl;
|
|
|
|
// Testing RightKidney Material
|
|
G4String RightKidneyMat = logicRightKidney->GetMaterial()->GetName();
|
|
G4cout << "Material of RightKidney = " << RightKidneyMat << G4endl;
|
|
|
|
// Testing Density
|
|
G4double RightKidneyDensity = logicRightKidney->GetMaterial()->GetDensity();
|
|
G4cout << "Density of Material = " << RightKidneyDensity * cm3 / g << " g/cm^3" << G4endl;
|
|
|
|
// Testing Mass
|
|
G4double RightKidneyMass = (RightKidneyVol)*RightKidneyDensity;
|
|
G4cout << "Mass of RightKidney = " << RightKidneyMass / gram << " g" << G4endl;
|
|
|
|
return physRightKidney;
|
|
}
|