File indexing completed on 2025-10-13 08:27:57
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
0033
0034
0035
0036
0037
0038
0039 #include "DetectorConstruction.hh"
0040
0041 #include "ScoreBasicMoleculeCounts.hh"
0042 #include "ScoreBasicReactionCounts.hh"
0043
0044 #include "G4Box.hh"
0045 #include "G4Ellipsoid.hh"
0046 #include "G4LogicalVolume.hh"
0047 #include "G4MultiFunctionalDetector.hh"
0048 #include "G4NistManager.hh"
0049 #include "G4Orb.hh"
0050 #include "G4PVPlacement.hh"
0051 #include "G4PhysicalConstants.hh"
0052 #include "G4PhysicalVolumeStore.hh"
0053 #include "G4SDManager.hh"
0054 #include "G4SystemOfUnits.hh"
0055 #include "G4VPhysicalVolume.hh"
0056 #include "G4VPrimitiveScorer.hh"
0057 #include "Randomize.hh"
0058
0059 G4VPhysicalVolume* DetectorConstruction::Construct()
0060 {
0061
0062 auto man = G4NistManager::Instance();
0063 auto water = man->FindOrBuildMaterial("G4_WATER");
0064
0065
0066
0067
0068 const G4double worldXYZ = 1 * m;
0069
0070 auto solidWorld = new G4Box("World", 0.5 * worldXYZ, 0.5 * worldXYZ, 0.5 * worldXYZ);
0071 auto lvWorld = new G4LogicalVolume(solidWorld, water, "World");
0072 auto pvWorld =
0073 new G4PVPlacement(nullptr, G4ThreeVector(), lvWorld, "World", nullptr, false, 0, true);
0074
0075
0076
0077 ConstructCell(pvWorld);
0078
0079
0080 return pvWorld;
0081 }
0082
0083 void DetectorConstruction::ConstructCell(G4VPhysicalVolume* pvWorld)
0084 {
0085 const G4double cellRadius = 7 * um;
0086 const G4double nucleusRadius = 4 * um;
0087
0088 const G4int nMitochondria = 100;
0089 const G4double mitoA = 0.55 * micrometer;
0090 const G4double mitoB = 0.25 * micrometer;
0091 const G4double mitoC = 0.90 * micrometer;
0092
0093 auto solidCell = new G4Orb("Cell", cellRadius);
0094 auto lvCell = new G4LogicalVolume(
0095 solidCell, G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER"), "Cell");
0096 auto pvCell =
0097 new G4PVPlacement(nullptr, G4ThreeVector(), "Cell", lvCell, pvWorld, false, 0, true);
0098
0099 auto solidNucleus = new G4Orb("Nucleus", nucleusRadius);
0100 auto lvNucleus = new G4LogicalVolume(
0101 solidNucleus, G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER"), "Nucleus");
0102
0103 new G4PVPlacement(nullptr, G4ThreeVector(), "Nucleus", lvNucleus, pvCell, false, 0, true);
0104
0105 auto solidMito = new G4Ellipsoid("Mitochondria", mitoA, mitoB, mitoC);
0106 auto lvMito = new G4LogicalVolume(
0107 solidMito, G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER"), "Mitochondria");
0108
0109 for (auto i = 0; i < nMitochondria; ++i) {
0110 G4bool overlap = true;
0111 do {
0112 auto u = twopi * G4UniformRand();
0113 auto v = std::acos(2 * G4UniformRand() - 1);
0114 auto dr = G4UniformRand() * cellRadius;
0115 auto x = dr * std::cos(u) * std::sin(v);
0116 auto y = dr * std::sin(u) * std::sin(v);
0117 auto z = dr * std::cos(v);
0118 auto pos = G4ThreeVector(x, y, z);
0119
0120 auto phi = G4UniformRand() * 2 * pi;
0121 auto psi = G4UniformRand() * 2 * pi;
0122 auto rot = new G4RotationMatrix();
0123 rot->rotateX(psi);
0124 rot->rotateY(phi);
0125
0126 auto pvMito = new G4PVPlacement(rot, pos, "Mitochondria", lvMito, pvCell, false, i, false);
0127
0128 overlap = pvMito->CheckOverlaps(1000, 0, false);
0129 if (overlap) {
0130 G4PhysicalVolumeStore::DeRegister(pvMito);
0131 }
0132 } while (overlap);
0133 }
0134 }
0135
0136 void DetectorConstruction::ConstructSDandField()
0137 {
0138 G4SDManager::GetSDMpointer()->SetVerboseLevel(1);
0139
0140 auto mfDetector = new G4MultiFunctionalDetector("mfDetector");
0141
0142 G4VPrimitiveScorer* primitivSpecies = new ScoreBasicMoleculeCounts("BasicMoleculeCounts", 1, "BasicCounter");
0143 mfDetector->RegisterPrimitive(primitivSpecies);
0144 primitivSpecies = new ScoreBasicMoleculeCounts("BasicCounter_VariablePrecision", 1,
0145 "BasicCounter_VariablePrecision");
0146 mfDetector->RegisterPrimitive(primitivSpecies);
0147 primitivSpecies = new ScoreBasicReactionCounts("BasicReactionCounts", 1, "Reactions");
0148 mfDetector->RegisterPrimitive(primitivSpecies);
0149
0150 G4SDManager::GetSDMpointer()->AddNewDetector(mfDetector);
0151 SetSensitiveDetector("Cell", mfDetector);
0152 }
0153