G4-DESCSS/src/G4MIRDSpleen.cpp

70 lines
2.5 KiB
C++

#include "G4HumanPhantomColour.h"
#include "G4HumanPhantomMaterial.h"
#include "G4MIRDSpleen.h"
#include "G4Ellipsoid.hh"
#include "G4LogicalVolume.hh"
#include "G4Material.hh"
#include "G4PVPlacement.hh"
#include "G4RotationMatrix.hh"
#include "G4SDManager.hh"
#include "G4SystemOfUnits.hh"
#include "G4ThreeVector.hh"
#include "G4VPhysicalVolume.hh"
#include "G4VisAttributes.hh"
#include "globals.hh"
G4MIRDSpleen::G4MIRDSpleen() {}
G4MIRDSpleen::~G4MIRDSpleen() {}
G4VPhysicalVolume* G4MIRDSpleen::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 = 3.5 * cm;
G4double by = 2. * cm;
G4double cz = 6. * cm;
G4Ellipsoid* spleen = new G4Ellipsoid("spleen", ax, by, cz);
G4LogicalVolume* logicSpleen = new G4LogicalVolume(spleen, soft, "logical" + volumeName, 0, 0, 0);
// Define rotation and position here!
G4VPhysicalVolume* physSpleen =
new G4PVPlacement(0, G4ThreeVector(11. * cm, 3. * cm, 2. * cm), // ztrans = half trunk lenght - z0
"physicalSpleen", logicSpleen, mother, false, 0, true);
// Visualization Attributes
// G4VisAttributes* SpleenVisAtt = new G4VisAttributes(G4Colour(0.41,0.41,0.41));
G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour();
G4Colour colour = colourPointer->GetColour(colourName);
G4VisAttributes* SpleenVisAtt = new G4VisAttributes(colour);
SpleenVisAtt->SetForceSolid(wireFrame);
logicSpleen->SetVisAttributes(SpleenVisAtt);
G4cout << "Spleen created !!!!!!" << G4endl;
// Testing Spleen Volume
G4double SpleenVol = logicSpleen->GetSolid()->GetCubicVolume();
G4cout << "Volume of Spleen = " << SpleenVol / cm3 << " cm^3" << G4endl;
// Testing Spleen Material
G4String SpleenMat = logicSpleen->GetMaterial()->GetName();
G4cout << "Material of Spleen = " << SpleenMat << G4endl;
// Testing Density
G4double SpleenDensity = logicSpleen->GetMaterial()->GetDensity();
G4cout << "Density of Material = " << SpleenDensity * cm3 / g << " g/cm^3" << G4endl;
// Testing Mass
G4double SpleenMass = (SpleenVol)*SpleenDensity;
G4cout << "Mass of Spleen = " << SpleenMass / gram << " g" << G4endl;
return physSpleen;
}