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
0031 #include "G4MIRDLeftKidney.hh"
0032 #include "globals.hh"
0033 #include "G4SystemOfUnits.hh"
0034 #include "G4SDManager.hh"
0035 #include "G4VisAttributes.hh"
0036 #include "G4HumanPhantomMaterial.hh"
0037 #include "G4SDManager.hh"
0038 #include "G4PVPlacement.hh"
0039 #include "G4SubtractionSolid.hh"
0040 #include "G4Ellipsoid.hh"
0041 #include "G4ThreeVector.hh"
0042 #include "G4VPhysicalVolume.hh"
0043 #include "G4RotationMatrix.hh"
0044 #include "G4Material.hh"
0045 #include "G4EllipticalTube.hh"
0046 #include "G4Box.hh"
0047 #include "G4UnionSolid.hh"
0048 #include "G4HumanPhantomColour.hh"
0049
0050 G4VPhysicalVolume* G4MIRDLeftKidney::Construct(const G4String& volumeName,G4VPhysicalVolume* mother,
0051 const G4String& colourName, 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 ax= 4.5 *cm;
0060 G4double by= 1.5 *cm;
0061 G4double cz= 5.5 *cm;
0062
0063 auto* oneLeftKidney = new G4Ellipsoid("OneLeftKidney",ax, by, cz);
0064
0065 G4double xx = 6. * cm;
0066 G4double yy = 12.00*cm;
0067 G4double zz = 12.00*cm;
0068 auto* subtrLeftKidney = new G4Box("SubtrLeftKidney",xx/2., yy/2., zz/2.);
0069
0070 auto* leftKidney = new G4SubtractionSolid("LeftKidney",
0071 oneLeftKidney,
0072 subtrLeftKidney,
0073 nullptr,
0074 G4ThreeVector(-6. *cm,
0075 0.0 *cm,
0076 0.0 * cm));
0077
0078 auto* logicLeftKidney = new G4LogicalVolume(leftKidney,
0079 soft,
0080 "logical" + volumeName,
0081 nullptr, nullptr, nullptr);
0082
0083 G4VPhysicalVolume* physLeftKidney = new G4PVPlacement(nullptr,G4ThreeVector(6.*cm,
0084 6. *cm,
0085 -2.50 *cm),
0086 "physicalLeftKidney", logicLeftKidney,
0087 mother,
0088 false,
0089 0, true);
0090
0091
0092
0093
0094
0095 auto* colourPointer = new G4HumanPhantomColour();
0096 G4Colour colour = colourPointer -> GetColour(colourName);
0097 auto* LeftKidneyVisAtt = new G4VisAttributes(colour);
0098 LeftKidneyVisAtt->SetForceSolid(wireFrame);
0099 logicLeftKidney->SetVisAttributes(LeftKidneyVisAtt);
0100
0101 G4cout << "Left LeftKidney created !!!!!!" << G4endl;
0102
0103
0104 G4double LeftKidneyVol = logicLeftKidney->GetSolid()->GetCubicVolume();
0105 G4cout << "Volume of LeftKidney = " << LeftKidneyVol/cm3 << " cm^3" << G4endl;
0106
0107
0108 G4String LeftKidneyMat = logicLeftKidney->GetMaterial()->GetName();
0109 G4cout << "Material of LeftKidney = " << LeftKidneyMat << G4endl;
0110
0111
0112 G4double LeftKidneyDensity = logicLeftKidney->GetMaterial()->GetDensity();
0113 G4cout << "Density of Material = " << LeftKidneyDensity*cm3/g << " g/cm^3" << G4endl;
0114
0115
0116 G4double LeftKidneyMass = (LeftKidneyVol)*LeftKidneyDensity;
0117 G4cout << "Mass of LeftKidney = " << LeftKidneyMass/gram << " g" << G4endl;
0118
0119 return physLeftKidney;
0120 }