source: trunk/documents/UserDoc/DocBookUsersGuides/ForApplicationDeveloper/xml/FAQ/geometry.xml @ 1013

Last change on this file since 1013 was 904, checked in by garnier, 16 years ago

ajout de la doc

File size: 4.9 KB
Line 
1<!-- ******************************************************** -->
2<!--                                                          -->
3<!--  [History]                                               -->
4<!--    1st version created:  Katsuya Amako, Dec-2006         -->
5<!--                                                          -->
6<!-- ******************************************************** -->
7
8
9<!-- ********************* QandA Section ******************** -->
10<section id="qanda.Geometry">
11<title>
12Geometry
13</title>
14
15<qandaset defaultlabel="qanda">
16
17<!-- ******* QandA Entry ******** -->
18<qandaentry  id="qanda.Geometry.GenPnt">
19
20<question><para>
21  I have a generic point and I would like to know in which physical
22  volume I'm located in my detector geometry.
23</para></question>
24
25<answer><para>
26  The best way of doing this is by invoking the G4Navigator. First
27  get a pointer of the navigator through the G4TransportationManager,
28  and then locate the point. i.e.
29  <informalexample><programlisting>
30     #include "G4TransportationManager.hh"
31     #include "G4Navigator.hh"
32     G4ThreeVector myPoint = ....;
33     G4Navigator* theNavigator = G4TransportationManager::GetTransportationManager()
34                                 -&gt;GetNavigatorForTracking();
35     G4VPhysicalVolume* myVolume = theNavigator-&gt;LocateGlobalPointAndSetup(myPoint);
36  </programlisting></informalexample>
37
38  <note>
39  <title>Note</title>
40  <para>
41    by using the navigator for tracking as shown above, the actual
42    particle gets also -relocated- in the specified position. Therefore,
43    if this information is needed during tracking time, in order to
44    avoid affecting tracking, you should either use an alternative
45    G4Navigator object (which you then assign to your world-volume),
46    or you access the information through the track or touchable as
47    specified in the FAQ for tracking and steps.
48  </para>
49  </note> 
50</para></answer>
51
52</qandaentry>
53
54<!-- ******* QandA Entry ******** -->
55<qandaentry id="qanda.Geometry.DghtVol">
56
57<question><para>
58  How can I access the daughter volumes of a specific physical volume?
59</para></question>
60
61<answer><para>
62  Through the associated logical volume.
63  <informalexample><programlisting>
64      G4VPhysicalVolume* myPVolume = ....;
65      G4LogicalVolume* myLVolume = myPVolume-&gt;GetLogicalVolume();
66      for (G4int i=0; iGetNoDaughters(); i++) 
67        myPVolume = myLVolume-&gt;GetDaughter(i);
68  </programlisting></informalexample>
69</para></answer>
70
71</qandaentry>
72
73<!-- ******* QandA Entry ******** -->
74<qandaentry id="qanda.Geometry.CpyNum">
75
76<question><para>
77  How can I identify the exact copy-number of a specific physical volume
78  in my mass geometry? I tried with GetCopyNo()  from my physical volume
79  pointer, but it doesn't seem to work!
80</para></question>
81
82<answer><para>
83  The correct way to identify -uniquely- a physical volume in your mass
84  geometry is by using the touchables (see also section 4.1.5 of the
85  User's Guide for Application Developers), as follows:
86  <informalexample><programlisting>
87      G4Step* aStep = ..;
88      G4StepPoint* preStepPoint = aStep-&gt;GetPreStepPoint();
89      G4TouchableHandle theTouchable = preStepPoint-&gt;GetTouchableHandle();
90      G4int copyNo = theTouchable-&gt;GetCopyNumber();
91      G4int motherCopyNo = theTouchable-&gt;GetCopyNumber(1);
92  </programlisting></informalexample>
93
94  where <literal>Copy</literal> here stays for any duplicated instance of
95  a physical volume, either if it is a <literal>G4PVPlacement</literal> 
96  (multiple placements of the same logical volume) or a
97  <literal>G4PVReplica/G4PVParameterised</literal>.
98  The method <literal>GetCopyNo()</literal> is meant to return only the
99  serial number of placements not duplicated in the geometry tree.
100</para></answer>
101
102</qandaentry>
103
104<!-- ******* QandA Entry ******** -->
105<qandaentry id="qanda.Geometry.ConvGtoL">
106
107<question><para>
108  How can I determine the exact position in global coordinates in my mass
109  geometry during tracking and how can I convert it to coordinates local
110  to the current volume ?
111</para></question>
112
113<answer><para>
114  You need again to do it through the touchables (see also section 4.1.5
115  of the User's Guide for Application Developers), as follows:
116  <informalexample><programlisting>
117      G4Step* aStep = ..;
118      G4StepPoint* preStepPoint = aStep-&gt;GetPreStepPoint();
119      G4TouchableHandle theTouchable = preStepPoint-&gt;GetTouchableHandle();
120      G4ThreeVector worldPosition = preStepPoint-&gt;GetPosition();
121      G4ThreeVector localPosition = theTouchable-&gt;GetHistory()-&gt;
122                                    GetTopTransform().TransformPoint(worldPosition);
123  </programlisting></informalexample>
124
125  where <literal>worldPosition</literal> here stays for the position related
126  to the world volume, while <literal>localPosition</literal> refers to the
127  coordinates local to the volume where the particle is currently placed.
128</para></answer>
129
130</qandaentry>
131
132
133</qandaset>
134</section>
Note: See TracBrowser for help on using the repository browser.