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 "G4MIRDRightKidney.hh"
0033
0034 #include "globals.hh"
0035 #include "G4SystemOfUnits.hh"
0036 #include "G4SDManager.hh"
0037 #include "G4VisAttributes.hh"
0038 #include "G4HumanPhantomMaterial.hh"
0039 #include "G4SDManager.hh"
0040 #include "G4PVPlacement.hh"
0041 #include "G4SubtractionSolid.hh"
0042 #include "G4Ellipsoid.hh"
0043 #include "G4ThreeVector.hh"
0044 #include "G4VPhysicalVolume.hh"
0045 #include "G4RotationMatrix.hh"
0046 #include "G4Material.hh"
0047 #include "G4EllipticalTube.hh"
0048 #include "G4Box.hh"
0049 #include "G4UnionSolid.hh"
0050 #include "G4HumanPhantomColour.hh"
0051
0052 G4VPhysicalVolume* G4MIRDRightKidney::Construct(const G4String& volumeName,
0053 G4VPhysicalVolume* mother,
0054 const G4String& colourName,
0055 G4bool wireFrame,G4bool )
0056 {
0057 G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
0058
0059 auto* material = new G4HumanPhantomMaterial();
0060 auto* soft = material -> GetMaterial("soft_tissue");
0061 delete material;
0062
0063 G4double ax= 4.5 *cm;
0064 G4double by= 1.5 *cm;
0065 G4double cz= 5.5 *cm;
0066
0067 auto* oneRightKidney = new G4Ellipsoid("OneRightKidney",ax, by, cz);
0068
0069 G4double xx = 6. * cm;
0070 G4double yy = 12.00*cm;
0071 G4double zz = 12.00*cm;
0072 G4VSolid* subtrRightKidney = new G4Box("SubtrRightKidney",xx/2., yy/2., zz/2.);
0073
0074 auto* kidney = new G4SubtractionSolid("RightKidney",
0075 oneRightKidney,
0076 subtrRightKidney,
0077 nullptr,
0078 G4ThreeVector(6. *cm,
0079 0.0 *cm,
0080 0.0 * cm));
0081
0082 auto* logicRightKidney = new G4LogicalVolume(kidney,
0083 soft,
0084 "logical" + volumeName,
0085 nullptr, nullptr, nullptr);
0086
0087 G4VPhysicalVolume* physRightKidney = new G4PVPlacement( nullptr, G4ThreeVector(-6.*cm,
0088 6. *cm,
0089 -2.50 *cm),
0090 "physicalRightKidney", logicRightKidney,
0091 mother,
0092 false,
0093 0, true);
0094
0095
0096
0097 auto* colourPointer = new G4HumanPhantomColour();
0098 G4Colour colour = colourPointer -> GetColour(colourName);
0099 auto* RightKidneyVisAtt = new G4VisAttributes(colour);
0100 RightKidneyVisAtt->SetForceSolid(wireFrame);
0101 logicRightKidney->SetVisAttributes(RightKidneyVisAtt);
0102
0103 G4cout << "RightKidney created !!!!!!" << G4endl;
0104
0105
0106 G4double RightKidneyVol = logicRightKidney->GetSolid()->GetCubicVolume();
0107 G4cout << "Volume of RightKidney = " << RightKidneyVol/cm3 << " cm^3" << G4endl;
0108
0109
0110 G4String RightKidneyMat = logicRightKidney->GetMaterial()->GetName();
0111 G4cout << "Material of RightKidney = " << RightKidneyMat << G4endl;
0112
0113
0114 G4double RightKidneyDensity = logicRightKidney->GetMaterial()->GetDensity();
0115 G4cout << "Density of Material = " << RightKidneyDensity*cm3/g << " g/cm^3" << G4endl;
0116
0117
0118 G4double RightKidneyMass = (RightKidneyVol)*RightKidneyDensity;
0119 G4cout << "Mass of RightKidney = " << RightKidneyMass/gram << " g" << G4endl;
0120
0121 return physRightKidney;
0122 }