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 "G4MIRDLeftLung.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 "G4SubtractionSolid.hh"
0047 #include "G4Box.hh"
0048 #include "G4HumanPhantomColour.hh"
0049
0050 G4VPhysicalVolume* G4MIRDLeftLung::Construct(const G4String& volumeName,G4VPhysicalVolume* mother,
0051 const G4String& colourName, G4bool wireFrame,G4bool)
0052 {
0053
0054 G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
0055
0056 auto* material = new G4HumanPhantomMaterial();
0057 auto* lung_material = material -> GetMaterial("lung_material");
0058 delete material;
0059
0060 G4double ax = 5. *cm;
0061 G4double by = 7.5 *cm;
0062 G4double cz = 24.*cm;
0063 G4double zcut1 = 0.0 *cm;
0064 G4double zcut2=24. *cm;
0065
0066 auto* oneLung = new G4Ellipsoid("OneLung",ax, by, cz, zcut1,zcut2);
0067
0068 ax= 5.*cm;
0069 by= 7.5*cm;
0070 cz= 24.*cm;
0071
0072
0073 auto* subtrLung = new G4Ellipsoid("subtrLung",ax, by, cz);
0074
0075
0076
0077 G4double dx = 5.5* cm;
0078 G4double dy = 8.5 * cm;
0079 G4double dz = 24. * cm;
0080
0081 auto* box = new G4Box("Box", dx, dy, dz);
0082
0083
0084
0085 auto* section2 = new G4SubtractionSolid("BoxSub2", subtrLung, box, nullptr, G4ThreeVector(0.*cm, 8.5* cm, 0.*cm));
0086
0087
0088
0089
0090
0091 auto* lung2 = new G4SubtractionSolid("Lung2", oneLung,
0092 section2,
0093 nullptr, G4ThreeVector(-6.*cm,0*cm,0.0*cm));
0094
0095
0096
0097
0098
0099
0100
0101
0102
0103 auto* logicLeftLung = new G4LogicalVolume(lung2,lung_material,
0104 "logical" + volumeName, nullptr, nullptr, nullptr);
0105
0106
0107 G4VPhysicalVolume* physLeftLung = new G4PVPlacement(nullptr,G4ThreeVector(8.50 *cm, 0.0*cm, 8.5*cm),
0108 "physicalLeftLung",
0109 logicLeftLung,
0110 mother,
0111 false,
0112 0, true);
0113
0114
0115
0116
0117 auto* colourPointer = new G4HumanPhantomColour();
0118 G4Colour colour = colourPointer -> GetColour(colourName);
0119 auto* LeftLungVisAtt = new G4VisAttributes(colour);
0120 LeftLungVisAtt->SetForceSolid(wireFrame);
0121 logicLeftLung->SetVisAttributes(LeftLungVisAtt);
0122
0123 G4cout << "LeftLung created !!!!!!" << G4endl;
0124
0125
0126 G4double LeftLungVol = logicLeftLung->GetSolid()->GetCubicVolume();
0127
0128 G4cout << "Volume of LeftLung = " << (LeftLungVol)/cm3 << " cm^3" << G4endl;
0129
0130
0131 G4String LeftLungMat = logicLeftLung->GetMaterial()->GetName();
0132 G4cout << "Material of LeftLung = " << LeftLungMat << G4endl;
0133
0134
0135 G4double LeftLungDensity = logicLeftLung->GetMaterial()->GetDensity();
0136 G4cout << "Density of Material = " << LeftLungDensity*cm3/g << " g/cm^3" << G4endl;
0137
0138
0139 G4double LeftLungMass = (LeftLungVol)*LeftLungDensity;
0140 G4cout << "Mass of LeftLung = " << LeftLungMass/gram << " g" << G4endl;
0141
0142 return physLeftLung;
0143 }