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 "G4MIRDMiddleLowerSpine.hh"
0031
0032 #include "globals.hh"
0033 #include "G4SystemOfUnits.hh"
0034 #include "G4SDManager.hh"
0035 #include "G4VisAttributes.hh"
0036 #include "G4HumanPhantomMaterial.hh"
0037 #include "G4EllipticalTube.hh"
0038 #include "G4RotationMatrix.hh"
0039 #include "G4ThreeVector.hh"
0040 #include "G4VPhysicalVolume.hh"
0041 #include "G4PVPlacement.hh"
0042 #include "G4UnionSolid.hh"
0043 #include "G4HumanPhantomColour.hh"
0044
0045 G4VPhysicalVolume* G4MIRDMiddleLowerSpine::Construct(const G4String& volumeName,
0046 G4VPhysicalVolume* mother,
0047 const G4String& colourName
0048 , G4bool wireFrame, G4bool )
0049 {
0050 auto* material = new G4HumanPhantomMaterial();
0051
0052 G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
0053
0054 auto* skeleton = material -> GetMaterial("skeleton");
0055
0056 delete material;
0057
0058 G4double dx = 2. *cm;
0059 G4double dy = 2.5 *cm;
0060 G4double dz = 24. *cm;
0061
0062 auto* middleLowerSpine = new G4EllipticalTube("MiddleLowerSpine",dx, dy, dz);
0063
0064 auto* logicMiddleLowerSpine = new G4LogicalVolume( middleLowerSpine, skeleton,
0065 "logical" + volumeName,
0066 nullptr, nullptr, nullptr);
0067
0068 G4VPhysicalVolume* physMiddleLowerSpine = new G4PVPlacement(nullptr,G4ThreeVector(0.0 *cm, 5.5 * cm,11. * cm),
0069 "physicalMiddleLowerSpine",
0070 logicMiddleLowerSpine,
0071 mother,
0072 false,
0073 0, true);
0074
0075
0076
0077
0078
0079 auto* colourPointer = new G4HumanPhantomColour();
0080 G4Colour colour = colourPointer -> GetColour(colourName);
0081 auto* MiddleLowerSpineVisAtt = new G4VisAttributes(colour);
0082 MiddleLowerSpineVisAtt->SetForceSolid(wireFrame);
0083 logicMiddleLowerSpine->SetVisAttributes(MiddleLowerSpineVisAtt);
0084
0085 G4cout << "MiddleLowerSpine created !!!!!!" << G4endl;
0086
0087
0088 G4double MiddleLowerSpineVol = logicMiddleLowerSpine->GetSolid()->GetCubicVolume();
0089 G4cout << "Volume of MiddleLowerSpine = " << MiddleLowerSpineVol/cm3 << " cm^3" << G4endl;
0090
0091
0092 G4String MiddleLowerSpineMat = logicMiddleLowerSpine->GetMaterial()->GetName();
0093 G4cout << "Material of MiddleLowerSpine = " << MiddleLowerSpineMat << G4endl;
0094
0095
0096 G4double MiddleLowerSpineDensity = logicMiddleLowerSpine->GetMaterial()->GetDensity();
0097 G4cout << "Density of Material = " << MiddleLowerSpineDensity*cm3/g << " g/cm^3" << G4endl;
0098
0099
0100 G4double MiddleLowerSpineMass = (MiddleLowerSpineVol)*MiddleLowerSpineDensity;
0101 G4cout << "Mass of MiddleLowerSpine = " << MiddleLowerSpineMass/gram << " g" << G4endl;
0102
0103 return physMiddleLowerSpine;
0104 }