File indexing completed on 2025-02-23 09:21:19
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 #include "F06DetectorConstruction.hh"
0035
0036 #include "G4Box.hh"
0037 #include "G4Colour.hh"
0038 #include "G4FieldManager.hh"
0039 #include "G4GeometryManager.hh"
0040 #include "G4LogicalVolume.hh"
0041 #include "G4LogicalVolumeStore.hh"
0042 #include "G4Material.hh"
0043 #include "G4NistManager.hh"
0044 #include "G4PVPlacement.hh"
0045 #include "G4PhysicalVolumeStore.hh"
0046 #include "G4RepleteEofM.hh"
0047 #include "G4SolidStore.hh"
0048 #include "G4SystemOfUnits.hh"
0049 #include "G4TransportationManager.hh"
0050 #include "G4UniformGravityField.hh"
0051 #include "G4UserLimits.hh"
0052 #include "G4VisAttributes.hh"
0053
0054
0055 #include "G4ChordFinder.hh"
0056 #include "G4ClassicalRK4.hh"
0057 #include "G4IntegrationDriver.hh"
0058 #include "G4MagIntegratorStepper.hh"
0059 #include "G4PropagatorInField.hh"
0060
0061
0062
0063 F06DetectorConstruction::F06DetectorConstruction() : fVacuum(nullptr)
0064 {
0065
0066 DefineMaterials();
0067 }
0068
0069
0070
0071 F06DetectorConstruction::~F06DetectorConstruction()
0072 {
0073 delete fField;
0074 }
0075
0076
0077
0078 void F06DetectorConstruction::DefineMaterials()
0079 {
0080 G4NistManager* nistMan = G4NistManager::Instance();
0081
0082 fVacuum = nistMan->FindOrBuildMaterial("G4_Galactic");
0083
0084 G4cout << *(G4Material::GetMaterialTable()) << G4endl;
0085 }
0086
0087
0088
0089 G4VPhysicalVolume* F06DetectorConstruction::Construct()
0090 {
0091
0092
0093
0094
0095 G4double expHall_x = 1.0 * m;
0096 G4double expHall_y = 1.0 * m;
0097 G4double expHall_z = 1.0 * m;
0098
0099 auto solidWorld = new G4Box("World",
0100 expHall_x, expHall_y, expHall_z);
0101
0102 auto logicWorld = new G4LogicalVolume(solidWorld,
0103 fVacuum,
0104 "World");
0105
0106 auto physiWorld = new G4PVPlacement(nullptr,
0107 G4ThreeVector(),
0108 logicWorld,
0109 "World",
0110 nullptr,
0111 false,
0112 0);
0113
0114 G4double maxStep = 1.0 * mm;
0115 G4double maxTime = 41. * s;
0116
0117 auto stepLimit = new G4UserLimits(maxStep, DBL_MAX, maxTime);
0118
0119 logicWorld->SetUserLimits(stepLimit);
0120
0121
0122
0123
0124
0125
0126
0127
0128
0129 return physiWorld;
0130 }
0131
0132
0133
0134 G4ThreadLocal G4UniformGravityField* F06DetectorConstruction::fField = nullptr;
0135
0136 void F06DetectorConstruction::ConstructSDandField()
0137 {
0138 using StepperType = G4ClassicalRK4;
0139
0140 if (!fField) {
0141 fField = new G4UniformGravityField();
0142
0143 auto equation = new G4RepleteEofM(fField);
0144
0145
0146 G4TransportationManager* transportMgr = G4TransportationManager::GetTransportationManager();
0147
0148 G4FieldManager* fieldManager = transportMgr->GetFieldManager();
0149 fieldManager->SetDetectorField(fField);
0150
0151 const int nVar = 8;
0152 auto stepper = new StepperType(equation, nVar);
0153
0154 G4double minStep = 0.01 * mm;
0155 G4ChordFinder* chordFinder = nullptr;
0156 if (stepper) {
0157 auto intgrDriver =
0158 new G4IntegrationDriver<StepperType>(minStep, stepper, stepper->GetNumberOfVariables());
0159 if (intgrDriver) {
0160 chordFinder = new G4ChordFinder(intgrDriver);
0161 }
0162 }
0163
0164
0165
0166
0167
0168 G4double deltaChord = 3.0 * mm;
0169 chordFinder->SetDeltaChord(deltaChord);
0170
0171 G4double deltaIntersection = 0.1 * mm;
0172 fieldManager->SetDeltaIntersection(deltaIntersection);
0173
0174
0175
0176 G4double deltaOneStep = 0.01 * mm;
0177 fieldManager->SetAccuraciesWithDeltaOneStep(deltaOneStep);
0178
0179 G4double epsMax = 1.0e-4;
0180 G4double epsMin = 2.5e-7;
0181 fieldManager->SetMinimumEpsilonStep(epsMin);
0182 fieldManager->SetMaximumEpsilonStep(epsMax);
0183
0184
0185
0186 fieldManager->SetChordFinder(chordFinder);
0187 }
0188 }
0189
0190