File indexing completed on 2025-10-13 08:26:57
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 #include "DetectorConstruction.hh"
0032
0033 #include "G4RunManager.hh"
0034 #include "G4NistManager.hh"
0035 #include "G4Box.hh"
0036 #include "G4LogicalVolume.hh"
0037 #include "G4RegionStore.hh"
0038 #include "G4VisAttributes.hh"
0039 #include "PrimaryGeneratorAction.hh"
0040 #include <CLHEP/Units/SystemOfUnits.h>
0041
0042
0043
0044 DetectorConstruction::DetectorConstruction()
0045 {}
0046
0047
0048
0049 G4VPhysicalVolume* DetectorConstruction::Construct()
0050 {
0051
0052 G4bool checkOverlaps = true;
0053
0054
0055 G4NistManager* nist = G4NistManager::Instance();
0056 G4Material* world_mat = nist->FindOrBuildMaterial("G4_Galactic");
0057 G4Material* silicon = nist->FindOrBuildMaterial("G4_Si");
0058
0059
0060 G4Box* solidWorld = new G4Box("World", 0.2*CLHEP::m, 0.2*CLHEP::m, 10.*CLHEP::m);
0061 G4LogicalVolume* logicWorld = new G4LogicalVolume(solidWorld, world_mat, "World");
0062 G4VPhysicalVolume* physWorld = new G4PVPlacement
0063 (0,
0064 G4ThreeVector(),
0065 logicWorld,
0066 "World",
0067 0,
0068 false,
0069 0,
0070 checkOverlaps);
0071 logicWorld->SetVisAttributes(G4VisAttributes::GetInvisible());
0072
0073
0074
0075
0076
0077
0078
0079
0080
0081
0082
0083
0084 fCrystalMaterial = nist->FindOrBuildMaterial("G4_Si");
0085
0086
0087
0088 G4double angleX = 0.*1e-6;
0089 G4RotationMatrix* crystalRotationMatrix = new G4RotationMatrix;
0090 crystalRotationMatrix->rotateY(-angleX);
0091
0092
0093 fBendingAngle = 0.905*CLHEP::mrad;
0094
0095
0096 G4ThreeVector crystalSize = G4ThreeVector(20.*CLHEP::mm,
0097 20.*CLHEP::mm,
0098 0.0305*CLHEP::mm);
0099
0100
0101 G4ThreeVector posCrystal = G4ThreeVector(0., 0., crystalSize.z()/2.);
0102
0103
0104 G4Box* solidCrystal = new G4Box("Crystal",
0105 crystalSize.x()/2,
0106 crystalSize.y()/2,
0107 crystalSize.z()/2.);
0108
0109 fLogicCrystal = new G4LogicalVolume(solidCrystal,
0110 fCrystalMaterial,
0111 "Crystal");
0112 new G4PVPlacement(crystalRotationMatrix,
0113 posCrystal,
0114 fLogicCrystal,
0115 "Crystal",
0116 logicWorld,
0117 false,
0118 0,
0119 checkOverlaps);
0120
0121
0122 G4Region* regionCh = new G4Region("Crystal");
0123 regionCh->AddRootLogicalVolume(fLogicCrystal);
0124
0125
0126 G4VisAttributes* crystalVisAttribute =
0127 new G4VisAttributes(G4Colour(1., 0., 0.));
0128 crystalVisAttribute->SetForceSolid(true);
0129 fLogicCrystal->SetVisAttributes(crystalVisAttribute);
0130
0131
0132 G4cout << "Crystal size: " << crystalSize.x()/CLHEP::mm
0133 << " " << crystalSize.y()/CLHEP::mm
0134 << " " << crystalSize.z()/CLHEP::mm << " mm3" << G4endl;
0135 G4cout << "Crystal bending angle: " << fBendingAngle << " rad" << G4endl;
0136 G4cout << "Crystal angleX: " << angleX << " rad" << G4endl;
0137
0138
0139
0140 G4ThreeVector posDetector = G4ThreeVector(0, 0, 5973*CLHEP::mm);
0141
0142
0143 G4Box* detector = new G4Box("Detector",
0144 10*CLHEP::cm/2,
0145 10*CLHEP::cm/2,
0146 0.3*CLHEP::mm/2);
0147
0148 G4LogicalVolume* logicDetector = new G4LogicalVolume(detector,
0149 silicon,
0150 "Detector");
0151 new G4PVPlacement(0,
0152 posDetector,
0153 logicDetector,
0154 "Detector",
0155 logicWorld,
0156 false,
0157 0,
0158 checkOverlaps);
0159
0160
0161 G4VisAttributes* detectorVisAttribute =
0162 new G4VisAttributes(G4Colour(0., 0., 1));
0163 detectorVisAttribute->SetForceSolid(true);
0164 logicDetector->SetVisAttributes(detectorVisAttribute);
0165
0166
0167 return physWorld;
0168 }
0169
0170
0171
0172 void DetectorConstruction::ConstructSDandField()
0173 {
0174
0175
0176 G4RegionStore* regionStore = G4RegionStore::GetInstance();
0177 G4Region* regionCh = regionStore->GetRegion("Crystal");
0178
0179
0180 G4ChannelingFastSimModel* channelingModel =
0181 new G4ChannelingFastSimModel("ChannelingModel", regionCh);
0182
0183
0184
0185 G4String lattice = "(111)";
0186
0187
0188 channelingModel->Input(fCrystalMaterial, lattice);
0189
0190 channelingModel->GetCrystalData()->SetBendingAngle(fBendingAngle,fLogicCrystal);
0191
0192
0193
0194
0195
0196 G4bool activateRadiationModel = true;
0197 if (activateRadiationModel)
0198 {
0199 channelingModel->RadiationModelActivate();
0200 G4cout << "Radiation model activated" << G4endl;
0201 }
0202 }
0203
0204
0205