File indexing completed on 2025-02-23 09:20:12
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 "G4MIRDRightOvary.hh"
0033
0034 #include "globals.hh"
0035 #include "G4SystemOfUnits.hh"
0036 #include "G4SDManager.hh"
0037 #include "G4VisAttributes.hh"
0038 #include "G4Ellipsoid.hh"
0039 #include "G4ThreeVector.hh"
0040 #include "G4VPhysicalVolume.hh"
0041 #include "G4RotationMatrix.hh"
0042 #include "G4Material.hh"
0043 #include "G4LogicalVolume.hh"
0044 #include "G4HumanPhantomMaterial.hh"
0045 #include "G4VPhysicalVolume.hh"
0046 #include "G4PVPlacement.hh"
0047 #include "G4UnionSolid.hh"
0048 #include "G4HumanPhantomColour.hh"
0049
0050 G4VPhysicalVolume* G4MIRDRightOvary::Construct(const G4String& volumeName, G4VPhysicalVolume* mother,
0051 const G4String& colourName, G4bool wireFrame,G4bool)
0052 {
0053 G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
0054
0055 auto* material = new G4HumanPhantomMaterial();
0056 auto* soft = material -> GetMaterial("soft_tissue");
0057 delete material;
0058
0059 G4double ax= 1. *cm;
0060 G4double by= 0.5*cm;
0061 G4double cz= 2.*cm;
0062
0063 auto* OneOvary = new G4Ellipsoid("OneOvary",
0064 ax, by, cz);
0065
0066 auto* logicRightOvary = new G4LogicalVolume(OneOvary,
0067 soft,
0068 "logical" + volumeName,
0069 nullptr, nullptr, nullptr);
0070
0071
0072 G4VPhysicalVolume* physRightOvary = new G4PVPlacement(nullptr,
0073 G4ThreeVector(6. *cm,0.5*cm, -20*cm),
0074 "physicalRightOvary",
0075 logicRightOvary,
0076 mother,
0077 false,
0078 0, true);
0079
0080
0081
0082
0083 auto* colourPointer = new G4HumanPhantomColour();
0084 G4Colour colour = colourPointer -> GetColour(colourName);
0085 auto* RightOvaryVisAtt = new G4VisAttributes(colour);
0086 RightOvaryVisAtt->SetForceSolid(wireFrame);
0087 logicRightOvary->SetVisAttributes(RightOvaryVisAtt);
0088
0089 G4cout << "RightOvary created !!!!!!" << G4endl;
0090
0091
0092 G4double RightOvaryVol = logicRightOvary->GetSolid()->GetCubicVolume();
0093 G4cout << "Volume of RightOvary = " << RightOvaryVol/cm3 << " cm^3" << G4endl;
0094
0095
0096 G4String RightOvaryMat = logicRightOvary->GetMaterial()->GetName();
0097 G4cout << "Material of RightOvary = " << RightOvaryMat << G4endl;
0098
0099
0100 G4double RightOvaryDensity = logicRightOvary->GetMaterial()->GetDensity();
0101 G4cout << "Density of Material = " << RightOvaryDensity*cm3/g << " g/cm^3" << G4endl;
0102
0103
0104 G4double RightOvaryMass = (RightOvaryVol)*RightOvaryDensity;
0105 G4cout << "Mass of RightOvary = " << RightOvaryMass/gram << " g" << G4endl;
0106
0107 return physRightOvary;
0108 }