Geant4 Cross Reference |
1 // 2 // ******************************************************************** 3 // * License and Disclaimer * 4 // * * 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. * 10 // * * 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitation of liability. * 17 // * * 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************************************** 25 // 26 /// \file DetectorConstruction.cc 27 /// \brief Implementation of the DetectorConstruction class 28 // 29 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 30 31 #include "DetectorConstruction.hh" 32 33 #include "G4RunManager.hh" 34 #include "G4NistManager.hh" 35 #include "G4Box.hh" 36 #include "G4LogicalVolume.hh" 37 #include "G4RegionStore.hh" 38 #include "G4VisAttributes.hh" 39 #include "PrimaryGeneratorAction.hh" 40 #include <CLHEP/Units/SystemOfUnits.h> 41 42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 43 44 G4VPhysicalVolume* DetectorConstruction::Construct() 45 { 46 //Check overlap option 47 G4bool checkOverlaps = true; 48 49 //Materials 50 G4NistManager* nist = G4NistManager::Instance(); 51 G4Material* world_mat = nist->FindOrBuildMaterial("G4_Galactic"); 52 G4Material* silicon = nist->FindOrBuildMaterial("G4_Si"); 53 54 //World 55 G4Box* solidWorld = new G4Box("World", 0.2*CLHEP::m, 0.2*CLHEP::m, 0.3*CLHEP::m); 56 G4LogicalVolume* logicWorld = new G4LogicalVolume(solidWorld, world_mat, "World"); 57 G4VPhysicalVolume* physWorld = new G4PVPlacement 58 (0, // no rotation 59 G4ThreeVector(), // centre position 60 logicWorld, // its logical volume 61 "World", // its name 62 0, // its mother volume 63 false, // no boolean operation 64 0, // copy number 65 checkOverlaps); // overlaps checking 66 logicWorld->SetVisAttributes(G4VisAttributes::GetInvisible()); 67 68 69 // --------------- Crystal ------------------------------------ 70 71 //Select crystal material 72 fCrystalMaterial = nist->FindOrBuildMaterial("G4_W"); 73 74 //Setting crystal rotation angle (also the angle of crystal planes vs the beam) 75 //Crystal rotation angle (also the angle of crystal planes vs the beam) 76 G4double angleX = 0.*1e-6; //rad 77 G4RotationMatrix* crystalRotationMatrix = new G4RotationMatrix; 78 crystalRotationMatrix->rotateY(-angleX); 79 80 //setting crystal dimensions 81 /*at high energies the electromagnetic shower in 82 oriented tungsten should behave similarly to the 83 e.m. shower in several times thicker amorphous tungsten 84 (several times reduction of the effective radiation length)*/ 85 G4ThreeVector crystalSize = G4ThreeVector(10.*CLHEP::mm, 86 10.*CLHEP::mm, 87 0.1*CLHEP::mm); 88 89 //Setting crystal position 90 G4ThreeVector posCrystal = G4ThreeVector(0., 0., crystalSize.z()/2.); 91 92 //crystal volume 93 G4Box* solidCrystal = new G4Box("Crystal", 94 crystalSize.x()/2, 95 crystalSize.y()/2, 96 crystalSize.z()/2.); 97 98 fLogicCrystal = new G4LogicalVolume(solidCrystal, 99 fCrystalMaterial, 100 "Crystal"); 101 new G4PVPlacement(crystalRotationMatrix, 102 posCrystal, 103 fLogicCrystal, 104 "Crystal", 105 logicWorld, 106 false, 107 0, 108 checkOverlaps); 109 110 //crystal region (necessary for the FastSim model) 111 G4Region* regionCh = new G4Region("Crystal"); 112 regionCh->AddRootLogicalVolume(fLogicCrystal); 113 114 //visualization attributes 115 G4VisAttributes* crystalVisAttribute = 116 new G4VisAttributes(G4Colour(1., 0., 0.)); 117 crystalVisAttribute->SetForceSolid(true); 118 fLogicCrystal->SetVisAttributes(crystalVisAttribute); 119 120 //print crystal info 121 G4cout << "Crystal size: " << crystalSize.x()/CLHEP::mm 122 << " " << crystalSize.y()/CLHEP::mm 123 << " " << crystalSize.z()/CLHEP::mm << " mm3" << G4endl; 124 G4cout << "Crystal angleX: " << angleX << " rad" << G4endl; 125 126 // --------------- Detector ----------------------------------- 127 //Setting detector position 128 G4ThreeVector posDetector = G4ThreeVector(0, 0, 0.1*CLHEP::m); 129 130 //particle detector volume 131 G4Box* detector = new G4Box("Detector", 132 10*CLHEP::cm/2, 133 10*CLHEP::cm/2, 134 0.3*CLHEP::mm/2); 135 136 G4LogicalVolume* logicDetector = new G4LogicalVolume(detector, 137 silicon, 138 "Detector"); 139 new G4PVPlacement(0, 140 posDetector, 141 logicDetector, 142 "Detector", 143 logicWorld, 144 false, 145 0, 146 checkOverlaps); 147 148 //visualization attributes 149 G4VisAttributes* detectorVisAttribute = 150 new G4VisAttributes(G4Colour(0., 0., 1)); 151 detectorVisAttribute->SetForceSolid(true); 152 logicDetector->SetVisAttributes(detectorVisAttribute); 153 154 //always return the physical World 155 return physWorld; 156 } 157 158 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 159 160 void DetectorConstruction::ConstructSDandField() 161 { 162 // --------------- fast simulation ---------------------------- 163 //extract the region of the crystal from the store 164 G4RegionStore* regionStore = G4RegionStore::GetInstance(); 165 G4Region* regionCh = regionStore->GetRegion("Crystal"); 166 167 //create the channeling model for this region 168 G4ChannelingFastSimModel* channelingModel = 169 new G4ChannelingFastSimModel("ChannelingModel", regionCh); 170 171 //Crystal planes or axes considered 172 //Use brackets (...) for planes and <...> for axes 173 G4String lattice = "<111>"; 174 175 //activate the channeling model 176 channelingModel->Input(fCrystalMaterial, lattice); 177 178 /* 179 activate radiation model (do it only when you want to take into account the 180 radiation production in an oriented crystal; it reduces simulation speed.) 181 */ 182 G4bool activateRadiationModel = true; 183 if (activateRadiationModel) 184 { 185 channelingModel->RadiationModelActivate(); 186 G4cout << "Radiation model activated" << G4endl; 187 } 188 } 189 190 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 191 192