![]() |
|
|||
File indexing completed on 2025-02-23 09:20:46
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 /// \file GB06/src/GB06ParallelWorldForSlices.cc 0027 /// \brief Implementation of the GB06ParallelWorldForSlices class 0028 // 0029 // 0030 #include "GB06ParallelWorldForSlices.hh" 0031 0032 #include "GB06BOptrSplitAndKillByImportance.hh" 0033 0034 #include "G4Box.hh" 0035 #include "G4LogicalVolume.hh" 0036 #include "G4LogicalVolumeStore.hh" 0037 #include "G4PVPlacement.hh" 0038 #include "G4PVReplica.hh" 0039 #include "G4PhysicalVolumeStore.hh" 0040 #include "G4SystemOfUnits.hh" 0041 #include "G4ThreeVector.hh" 0042 0043 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0044 0045 GB06ParallelWorldForSlices::GB06ParallelWorldForSlices(G4String worldName) 0046 : G4VUserParallelWorld(worldName) 0047 { 0048 ; 0049 } 0050 0051 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0052 0053 GB06ParallelWorldForSlices::~GB06ParallelWorldForSlices() 0054 { 0055 ; 0056 } 0057 0058 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0059 0060 void GB06ParallelWorldForSlices::Construct() 0061 { 0062 // -- Inform about construction: 0063 // -- (fWorldName is a protected data member of the base parallel world class) 0064 G4cout << "Parallel World `" << fWorldName << "' constructed." << G4endl; 0065 0066 // ------------------------- 0067 // Build parallel geometry: 0068 // ------------------------- 0069 0070 // -- Obtain clone of mass geometry world from GetWorld() base class utility: 0071 G4VPhysicalVolume* physicalParallelWorld = GetWorld(); 0072 G4LogicalVolume* logicalParallelWorld = physicalParallelWorld->GetLogicalVolume(); 0073 0074 // -- We overlay a sliced geometry on top of the block of concrete in the mass geometry 0075 // -- (ie, in the detector construction class), using the same dimensions. 0076 // -- [Note that this is a choice : we can use different dimensions and shapes, creating 0077 // -- a new solid for that.] 0078 // -- For this we: 0079 // -- - 1) get back the solid used to create the concrete shield; 0080 // -- - 2) create a new logical volume of same shape than the shield and we place 0081 // -- inside the slices 0082 // -- - 3) place the sliced structure, using the placement of the physical volume of 0083 // -- the concrete shield 0084 // -- In all this construction, no materials are used, as only the volumes boundaries 0085 // -- are of interest. Note that the absence of materials is only possible in parallel 0086 // -- geometries. 0087 0088 // -- 1) get back the solid used to create the concrete shield: 0089 // ------------------------------------------------------ 0090 0091 // -- get back the logical volume of the shield, using its name: 0092 G4LogicalVolume* shieldLogical = G4LogicalVolumeStore::GetInstance()->GetVolume("shield.logical"); 0093 0094 // -- get back the solid, a G4box in this case. We cast the pointer to access later on 0095 // -- the G4Box class specific methods: 0096 G4Box* shieldSolid = (G4Box*)shieldLogical->GetSolid(); 0097 0098 // -- we now re-create a logical volume for the mother volume of the slices: 0099 G4LogicalVolume* motherForSlicesLogical = 0100 new G4LogicalVolume(shieldSolid, // its solid 0101 nullptr, // no material 0102 "motherForSlices.logical"); // its name 0103 0104 // -- 2) new logical volume of same shape than the shield and place inside the slices: 0105 // ----------------------------------------------------------------------------- 0106 0107 // -- We create now the slices; we choose 20 slices: 0108 const G4int nSlices(20); 0109 // -- the solid for slices: 0110 G4double halfSliceZ = shieldSolid->GetZHalfLength() / nSlices; 0111 G4Box* sliceSolid = new G4Box("slice.solid", shieldSolid->GetXHalfLength(), 0112 shieldSolid->GetYHalfLength(), halfSliceZ); 0113 0114 // -- the logical volume for slices: 0115 G4LogicalVolume* sliceLogical = new G4LogicalVolume(sliceSolid, // its solid 0116 nullptr, // no material 0117 "slice.logical"); // its name 0118 0119 // -- we use a replica, to place the 20 slices in one go, along the Z axis: 0120 new G4PVReplica("slice.physical", // its name 0121 sliceLogical, // its logical volume 0122 motherForSlicesLogical, // its mother volume 0123 kZAxis, // axis of replication 0124 nSlices, // number of replica 0125 2 * halfSliceZ); // width of replica 0126 0127 // -- 3) place the sliced structure, using the concrete shield placement: 0128 // ---------------------------------------------------------------- 0129 0130 // -- get back the physical volume of the shield, using its name: 0131 // -- (note that we know we have only one physical volume with this name. If we had 0132 // -- several, we should loop by ourselves on the store which is of 0133 // -- std::vector<G4VPhysicalVolume*> type.) 0134 G4VPhysicalVolume* shieldPhysical = 0135 G4PhysicalVolumeStore::GetInstance()->GetVolume("shield.physical"); 0136 0137 // -- get back the translation 0138 // -- (we don't try to get back the rotation, we know we used nullptr): 0139 G4ThreeVector translation = shieldPhysical->GetObjectTranslation(); 0140 0141 // -- finally, we place the sliced structure: 0142 new G4PVPlacement(nullptr, // no rotation 0143 translation, // translate as for the shield 0144 motherForSlicesLogical, // its logical volume 0145 "motherForSlices.physical", // its name 0146 logicalParallelWorld, // its mother volume 0147 false, // no boolean operation 0148 0); // copy number 0149 } 0150 0151 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0152 0153 void GB06ParallelWorldForSlices::ConstructSD() 0154 { 0155 // -- Create the biasing operator: 0156 auto biasingOperator = new GB06BOptrSplitAndKillByImportance("neutron", "parallelOptr"); 0157 // -- Tell it it is active for this parallel geometry, passing the world 0158 // -- volume of this geometry : 0159 biasingOperator->SetParallelWorld(GetWorld()); 0160 0161 // -- Attach to the logical volume where the biasing has to be applied: 0162 auto slice = G4LogicalVolumeStore::GetInstance()->GetVolume("slice.logical"); 0163 biasingOperator->AttachTo(slice); 0164 0165 // -- Create a simple "volume importance" map, linking replica numbers to importances: 0166 // -------------------------------------------------------------------------------- 0167 // -- we define the map as going from an importance to 2*importance when going from 0168 // -- a slice to the next one, in the Z direction. 0169 // -- Get back the replica of slices: 0170 G4PVReplica* slicePhysical = 0171 (G4PVReplica*)(G4PhysicalVolumeStore::GetInstance()->GetVolume("slice.physical")); 0172 G4int nReplica = slicePhysical->GetMultiplicity(); 0173 // -- We use and fill the map we defined in the biasing operator: 0174 G4int importance = 1; 0175 for (G4int iReplica = 0; iReplica < nReplica; iReplica++) { 0176 (biasingOperator->GetImportanceMap())[iReplica] = importance; 0177 importance *= 2; 0178 } 0179 }
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |
![]() ![]() |