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
0032 #include "G4MIRDRightClavicle.hh"
0033
0034 #include "globals.hh"
0035 #include "G4SystemOfUnits.hh"
0036 #include "G4SDManager.hh"
0037 #include "G4Cons.hh"
0038 #include "G4VisAttributes.hh"
0039 #include "G4HumanPhantomMaterial.hh"
0040 #include "G4EllipticalTube.hh"
0041 #include "G4ThreeVector.hh"
0042 #include "G4VPhysicalVolume.hh"
0043 #include "G4RotationMatrix.hh"
0044 #include "G4LogicalVolume.hh"
0045 #include "G4PVPlacement.hh"
0046 #include "G4SubtractionSolid.hh"
0047 #include "G4HumanPhantomColour.hh"
0048 #include "G4Box.hh"
0049 #include "G4Torus.hh"
0050
0051 G4VPhysicalVolume* G4MIRDRightClavicle::Construct(const G4String& volumeName, G4VPhysicalVolume* mother,
0052 const G4String& colourName, G4bool wireFrame,G4bool)
0053 {
0054 G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
0055
0056 auto* material = new G4HumanPhantomMaterial();
0057 auto* skeleton = material -> GetMaterial("skeleton");
0058
0059 G4double rMin = 0*cm;
0060 G4double rMax = 0.7883*cm;
0061 G4double rTor = 10*cm;
0062 G4double pSPhi = 201.75*degree;
0063 G4double pDPhi = 0.7*rad;
0064
0065 auto* clavicle = new G4Torus("Clavicle",rMin,rMax,rTor,pSPhi,pDPhi);
0066
0067 auto* logicRightClavicle = new G4LogicalVolume(clavicle,
0068 skeleton,
0069 "logical" + volumeName,
0070 nullptr, nullptr, nullptr);
0071
0072 G4VPhysicalVolume* physRightClavicle = new G4PVPlacement(nullptr,
0073 G4ThreeVector(0.*cm,2.*cm,33.25*cm),
0074 "physicalRightClavicle",
0075 logicRightClavicle,
0076 mother,
0077 false,
0078 0, true);
0079
0080
0081 auto* colourPointer = new G4HumanPhantomColour();
0082 G4Colour colour = colourPointer -> GetColour(colourName);
0083 auto* RightClavicleVisAtt = new G4VisAttributes(colour);
0084 RightClavicleVisAtt->SetForceSolid(wireFrame);
0085 logicRightClavicle->SetVisAttributes(RightClavicleVisAtt);
0086 G4cout << "RightClavicle created !!!!!!" << G4endl;
0087
0088
0089 G4double RightClavicleVol = logicRightClavicle->GetSolid()->GetCubicVolume();
0090 G4cout << "Volume of RightClavicle = " << RightClavicleVol/cm3 << " cm^3" << G4endl;
0091
0092
0093 G4String RightClavicleMat = logicRightClavicle->GetMaterial()->GetName();
0094 G4cout << "Material of RightClavicle = " << RightClavicleMat << G4endl;
0095
0096
0097 G4double RightClavicleDensity = logicRightClavicle->GetMaterial()->GetDensity();
0098 G4cout << "Density of Material = " << RightClavicleDensity*cm3/g << " g/cm^3" << G4endl;
0099
0100
0101 G4double RightClavicleMass = (RightClavicleVol)*RightClavicleDensity;
0102 G4cout << "Mass of RightClavicle = " << RightClavicleMass/gram << " g" << G4endl;
0103
0104 return physRightClavicle;
0105 }