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 #include "G4MIRDRightArmBone.hh"
0031 #include "globals.hh"
0032 #include "G4SystemOfUnits.hh"
0033 #include "G4SDManager.hh"
0034 #include "G4VisAttributes.hh"
0035 #include "G4EllipticalTube.hh"
0036 #include "G4RotationMatrix.hh"
0037 #include "G4ThreeVector.hh"
0038 #include "G4VPhysicalVolume.hh"
0039 #include "G4PVPlacement.hh"
0040 #include "G4UnionSolid.hh"
0041 #include "G4EllipticalCone.hh"
0042 #include "G4HumanPhantomColour.hh"
0043 #include "G4HumanPhantomMaterial.hh"
0044
0045 G4VPhysicalVolume* G4MIRDRightArmBone::Construct(const G4String& volumeName,G4VPhysicalVolume* mother,
0046 const G4String& colourName,G4bool wireFrame, G4bool)
0047 {
0048 auto* material = new G4HumanPhantomMaterial();
0049
0050 G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
0051
0052 auto* skeleton = material -> GetMaterial("skeleton");
0053
0054 delete material;
0055
0056 G4double dx = 1.4 * cm;
0057 G4double dy = 2.7 * cm;
0058
0059
0060 auto* rightArm = new G4EllipticalTube("OneRightArmBone",dx,dy,34.5 *cm);
0061
0062 auto* logicRightArmBone = new G4LogicalVolume(rightArm,
0063 skeleton,
0064 "logical" + volumeName,
0065 nullptr, nullptr, nullptr);
0066
0067 auto* rm = new G4RotationMatrix();
0068 rm->rotateX(180.*degree);
0069 G4VPhysicalVolume* physRightArmBone = new G4PVPlacement(rm,
0070 G4ThreeVector(-18.4 * cm, 0.0, -0.5*cm),
0071
0072 "physicalRightArmBone",
0073 logicRightArmBone,
0074 mother,
0075 false,0,true);
0076
0077
0078
0079 auto* colourPointer = new G4HumanPhantomColour();
0080 G4Colour colour = colourPointer -> GetColour(colourName);
0081 auto* RightArmBoneVisAtt = new G4VisAttributes(colour);
0082
0083 RightArmBoneVisAtt->SetForceSolid(wireFrame);
0084 logicRightArmBone->SetVisAttributes(RightArmBoneVisAtt);
0085
0086 G4cout << "RightArmBone created !!!!!!" << G4endl;
0087
0088
0089 G4double RightArmBoneVol = logicRightArmBone->GetSolid()->GetCubicVolume();
0090 G4cout << "Volume of RightArmBone = " << RightArmBoneVol/cm3 << " cm^3" << G4endl;
0091
0092
0093 G4String RightArmBoneMat = logicRightArmBone->GetMaterial()->GetName();
0094 G4cout << "Material of RightArmBone = " << RightArmBoneMat << G4endl;
0095
0096
0097 G4double RightArmBoneDensity = logicRightArmBone->GetMaterial()->GetDensity();
0098 G4cout << "Density of Material = " << RightArmBoneDensity*cm3/g << " g/cm^3" << G4endl;
0099
0100
0101 G4double RightArmBoneMass = (RightArmBoneVol)*RightArmBoneDensity;
0102 G4cout << "Mass of RightArmBone = " << RightArmBoneMass/gram << " g" << G4endl;
0103
0104 return physRightArmBone;
0105 }