Loading...
Searching...
No Matches
Public Member Functions | Private Member Functions | List of all members
GB06ParallelWorldForSlices Class Reference

#include <Doxymodules_biasing.h>

Inheritance diagram for GB06ParallelWorldForSlices:
G4VUserParallelWorld

Public Member Functions

 GB06ParallelWorldForSlices (G4String worldName)
 
 ~GB06ParallelWorldForSlices ()
 

Private Member Functions

virtual void Construct ()
 
virtual void ConstructSD ()
 

Detailed Description

Definition at line 141 of file Doxymodules_biasing.h.

Constructor & Destructor Documentation

◆ GB06ParallelWorldForSlices()

GB06ParallelWorldForSlices::GB06ParallelWorldForSlices ( G4String  worldName)

Definition at line 45 of file GB06ParallelWorldForSlices.cc.

◆ ~GB06ParallelWorldForSlices()

GB06ParallelWorldForSlices::~GB06ParallelWorldForSlices ( )

Definition at line 51 of file GB06ParallelWorldForSlices.cc.

52{;}

Member Function Documentation

◆ Construct()

void GB06ParallelWorldForSlices::Construct ( )
privatevirtual

Definition at line 56 of file GB06ParallelWorldForSlices.cc.

57{
58 // -- Inform about construction:
59 // -- (fWorldName is a protected data member of the base parallel world class)
60 G4cout << "Parallel World `" << fWorldName << "' constructed." << G4endl;
61
62 // -------------------------
63 // Build parallel geometry:
64 // -------------------------
65
66 // -- Obtain clone of mass geometry world from GetWorld() base class utility:
67 G4VPhysicalVolume* physicalParallelWorld = GetWorld();
68 G4LogicalVolume* logicalParallelWorld = physicalParallelWorld->GetLogicalVolume();
69
70
71 // -- We overlay a sliced geometry on top of the block of concrete in the mass geometry
72 // -- (ie, in the detector construction class), using the same dimensions.
73 // -- [Note that this is a choice : we can use different dimensions and shapes, creating
74 // -- a new solid for that.]
75 // -- For this we:
76 // -- - 1) get back the solid used to create the concrete shield;
77 // -- - 2) create a new logical volume of same shape than the shield and we place
78 // -- inside the slices
79 // -- - 3) place the sliced structure, using the placement of the physical volume of
80 // -- the concrete shield
81 // -- In all this construction, no materials are used, as only the volumes boundaries
82 // -- are of interest. Note that the absence of materials is only possible in parallel
83 // -- geometries.
84
85
86 // -- 1) get back the solid used to create the concrete shield:
87 // ------------------------------------------------------
88
89 // -- get back the logical volume of the shield, using its name:
90 G4LogicalVolume* shieldLogical =
91 G4LogicalVolumeStore::GetInstance()->GetVolume("shield.logical");
92
93 // -- get back the solid, a G4box in this case. We cast the pointer to access later on
94 // -- the G4Box class specific methods:
95 G4Box* shieldSolid = (G4Box*) shieldLogical->GetSolid();
96
97 // -- we now re-create a logical volume for the mother volume of the slices:
98 G4LogicalVolume* motherForSlicesLogical =
99 new G4LogicalVolume(shieldSolid, // its solid
100 nullptr, // no material
101 "motherForSlices.logical"); // its name
102
103
104
105 // -- 2) new logical volume of same shape than the shield and place inside the slices:
106 // -----------------------------------------------------------------------------
107
108 // -- We create now the slices; we choose 20 slices:
109 const G4int nSlices(20);
110 // -- the solid for slices:
111 G4double halfSliceZ = shieldSolid->GetZHalfLength() / nSlices;
112 G4Box* sliceSolid = new G4Box("slice.solid",
113 shieldSolid->GetXHalfLength(),
114 shieldSolid->GetYHalfLength(),
115 halfSliceZ );
116
117 // -- the logical volume for slices:
118 G4LogicalVolume* sliceLogical = new G4LogicalVolume(sliceSolid, // its solid
119 nullptr, // no material
120 "slice.logical"); // its name
121
122 // -- we use a replica, to place the 20 slices in one go, along the Z axis:
123 new G4PVReplica( "slice.physical", // its name
124 sliceLogical, // its logical volume
125 motherForSlicesLogical, // its mother volume
126 kZAxis, // axis of replication
127 nSlices, // number of replica
128 2*halfSliceZ); // width of replica
129
130
131 // -- 3) place the sliced structure, using the concrete shield placement:
132 // ----------------------------------------------------------------
133
134 // -- get back the physical volume of the shield, using its name:
135 // -- (note that we know we have only one physical volume with this name. If we had
136 // -- several, we should loop by ourselves on the store which is of
137 // -- std::vector<G4VPhysicalVolume*> type.)
139 shieldPhysical = G4PhysicalVolumeStore::GetInstance()->GetVolume("shield.physical");
140
141 // -- get back the translation
142 // -- (we don't try to get back the rotation, we know we used nullptr):
143 G4ThreeVector translation = shieldPhysical->GetObjectTranslation();
144
145 // -- finally, we place the sliced structure:
146 new G4PVPlacement( nullptr, // no rotation
147 translation, // translate as for the shield
148 motherForSlicesLogical, // its logical volume
149 "motherForSlices.physical", // its name
150 logicalParallelWorld, // its mother volume
151 false, // no boolean operation
152 0); // copy number
153
154}

◆ ConstructSD()

void GB06ParallelWorldForSlices::ConstructSD ( )
privatevirtual

Definition at line 158 of file GB06ParallelWorldForSlices.cc.

159{
160 // -- Create the biasing operator:
161 auto biasingOperator = new GB06BOptrSplitAndKillByImportance("neutron","parallelOptr");
162 // -- Tell it it is active for this parallel geometry, passing the world
163 // -- volume of this geometry :
164 biasingOperator->SetParallelWorld( GetWorld() );
165
166 // -- Attach to the logical volume where the biasing has to be applied:
167 auto slice = G4LogicalVolumeStore::GetInstance()->GetVolume("slice.logical");
168 biasingOperator->AttachTo(slice);
169
170 // -- Create a simple "volume importance" map, linking replica numbers to importances:
171 // --------------------------------------------------------------------------------
172 // -- we define the map as going from an importance to 2*importance when going from
173 // -- a slice to the next one, in the Z direction.
174 // -- Get back the replica of slices:
175 G4PVReplica* slicePhysical =
176 (G4PVReplica*)(G4PhysicalVolumeStore::GetInstance()->GetVolume("slice.physical"));
177 G4int nReplica = slicePhysical->GetMultiplicity();
178 // -- We use and fill the map we defined in the biasing operator:
179 G4int importance = 1;
180 for ( G4int iReplica = 0 ; iReplica < nReplica ; iReplica++ )
181 {
182 (biasingOperator->GetImportanceMap())[ iReplica ] = importance;
183 importance *= 2;
184 }
185
186}

The documentation for this class was generated from the following files:

Applications | User Support | Publications | Collaboration