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
0032 #include "G4MIRDUpperLargeIntestine.hh"
0033
0034 #include "globals.hh"
0035 #include "G4SystemOfUnits.hh"
0036 #include "G4SDManager.hh"
0037 #include "G4VisAttributes.hh"
0038 #include "G4EllipticalTube.hh"
0039 #include "G4UnionSolid.hh"
0040 #include "G4RotationMatrix.hh"
0041 #include "G4ThreeVector.hh"
0042 #include "G4VPhysicalVolume.hh"
0043 #include "G4PVPlacement.hh"
0044 #include "G4LogicalVolume.hh"
0045 #include "G4HumanPhantomMaterial.hh"
0046 #include "G4HumanPhantomColour.hh"
0047
0048 G4VPhysicalVolume* G4MIRDUpperLargeIntestine::Construct(const G4String& volumeName,
0049 G4VPhysicalVolume* mother,
0050 const G4String& colourName
0051 , G4bool wireFrame,G4bool)
0052 {
0053 G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
0054
0055 auto* material = new G4HumanPhantomMaterial();
0056 auto* soft = material -> GetMaterial("soft_tissue");
0057 delete material;
0058
0059 G4double dx = 2.5 * cm;
0060 G4double dy = 2.5* cm;
0061 G4double dz = 4.775 * cm;
0062
0063 auto* AscendingColonUpperLargeIntestine = new G4EllipticalTube("AscendingColon",dx, dy, dz);
0064
0065 dx = 2.5 * cm;
0066 dy = 1.5 *cm;
0067 dz = 10.5* cm;
0068
0069 auto* TraverseColonUpperLargeIntestine = new G4EllipticalTube("TraverseColon",dx, dy, dz);
0070
0071 auto* relative_rm = new G4RotationMatrix();
0072 relative_rm -> rotateX(90. * degree);
0073 relative_rm -> rotateZ(0. * degree);
0074 relative_rm -> rotateY(90. * degree);
0075 auto* upperLargeIntestine = new G4UnionSolid("UpperLargeIntestine",
0076 AscendingColonUpperLargeIntestine,
0077 TraverseColonUpperLargeIntestine,
0078 relative_rm,
0079 G4ThreeVector(8.0 *cm, 0.0,6.275 * cm));
0080
0081
0082 auto* logicUpperLargeIntestine = new G4LogicalVolume(upperLargeIntestine, soft,
0083 "logical" + volumeName,
0084 nullptr, nullptr, nullptr);
0085
0086 G4VPhysicalVolume* physUpperLargeIntestine = new G4PVPlacement(nullptr,
0087 G4ThreeVector(-8.0 * cm, -2.36 *cm,-15.775 *cm),
0088 "physicalUpperLargeIntestine",
0089 logicUpperLargeIntestine,
0090 mother,
0091 false,
0092 0, true);
0093
0094
0095
0096
0097 auto* colourPointer = new G4HumanPhantomColour();
0098 G4Colour colour = colourPointer -> GetColour(colourName);
0099 auto* UpperLargeIntestineVisAtt = new G4VisAttributes(colour);
0100 UpperLargeIntestineVisAtt->SetForceSolid(wireFrame);
0101 logicUpperLargeIntestine->SetVisAttributes(UpperLargeIntestineVisAtt);
0102
0103 G4cout << "UpperLargeIntestine created !!!!!!" << G4endl;
0104
0105
0106 G4double UpperLargeIntestineVol = logicUpperLargeIntestine->GetSolid()->GetCubicVolume();
0107 G4cout << "Volume of UpperLargeIntestine = " << UpperLargeIntestineVol/cm3 << " cm^3" << G4endl;
0108
0109
0110 G4String UpperLargeIntestineMat = logicUpperLargeIntestine->GetMaterial()->GetName();
0111 G4cout << "Material of UpperLargeIntestine = " << UpperLargeIntestineMat << G4endl;
0112
0113
0114 G4double UpperLargeIntestineDensity = logicUpperLargeIntestine->GetMaterial()->GetDensity();
0115 G4cout << "Density of Material = " << UpperLargeIntestineDensity*cm3/g << " g/cm^3" << G4endl;
0116
0117
0118 G4double UpperLargeIntestineMass = (UpperLargeIntestineVol)*UpperLargeIntestineDensity;
0119 G4cout << "Mass of UpperLargeIntestine = " << UpperLargeIntestineMass/gram << " g" << G4endl;
0120
0121 return physUpperLargeIntestine;
0122 }