File indexing completed on 2025-02-23 09:22:30
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
0040
0041
0042
0043
0044
0045
0046
0047
0048
0049
0050
0051
0052
0053
0054 #include "TSDetectorConstruction.hh"
0055
0056 #include "G4Box.hh"
0057 #include "G4Colour.hh"
0058 #include "G4LogicalVolume.hh"
0059 #include "G4Material.hh"
0060 #include "G4MultiFunctionalDetector.hh"
0061 #include "G4NistManager.hh"
0062 #include "G4PSEnergyDeposit.hh"
0063 #include "G4PSNofStep.hh"
0064 #include "G4PVPlacement.hh"
0065 #include "G4RunManager.hh"
0066 #include "G4SDManager.hh"
0067 #include "G4UnitsTable.hh"
0068 #include "G4UserLimits.hh"
0069 #include "G4VPhysicalVolume.hh"
0070 #include "G4VisAttributes.hh"
0071
0072 using namespace CLHEP;
0073
0074
0075
0076 TSDetectorConstruction* TSDetectorConstruction::fgInstance = 0;
0077
0078
0079
0080 TSDetectorConstruction* TSDetectorConstruction::Instance()
0081 {
0082 return fgInstance;
0083 }
0084
0085
0086
0087 TSDetectorConstruction::TSDetectorConstruction()
0088 : fWorldPhys(0),
0089 fWorldMaterialName("G4_Galactic"),
0090 fTargetMaterialName("G4_B"),
0091 fCasingMaterialName("G4_WATER"),
0092 fWorldDim(G4ThreeVector(0.5 * m, 0.5 * m, 0.5 * m)),
0093 fTargetDim(G4ThreeVector(0.5 * m, 0.5 * m, 0.5 * m)),
0094 fTargetSections(G4ThreeVector(5, 5, 5)),
0095 fMfdName("Target_MFD")
0096 {
0097 fgInstance = this;
0098 }
0099
0100
0101
0102 TSDetectorConstruction::~TSDetectorConstruction()
0103 {
0104 fgInstance = 0;
0105 }
0106
0107
0108
0109 G4VPhysicalVolume* TSDetectorConstruction::Construct()
0110 {
0111 return ConstructWorld(ConstructMaterials());
0112 }
0113
0114
0115
0116 TSDetectorConstruction::MaterialCollection_t TSDetectorConstruction::ConstructMaterials()
0117 {
0118 MaterialCollection_t materials;
0119 G4NistManager* nist = G4NistManager::Instance();
0120
0121 materials["World"] = nist->FindOrBuildMaterial(fWorldMaterialName);
0122 materials["Target"] = nist->FindOrBuildMaterial(fTargetMaterialName);
0123 materials["Casing"] = nist->FindOrBuildMaterial(fCasingMaterialName);
0124
0125 return materials;
0126 }
0127
0128
0129
0130 G4VPhysicalVolume* TSDetectorConstruction::ConstructWorld(const MaterialCollection_t& materials)
0131 {
0132 G4UserLimits* steplimit = new G4UserLimits(0.1 * (fTargetDim.z() / fTargetSections.z()));
0133 G4bool check_overlap = false;
0134
0135 G4Box* world_solid =
0136 new G4Box("World", 0.5 * fWorldDim.x(), 0.5 * fWorldDim.y(), 0.5 * fWorldDim.z());
0137 G4LogicalVolume* world_log =
0138 new G4LogicalVolume(world_solid, materials.find("World")->second, "World");
0139 fWorldPhys =
0140 new G4PVPlacement(0, G4ThreeVector(0.), "World", world_log, 0, false, 0, check_overlap);
0141
0142 G4int nz = fTargetSections.z();
0143 G4int ny = fTargetSections.y();
0144 G4int nx = fTargetSections.x();
0145
0146
0147 G4double sx = fTargetDim.x() / fTargetSections.x();
0148 G4double sy = fTargetDim.y() / fTargetSections.y();
0149 G4double sz = fTargetDim.z() / fTargetSections.z();
0150
0151
0152
0153
0154
0155
0156
0157 G4VisAttributes* red = new G4VisAttributes(G4Color(1., 0., 0., 1.0));
0158 G4VisAttributes* green = new G4VisAttributes(G4Color(0., 1., 0., 0.25));
0159 G4VisAttributes* blue = new G4VisAttributes(G4Color(0., 0., 1., 0.1));
0160 G4VisAttributes* white = new G4VisAttributes(G4Color(1., 1., 1., 1.));
0161
0162 white->SetVisibility(true);
0163 red->SetVisibility(true);
0164 green->SetVisibility(true);
0165 blue->SetVisibility(true);
0166
0167 white->SetForceWireframe(true);
0168 red->SetForceSolid(true);
0169 green->SetForceSolid(true);
0170 blue->SetForceSolid(true);
0171
0172 world_log->SetVisAttributes(white);
0173
0174 for (G4int k = 0; k < nz; ++k)
0175 for (G4int j = 0; j < ny; ++j)
0176 for (G4int i = 0; i < nx; ++i) {
0177
0178 G4double dx = 0.5 * sx + static_cast<G4double>(i) * sx - 0.5 * fWorldDim.x();
0179 G4double dy = 0.5 * sy + static_cast<G4double>(j) * sy - 0.5 * fWorldDim.y();
0180 G4double dz = 0.5 * sz + static_cast<G4double>(k) * sz - 0.5 * fWorldDim.z();
0181 G4ThreeVector td = G4ThreeVector(dx, dy, -dz);
0182
0183 std::stringstream ss_name;
0184 ss_name << "Target_" << i << "_" << j << "_" << k;
0185
0186 G4Box* target_solid = new G4Box(ss_name.str(), 0.5 * sx, 0.5 * sy, 0.5 * sz);
0187
0188 G4Material* target_material = 0;
0189 G4bool is_casing = true;
0190
0191 if (j == 0 || j + 1 == ny || i == 0 || i + 1 == nx || (nz > 1 && (k == 0 || k + 1 == nz)))
0192 target_material = materials.find("Casing")->second;
0193 else {
0194 target_material = materials.find("Target")->second;
0195 is_casing = false;
0196 }
0197
0198 G4LogicalVolume* target_log =
0199 new G4LogicalVolume(target_solid, target_material, ss_name.str());
0200
0201 target_log->SetUserLimits(steplimit);
0202
0203 new G4PVPlacement(0, td, ss_name.str(), target_log, fWorldPhys, true,
0204 k * nx * ny + j * nx + i, check_overlap);
0205
0206 fScoringVolumes.insert(target_log);
0207
0208 if (is_casing)
0209 target_log->SetVisAttributes(blue);
0210 else {
0211
0212 G4bool even_z = (k % 2 == 0) ? true : false;
0213 G4bool even_y = (j % 2 == 0) ? true : false;
0214 G4bool even_x = (i % 2 == 0) ? true : false;
0215
0216 G4VisAttributes* theColor = nullptr;
0217
0218 if ((even_z)) {
0219 if ((even_y && even_x) || (!even_y && !even_x))
0220 theColor = red;
0221 else
0222 theColor = green;
0223 }
0224 else
0225 {
0226 if ((!even_y && even_x) || (even_y && !even_x))
0227 theColor = red;
0228 else
0229 theColor = green;
0230 }
0231
0232 target_log->SetVisAttributes(theColor);
0233 }
0234 }
0235
0236 return fWorldPhys;
0237 }
0238
0239
0240
0241 void TSDetectorConstruction::ConstructSDandField()
0242 {
0243
0244
0245
0246
0247 G4MultiFunctionalDetector* MFDet = new G4MultiFunctionalDetector(fMfdName);
0248 G4SDManager::GetSDMpointer()->AddNewDetector(MFDet);
0249 G4VPrimitiveScorer* edep = new G4PSEnergyDeposit("EnergyDeposit");
0250 MFDet->RegisterPrimitive(edep);
0251 G4VPrimitiveScorer* nstep = new G4PSNofStep("NumberOfSteps");
0252 MFDet->RegisterPrimitive(nstep);
0253
0254
0255 for (auto ite : fScoringVolumes) {
0256 SetSensitiveDetector(ite, MFDet);
0257 }
0258 }
0259
0260