File indexing completed on 2025-02-23 09:20:12
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031 #include "G4MIRDSpleen.hh"
0032
0033 #include "globals.hh"
0034 #include "G4SystemOfUnits.hh"
0035 #include "G4SDManager.hh"
0036 #include "G4VisAttributes.hh"
0037 #include "G4ThreeVector.hh"
0038 #include "G4VPhysicalVolume.hh"
0039 #include "G4RotationMatrix.hh"
0040 #include "G4Material.hh"
0041 #include "G4LogicalVolume.hh"
0042 #include "G4HumanPhantomMaterial.hh"
0043 #include "G4VPhysicalVolume.hh"
0044 #include "G4PVPlacement.hh"
0045 #include "G4Ellipsoid.hh"
0046 #include "G4HumanPhantomColour.hh"
0047
0048 G4VPhysicalVolume* G4MIRDSpleen::Construct(const G4String& volumeName,G4VPhysicalVolume* mother,
0049 const G4String& colourName, G4bool wireFrame, G4bool)
0050 {
0051 G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
0052
0053 auto* material = new G4HumanPhantomMaterial();
0054 auto* soft = material -> GetMaterial("soft_tissue");
0055 delete material;
0056
0057 G4double ax= 3.5 *cm;
0058 G4double by= 2. *cm;
0059 G4double cz= 6. * cm;
0060
0061 auto* spleen = new G4Ellipsoid("spleen", ax, by, cz);
0062
0063
0064 auto* logicSpleen = new G4LogicalVolume(spleen, soft,
0065 "logical" + volumeName,
0066 nullptr, nullptr, nullptr);
0067
0068
0069 G4VPhysicalVolume* physSpleen = new G4PVPlacement(nullptr,
0070 G4ThreeVector(11. *cm, 3. *cm, 2.*cm),
0071 "physicalSpleen",
0072 logicSpleen,
0073 mother,
0074 false,
0075 0, true);
0076
0077
0078
0079 auto* colourPointer = new G4HumanPhantomColour();
0080 G4Colour colour = colourPointer -> GetColour(colourName);
0081 auto* SpleenVisAtt = new G4VisAttributes(colour);
0082 SpleenVisAtt->SetForceSolid(wireFrame);
0083 logicSpleen->SetVisAttributes(SpleenVisAtt);
0084
0085 G4cout << "Spleen created !!!!!!" << G4endl;
0086
0087
0088 G4double SpleenVol = logicSpleen->GetSolid()->GetCubicVolume();
0089 G4cout << "Volume of Spleen = " << SpleenVol/cm3 << " cm^3" << G4endl;
0090
0091
0092 G4String SpleenMat = logicSpleen->GetMaterial()->GetName();
0093 G4cout << "Material of Spleen = " << SpleenMat << G4endl;
0094
0095
0096 G4double SpleenDensity = logicSpleen->GetMaterial()->GetDensity();
0097 G4cout << "Density of Material = " << SpleenDensity*cm3/g << " g/cm^3" << G4endl;
0098
0099
0100 G4double SpleenMass = (SpleenVol)*SpleenDensity;
0101 G4cout << "Mass of Spleen = " << SpleenMass/gram << " g" << G4endl;
0102
0103 return physSpleen;
0104 }