1 | <!-- ******************************************************** --> |
---|
2 | <!-- --> |
---|
3 | <!-- [History] --> |
---|
4 | <!-- Changed by: Gabriele Cosmo, 18-Apr-2005 --> |
---|
5 | <!-- Converted to DocBook: Katsuya Amako, Aug-2006 --> |
---|
6 | <!-- --> |
---|
7 | <!-- ******************************************************** --> |
---|
8 | |
---|
9 | |
---|
10 | <!-- ******************* Section (Level#2) ****************** --> |
---|
11 | <sect2 id="sect.Geom.Navig"> |
---|
12 | <title> |
---|
13 | The Geometry Navigator |
---|
14 | </title> |
---|
15 | |
---|
16 | <para> |
---|
17 | Navigation through the geometry at tracking time is implemented |
---|
18 | by the class <literal>G4Navigator</literal>. The navigator is used to locate |
---|
19 | points in the geometry and compute distances to geometry |
---|
20 | boundaries. At tracking time, the navigator is intended to be the |
---|
21 | only point of interaction with tracking. |
---|
22 | </para> |
---|
23 | |
---|
24 | <para> |
---|
25 | Internally, the G4Navigator has several private helper/utility |
---|
26 | classes: |
---|
27 | |
---|
28 | <itemizedlist spacing="compact"> |
---|
29 | <listitem><para> |
---|
30 | <emphasis role="bold">G4NavigationHistory</emphasis> - stores the compounded |
---|
31 | transformations, replication/parameterisation information, and |
---|
32 | volume pointers at each level of the hierarchy to the current |
---|
33 | location. The volume types at each level are also stored - whether |
---|
34 | normal (placement), replicated or parameterised. |
---|
35 | </para></listitem> |
---|
36 | <listitem><para> |
---|
37 | <emphasis role="bold">G4NormalNavigation</emphasis> - provides location & |
---|
38 | distance computation functions for geometries containing 'placement' |
---|
39 | volumes, with no voxels. |
---|
40 | </para></listitem> |
---|
41 | <listitem><para> |
---|
42 | <emphasis role="bold">G4VoxelNavigation</emphasis> - provides location and distance |
---|
43 | computation functions for geometries containing 'placement' |
---|
44 | physical volumes with voxels. Internally a stack of voxel |
---|
45 | information is maintained. Private functions allow for isotropic |
---|
46 | distance computation to voxel boundaries and for computation of the |
---|
47 | 'next voxel' in a specified direction. |
---|
48 | </para></listitem> |
---|
49 | <listitem><para> |
---|
50 | <emphasis role="bold">G4ParameterisedNavigation</emphasis> - provides location and |
---|
51 | distance computation functions for geometries containing |
---|
52 | parameterised volumes with voxels. Voxel information is maintained |
---|
53 | similarly to <literal>G4VoxelNavigation</literal>, but computation can also |
---|
54 | be simpler by adopting voxels to be one level deep only |
---|
55 | (<emphasis>unrefined</emphasis>, or 1D optimisation) |
---|
56 | </para></listitem> |
---|
57 | <listitem><para> |
---|
58 | <emphasis role="bold">G4ReplicaNavigation</emphasis> - provides location and distance |
---|
59 | computation functions for replicated volumes. |
---|
60 | </para></listitem> |
---|
61 | </itemizedlist> |
---|
62 | </para> |
---|
63 | |
---|
64 | <para> |
---|
65 | In addition, the navigator maintains a set of flags for |
---|
66 | exiting/entry optimisation. A navigator is not a singleton class; |
---|
67 | this is mainly to allow a design extension in future (e.g |
---|
68 | geometrical event biasing). |
---|
69 | </para> |
---|
70 | |
---|
71 | |
---|
72 | <!-- ******************* Section (Level#3) ****************** --> |
---|
73 | <sect3 id="sect.Geom.Navig.NavigTrack"> |
---|
74 | <title> |
---|
75 | Navigation and Tracking |
---|
76 | </title> |
---|
77 | |
---|
78 | <para> |
---|
79 | The main functions required for tracking in the geometry are |
---|
80 | described below. Additional functions are provided to return the |
---|
81 | net transformation of volumes and for the creation of touchables. |
---|
82 | None of the functions implicitly requires that the geometry be |
---|
83 | described hierarchically. |
---|
84 | |
---|
85 | <itemizedlist spacing="compact"> |
---|
86 | <listitem><para> |
---|
87 | <emphasis role="bold">SetWorldVolume()</emphasis> |
---|
88 | <para> |
---|
89 | Sets the first volume in the hierarchy. It must be unrotated and |
---|
90 | untranslated from the origin. |
---|
91 | </para> |
---|
92 | </para></listitem> |
---|
93 | <listitem><para> |
---|
94 | <emphasis role="bold">LocateGlobalPointAndSetup()</emphasis> |
---|
95 | <para> |
---|
96 | Locates the volume containing the specified global point. This |
---|
97 | involves a traverse of the hierarchy, requiring the computation of |
---|
98 | compound transformations, testing replicated and parameterised |
---|
99 | volumes (etc). To improve efficiency this search may be performed |
---|
100 | relative to the last, and this is the recommended way of calling |
---|
101 | the function. A 'relative' search may be used for the first call of |
---|
102 | the function which will result in the search defaulting to a search |
---|
103 | from the root node of the hierarchy. Searches may also be performed |
---|
104 | using a <literal>G4TouchableHistory</literal>. |
---|
105 | </para> |
---|
106 | </para></listitem> |
---|
107 | <listitem><para> |
---|
108 | <emphasis role="bold">LocateGlobalPointAndUpdateTouchableHandle()</emphasis> |
---|
109 | <para> |
---|
110 | First, search the geometrical hierarchy like the above method |
---|
111 | <literal>LocateGlobalPointAndSetup()</literal>. Then use the volume found and |
---|
112 | its navigation history to update the touchable. |
---|
113 | </para> |
---|
114 | </para></listitem> |
---|
115 | <listitem><para> |
---|
116 | <emphasis role="bold">ComputeStep()</emphasis> |
---|
117 | <para> |
---|
118 | Computes the distance to the next boundary intersected along the |
---|
119 | specified unit direction from a specified point. The point must be |
---|
120 | have been located prior to calling <literal>ComputeStep()</literal>. |
---|
121 | </para> |
---|
122 | <para> |
---|
123 | When calling <literal>ComputeStep()</literal>, a proposed physics step is |
---|
124 | passed. If it can be determined that the first intersection lies at |
---|
125 | or beyond that distance then <literal>kInfinity</literal> is returned. In any |
---|
126 | case, if the returned step is greater than the physics step, the |
---|
127 | physics step must be taken. |
---|
128 | </para> |
---|
129 | </para></listitem> |
---|
130 | <listitem><para> |
---|
131 | <emphasis role="bold">SetGeometricallyLimitedStep()</emphasis> |
---|
132 | <para> |
---|
133 | Informs the navigator that the last computed step was taken in its |
---|
134 | entirety. This enables entering/exiting optimisation, and should be |
---|
135 | called prior to calling <literal>LocateGlobalPointAndSetup()</literal>. |
---|
136 | </para> |
---|
137 | </para></listitem> |
---|
138 | <listitem><para> |
---|
139 | <emphasis role="bold">CreateTouchableHistory()</emphasis> |
---|
140 | <para> |
---|
141 | Creates a <literal>G4TouchableHistory</literal> object, for which the caller |
---|
142 | has deletion responsibility. The 'touchable' volume is the volume |
---|
143 | returned by the last Locate operation. The object includes a copy |
---|
144 | of the current NavigationHistory, enabling the efficient relocation |
---|
145 | of points in/close to the current volume in the hierarchy. |
---|
146 | </para> |
---|
147 | </para></listitem> |
---|
148 | </itemizedlist> |
---|
149 | </para> |
---|
150 | |
---|
151 | <para> |
---|
152 | As stated previously, the navigator makes use of utility classes to |
---|
153 | perform location and step computation functions. The different |
---|
154 | navigation utilities manipulate the <literal>G4NavigationHistory</literal> |
---|
155 | object. |
---|
156 | </para> |
---|
157 | |
---|
158 | <para> |
---|
159 | In <literal>LocateGlobalPointAndSetup()</literal> the process of locating a |
---|
160 | point breaks down into three main stages - optimisation, |
---|
161 | determination that the point is contained with a subtree (mother |
---|
162 | and daughters), and determination of the actual containing |
---|
163 | daughter. The latter two can be thought of as scanning first 'up' |
---|
164 | the hierarchy until a volume that is guaranteed to contain the |
---|
165 | point is found, and then scanning 'down' until the actual volume |
---|
166 | that contains the point is found. |
---|
167 | </para> |
---|
168 | |
---|
169 | <para> |
---|
170 | In <literal>ComputeStep()</literal> three types of computation are treated |
---|
171 | depending on the current containing volume: |
---|
172 | |
---|
173 | <itemizedlist spacing="compact"> |
---|
174 | <listitem><para> |
---|
175 | The volume contains normal (placement) daughters (or none) |
---|
176 | </para></listitem> |
---|
177 | <listitem><para> |
---|
178 | The volume contains a single parameterised volume object, |
---|
179 | representing many volumes |
---|
180 | </para></listitem> |
---|
181 | <listitem><para> |
---|
182 | The volume is a replica and contains normal (placement) daughters |
---|
183 | </para></listitem> |
---|
184 | </itemizedlist> |
---|
185 | </para> |
---|
186 | |
---|
187 | </sect3> |
---|
188 | |
---|
189 | <!-- ******************* Section (Level#3) ****************** --> |
---|
190 | <sect3 id="sect.Geom.Navig.LocPnt"> |
---|
191 | <title> |
---|
192 | Using the navigator to locate points |
---|
193 | </title> |
---|
194 | |
---|
195 | <para> |
---|
196 | More than one navigator objects can be created inside an application; these |
---|
197 | navigators can act independently for different purposes. The main navigator |
---|
198 | which is "<emphasis>activated</emphasis> automatically at the startup of a |
---|
199 | simulation program is the navigator used for the |
---|
200 | <emphasis>tracking</emphasis> and attached the world volume |
---|
201 | of the main tracking (or <emphasis>mass</emphasis>) geometry. |
---|
202 | </para> |
---|
203 | |
---|
204 | <para> |
---|
205 | The navigator for tracking can be retrieved at any state of the application |
---|
206 | by messagging the <literal>G4TransportationManager</literal>: |
---|
207 | |
---|
208 | <informalexample> |
---|
209 | <programlisting> |
---|
210 | G4Navigator* tracking_navigator = |
---|
211 | G4TransportationManager::GetInstance()->GetNavigatorForTracking(); |
---|
212 | </programlisting> |
---|
213 | </informalexample> |
---|
214 | |
---|
215 | The navigator for tracking also retains all the information of the current |
---|
216 | history of volumes transversed at a precise moment of the tracking during a |
---|
217 | run. Therefore, if the navigator for tracking is used during tracking for |
---|
218 | locating a generic point in the tree of volumes, the actual particle gets |
---|
219 | also -relocated- in the specified position and tracking will be of course |
---|
220 | affected ! |
---|
221 | </para> |
---|
222 | |
---|
223 | <para> |
---|
224 | In order to avoid the problem above and provide information about location |
---|
225 | of a point without affecting the tracking, it is suggested to either use an |
---|
226 | alternative <literal>G4Navigator</literal> object (which can then be assigned |
---|
227 | to the world-volume), or access the information through the step. |
---|
228 | </para> |
---|
229 | |
---|
230 | <!-- ******* Bridgehead ******* --> |
---|
231 | <bridgehead renderas='sect4'> |
---|
232 | Using the 'step' to retrieve geometrical information |
---|
233 | </bridgehead> |
---|
234 | |
---|
235 | <para> |
---|
236 | During the tracking run, geometrical information can be retrieved through |
---|
237 | the touchable handle associated to the current step. For example, to identify |
---|
238 | the exact copy-number of a specific physical volume in the mass geometry, |
---|
239 | one should do the following: |
---|
240 | |
---|
241 | <informalexample> |
---|
242 | <programlisting> |
---|
243 | // Given the pointer to the step object ... |
---|
244 | // |
---|
245 | G4Step* aStep = ..; |
---|
246 | |
---|
247 | // ... retrieve the 'pre-step' point |
---|
248 | // |
---|
249 | G4StepPoint* preStepPoint = aStep->GetPreStepPoint(); |
---|
250 | |
---|
251 | // ... retrieve a touchable handle and access to the information |
---|
252 | // |
---|
253 | G4TouchableHandle theTouchable = preStepPoint->GetTouchableHandle(); |
---|
254 | G4int copyNo = theTouchable->GetCopyNumber(); |
---|
255 | G4int motherCopyNo = theTouchable->GetCopyNumber(1); |
---|
256 | </programlisting> |
---|
257 | </informalexample> |
---|
258 | </para> |
---|
259 | |
---|
260 | <para> |
---|
261 | To determine the exact position in global coordinates in the mass geometry |
---|
262 | and convert to local coordinates (local to the current volume): |
---|
263 | |
---|
264 | <informalexample> |
---|
265 | <programlisting> |
---|
266 | G4ThreeVector worldPosition = preStepPoint->GetPosition(); |
---|
267 | G4ThreeVector localPosition = theTouchable->GetHistory()-> |
---|
268 | GetTopTransform().TransformPoint(worldPosition); |
---|
269 | </programlisting> |
---|
270 | </informalexample> |
---|
271 | </para> |
---|
272 | |
---|
273 | <!-- ******* Bridgehead ******* --> |
---|
274 | <bridgehead renderas='sect4'> |
---|
275 | Using an alternative navigator to locate points |
---|
276 | </bridgehead> |
---|
277 | |
---|
278 | <para> |
---|
279 | In order to know (when in the <literal>idle</literal> state of the |
---|
280 | application) in which physical volume a given point is located |
---|
281 | in the detector geometry, it is necessary to create an alternative |
---|
282 | navigator object first and assign it to the world volume: |
---|
283 | |
---|
284 | <informalexample> |
---|
285 | <programlisting> |
---|
286 | G4Navigator* aNavigator = new G4Navigator(); |
---|
287 | aNavigator->SetWorldVolume(worldVolumePointer); |
---|
288 | </programlisting> |
---|
289 | </informalexample> |
---|
290 | </para> |
---|
291 | |
---|
292 | <para> |
---|
293 | Then, locate the point <literal>myPoint</literal> (defined in global coordinates), |
---|
294 | retrieve a <emphasis>touchable handle</emphasis> and do whatever you need with it: |
---|
295 | |
---|
296 | <informalexample> |
---|
297 | <programlisting> |
---|
298 | aNavigator->LocateGlobalPointAndSetup(myPoint); |
---|
299 | G4TouchableHistoryHandle aTouchable = |
---|
300 | aNavigator->CreateTouchableHistoryHandle(); |
---|
301 | |
---|
302 | // Do whatever you need with it ... |
---|
303 | // ... convert point in local coordinates (local to the current volume) |
---|
304 | // |
---|
305 | G4ThreeVector localPosition = aTouchable->GetHistory()-> |
---|
306 | GetTopTransform().TransformPoint(myPoint); |
---|
307 | |
---|
308 | // ... convert back to global coordinates system |
---|
309 | G4ThreeVector globalPosition = aTouchable->GetHistory()-> |
---|
310 | GetTopTransform().Inverse().TransformPoint(localPosition); |
---|
311 | </programlisting> |
---|
312 | </informalexample> |
---|
313 | </para> |
---|
314 | |
---|
315 | <para> |
---|
316 | If outside of the tracking run and given a generic local position (local to a |
---|
317 | given volume in the geometry tree), it is -not- possible to determine a priori |
---|
318 | its global position and convert it to the global coordinates system. |
---|
319 | The reason for this is rather simple, nobody can guarantee that the given |
---|
320 | (local) point is located in the right -copy- of the physical volume ! |
---|
321 | In order to retrieve this information, some extra knowledge related to the |
---|
322 | absolute position of the physical volume is required first, i.e. one should |
---|
323 | first determine a global point belonging to that volume, eventually making |
---|
324 | a dedicated scan of the geometry tree through a dedicated |
---|
325 | <literal>G4Navigator</literal> object and then apply the method above after |
---|
326 | having created the touchable for it. |
---|
327 | </para> |
---|
328 | |
---|
329 | </sect3> |
---|
330 | |
---|
331 | <!-- ******************* Section (Level#3) ****************** --> |
---|
332 | <sect3 id="sect.Geom.Navig.ParaGeom"> |
---|
333 | <title> |
---|
334 | Navigation in parallel geometries |
---|
335 | </title> |
---|
336 | |
---|
337 | <para> |
---|
338 | Since release 8.2 of Geant4, it is possible to define geometry trees which |
---|
339 | are <literal>parallel</literal> to the tracking geometry and having them |
---|
340 | assigned to navigator objects that transparently communicate in sync with |
---|
341 | the normal tracking geometry. |
---|
342 | </para> |
---|
343 | |
---|
344 | <para> |
---|
345 | Parallel geometries can be defined for several uses (fast shower |
---|
346 | parameterisation, geometrical biasing, particle scoring, readout |
---|
347 | geometries, etc ...) and can <emphasis>overlap</emphasis> with the mass |
---|
348 | geometry defined for the tracking. The <literal>parallel</literal> |
---|
349 | transportation will be activated only after the registration of the |
---|
350 | parallel geometry in the detector description |
---|
351 | setup; see Section <xref linkend="sect.ParaGeom" /> for how to define a parallel |
---|
352 | geometry and register it to the run-manager. |
---|
353 | </para> |
---|
354 | |
---|
355 | <para> |
---|
356 | The <literal>G4TransportationManager</literal> provides all the utilities to verify, |
---|
357 | retrieve and activate the navigators associated to the various parallel |
---|
358 | geometries defined. |
---|
359 | </para> |
---|
360 | |
---|
361 | </sect3> |
---|
362 | |
---|
363 | <!-- ******************* Section (Level#3) ****************** --> |
---|
364 | <sect3 id="sect.Geom.Navig.PhantomGeom"> |
---|
365 | <title> |
---|
366 | Fast navigation in regular patterned geometries and phantoms |
---|
367 | </title> |
---|
368 | |
---|
369 | <para> |
---|
370 | Since release 9.1 of Geant4, a specialised navigation algorithm has |
---|
371 | been introduced to allow for optimal memory use and extremely efficient |
---|
372 | navigation in geometries represented by a regular pattern of volumes |
---|
373 | and particularly three-dimensional grids of boxes. A typical application |
---|
374 | of this kind is the case of DICOM phantoms for medical physics studies. |
---|
375 | </para> |
---|
376 | |
---|
377 | <para> |
---|
378 | The class <literal>G4RegularNavigation</literal> is used and automatically |
---|
379 | activated when such geometries are defined. It is required to the user to |
---|
380 | implement a parameterisation of the kind <literal>G4PhantomParameterisation</literal> |
---|
381 | and place the parameterised volume containing it in a container volume, so that |
---|
382 | all cells in the three-dimensional grid (<emphasis>voxels</emphasis>) completely |
---|
383 | fill the container volume. |
---|
384 | This way the location of a point inside a voxel can be done in a fast way, |
---|
385 | transforming the position to the coordinate system of the container volume |
---|
386 | and doing a simple calculation of the kind: |
---|
387 | |
---|
388 | <informalexample> |
---|
389 | <programlisting> |
---|
390 | copyNo_x = (localPoint.x()+fVoxelHalfX*fNoVoxelX)/(fVoxelHalfX*2.) |
---|
391 | </programlisting> |
---|
392 | </informalexample> |
---|
393 | |
---|
394 | where <literal>fVoxelHalfX</literal> is the half dimension of the voxel along |
---|
395 | <literal>X</literal> and <literal>fNoVoxelX</literal> is the number of voxels |
---|
396 | in the <literal>X</literal> dimension. |
---|
397 | Voxel <literal>0</literal> will be the one closest to the corner |
---|
398 | <literal>(fVoxelHalfX*fNoVoxelX, fVoxelHalfY*fNoVoxelY, fVoxelHalfZ*fNoVoxelZ)</literal>. |
---|
399 | </para> |
---|
400 | |
---|
401 | <para> |
---|
402 | Having the voxels filling completely the container volume allows to avoid |
---|
403 | the lengthy computation of <literal>ComputeStep()</literal> and |
---|
404 | <literal>ComputeSafety</literal> methods required in the traditional |
---|
405 | navigation algorithm. |
---|
406 | In this case, when a track is inside the parent volume, it has always to |
---|
407 | be inside one of the voxels and it will be only necessary to calculate |
---|
408 | the distance to the walls of the current voxel. |
---|
409 | </para> |
---|
410 | |
---|
411 | <!-- ******* Bridgehead ******* --> |
---|
412 | <bridgehead renderas='sect4'> |
---|
413 | Skipping borders of voxels with same material |
---|
414 | </bridgehead> |
---|
415 | |
---|
416 | <para> |
---|
417 | Another speed optimisation can be provided by skipping the frontiers |
---|
418 | of two voxels which the same material assigned, so that bigger steps |
---|
419 | can be done. This optimisation may be not very useful when the number of |
---|
420 | materials is very big (in which case the probability of having contiguous |
---|
421 | voxels with same material is reduced), or when the physical step is small |
---|
422 | compared to the voxel dimensions (very often the case of electrons). |
---|
423 | The optimisation can be switched off in such cases, by invoking the |
---|
424 | following method with argument <literal>skip = 0</literal>: |
---|
425 | |
---|
426 | <informalexample> |
---|
427 | <programlisting> |
---|
428 | G4RegularParameterisation::SetSkipEqualMaterials( G4bool skip ); |
---|
429 | </programlisting> |
---|
430 | </informalexample> |
---|
431 | </para> |
---|
432 | |
---|
433 | <!-- ******* Bridgehead ******* --> |
---|
434 | <bridgehead renderas='sect4'> |
---|
435 | Example |
---|
436 | </bridgehead> |
---|
437 | |
---|
438 | <para> |
---|
439 | To use the specialised navigation, it is required to first create an object |
---|
440 | of type <literal>G4PhantomParameterisation</literal>: |
---|
441 | |
---|
442 | <informalexample> |
---|
443 | <programlisting> |
---|
444 | G4PhantomParameterisation* param = new G4PhantomParameterisation(); |
---|
445 | </programlisting> |
---|
446 | </informalexample> |
---|
447 | |
---|
448 | Then, fill it with the all the necessary data: |
---|
449 | |
---|
450 | <informalexample> |
---|
451 | <programlisting> |
---|
452 | // Voxel dimensions in the three dimensions |
---|
453 | // |
---|
454 | G4double halfX = ...; |
---|
455 | G4double halfY = ...; |
---|
456 | G4double halfZ = ...; |
---|
457 | param->SetVoxelDimensions( halfX, halfY, halfZ ); |
---|
458 | |
---|
459 | // Number of voxels in the three dimensions |
---|
460 | // |
---|
461 | G4int nVoxelX = ...; |
---|
462 | G4int nVoxelY = ...; |
---|
463 | G4int nVoxelZ = ...; |
---|
464 | param->SetNoVoxel( nVoxelX, nVoxelY, nVoxelZ ); |
---|
465 | |
---|
466 | // Vector of materials of the voxels |
---|
467 | // |
---|
468 | std::vector < G4Material* > theMaterials; |
---|
469 | theMaterials.push_back( new G4Material( ... |
---|
470 | theMaterials.push_back( new G4Material( ... |
---|
471 | param->SetMaterials( theMaterials ); |
---|
472 | |
---|
473 | // List of material indices |
---|
474 | // For each voxel it is a number that correspond to the index of its |
---|
475 | // material in the vector of materials defined above; |
---|
476 | // |
---|
477 | size_t* mateIDs = new size_t[nVoxelX*nVoxelY*nVoxelZ]; |
---|
478 | mateIDs[0] = n0; |
---|
479 | mateIDs[1] = n1; |
---|
480 | ... |
---|
481 | param->SetMaterialIndices( mateIDs ); |
---|
482 | </programlisting> |
---|
483 | </informalexample> |
---|
484 | |
---|
485 | Then, define the volume that contains all the voxels: |
---|
486 | |
---|
487 | <informalexample> |
---|
488 | <programlisting> |
---|
489 | G4Box* cont_solid = new G4Box("PhantomContainer",nVoxelX*halfX.,nVoxelY*halfY.,nVoxelZ*halfZ); |
---|
490 | G4LogicalVolume* cont_logic = |
---|
491 | new G4LogicalVolume( cont_solid, |
---|
492 | matePatient, // material is not relevant here... |
---|
493 | "PhantomContainer", |
---|
494 | 0, 0, 0 ); |
---|
495 | G4VPhysicalVolume * cont_phys = |
---|
496 | new G4PVPlacement(rotm, // rotation |
---|
497 | pos, // translation |
---|
498 | cont_logic, // logical volume |
---|
499 | "PhantomContainer", // name |
---|
500 | world_logic, // mother volume |
---|
501 | false, // No op. bool. |
---|
502 | 1); // Copy number |
---|
503 | |
---|
504 | </programlisting> |
---|
505 | </informalexample> |
---|
506 | |
---|
507 | The physical volume should be assigned as the container volume of the |
---|
508 | parameterisation: |
---|
509 | |
---|
510 | <informalexample> |
---|
511 | <programlisting> |
---|
512 | param->BuildContainerSolid(cont_phys); |
---|
513 | |
---|
514 | // Assure that the voxels are completely filling the container volume |
---|
515 | // |
---|
516 | param->CheckVoxelsFillContainer( cont_solid->GetXHalfLength(), |
---|
517 | cont_solid->GetyHalfLength(), |
---|
518 | cont_solid->GetzHalfLength() ); |
---|
519 | |
---|
520 | // The parameterised volume which uses this parameterisation is placed |
---|
521 | // in the container logical volume |
---|
522 | // |
---|
523 | G4PVParameterised * patient_phys = |
---|
524 | new G4PVParameterised("Patient", // name |
---|
525 | patient_logic, // logical volume |
---|
526 | cont_logic, // mother volume |
---|
527 | kXAxis, // optimisation hint |
---|
528 | nVoxelX*nVoxelY*nVoxelZ, // number of voxels |
---|
529 | param); // parameterisation |
---|
530 | |
---|
531 | // Indicate that this physical volume is having a regular structure |
---|
532 | // |
---|
533 | patient_phys->SetRegularStructureId(1); |
---|
534 | </programlisting> |
---|
535 | </informalexample> |
---|
536 | </para> |
---|
537 | |
---|
538 | An example showing the application of the optimised navigation algorithm |
---|
539 | for phantoms geometries is available in |
---|
540 | <literal>examples/extended/medical/DICOM</literal>. It implements a real |
---|
541 | application for reading <literal>DICOM</literal> images and convert them to |
---|
542 | Geant4 geometries with defined materials and densities, allowing for |
---|
543 | different implementation solutions to be chosen (non optimised, classical |
---|
544 | 3D optimisation, nested parameterisations and use of |
---|
545 | <literal>G4PhantomParameterisation</literal>). |
---|
546 | |
---|
547 | </sect3> |
---|
548 | |
---|
549 | <!-- ******************* Section (Level#3) ****************** --> |
---|
550 | <sect3 id="sect.Geom.Navig.RunTime"> |
---|
551 | <title> |
---|
552 | Run-time commands |
---|
553 | </title> |
---|
554 | |
---|
555 | <para> |
---|
556 | When running in <emphasis>verbose</emphasis> mode (i.e. the default, |
---|
557 | <literal>G4VERBOSE</literal> set while installing the Geant4 kernel |
---|
558 | libraries), the navigator provides a few commands to control its |
---|
559 | behavior. It is possible to select different verbosity levels (up |
---|
560 | to 5), with the command: |
---|
561 | |
---|
562 | <informalexample> |
---|
563 | <programlisting> |
---|
564 | geometry/navigator/verbose [verbose_level] |
---|
565 | </programlisting> |
---|
566 | </informalexample> |
---|
567 | |
---|
568 | or to force the navigator to run in <emphasis>check</emphasis> mode: |
---|
569 | |
---|
570 | <informalexample> |
---|
571 | <programlisting> |
---|
572 | geometry/navigator/check_mode [true/false] |
---|
573 | </programlisting> |
---|
574 | </informalexample> |
---|
575 | </para> |
---|
576 | |
---|
577 | <para> |
---|
578 | The latter will force more strict and less tolerant checks in |
---|
579 | step/safety computation to verify the correctness of the solids' |
---|
580 | response in the geometry. |
---|
581 | </para> |
---|
582 | |
---|
583 | <para> |
---|
584 | By combining <emphasis>check_mode</emphasis> with verbosity level-1, additional |
---|
585 | verbosity checks on the response from the solids can be activated. |
---|
586 | </para> |
---|
587 | |
---|
588 | </sect3> |
---|
589 | |
---|
590 | <!-- ******************* Section (Level#3) ****************** --> |
---|
591 | <sect3 id="sect.Geom.Tolerance"> |
---|
592 | <title> |
---|
593 | Setting Geometry Tolerance to be relative |
---|
594 | </title> |
---|
595 | |
---|
596 | <para> |
---|
597 | The tolerance value defining the accuracy of tracking on the surfaces |
---|
598 | is by default set to a reasonably small value of <emphasis>10E-9 mm</emphasis>. |
---|
599 | Such accuracy may be however redundant for use on simulation of detectors |
---|
600 | of big size or macroscopic dimensions. Since release 9.0, it is possible to |
---|
601 | specify the surface tolerance to be relative to the extent of the world volume |
---|
602 | defined for containing the geometry setup. |
---|
603 | </para> |
---|
604 | |
---|
605 | <para> |
---|
606 | The class <literal>G4GeometryManager</literal> can be used to activate |
---|
607 | the computation of the surface tolerance to be relative to the geometry |
---|
608 | setup which has been defined. It can be done this way: |
---|
609 | |
---|
610 | <informalexample> |
---|
611 | <programlisting> |
---|
612 | G4GeometryManager::GetInstance()->SetWorldMaximumExtent(WorldExtent); |
---|
613 | </programlisting> |
---|
614 | </informalexample> |
---|
615 | |
---|
616 | where, <literal>WorldExtent</literal> is the actual maximum extent of the |
---|
617 | world volume used for placing the whole geometry setup. |
---|
618 | </para> |
---|
619 | |
---|
620 | <para> |
---|
621 | Such call to <literal>G4GeometryManager</literal> must be done |
---|
622 | <emphasis role="bold">before</emphasis> defining any geometrical component of |
---|
623 | the setup (solid shape or volume), and can be done only |
---|
624 | <emphasis role="bold">once</emphasis> ! |
---|
625 | </para> |
---|
626 | |
---|
627 | <para> |
---|
628 | The class <literal>G4GeometryTolerance</literal> is to be used for retrieving the |
---|
629 | actual values defined for tolerances, surface (Cartesian), angular or radial |
---|
630 | respectively: |
---|
631 | |
---|
632 | <informalexample> |
---|
633 | <programlisting> |
---|
634 | G4GeometryTolerance::GetInstance()->GetSurfaceTolerance(); |
---|
635 | G4GeometryTolerance::GetInstance()->GetAngularTolerance(); |
---|
636 | G4GeometryTolerance::GetInstance()->GetRadialTolerance(); |
---|
637 | </programlisting> |
---|
638 | </informalexample> |
---|
639 | |
---|
640 | </para> |
---|
641 | |
---|
642 | </sect3> |
---|
643 | |
---|
644 | </sect2> |
---|