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 // 27 /// \file GB05DetectorConstruction.cc 28 /// \brief Implementation of the GB05DetectorConstruction class 29 30 #include "GB05DetectorConstruction.hh" 31 32 #include "GB05BOptrSplitAndKillByCrossSection.hh" 33 #include "GB05SD.hh" 34 35 #include "G4Box.hh" 36 #include "G4LogicalVolume.hh" 37 #include "G4LogicalVolumeStore.hh" 38 #include "G4Material.hh" 39 #include "G4NistManager.hh" 40 #include "G4PVPlacement.hh" 41 #include "G4SDManager.hh" 42 #include "G4SystemOfUnits.hh" 43 44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 45 46 GB05DetectorConstruction::GB05DetectorConstruction() : G4VUserDetectorConstruction() {} 47 48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 49 50 GB05DetectorConstruction::~GB05DetectorConstruction() {} 51 52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 53 54 G4VPhysicalVolume* GB05DetectorConstruction::Construct() 55 { 56 // -- Collect a few materials from NIST database: 57 auto nistManager = G4NistManager::Instance(); 58 G4Material* worldMaterial = nistManager->FindOrBuildMaterial("G4_Galactic"); 59 G4Material* defaultMaterial = nistManager->FindOrBuildMaterial("G4_CONCRETE"); 60 61 G4VSolid* solidWorld = new G4Box("World", 10 * m, 10 * m, 10 * m); 62 63 G4LogicalVolume* logicWorld = new G4LogicalVolume(solidWorld, // its solid 64 worldMaterial, // its material 65 "World"); // its name 66 67 G4PVPlacement* physiWorld = new G4PVPlacement(0, // no rotation 68 G4ThreeVector(), // at (0,0,0) 69 logicWorld, // its logical volume 70 "World", // its name 71 0, // its mother volume 72 false, // no boolean operation 73 0); // copy number 74 75 // ----------------------------------- 76 // -- volume where biasing is applied: 77 // ----------------------------------- 78 G4double halfXY = 5.0 * m; 79 G4double halfZ = 1.0 * m; 80 G4VSolid* solidShield = new G4Box("shield.solid", halfXY, halfXY, halfZ); 81 82 G4LogicalVolume* logicShield = new G4LogicalVolume(solidShield, // its solid 83 defaultMaterial, // its material 84 "shield.logical"); // its name 85 86 new G4PVPlacement(0, // no rotation 87 G4ThreeVector(0, 0, halfZ), // volume entrance at (0,0,0) 88 logicShield, // its logical volume 89 "shield.phys", // its name 90 logicWorld, // its mother volume 91 false, // no boolean operation 92 0); // copy number 93 94 // -------------------------------------------------------------------- 95 // -- dummy volume to print information on particle exiting the shield: 96 // -------------------------------------------------------------------- 97 G4double halfz = 1 * cm; 98 G4VSolid* solidMeasurement = new G4Box("meas.solid", halfXY, halfXY, halfz); 99 100 G4LogicalVolume* logicMeasurement = new G4LogicalVolume(solidMeasurement, // its solid 101 worldMaterial, // its material 102 "meas.logical"); // its name 103 104 new G4PVPlacement(0, // no rotation 105 G4ThreeVector(0, 0, 2 * halfZ + halfz), // volume entrance at (0,0,0) 106 logicMeasurement, // its logical volume 107 "meas.phys", // its name 108 logicWorld, // its mother volume 109 false, // no boolean operation 110 0); // copy number 111 112 return physiWorld; 113 } 114 115 void GB05DetectorConstruction::ConstructSDandField() 116 { 117 // -- Fetch volume for biasing: 118 G4LogicalVolume* logicShield = G4LogicalVolumeStore::GetInstance()->GetVolume("shield.logical"); 119 120 // ------------------------------------------------------------- 121 // -- operator creation, configuration and attachment to volume: 122 // ------------------------------------------------------------- 123 GB05BOptrSplitAndKillByCrossSection* biasingOperator = 124 new GB05BOptrSplitAndKillByCrossSection("neutron"); 125 // -- Now, we declare to our biasing operator all the processes we 126 // -- their disapperance effect on neutrons to be counterbalanced 127 // -- by the splitting by cross-section : 128 biasingOperator->AddProcessToEquipoise("Decay"); 129 biasingOperator->AddProcessToEquipoise("nCapture"); 130 biasingOperator->AddProcessToEquipoise("neutronInelastic"); 131 132 biasingOperator->AttachTo(logicShield); 133 134 G4cout << " Attaching biasing operator " << biasingOperator->GetName() << " to logical volume " 135 << biasingOperator->GetName() << G4endl; 136 137 // ------------------------------------------------------------------------------------ 138 // -- Attach a sensitive detector to print information on particles exiting the shield: 139 // ------------------------------------------------------------------------------------ 140 // -- create and register sensitive detector code module: 141 G4SDManager* SDman = G4SDManager::GetSDMpointer(); 142 G4VSensitiveDetector* sd = new GB05SD("Scorer"); 143 SDman->AddNewDetector(sd); 144 // -- Fetch volume for sensitivity and attach sensitive module to it: 145 G4LogicalVolume* logicSD = G4LogicalVolumeStore::GetInstance()->GetVolume("meas.logical"); 146 logicSD->SetSensitiveDetector(sd); 147 } 148