File indexing completed on 2026-04-17 07:52:13
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 #ifndef DetectorConstruction_h
0034 #define DetectorConstruction_h 1
0035
0036 #include "G4VUserDetectorConstruction.hh"
0037 #include "G4ios.hh"
0038 #include "globals.hh"
0039 #include "G4SystemOfUnits.hh"
0040 #include <vector>
0041
0042 #include "G4Region.hh"
0043 #include "G4PVPlacement.hh"
0044
0045 #include "DetectorConstructionMessenger.hh"
0046 #include "G4ChannelingFastSimModel.hh"
0047
0048 #define NSpheresMax 10000
0049
0050 class G4VPhysicalVolume;
0051 class G4LogicalVolume;
0052
0053
0054
0055
0056
0057 class DetectorConstruction : public G4VUserDetectorConstruction
0058 {
0059 public:
0060 DetectorConstruction();
0061 ~DetectorConstruction() override = default;
0062
0063 G4VPhysicalVolume* Construct() override;
0064 void ConstructSDandField() override;
0065
0066
0067 std::vector<G4LogicalVolume*> GetScoringVolume() const {
0068 return fScoringVolume;}
0069
0070
0071 void SetHybridSource(G4bool val) {fHybridSource = val;}
0072
0073
0074 void SetCrystalMaterial(G4String val) {fCrystalMaterialStr = val;}
0075 void SetCrystalSize(G4ThreeVector val) {fCrystalSize = val;}
0076 void SetCrystalBendingAngle(G4double val) {fBendingAngle = val;}
0077 void SetCrystalLattice(G4String val) {fLattice = val;}
0078 void SetCrystalAngleX(G4double val) {fAngleX = val;}
0079 void SetCrystalAngleY(G4double val) {fAngleY = val;}
0080 G4double GetCrystalZ() const {return fCrystalZ;}
0081 void SetRadiationModel(G4bool val) {fActivateRadiationModel = val;}
0082 void SetOCeffects(G4bool val) {fActivateOCeffects = val;}
0083 G4bool GetOCeffects() const {return fActivateOCeffects;}
0084 G4LogicalVolume* GetCrystalVolume() const {return fCrystalLogic;}
0085 void SetPotentialPath(const G4String path){fPotentialPath = path;}
0086
0087
0088 void SetRadiatorConverterSepDistance(G4double val) {
0089 fRadiatorConverterSepDistance = val;}
0090 G4double GetRadiatorConverterSepDistance() const {
0091 return fRadiatorConverterSepDistance;}
0092 void SetConverterSize(G4ThreeVector val) {fConverterSize = val;}
0093 void SetConverterMaterial(G4String val) {fConverterMaterialStr = val;}
0094 void SetGranularConverter(G4bool val) {fGranularConverter = val;}
0095 void SetSphereRadius(G4double val) {fSphereRadius = val;}
0096 G4int GetNSpheres() const {return fNSpheres;}
0097 G4LogicalVolume* GetConverterVolume() const {return fConverterLogic;}
0098
0099
0100 void SetMagneticField(G4bool val) {fSetMagneticField = val;}
0101 void SetFieldValue(G4double val) {fFieldValue = val;}
0102 void SetFieldRegionLength(G4double val) {fFieldRegionLength = val;}
0103
0104
0105 void SetCollimator(G4bool val) {fSetCollimator = val;}
0106 void SetCollimatorHole(G4String val) {fCollimatorHole = val;}
0107 void SetCollimatorAperture(G4double val) {fCollimatorAperture = val;}
0108 void SetCollimatorThickness(G4double val) {fCollimatorThickness = val;}
0109 void SetCollimatorSide(G4double val) {fCollimatorSide = val;}
0110 void SetRadiatorCollimatorSepDistance(G4double val) {
0111 fRadiatorCollimatorSepDistance = val;}
0112 G4double GetRadiatorCollimatorSepDistance() const {
0113 return fRadiatorCollimatorSepDistance;}
0114
0115
0116 void SetVirtualDetectorSize(G4ThreeVector val) {fVirtualDetectorSize = val;}
0117 std::vector<G4ThreeVector> GetVirtualDetectorPositionVector() const {
0118 return fVirtualDetectorPositionVector;}
0119
0120
0121 void SetScoringCrystalExit(G4bool bval) {fScoringCrystalExit = bval;}
0122 G4bool GetScoringCrystalExit() const {return fScoringCrystalExit;}
0123
0124 protected:
0125 std::vector<G4LogicalVolume*> fScoringVolume;
0126
0127 private:
0128 DetectorConstructionMessenger* fMessenger;
0129
0130 G4bool fHybridSource = true;
0131
0132 G4Region* fCrystalRegion{nullptr};
0133 G4LogicalVolume* fCrystalLogic{nullptr};
0134 G4String fCrystalMaterialStr = "W";
0135 G4Material* fCrystalMaterial{nullptr};
0136 G4ThreeVector fCrystalSize = G4ThreeVector(7.*mm, 7.*mm, 2.*mm);
0137 G4double fBendingAngle = 0.e-6;
0138 G4String fLattice = "<111>";
0139 G4double fAngleX = 0.e-6;
0140 G4double fAngleY = 0.e-6;
0141 G4double fCrystalZ = 0.;
0142 G4bool fActivateRadiationModel = true;
0143 G4bool fActivateOCeffects = true;
0144 G4String fPotentialPath = "";
0145
0146 G4double fRadiatorConverterSepDistance = 60.*cm;
0147 G4ThreeVector fConverterSize = G4ThreeVector(199.75*mm, 199.75*mm, 11.6*mm);
0148 G4double fConverterZ = 0.;
0149 G4LogicalVolume* fConverterLogic{nullptr};
0150 G4bool fGranularConverter = false;
0151 G4String fConverterMaterialStr = "W";
0152 G4Material* fConverterMaterial{nullptr};
0153 G4double fSphereRadius = 1.1*mm;
0154 G4LogicalVolume* fSphereLogic[NSpheresMax];
0155 G4int fNSpheres = 0;
0156 G4bool fConverter = true;
0157
0158 G4bool fSetMagneticField = false;
0159 G4double fFieldValue = 100.*tesla;
0160 G4double fFieldRegionLength = 90.*cm;
0161 G4LogicalVolume* fMFlogic{nullptr};
0162
0163 G4bool fSetCollimator = false;
0164 G4double fCollimatorAperture = 2.*mm;
0165 G4String fCollimatorHole = "squared";
0166 G4double fCollimatorThickness = 50.*cm;
0167 G4double fCollimatorSide = 2.5*m;
0168 G4double fRadiatorCollimatorSepDistance = 5.*cm;
0169 G4LogicalVolume* fCollimatorLogic{nullptr};
0170
0171 G4ThreeVector fVirtualDetectorSize = G4ThreeVector(40.*cm, 40.*cm, 0.01*mm);
0172 std::vector<G4ThreeVector> fVirtualDetectorPositionVector;
0173 G4LogicalVolume* fVirtualDetectorLogic0{nullptr};
0174 G4LogicalVolume* fVirtualDetectorLogic1{nullptr};
0175 G4LogicalVolume* fVirtualDetectorLogic2{nullptr};
0176
0177 G4bool fScoringCrystalExit = false;
0178 };
0179
0180
0181
0182 #endif