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 #include "G4MIRDLeftLeg.hh"
0030
0031 #include "globals.hh"
0032 #include "G4SystemOfUnits.hh"
0033 #include "G4SDManager.hh"
0034 #include "G4Cons.hh"
0035 #include "G4VisAttributes.hh"
0036 #include "G4HumanPhantomMaterial.hh"
0037 #include "G4EllipticalTube.hh"
0038 #include "G4ThreeVector.hh"
0039 #include "G4VPhysicalVolume.hh"
0040 #include "G4RotationMatrix.hh"
0041 #include "G4LogicalVolume.hh"
0042 #include "G4PVPlacement.hh"
0043 #include "G4UnionSolid.hh"
0044 #include "G4HumanPhantomColour.hh"
0045
0046 G4VPhysicalVolume* G4MIRDLeftLeg::Construct(const G4String& volumeName, G4VPhysicalVolume* mother,
0047 const G4String& colourName, G4bool wireFrame,G4bool)
0048 {
0049
0050 G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
0051
0052 auto* material = new G4HumanPhantomMaterial();
0053 auto* soft = material -> GetMaterial("soft_tissue");
0054
0055 G4double rmin1 = 0.* cm;
0056 G4double rmin2 = 0.* cm;
0057 G4double dz= 80.0 * cm;
0058 G4double rmax1= 2.0 * cm;
0059 G4double rmax2= 10. * cm;
0060 G4double startphi= 0.* degree;
0061 G4double deltaphi= 360. * degree;
0062
0063 auto* leg1 = new G4Cons("Leg1",
0064 rmin1, rmax1,
0065 rmin2, rmax2, dz/2.,
0066 startphi, deltaphi);
0067
0068 auto* logicLeftLeg = new G4LogicalVolume(leg1,
0069 soft,
0070 "logical" + volumeName,
0071 nullptr, nullptr, nullptr);
0072
0073 auto* rm = new G4RotationMatrix();
0074 rm->rotateX(180.*degree);
0075 rm->rotateY(180.*degree);
0076 G4VPhysicalVolume* physLeftLeg = new G4PVPlacement(rm,
0077
0078 G4ThreeVector(10. * cm, 0. * cm, -40. *cm),
0079 "physicalLeftLeg",
0080 logicLeftLeg,
0081 mother,
0082 false,
0083 0, true);
0084
0085
0086
0087
0088 auto* colourPointer = new G4HumanPhantomColour();
0089 G4Colour colour = colourPointer -> GetColour(colourName);
0090 auto* LeftLegVisAtt = new G4VisAttributes(colour);
0091 LeftLegVisAtt->SetForceSolid(wireFrame);
0092 logicLeftLeg->SetVisAttributes(LeftLegVisAtt);
0093
0094 G4cout << "LeftLeg created !!!!!!" << G4endl;
0095
0096
0097 G4double LeftLegVol = logicLeftLeg->GetSolid()->GetCubicVolume();
0098 G4cout << "Volume of LeftLeg = " << LeftLegVol/cm3 << " cm^3" << G4endl;
0099
0100
0101 G4String LeftLegMat = logicLeftLeg->GetMaterial()->GetName();
0102 G4cout << "Material of LeftLeg = " << LeftLegMat << G4endl;
0103
0104
0105 G4double LeftLegDensity = logicLeftLeg->GetMaterial()->GetDensity();
0106 G4cout << "Density of Material = " << LeftLegDensity*cm3/g << " g/cm^3" << G4endl;
0107
0108
0109 G4double LeftLegMass = (LeftLegVol)*LeftLegDensity;
0110 G4cout << "Mass of LeftLeg = " << LeftLegMass/gram << " g" << G4endl;
0111
0112 return physLeftLeg;
0113 }