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 "G4MIRDPelvis.hh"
0032
0033 #include "globals.hh"
0034 #include "G4SystemOfUnits.hh"
0035 #include "G4SDManager.hh"
0036 #include "G4VisAttributes.hh"
0037 #include "G4HumanPhantomMaterial.hh"
0038 #include "G4EllipticalTube.hh"
0039 #include "G4RotationMatrix.hh"
0040 #include "G4ThreeVector.hh"
0041 #include "G4VPhysicalVolume.hh"
0042 #include "G4PVPlacement.hh"
0043 #include "G4SubtractionSolid.hh"
0044 #include "G4Box.hh"
0045 #include "G4VSolid.hh"
0046 #include "G4LogicalVolume.hh"
0047 #include "G4HumanPhantomColour.hh"
0048
0049 G4VPhysicalVolume* G4MIRDPelvis::Construct(const G4String& volumeName,G4VPhysicalVolume* mother,
0050 const G4String& colourName, G4bool wireFrame,G4bool)
0051 {
0052 auto* material = new G4HumanPhantomMaterial();
0053
0054 G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
0055
0056 auto* skeleton = material -> GetMaterial("skeleton");
0057
0058 delete material;
0059
0060
0061
0062
0063
0064
0065
0066
0067
0068
0069
0070
0071 G4double dx= 12. *cm;
0072 G4double dy= 12. * cm;
0073 G4double dz= 11. * cm;
0074
0075 G4VSolid* outPelvis = new G4EllipticalTube("OutPelvis",dx, dy, dz);
0076
0077 dx = 11.3 * cm;
0078 dy = 11.3* cm;
0079 dz = 12.0 *cm;
0080
0081 G4VSolid* inPelvis = new G4EllipticalTube("InPelvis",dx, dy, dz);
0082
0083 G4double x = 28. * cm;
0084 G4double y = 28. * cm;
0085 G4double z = 24. *cm;
0086
0087 auto* subPelvis = new G4Box("SubtrPelvis", x/2., y/2., z/2.);
0088
0089 auto* firstPelvis = new G4SubtractionSolid("FirstPelvis",
0090 outPelvis,
0091 inPelvis, nullptr, G4ThreeVector(0.*cm, -0.8 *cm, 0. * cm));
0092
0093 auto* secondPelvis = new G4SubtractionSolid("SecondPelvis",
0094 firstPelvis,
0095 subPelvis, nullptr,
0096 G4ThreeVector(0.0,
0097 -14. * cm, 0.*cm));
0098
0099 auto* pelvis = new G4SubtractionSolid("Pelvis", secondPelvis, subPelvis,
0100 nullptr,
0101 G4ThreeVector(0.0,
0102 22. * cm, -9. *cm));
0103
0104
0105 auto* logicPelvis = new G4LogicalVolume(pelvis, skeleton,
0106 "logical" + volumeName, nullptr, nullptr, nullptr);
0107
0108
0109 G4VPhysicalVolume* physPelvis = new G4PVPlacement(nullptr,G4ThreeVector(0.0, -3. * cm,-24. * cm),
0110
0111 "physicalPelvis",
0112 logicPelvis,
0113 mother,
0114 false,
0115 0, true);
0116
0117
0118
0119
0120 auto* colourPointer = new G4HumanPhantomColour();
0121 G4Colour colour = colourPointer -> GetColour(colourName);
0122 auto* PelvisVisAtt = new G4VisAttributes(colour);
0123 PelvisVisAtt->SetForceSolid(wireFrame);
0124 logicPelvis->SetVisAttributes(PelvisVisAtt);
0125
0126 G4cout << "Pelvis created !!!!!!" << G4endl;
0127
0128
0129 G4double PelvisVol = logicPelvis->GetSolid()->GetCubicVolume();
0130 G4cout << "Volume of Pelvis = " << PelvisVol/cm3 << " cm^3" << G4endl;
0131
0132
0133 G4String PelvisMat = logicPelvis->GetMaterial()->GetName();
0134 G4cout << "Material of Pelvis = " << PelvisMat << G4endl;
0135
0136
0137 G4double PelvisDensity = logicPelvis->GetMaterial()->GetDensity();
0138 G4cout << "Density of Material = " << PelvisDensity*cm3/g << " g/cm^3" << G4endl;
0139
0140
0141 G4double PelvisMass = (PelvisVol)*PelvisDensity;
0142 G4cout << "Mass of Pelvis = " << PelvisMass/gram << " g" << G4endl;
0143
0144 return physPelvis;
0145 }