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 #include "G4MIRDLeftScapula.hh"
0031
0032 #include "globals.hh"
0033 #include "G4SystemOfUnits.hh"
0034 #include "G4SDManager.hh"
0035 #include "G4Cons.hh"
0036 #include "G4VisAttributes.hh"
0037 #include "G4HumanPhantomMaterial.hh"
0038 #include "G4EllipticalTube.hh"
0039 #include "G4ThreeVector.hh"
0040 #include "G4VPhysicalVolume.hh"
0041 #include "G4RotationMatrix.hh"
0042 #include "G4LogicalVolume.hh"
0043 #include "G4PVPlacement.hh"
0044 #include "G4SubtractionSolid.hh"
0045 #include "G4HumanPhantomColour.hh"
0046 #include "G4Box.hh"
0047
0048 G4VPhysicalVolume* G4MIRDLeftScapula::Construct(const G4String& volumeName, G4VPhysicalVolume* mother,
0049 const G4String& colourName, G4bool wireFrame,G4bool)
0050 {
0051
0052 G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
0053
0054 auto* material = new G4HumanPhantomMaterial();
0055 auto* skeleton = material -> GetMaterial("skeleton");
0056
0057 G4double ax_in = 17.* cm;
0058 G4double by_in = 9.8* cm;
0059 G4double ax_out = 19.*cm;
0060 G4double by_out = 9.8*cm;
0061 G4double dz= 16.4* cm;
0062
0063 auto* inner_scapula = new G4EllipticalTube("ScapulaIn", ax_in, by_in, (dz+ 1.*cm)/2);
0064 auto* outer_scapula = new G4EllipticalTube("ScapulaOut", ax_out, by_out, dz/2);
0065
0066 auto* subtraction = new G4Box("subtraction",ax_out, ax_out, ax_out);
0067
0068 G4double xx = ax_out * 0.242 ;
0069 G4double yy = - ax_out * 0.97;
0070
0071 auto* rm = new G4RotationMatrix();
0072 rm -> rotateZ(-14.* degree);
0073
0074 auto* scapula_first = new G4SubtractionSolid("Scapula_first",
0075 outer_scapula,
0076 subtraction,
0077 rm,
0078 G4ThreeVector(xx, yy, 0. *cm));
0079
0080 G4double xx2 = -ax_out * 0.62470 ;
0081 G4double yy2 = ax_out * 0.78087;
0082
0083 auto* rm2 = new G4RotationMatrix();
0084 rm2 -> rotateZ(-38.6598* degree);
0085
0086 auto* scapula_bone = new G4SubtractionSolid("Scapula",
0087 scapula_first,
0088 subtraction,
0089 rm2,
0090 G4ThreeVector(xx2, yy2, 0. *cm));
0091
0092 auto* scapula = new G4SubtractionSolid("Scapula",
0093 scapula_bone,
0094 inner_scapula);
0095
0096 auto* logicLeftScapula = new G4LogicalVolume(scapula,
0097 skeleton,
0098 "logical" + volumeName,
0099 nullptr, nullptr, nullptr);
0100
0101 G4VPhysicalVolume* physLeftScapula = new G4PVPlacement(nullptr,
0102 G4ThreeVector(0. * cm, 0. * cm, 24.1 *cm),
0103 "physicalLeftScapula",
0104 logicLeftScapula,
0105 mother,
0106 false,
0107 0, true);
0108
0109
0110
0111 auto* colourPointer = new G4HumanPhantomColour();
0112 G4Colour colour = colourPointer -> GetColour(colourName);
0113 auto* LeftScapulaVisAtt = new G4VisAttributes(colour);
0114 LeftScapulaVisAtt->SetForceSolid(wireFrame);
0115 logicLeftScapula->SetVisAttributes(LeftScapulaVisAtt);
0116
0117 G4cout << "LeftScapula created !!!!!!" << G4endl;
0118
0119
0120 G4double LeftScapulaVol = logicLeftScapula->GetSolid()->GetCubicVolume();
0121 G4cout << "Volume of LeftScapula = " << LeftScapulaVol/cm3 << " cm^3" << G4endl;
0122
0123
0124 G4String LeftScapulaMat = logicLeftScapula->GetMaterial()->GetName();
0125 G4cout << "Material of LeftScapula = " << LeftScapulaMat << G4endl;
0126
0127
0128 G4double LeftScapulaDensity = logicLeftScapula->GetMaterial()->GetDensity();
0129 G4cout << "Density of Material = " << LeftScapulaDensity*cm3/g << " g/cm^3" << G4endl;
0130
0131
0132 G4double LeftScapulaMass = (LeftScapulaVol)*LeftScapulaDensity;
0133 G4cout << "Mass of LeftScapula = " << LeftScapulaMass/gram << " g" << G4endl;
0134
0135 return physLeftScapula;
0136 }