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 "G4MIRDLeftBreast.hh"
0031 #include "globals.hh"
0032 #include "G4SystemOfUnits.hh"
0033 #include "G4SDManager.hh"
0034 #include "G4VisAttributes.hh"
0035 #include "G4Ellipsoid.hh"
0036 #include "G4ThreeVector.hh"
0037 #include "G4VPhysicalVolume.hh"
0038 #include "G4RotationMatrix.hh"
0039 #include "G4Material.hh"
0040 #include "G4LogicalVolume.hh"
0041 #include "G4HumanPhantomMaterial.hh"
0042 #include "G4VPhysicalVolume.hh"
0043 #include "G4PVPlacement.hh"
0044 #include "G4EllipticalTube.hh"
0045 #include "G4UnionSolid.hh"
0046 #include "G4SubtractionSolid.hh"
0047 #include "G4HumanPhantomColour.hh"
0048
0049 G4VPhysicalVolume* G4MIRDLeftBreast::Construct(const G4String& volumeName,G4VPhysicalVolume* mother,
0050 const G4String& colourName, G4bool wireFrame, G4bool)
0051 {
0052 G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
0053
0054 auto* material = new G4HumanPhantomMaterial();
0055 auto* soft = material -> GetMaterial("soft_tissue");
0056 delete material;
0057
0058 G4double ax= 4.95* cm;
0059 G4double by= 4.35* cm;
0060 G4double cz= 4.15*cm;
0061
0062 auto* oneLeftBreast = new G4Ellipsoid("OneLeftBreast",
0063 ax, by, cz);
0064
0065 G4double dx= 20.* cm;
0066 G4double dy= 10.* cm;
0067 G4double dz= 35.*cm;
0068
0069 auto* Trunk = new G4EllipticalTube("Trunk",dx, dy, dz );
0070
0071
0072
0073 auto* rm_relative = new G4RotationMatrix();
0074 rm_relative -> rotateX(90. * degree);
0075
0076 auto* breast = new G4SubtractionSolid("LeftBreast",
0077 oneLeftBreast,
0078 Trunk,
0079 rm_relative,
0080 G4ThreeVector(-10.*cm,
0081 0.0*cm,
0082 -8.66*cm));
0083
0084
0085 auto* logicLeftBreast = new G4LogicalVolume(breast, soft,
0086 "logical" + volumeName, nullptr, nullptr, nullptr);
0087
0088
0089
0090 auto* rm = new G4RotationMatrix();
0091 rm->rotateX(90.*degree);
0092 rm->rotateY(0.*degree);
0093 rm->rotateZ(16.*degree);
0094 G4VPhysicalVolume* physLeftBreast = new G4PVPlacement(rm,G4ThreeVector(10*cm, 9.1 *cm, 52.* cm),
0095 "physicalLeftBreast",
0096 logicLeftBreast,
0097 mother,
0098 false,
0099 0, true);
0100
0101
0102
0103 auto* colourPointer = new G4HumanPhantomColour();
0104 G4Colour colour = colourPointer -> GetColour(colourName);
0105 auto* LeftBreastVisAtt = new G4VisAttributes(colour);
0106 LeftBreastVisAtt->SetForceSolid(wireFrame);
0107
0108 logicLeftBreast->SetVisAttributes(LeftBreastVisAtt);
0109
0110 G4cout << "LeftBreast created !!!!!!" << G4endl;
0111
0112
0113 G4double LeftBreastVol = logicLeftBreast->GetSolid()->GetCubicVolume();
0114 G4cout << "Volume of LeftBreast = " << LeftBreastVol/cm3 << " cm^3" << G4endl;
0115
0116
0117 G4String LeftBreastMat = logicLeftBreast->GetMaterial()->GetName();
0118 G4cout << "Material of LeftBreast = " << LeftBreastMat << G4endl;
0119
0120
0121 G4double LeftBreastDensity = logicLeftBreast->GetMaterial()->GetDensity();
0122 G4cout << "Density of Material = " << LeftBreastDensity*cm3/g << " g/cm^3" << G4endl;
0123
0124
0125 G4double LeftBreastMass = (LeftBreastVol)*LeftBreastDensity;
0126 G4cout << "Mass of LeftBreast = " << LeftBreastMass/gram << " g" << G4endl;
0127
0128 return physLeftBreast;
0129 }