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 "G4MIRDLowerLargeIntestine.hh"
0030
0031 #include "globals.hh"
0032 #include "G4SystemOfUnits.hh"
0033 #include "G4SDManager.hh"
0034 #include "G4VisAttributes.hh"
0035 #include "G4EllipticalTube.hh"
0036 #include "G4UnionSolid.hh"
0037 #include "G4RotationMatrix.hh"
0038 #include "G4ThreeVector.hh"
0039 #include "G4VPhysicalVolume.hh"
0040 #include "G4PVPlacement.hh"
0041 #include "G4LogicalVolume.hh"
0042 #include "G4Torus.hh"
0043 #include "G4HumanPhantomMaterial.hh"
0044 #include "G4HumanPhantomColour.hh"
0045
0046 G4VPhysicalVolume* G4MIRDLowerLargeIntestine::Construct(const G4String& volumeName,
0047 G4VPhysicalVolume* mother,
0048 const G4String& colourName, G4bool wireFrame,G4bool)
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 delete material;
0055
0056 G4double dx = 1.88 * cm;
0057 G4double dy = 2.13 *cm;
0058 G4double dz = 7.64 *cm;
0059
0060 auto* DescendingColonLowerLargeIntestine = new G4EllipticalTube("DiscendingColon",dx, dy, dz);
0061
0062 G4double rmin= 0.0 *cm;
0063 G4double rmax = 1.88 * cm;
0064 G4double rtor= 5.72*cm;
0065 G4double startphi= 0. * degree;
0066 G4double deltaphi= 90. * degree;
0067
0068 auto* SigmoidColonUpLowerLargeIntestine = new G4Torus("SigmoidColonUpLowerLargeIntestine",
0069 rmin, rmax,rtor,
0070 startphi, deltaphi);
0071
0072 rtor = 3. * cm;
0073 G4VSolid* SigmoidColonDownLowerLargeIntestine = new G4Torus("SigmoidColonDownLowerLargeIntestine",
0074 rmin, rmax,
0075 rtor,startphi,deltaphi);
0076
0077 auto* relative_rm = new G4RotationMatrix();
0078 relative_rm -> rotateY(180. * degree);
0079 relative_rm -> rotateZ(90. * degree);
0080
0081 auto* SigmoidColonLowerLargeIntestine = new G4UnionSolid( "SigmoidColonLowerLargeIntestine",
0082 SigmoidColonUpLowerLargeIntestine,
0083 SigmoidColonDownLowerLargeIntestine,
0084 relative_rm,
0085 G4ThreeVector(0.0,8.72*cm,0.0));
0086
0087
0088 auto* relative_rm_2 = new G4RotationMatrix();
0089 relative_rm_2 -> rotateX(90. * degree);
0090
0091 auto* LowerLargeIntestine = new G4UnionSolid( "LowerLargeIntestine",
0092 DescendingColonLowerLargeIntestine,
0093 SigmoidColonLowerLargeIntestine,
0094 relative_rm_2,
0095 G4ThreeVector(-5.72*cm,0.0*cm, -7.64*cm)
0096 );
0097
0098
0099 auto* logicLowerLargeIntestine = new G4LogicalVolume( LowerLargeIntestine, soft,
0100 "logical" + volumeName,
0101 nullptr, nullptr, nullptr);
0102
0103 G4VPhysicalVolume* physLowerLargeIntestine = new G4PVPlacement(nullptr,
0104 G4ThreeVector(8.72*cm, -2.36*cm,-18.64 *cm),
0105 "physicalLowerLargeIntestine",
0106 logicLowerLargeIntestine,
0107 mother,
0108 false,
0109 0, true);
0110
0111
0112
0113
0114 auto* colourPointer = new G4HumanPhantomColour();
0115 G4Colour colour = colourPointer -> GetColour(colourName);
0116 auto* LowerLargeIntestineVisAtt = new G4VisAttributes(colour);
0117 LowerLargeIntestineVisAtt->SetForceSolid(wireFrame);
0118 logicLowerLargeIntestine->SetVisAttributes(LowerLargeIntestineVisAtt);
0119
0120 G4cout << "LowerLargeIntestine created !!!!!!" << G4endl;
0121
0122
0123 G4double LowerLargeIntestineVol = logicLowerLargeIntestine->GetSolid()->GetCubicVolume();
0124 G4cout << "Volume of LowerLargeIntestine = " << LowerLargeIntestineVol/cm3 << " cm^3" << G4endl;
0125
0126
0127 G4String LowerLargeIntestineMat = logicLowerLargeIntestine->GetMaterial()->GetName();
0128 G4cout << "Material of LowerLargeIntestine = " << LowerLargeIntestineMat << G4endl;
0129
0130
0131 G4double LowerLargeIntestineDensity = logicLowerLargeIntestine->GetMaterial()->GetDensity();
0132 G4cout << "Density of Material = " << LowerLargeIntestineDensity*cm3/g << " g/cm^3" << G4endl;
0133
0134
0135 G4double LowerLargeIntestineMass = (LowerLargeIntestineVol)*LowerLargeIntestineDensity;
0136 G4cout << "Mass of LowerLargeIntestine = " << LowerLargeIntestineMass/gram << " g" << G4endl;
0137
0138 return physLowerLargeIntestine;
0139 }