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