File indexing completed on 2026-04-29 07:39:22
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 #include "DetectorConstruction.hh"
0030
0031 #include "DetectorMessenger.hh"
0032 #include "PrimaryGeneratorAction.hh"
0033
0034 #include "G4Box.hh"
0035 #include "G4GeometryManager.hh"
0036 #include "G4LogicalVolume.hh"
0037 #include "G4LogicalVolumeStore.hh"
0038 #include "G4Material.hh"
0039 #include "G4NistManager.hh"
0040 #include "G4PVPlacement.hh"
0041 #include "G4PhysicalConstants.hh"
0042 #include "G4PhysicalVolumeStore.hh"
0043 #include "G4RunManager.hh"
0044 #include "G4SolidStore.hh"
0045 #include "G4SystemOfUnits.hh"
0046 #include "G4ThreeVector.hh"
0047 #include "G4Tubs.hh"
0048 #include "globals.hh"
0049
0050
0051
0052 DetectorConstruction::DetectorConstruction()
0053 : fMaterial(nullptr),
0054 fExperimentalHall_log(nullptr),
0055 fExperimentalHall_phys(nullptr),
0056 fLogicLayer(nullptr),
0057 fPhysiLayer(nullptr),
0058 fLogicScoringUpDown(nullptr),
0059 fPhysiScoringUpstream(nullptr),
0060 fPhysiScoringDownstream(nullptr),
0061 fLogicScoringSide(nullptr),
0062 fPhysiScoringSide(nullptr),
0063 fDetectorMessenger(nullptr),
0064 fThickness(2.0 * CLHEP::m),
0065 fDiameter(2.0 * CLHEP::m)
0066 {
0067 fMaterial = G4NistManager::Instance()->FindOrBuildMaterial("G4_Fe");
0068
0069 fDetectorMessenger = new DetectorMessenger(this);
0070 }
0071
0072
0073
0074 DetectorConstruction::~DetectorConstruction()
0075 {
0076 delete fDetectorMessenger;
0077 }
0078
0079
0080
0081 G4VPhysicalVolume* DetectorConstruction::Construct()
0082 {
0083 return ConstructLayer();
0084 }
0085
0086
0087
0088 G4VPhysicalVolume* DetectorConstruction::ConstructLayer()
0089 {
0090
0091 G4GeometryManager::GetInstance()->OpenGeometry();
0092 G4PhysicalVolumeStore::GetInstance()->Clean();
0093 G4LogicalVolumeStore::GetInstance()->Clean();
0094 G4SolidStore::GetInstance()->Clean();
0095
0096
0097
0098
0099
0100
0101 G4double expHall_x = 0.6 * fDiameter;
0102
0103 G4double expHall_y = 0.6 * fDiameter;
0104
0105 G4double expHall_z = 0.6 * fThickness;
0106
0107
0108 G4Material* vacuum = G4NistManager::Instance()->FindOrBuildMaterial("G4_Galactic");
0109
0110
0111 G4Box* experimentalHall_box = new G4Box("expHall_box", expHall_x, expHall_y, expHall_z);
0112 fExperimentalHall_log = new G4LogicalVolume(experimentalHall_box,
0113 vacuum,
0114 "expHall_log",
0115 0,
0116 0,
0117 0);
0118 fExperimentalHall_phys = new G4PVPlacement(0,
0119 G4ThreeVector(),
0120 "expHall",
0121 fExperimentalHall_log,
0122 0,
0123 false,
0124 0);
0125
0126
0127 G4Tubs* solidLayer = new G4Tubs("solidLayer",
0128 0.0,
0129 0.5 * fDiameter,
0130 0.5 * fThickness,
0131 0.0,
0132 2.0 * pi);
0133 fLogicLayer = new G4LogicalVolume(solidLayer,
0134 fMaterial,
0135 "logicLayer",
0136 0,
0137 0,
0138 0);
0139 fPhysiLayer = new G4PVPlacement(0,
0140 G4ThreeVector(),
0141 "physiLayer",
0142 fLogicLayer,
0143 fExperimentalHall_phys,
0144 false,
0145 0);
0146
0147
0148
0149
0150 G4Tubs* solidScoringUpDown = new G4Tubs("solidScoringUpDown",
0151 0.0,
0152 0.5 * fDiameter,
0153 0.5 * fScoringThickness,
0154 0.0,
0155 2.0 * pi);
0156 fLogicScoringUpDown = new G4LogicalVolume(solidScoringUpDown,
0157 vacuum,
0158 "logicScoringUpDown",
0159 0,
0160 0,
0161 0);
0162 G4double zScoringUpDown = 0.5 * (fThickness + fScoringThickness);
0163 fPhysiScoringUpstream = new G4PVPlacement(0,
0164 G4ThreeVector(0.0, 0.0, -zScoringUpDown),
0165
0166 "physiScoringUpstream",
0167 fLogicScoringUpDown,
0168 fExperimentalHall_phys,
0169 false,
0170 0);
0171 fPhysiScoringDownstream = new G4PVPlacement(0,
0172 G4ThreeVector(0.0, 0.0, zScoringUpDown),
0173
0174 "physiScoringDownstream",
0175 fLogicScoringUpDown,
0176 fExperimentalHall_phys,
0177 false,
0178 0);
0179
0180 G4Tubs* solidScoringSide = new G4Tubs("solidScoringSide",
0181 0.5 * fDiameter,
0182 0.5 * fDiameter + fScoringThickness,
0183 0.5 * fThickness,
0184 0.0,
0185 2.0 * pi);
0186 fLogicScoringSide = new G4LogicalVolume(solidScoringSide,
0187 vacuum,
0188 "logicScoringSide",
0189 0,
0190 0,
0191 0);
0192 fPhysiScoringSide = new G4PVPlacement(0,
0193 G4ThreeVector(0.0, 0.0, 0.0),
0194 "physiScoringSide",
0195 fLogicScoringSide,
0196 fExperimentalHall_phys,
0197 false,
0198 0);
0199
0200 G4cout << G4endl << "DetectorConstruction::ConstructLayer() : " << G4endl
0201 << "\t World (box) size: " << G4endl << "\t \t x : -/+ " << expHall_x << " mm ;"
0202 << "\t y : -/+ " << expHall_y << " mm ;"
0203 << "\t z : -/+ " << expHall_z << " mm ;" << G4endl
0204 << "\t Target layer (cylinder) size: " << G4endl << "\t \t x : -/+ " << 0.5 * fDiameter
0205 << " mm ;"
0206 << "\t y : -/+ " << 0.5 * fDiameter << " mm ;"
0207 << "\t z : -/+ " << 0.5 * fThickness << " mm ;" << G4endl << G4endl << G4endl;
0208
0209 return fExperimentalHall_phys;
0210 }
0211
0212
0213
0214 void DetectorConstruction::SetMaterial(const G4String name)
0215 {
0216 fMaterial = G4NistManager::Instance()->FindOrBuildMaterial(name);
0217 if (!fMaterial) {
0218 G4cout << G4endl << G4endl << "WARNING: the name of the material has not been recognized!"
0219 << G4endl << " ===> the default * G4_Fe * will be used." << G4endl << G4endl;
0220 fMaterial = G4NistManager::Instance()->FindOrBuildMaterial("G4_Fe");
0221 }
0222 if (fLogicLayer) fLogicLayer->SetMaterial(fMaterial);
0223 }
0224
0225
0226
0227 void DetectorConstruction::UpdateGeometry()
0228 {
0229 G4RunManager::GetRunManager()->ReinitializeGeometry();
0230 PrintParameters();
0231
0232 const PrimaryGeneratorAction* pPrimaryAction = dynamic_cast<const PrimaryGeneratorAction*>(
0233 G4RunManager::GetRunManager()->GetUserPrimaryGeneratorAction());
0234 if (pPrimaryAction) pPrimaryAction->SetGunPosition();
0235 }
0236
0237
0238
0239 void DetectorConstruction::PrintParameters()
0240 {
0241 G4cout << G4endl << G4endl << " ------ DetectorConstruction::PrintParameters() ------ " << G4endl
0242 << " Material = " << fMaterial->GetName() << G4endl
0243 << " Thickness = " << fThickness << " mm" << G4endl
0244 << " Diameter = " << fDiameter << " mm" << G4endl
0245 << " ScoringThickness = " << fScoringThickness << " mm" << G4endl
0246 << " ------------------------------------------------------ " << G4endl << G4endl;
0247 }
0248
0249