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
0032 #include "G4MIRDRightAdrenal.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
0053 G4VPhysicalVolume* G4MIRDRightAdrenal::Construct(const G4String& volumeName,G4VPhysicalVolume* mother,
0054 const G4String& colourName, G4bool wireFrame, G4bool)
0055 {
0056 G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
0057
0058 auto* material = new G4HumanPhantomMaterial();
0059 auto* soft = material -> GetMaterial("soft_tissue");
0060 delete material;
0061
0062 G4double ax= 1.5 *cm;
0063 G4double by= 0.5 *cm;
0064 G4double cz= 5.0 *cm;
0065
0066 auto* rightAdrenal = new G4Ellipsoid("OneRightAdrenal",ax, by, cz, 0. *cm, cz);
0067
0068
0069 auto* logicRightAdrenal = new G4LogicalVolume(rightAdrenal,
0070 soft,
0071 "logical" + volumeName,
0072 nullptr, nullptr, nullptr);
0073
0074 G4VPhysicalVolume* physRightAdrenal = new G4PVPlacement(nullptr, G4ThreeVector(-4.5*cm,
0075 6.5 *cm,
0076 3. *cm),
0077 "physicalRightAdrenal", logicRightAdrenal,
0078 mother,
0079 false,
0080 0, true);
0081
0082
0083
0084 auto* colourPointer = new G4HumanPhantomColour();
0085 G4Colour colour = colourPointer -> GetColour(colourName);
0086 auto* RightAdrenalVisAtt = new G4VisAttributes(colour);
0087 RightAdrenalVisAtt->SetForceSolid(wireFrame);
0088 logicRightAdrenal->SetVisAttributes(RightAdrenalVisAtt);
0089
0090 G4cout << "Right RightAdrenal created !!!!!!" << G4endl;
0091
0092
0093 G4double RightAdrenalVol = logicRightAdrenal->GetSolid()->GetCubicVolume();
0094 G4cout << "Volume of RightAdrenal = " << RightAdrenalVol/cm3 << " cm^3" << G4endl;
0095
0096
0097 G4String RightAdrenalMat = logicRightAdrenal->GetMaterial()->GetName();
0098 G4cout << "Material of RightAdrenal = " << RightAdrenalMat << G4endl;
0099
0100
0101 G4double RightAdrenalDensity = logicRightAdrenal->GetMaterial()->GetDensity();
0102 G4cout << "Density of Material = " << RightAdrenalDensity*cm3/g << " g/cm^3" << G4endl;
0103
0104
0105 G4double RightAdrenalMass = (RightAdrenalVol)*RightAdrenalDensity;
0106 G4cout << "Mass of RightAdrenal = " << RightAdrenalMass/gram << " g" << G4endl;
0107
0108 return physRightAdrenal;
0109 }