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 "G4MIRDLeftAdrenal.hh"
0031
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* G4MIRDLeftAdrenal::Construct(const G4String& volumeName,G4VPhysicalVolume* mother,
0051 const G4String& colourName, G4bool wireFrame, G4bool)
0052 {
0053 G4cout << "Construct " << volumeName <<" with mother "<<mother->GetName()<<G4endl;
0054
0055 auto* material = new G4HumanPhantomMaterial();
0056 auto* soft = material -> GetMaterial("soft_tissue");
0057 delete material;
0058
0059 G4double ax= 1.5 *cm;
0060 G4double by= 0.5 *cm;
0061 G4double cz= 5.0 *cm;
0062
0063 auto* leftAdrenal = new G4Ellipsoid("OneLeftAdrenal",ax, by, cz, 0. *cm, cz);
0064
0065
0066 auto* logicLeftAdrenal = new G4LogicalVolume(leftAdrenal,
0067 soft,
0068 "logical" + volumeName,
0069 nullptr, nullptr, nullptr);
0070
0071 G4VPhysicalVolume* physLeftAdrenal = new G4PVPlacement(nullptr,G4ThreeVector(4.5*cm,
0072 6.5 *cm,
0073 3. *cm),
0074 "physicalLeftAdrenal", logicLeftAdrenal,
0075 mother,
0076 false,
0077 0, true);
0078
0079
0080
0081 auto* colourPointer = new G4HumanPhantomColour();
0082 G4Colour colour = colourPointer -> GetColour(colourName);
0083 auto* leftAdrenalVisAtt = new G4VisAttributes(colour);
0084 leftAdrenalVisAtt->SetForceSolid(wireFrame);
0085 logicLeftAdrenal->SetVisAttributes(leftAdrenalVisAtt);
0086
0087 G4cout << "Left LeftAdrenal created !!!!!!" << G4endl;
0088
0089
0090 G4double LeftAdrenalVol = logicLeftAdrenal->GetSolid()->GetCubicVolume();
0091 G4cout << "Volume of LeftAdrenal = " << LeftAdrenalVol/cm3 << " cm^3" << G4endl;
0092
0093
0094 G4String LeftAdrenalMat = logicLeftAdrenal->GetMaterial()->GetName();
0095 G4cout << "Material of LeftAdrenal = " << LeftAdrenalMat << G4endl;
0096
0097
0098 G4double LeftAdrenalDensity = logicLeftAdrenal->GetMaterial()->GetDensity();
0099 G4cout << "Density of Material = " << LeftAdrenalDensity*cm3/g << " g/cm^3" << G4endl;
0100
0101
0102 G4double LeftAdrenalMass = (LeftAdrenalVol)*LeftAdrenalDensity;
0103 G4cout << "Mass of LeftAdrenal = " << LeftAdrenalMass/gram << " g" << G4endl;
0104
0105
0106 return physLeftAdrenal;
0107 }