File indexing completed on 2025-01-18 09:16:53
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 "ICRP110PhantomVisAction.hh"
0032
0033 #include "G4VisManager.hh"
0034 #include "ICRP110PhantomConstruction.hh"
0035 #include "ICRP110PhantomPseudoScene.hh"
0036 #include "G4PhysicalVolumeModel.hh"
0037 #include "G4Polymarker.hh"
0038 #include "G4SystemOfUnits.hh"
0039 #include "Randomize.hh"
0040
0041 void ICRP110PhantomVisAction::Draw ()
0042 {
0043
0044 G4VVisManager* pVisManager = G4VVisManager::GetConcreteInstance();
0045 if (!pVisManager) return;
0046
0047 auto container = fpDetectorConstruction->GetPhantumContainer();
0048
0049 static G4bool first = true;
0050 if (first && container) {
0051 first = false;
0052
0053 G4ModelingParameters tmpMP;
0054 tmpMP.SetCulling(true);
0055 tmpMP.SetCullingInvisible(true);
0056 const G4bool useFullExtent = true;
0057 G4PhysicalVolumeModel tmpPVModel
0058 (container,
0059 G4PhysicalVolumeModel::UNLIMITED,
0060 G4Transform3D(),
0061 &tmpMP,
0062 useFullExtent);
0063 ICRP110PhantomPseudoScene tmpScene
0064 (&tmpPVModel
0065 ,fPositionByMaterial
0066 ,fColourByMaterial);
0067 tmpPVModel.DescribeYourselfTo(tmpScene);
0068
0069 for (const auto& entry: fColourByMaterial) {
0070 fSetOfMaterials.insert(entry.first);
0071 }
0072 for (const auto& material: fSetOfMaterials) {
0073 G4cout << "fSetOfMaterials: " << material->GetName() << G4endl;
0074 }
0075 }
0076
0077
0078 const G4double voxelHalfDimX = fpDetectorConstruction->GetVoxelHalfDimX()*mm;
0079 const G4double voxelHalfDimY = fpDetectorConstruction->GetVoxelHalfDimY()*mm;
0080 const G4double voxelHalfDimZ = fpDetectorConstruction->GetVoxelHalfDimZ()*mm;
0081
0082 static G4bool firstPrint = true;
0083 G4int nDotsTotal = 0;
0084 const G4int nDotsPerPositionMin = 1;
0085 const G4int nDotsPerPositionMax = 10;
0086 const G4double densityMin = 1.*g/cm3;
0087 const G4double densityMax = 3.*g/cm3;;
0088 for (const auto& material: fSetOfMaterials) {
0089 auto& colour = fColourByMaterial[material];
0090 G4int nDots = 0;
0091 G4Polymarker dots;
0092 dots.SetVisAttributes(G4Colour(colour));
0093 dots.SetMarkerType(G4Polymarker::dots);
0094 dots.SetSize(G4VMarker::screen,1.);
0095 dots.SetInfo(material->GetName());
0096 const G4double currentDensity = material->GetDensity();
0097 const auto range = fPositionByMaterial.equal_range(material);
0098 for (auto posByMat = range.first; posByMat != range.second; ++posByMat) {
0099 G4int nDotsPerPosition
0100 = nDotsPerPositionMin
0101 + (nDotsPerPositionMax-nDotsPerPositionMin)
0102 *(currentDensity-densityMin)/(densityMax-densityMin);
0103 nDotsPerPosition = std::max(nDotsPerPositionMin,nDotsPerPosition);
0104 nDotsPerPosition = std::min(nDotsPerPositionMax,nDotsPerPosition);
0105 for (G4int i = 0; i < nDotsPerPosition; ++i) {
0106 const G4double x = posByMat->second.getX() + (2.*G4UniformRand()-1.)*voxelHalfDimX;
0107 const G4double y = posByMat->second.getY() + (2.*G4UniformRand()-1.)*voxelHalfDimY;
0108 const G4double z = posByMat->second.getZ() + (2.*G4UniformRand()-1.)*voxelHalfDimZ;
0109 dots.push_back(G4ThreeVector(x,y,z));
0110 ++nDots;
0111 }
0112 }
0113 pVisManager->Draw(dots);
0114 if (firstPrint) {
0115 G4cout
0116 << "Number of dots for material " << material->GetName()
0117 << ", colour " << colour
0118 << ": " << nDots << G4endl;
0119 }
0120 nDotsTotal += nDots;
0121 }
0122 if (firstPrint) {
0123 G4cout << "Total number of dots: " << nDotsTotal << G4endl;
0124 firstPrint = false;
0125 }
0126 }