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 #include "G4MIRDRightTeste.hh"
0032 #include "globals.hh"
0033 #include "G4SystemOfUnits.hh"
0034 #include "G4SDManager.hh"
0035 #include "G4VisAttributes.hh"
0036 #include "G4Ellipsoid.hh"
0037 #include "G4ThreeVector.hh"
0038 #include "G4VPhysicalVolume.hh"
0039 #include "G4RotationMatrix.hh"
0040 #include "G4Material.hh"
0041 #include "G4LogicalVolume.hh"
0042 #include "G4HumanPhantomMaterial.hh"
0043 #include "G4VPhysicalVolume.hh"
0044 #include "G4PVPlacement.hh"
0045 #include "G4UnionSolid.hh"
0046 #include "G4HumanPhantomColour.hh"
0047
0048 G4VPhysicalVolume* G4MIRDRightTeste::Construct(const G4String& volumeName, G4VPhysicalVolume* mother,
0049 const G4String& colourName, G4bool wireFrame,G4bool)
0050 {
0051 G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
0052 auto* material = new G4HumanPhantomMaterial();
0053 auto* soft = material -> GetMaterial("soft_tissue");
0054 delete material;
0055
0056 G4double ax= 1.3*cm;
0057 G4double by= 1.5*cm;
0058 G4double cz= 2.3*cm;
0059
0060 auto* OneTeste = new G4Ellipsoid("OneTeste",
0061 ax, by, cz);
0062
0063 auto* logicRightTeste = new G4LogicalVolume(OneTeste,
0064 soft,
0065 "logical" + volumeName,
0066 nullptr, nullptr, nullptr);
0067
0068
0069 G4VPhysicalVolume* physRightTeste = new G4PVPlacement(nullptr,
0070 G4ThreeVector(-1.4*cm,3*cm, 0*cm),
0071 "physicalRightTeste",
0072 logicRightTeste,
0073 mother,
0074 false,
0075 0, true);
0076
0077
0078
0079
0080 auto* colourPointer = new G4HumanPhantomColour();
0081 G4Colour colour = colourPointer -> GetColour(colourName);
0082 auto* RightTesteVisAtt = new G4VisAttributes(colour);
0083 RightTesteVisAtt->SetForceSolid(wireFrame);
0084 logicRightTeste->SetVisAttributes(RightTesteVisAtt);
0085
0086 G4cout << "RightTeste created !!!!!!" << G4endl;
0087
0088
0089 G4double RightTesteVol = logicRightTeste->GetSolid()->GetCubicVolume();
0090 G4cout << "Volume of RightTeste = " << RightTesteVol/cm3 << " cm^3" << G4endl;
0091
0092
0093 G4String RightTesteMat = logicRightTeste->GetMaterial()->GetName();
0094 G4cout << "Material of RightTeste = " << RightTesteMat << G4endl;
0095
0096
0097 G4double RightTesteDensity = logicRightTeste->GetMaterial()->GetDensity();
0098 G4cout << "Density of Material = " << RightTesteDensity*cm3/g << " g/cm^3" << G4endl;
0099
0100
0101 G4double RightTesteMass = (RightTesteVol)*RightTesteDensity;
0102 G4cout << "Mass of RightTeste = " << RightTesteMass/gram << " g" << G4endl;
0103
0104 return physRightTeste;
0105 }