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 "G4MIRDLeftTeste.hh"
0031
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* G4MIRDLeftTeste::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
0053 auto* material = new G4HumanPhantomMaterial();
0054 auto* soft = material -> GetMaterial("soft_tissue");
0055 delete material;
0056
0057 G4double ax= 1.3*cm;
0058 G4double by= 1.5*cm;
0059 G4double cz= 2.3*cm;
0060
0061 auto* OneTeste = new G4Ellipsoid("OneTeste",
0062 ax, by, cz);
0063
0064
0065 auto* logicLeftTeste = new G4LogicalVolume(OneTeste,
0066 soft,
0067 "logical" + volumeName,
0068 nullptr, nullptr, nullptr);
0069
0070
0071 G4VPhysicalVolume* physLeftTeste = new G4PVPlacement(nullptr,
0072 G4ThreeVector(1.4 *cm,3.*cm, 0*cm),
0073 "physicalLeftTeste",
0074 logicLeftTeste,
0075 mother,
0076 false,
0077 0, true);
0078
0079
0080
0081
0082 auto* colourPointer = new G4HumanPhantomColour();
0083 G4Colour colour = colourPointer -> GetColour(colourName);
0084 delete colourPointer;
0085
0086 auto* LeftTesteVisAtt = new G4VisAttributes(colour);
0087 LeftTesteVisAtt->SetForceSolid(wireFrame);
0088 logicLeftTeste->SetVisAttributes(LeftTesteVisAtt);
0089
0090 G4cout << "LeftTeste created !!!!!!" << G4endl;
0091
0092
0093 G4double LeftTesteVol = logicLeftTeste->GetSolid()->GetCubicVolume();
0094 G4cout << "Volume of LeftTeste = " << LeftTesteVol/cm3 << " cm^3" << G4endl;
0095
0096
0097 G4String LeftTesteMat = logicLeftTeste->GetMaterial()->GetName();
0098 G4cout << "Material of LeftTeste = " << LeftTesteMat << G4endl;
0099
0100
0101 G4double LeftTesteDensity = logicLeftTeste->GetMaterial()->GetDensity();
0102 G4cout << "Density of Material = " << LeftTesteDensity*cm3/g << " g/cm^3" << G4endl;
0103
0104
0105 G4double LeftTesteMass = (LeftTesteVol)*LeftTesteDensity;
0106 G4cout << "Mass of LeftTeste = " << LeftTesteMass/gram << " g" << G4endl;
0107
0108 return physLeftTeste;
0109 }