File indexing completed on 2026-05-02 07:42:12
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 "DetectorALICE06.hh"
0032 #include "DetectorBari05.hh"
0033 #include "DetectorBarr90.hh"
0034 #include "DetectorHarris73.hh"
0035 #include "DetectorMessenger.hh"
0036 #include "DetectorSimpleALICE.hh"
0037 #include "DetectorWatase86.hh"
0038 #include "Materials.hh"
0039
0040 #include "G4LogicalVolume.hh"
0041 #include "G4LogicalVolumeStore.hh"
0042 #include "G4ProductionCuts.hh"
0043 #include "G4Region.hh"
0044 #include "G4RegionStore.hh"
0045 #include "G4SystemOfUnits.hh"
0046 #include "G4ios.hh"
0047
0048
0049
0050 DetectorConstruction::DetectorConstruction()
0051 : G4VUserDetectorConstruction(), fRadiatorDescription(0), fSetUp("simpleALICE")
0052 {
0053 fDetectorMessenger = new DetectorMessenger(this);
0054 }
0055
0056
0057
0058 DetectorConstruction::~DetectorConstruction()
0059 {
0060 delete fRadiatorDescription;
0061 delete fDetectorMessenger;
0062 delete Materials::GetInstance();
0063 }
0064
0065
0066
0067 G4VPhysicalVolume* DetectorConstruction::Construct()
0068 {
0069 G4VPhysicalVolume* world = nullptr;
0070 if (fSetUp == "simpleALICE") {
0071 DetectorSimpleALICE simpleALICE;
0072 world = simpleALICE.Construct();
0073 fRadiatorDescription = simpleALICE.GetRadiatorDescription();
0074 }
0075 else if (fSetUp == "alice06") {
0076 DetectorALICE06 alice06;
0077 world = alice06.Construct();
0078 fRadiatorDescription = alice06.GetRadiatorDescription();
0079 }
0080 else if (fSetUp == "bari05") {
0081 DetectorBari05 bari05;
0082 world = bari05.Construct();
0083 fRadiatorDescription = bari05.GetRadiatorDescription();
0084 }
0085 else if (fSetUp == "harris73") {
0086 DetectorHarris73 harris73;
0087 world = harris73.Construct();
0088 fRadiatorDescription = harris73.GetRadiatorDescription();
0089 }
0090 else if (fSetUp == "watase86") {
0091 DetectorWatase86 watase86;
0092 world = watase86.Construct();
0093 fRadiatorDescription = watase86.GetRadiatorDescription();
0094 }
0095 else if (fSetUp == "barr90") {
0096 DetectorBarr90 barr90;
0097 world = barr90.Construct();
0098 fRadiatorDescription = barr90.GetRadiatorDescription();
0099 }
0100 else {
0101 G4cout << "Experimental setup is unsupported. Check /XTRdetector/setup " << G4endl;
0102 G4cout << "Run default: simpleALICE" << G4endl;
0103
0104 DetectorSimpleALICE simpleALICE;
0105 world = simpleALICE.Construct();
0106 fRadiatorDescription = simpleALICE.GetRadiatorDescription();
0107 }
0108 G4RegionStore* regionStore = G4RegionStore::GetInstance();
0109 size_t nRegions = regionStore->size();
0110 G4double cut = 1. * CLHEP::mm;
0111 for (size_t i = 0; i < nRegions; ++i) {
0112 G4Region* reg = (*regionStore)[i];
0113 if (!reg->GetProductionCuts()) {
0114 G4ProductionCuts* cuts = new G4ProductionCuts();
0115 cuts->SetProductionCut(cut, "gamma");
0116 cuts->SetProductionCut(cut, "e-");
0117 cuts->SetProductionCut(cut, "e+");
0118 cuts->SetProductionCut(cut, "proton");
0119 reg->SetProductionCuts(cuts);
0120 }
0121 }
0122 return world;
0123 }
0124
0125
0126
0127 RadiatorDescription* DetectorConstruction::GetRadiatorDescription() const
0128 {
0129 if (!fRadiatorDescription) {
0130 G4cout << "RadiatorDescription is not defined" << G4endl;
0131 }
0132 return fRadiatorDescription;
0133 }
0134
0135
0136
0137 G4Material* DetectorConstruction::GetAbsorberMaterial() const
0138 {
0139 G4LogicalVolume* lv = G4LogicalVolumeStore::GetInstance()->GetVolume("Absorber");
0140
0141 if (!lv) {
0142 G4cerr << "Absorber logical volume is not defined." << G4endl;
0143 return 0;
0144 }
0145
0146 return lv->GetMaterial();
0147 }
0148
0149