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 #include "G4MIRDSmallIntestine.hh"
0032
0033 #include "globals.hh"
0034 #include "G4SystemOfUnits.hh"
0035 #include "G4SDManager.hh"
0036 #include "G4VisAttributes.hh"
0037 #include "G4EllipticalTube.hh"
0038 #include "G4UnionSolid.hh"
0039 #include "G4RotationMatrix.hh"
0040 #include "G4ThreeVector.hh"
0041 #include "G4VPhysicalVolume.hh"
0042 #include "G4PVPlacement.hh"
0043 #include "G4LogicalVolume.hh"
0044 #include "G4Torus.hh"
0045 #include "G4Tubs.hh"
0046 #include "G4Box.hh"
0047 #include "G4IntersectionSolid.hh"
0048 #include "G4SubtractionSolid.hh"
0049 #include "G4HumanPhantomMaterial.hh"
0050 #include "G4HumanPhantomColour.hh"
0051
0052 G4VPhysicalVolume* G4MIRDSmallIntestine::Construct(const G4String& volumeName,
0053 G4VPhysicalVolume* mother,
0054 const G4String& colourName, G4bool wireFrame,G4bool)
0055 {
0056 G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
0057
0058 auto* material = new G4HumanPhantomMaterial();
0059 auto* soft = material -> GetMaterial("soft_tissue");
0060 delete material;
0061
0062 G4double boxX = 11.*cm;
0063 G4double boxY = 3.53*cm;
0064 G4double boxZ = 5*cm;
0065
0066 auto* smallIntestineBox = new G4Box("smallIntestineBox",boxX,boxY,boxZ);
0067
0068 G4double tubsRmin = 0*cm;
0069 G4double tubsRmax = 11.*cm;
0070 G4double tubsZ = 5*cm;
0071 G4double tubsSphi = 0*degree;
0072 G4double tubsDphi = 360*degree;
0073
0074 auto* smallIntestineTubs = new G4Tubs("smallIntestineTubs",tubsRmin,tubsRmax,tubsZ,tubsSphi,tubsDphi);
0075
0076
0077 auto* filledSmallIntestine1 = new G4IntersectionSolid("filledSmallIntestine1",smallIntestineTubs,smallIntestineBox,
0078 nullptr,G4ThreeVector(0*cm,-1.33*cm, 0*cm));
0079
0080 auto* filledSmallIntestine = new G4IntersectionSolid("filledSmallIntestine",filledSmallIntestine1,smallIntestineTubs,
0081 nullptr,G4ThreeVector(0*cm,0.8*cm, 0*cm));
0082
0083 G4double dx = 2.50*cm;
0084 G4double dy = 2.50*cm;
0085 G4double dz = 4.775*cm;
0086
0087 auto* AscendingColonUpperLargeIntestine = new G4EllipticalTube("AscendingColon",dx, dy, dz);
0088
0089 dx = 2.50 * cm;
0090 dy = 1.50 *cm;
0091 dz = 10.50* cm;
0092
0093 auto* TraverseColonUpperLargeIntestine = new G4EllipticalTube("TraverseColon",dx, dy, dz);
0094
0095 auto* relative_rm = new G4RotationMatrix();
0096 relative_rm -> rotateX(90. * degree);
0097
0098 relative_rm -> rotateY(90. * degree);
0099 auto* upperLargeIntestine = new G4UnionSolid("UpperLargeIntestine",
0100 AscendingColonUpperLargeIntestine,
0101 TraverseColonUpperLargeIntestine,
0102 relative_rm,
0103 G4ThreeVector(-8.0 *cm, 0.0*cm,6.275* cm));
0104
0105 dx = 1.88 * cm;
0106 dy = 2.13 *cm;
0107 dz = 7.64 *cm;
0108
0109 auto* DescendingColonLowerLargeIntestine = new G4EllipticalTube("DiscendingColon",dx, dy, dz);
0110
0111 auto* upperlowerLargeIntestine = new G4UnionSolid("UpperLowerLargeIntestine",
0112 upperLargeIntestine,
0113 DescendingColonLowerLargeIntestine,
0114 nullptr,
0115 G4ThreeVector(-16.72*cm, 0.0*cm,-2.865* cm));
0116
0117
0118 auto* SmallIntestine = new G4SubtractionSolid("SmallIntestine",
0119 filledSmallIntestine,
0120 upperlowerLargeIntestine,
0121 nullptr,
0122 G4ThreeVector(8.0*cm,-0.3*cm,-2.775*cm));
0123
0124
0125 auto* logicSmallIntestine = new G4LogicalVolume( SmallIntestine,
0126 soft,
0127 "logical"+volumeName,
0128 nullptr, nullptr, nullptr);
0129 auto* rm = new G4RotationMatrix();
0130 rm->rotateX(180.*degree);
0131 rm->rotateY(180.*degree);
0132 G4VPhysicalVolume* physSmallIntestine = new G4PVPlacement(rm,
0133 G4ThreeVector(0*cm, -2.66*cm, -13*cm),
0134
0135
0136 "physical"+volumeName,
0137 logicSmallIntestine,
0138 mother,
0139 false,
0140 0, true);
0141
0142
0143
0144
0145 auto* colourPointer = new G4HumanPhantomColour();
0146 G4Colour colour = colourPointer -> GetColour(colourName);
0147 auto* SmallIntestineVisAtt = new G4VisAttributes(colour);
0148 SmallIntestineVisAtt->SetForceSolid(wireFrame);
0149 logicSmallIntestine->SetVisAttributes(SmallIntestineVisAtt);
0150
0151 G4cout << "SmallIntestine created !!!!!!" << G4endl;
0152
0153
0154 G4double SmallIntestineVol = logicSmallIntestine->GetSolid()->GetCubicVolume();
0155 G4cout << "Volume of SmallIntestine = " << SmallIntestineVol/cm3 << " cm^3" << G4endl;
0156
0157
0158 G4String SmallIntestineMat = logicSmallIntestine->GetMaterial()->GetName();
0159 G4cout << "Material of SmallIntestine = " << SmallIntestineMat << G4endl;
0160
0161
0162 G4double SmallIntestineDensity = logicSmallIntestine->GetMaterial()->GetDensity();
0163 G4cout << "Density of Material = " << SmallIntestineDensity*cm3/g << " g/cm^3" << G4endl;
0164
0165
0166 G4double SmallIntestineMass = (SmallIntestineVol)*SmallIntestineDensity;
0167 G4cout << "Mass of SmallIntestine = " << SmallIntestineMass/gram << " g" << G4endl;
0168
0169 return physSmallIntestine;
0170 }