source: trunk/documents/UserDoc/DocBookUsersGuides/ForApplicationDeveloper/xml/TrackingAndPhysics/trackErrorPropagation.xml @ 1208

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

ajout de la doc

File size: 19.0 KB
Line 
1<!-- ******************* Section (Level#1) ****************** -->
2<sect1 id="sect.G4Error">
3<title>
4Track Error Propagation
5</title>
6
7<para>
8The error propagation package serves to propagate one particle together with its error from a given trajectory state until a user-defined target is reached (a surface, a volume, a given track length,...).
9</para>
10
11<!-- ******************* Section (Level#2) ****************** -->
12<sect2 id="sect.G4Error.Physics">
13<title>
14Physics
15</title>
16
17<para>
18The error propagator package computes the average trajectory that a particle would follow. 
19This means that the physics list must have the following characteristics:
20</para>
21
22<itemizedlist>
23<listitem> No multiple scattering </listitem>
24<listitem> No random fluctuations for energy loss </listitem>
25<listitem> No creation of secondary tracks </listitem>
26<listitem> No hadronic processes </listitem>
27</itemizedlist>
28
29<para>
30It has also to be taken into account that when the propagation is done backwards (in the direction opposed to the one the original track traveled) the energy loss has to be changed into an energy gain.
31</para>
32
33<para>
34All this is done in the <literal>G4ErrorPhysicsList</literal> class, that is automatically set by <literal>G4ErrorPropagatorManager</literal> as the GEANT4 physics list.
35It sets <literal>G4ErrorEnergyLoss</literal> as unique electromagnetic process. This process uses the GEANT4 class <literal>G4EnergyLossForExtrapolator</literal> to compute the average energy loss for forwards or backwards propagation.
36To avoid getting too different energy loss calculation when the propagation is done forwards (when the energy at the beginning of the step is used) or backwards (when the energy at the end of the step is used, always smaller than at the beginning) <literal>G4ErrorEnergyLoss</literal> computes once the energy loss and then replaces the original energy loss by subtracting/adding half of this value (what is approximately the same as computing the energy loss with the energy at the middle of the step).
37In this way, a better calculation of the energy loss is obtained with a minimal impact on the total CPU time.
38</para>
39
40<para>
41The user may use his/her own physics list instead of <literal>G4ErrorPhysicsList</literal>. As it is not needed to define a physics list when running this package, the user may have not realized that somewhere else in his/her application it has been defined; therefore a warning will be sent to advert the user that he is using a physics list different to <literal>G4ErrorPhysicsList</literal>. If a new physics list is used, it should also initialize the <literal>G4ErrorMessenger</literal> with the classes that serve to limit the step:
42
43<informalexample>
44<programlisting>
45  G4ErrorEnergyLoss* eLossProcess = new G4ErrorEnergyLoss;
46  G4ErrorStepLengthLimitProcess* stepLengthLimitProcess = new G4ErrorStepLengthLimitProcess;
47  G4ErrorMagFieldLimitProcess* magFieldLimitProcess = new G4ErrorMagFieldLimitProcess;
48  new G4ErrorMessenger( stepLengthLimitProcess, magFieldLimitProcess, eLossProcess );
49</programlisting>
50</informalexample>
51
52</para>
53
54<para>
55To ease the use of this package in the reconstruction code, the physics list, whether <literal>G4ErrorPhysicsList</literal> or the user's one, will be automatically initialized before starting the track propagation if it has not been done by the user.
56</para>
57
58</sect2>
59
60<!-- ******************* Section (Level#2) ****************** -->
61<sect2 id="sect.G4Error.TrajState">
62<title>
63Trajectory state
64</title>
65
66<para>
67The user has to provide the particle trajectory state at the initial point.
68To do this it has to create an object of one of the children classes of <literal>G4ErrorTrajState</literal>, providing:
69</para>
70
71<itemizedlist>
72<listitem> Particle type </listitem>
73<listitem> Position </listitem>
74<listitem> Momentum </listitem>
75<listitem> Trajectory error matrix </listitem>
76</itemizedlist>
77
78<informalexample>
79<programlisting>
80 G4ErrorTrajState( const G4String&amp; partType,
81                   const G4Point3D&amp; pos,
82                   const G4Vector3D&amp; mom,
83                   const G4ErrorTrajErr&amp; errmat = G4ErrorTrajErr(5,0) );
84</programlisting>
85</informalexample>
86
87<para>
88A particle trajectory is characterized by five independent variables as a function of one parameter (e.g. the path length).
89Among the five variables, one is related to the curvature (to the absolute value of the momentum), two are related to the direction of the particle and the other two are related to the spatial location.
90</para>
91
92<para>
93There are two possible representations of these five parameters in the error propagator package: as a free trajectory state, class <literal>G4ErrorTrajStateFree</literal>, or as a trajectory state on a surface, class <literal>G4ErrorTrajStateonSurface</literal>.
94</para>
95
96
97<!-- ******************* Section (Level#3) ****************** -->
98<sect3 id="sect.G4Error.TrajState.FreTrajState">
99<title>
100Free trajectory state
101</title>
102
103In the free trajectory state representation the five trajectory parameters are
104
105<itemizedlist>
106<listitem> G4double fInvP </listitem> 
107<listitem> G4double fLambda </listitem> 
108<listitem> G4double fPhi </listitem> 
109<listitem> G4double fYPerp </listitem> 
110<listitem> G4double fZPerp </listitem> 
111</itemizedlist>
112
113<para>
114where <literal>fInvP</literal> is the inverse of the momentum.
115<literal>fLambda</literal> and <literal>fPhi</literal> are the dip and azimuthal angles related to the momentum components in the following way:
116</para>
117
118<literal>
119            p_x = p cos(lambda) cos(phi)
120            p_y = p cos(lambda) sin(phi)   
121            p_z = p sin(lambda)
122</literal>
123
124<para>
125that is, <literal>lambda = 90 - theta</literal>, where <literal>theta</literal> is the usual angle with respect to the Z axis.
126</para>
127
128<para>
129<literal>fYperp</literal> and <literal>fZperp</literal> are the coordinates of the trajectory in a local orthonormal reference frame with the X axis along the particle direction, the Y axis being parallel to the X-Y plane (obtained by the vectorial product of the global Z axis and the momentum).
130</para>
131</sect3>
132
133<!-- ******************* Section (Level#3) ****************** -->
134<sect3 id="sect.G4Error.TrajState.SurfaceTrajState">
135<title>
136Trajectory state on a surface
137</title>
138
139<para>
140In the trajectory state on a surface representation the five trajectory parameters are
141</para>
142
143<itemizedlist>
144<listitem> G4double fInvP </listitem> 
145<listitem> G4double fPV </listitem> 
146<listitem> G4double fPW </listitem> 
147<listitem> G4double fV </listitem> 
148<listitem> G4double fW </listitem> 
149</itemizedlist>
150
151<para>
152where <literal>fInvP</literal> is the inverse of the momentum;
153<literal>fPV</literal> and <literal>fPW</literal> are the momentum components in an orthonormal coordinate system with axis U, V and W;
154<literal>fV</literal> and <literal>fW</literal> are the position components on this coordinate system.
155</para>
156
157For this representation the user has to provide the plane where the parameters are calculated.
158This can be done by providing two vectors, V and W, contained in the plane:
159
160<informalexample>
161<programlisting>
162  G4ErrorSurfaceTrajState( const G4String&amp; partType,
163                           const G4Point3D&amp; pos,
164                           const G4Vector3D&amp; mom,
165                           const G4Vector3D&amp; vecV,
166                           const G4Vector3D&amp; vecW,
167                           const G4ErrorTrajErr&amp; errmat = G4ErrorTrajErr(5,0) );
168</programlisting>
169</informalexample>
170
171 or by providing a plane
172
173<informalexample>
174<programlisting>
175  G4ErrorSurfaceTrajState( const G4String&amp; partType,
176                           const G4Point3D&amp; pos,
177                           const G4Vector3D&amp; mom,
178                           const G4Plane3D&amp; plane,
179                           const G4ErrorTrajErr&amp; errmat = G4ErrorTrajErr(5,0) );
180</programlisting>
181</informalexample>
182
183In this second case the vector V is calculated as the vector in the plane perpendicular to the global vector X (if the plane normal is equal to X, Z is used instead) and W is calculated as the vector in the plane perpendicular to V.
184
185</sect3>
186</sect2>
187
188<!-- ******************* Section (Level#2) ****************** -->
189<sect2 id="sect.G4Error.Err">
190<title>
191Trajectory state error
192</title>
193
194<para>
195The 5X5 error matrix should also be provided at the creation of the trajectory state as a <literal>G4ErrorTrajErr</literal> object.
196If it is not provided a default object will be created filled with null values.
197</para>
198
199<para>
200Currently the <literal>G4ErrorTrajErr</literal> is a <literal>G4ErrorSymMatrix</literal>, a simplified version of <literal>CLHEP HepSymMatrix</literal>.
201</para>
202
203<para>
204The error matrix is given in units of GeV and cm. Therefore you should do the conversion if your code is using other units.
205</para>
206
207</sect2>
208
209<!-- ******************* Section (Level#2) ****************** -->
210<sect2 id="sect.G4Error.Target">
211<title>
212Targets
213</title>
214
215<para>
216The user has to define up to where the propagation must be done: the target.
217The target can be a surface <literal>G4ErrorSurfaceTarget</literal>, which is not part of the GEANT4 geometry.
218It can also be the surface of a GEANT4 volume  <literal>G4ErrorGeomVolumeTarget</literal>, so that the particle will be stopped when it enters this volume.
219Or it can be that the particle is stopped when a certain track length is reached, by implementing a <literal>G4ErrorTrackLengthTarget</literal>.
220</para>
221
222
223<!-- ******************* Section (Level#3) ****************** -->
224<sect3 id="sect.G4Error.Target.SurfaceTarget">
225<title>
226Surface target
227</title>
228
229When the user chooses a <literal>G4ErrorSurfaceTarget</literal> as target, the track is propagated until the surface is reached.
230This surface is not part of GEANT4 geometry, but usually traverses many GEANT4 volumes.
231The class <literal>G4ErrorNavigator</literal> takes care of the double navigation: for each step the step length is calculated as the minimum of the step length in the full geometry (up to a GEANT4 volume surface) and the distance to the user-defined surface.
232To do it, <literal>G4ErrorNavigator</literal> inherits from <literal>G4Navigator</literal> and overwrites the methods <literal>ComputeStep()</literal> and <literal>ComputeSafety()</literal>.
233
234Two types of surface are currently supported (more types could be easily implemented at user request): plane and cylindrical.
235
236<!-- ******************* Section (Level#4) ****************** -->
237<sect4 id="sect.G4Error.Target.SurfaceTarget.PlaneSurfaceTarget">
238<title>
239Plane surface target
240</title>
241
242<literal>G4ErrorPlaneSurfaceTarget</literal> implements an infinite plane surface.
243The surface can be given as the four coefficients of the plane equation <literal>ax+by+cz+d = 0</literal>:
244
245<informalexample>
246<programlisting>
247    G4ErrorPlaneSurfaceTarget(G4double a=0,
248                              G4double b=0,
249                              G4double c=0,
250                              G4double d=0);
251</programlisting>
252</informalexample>
253
254or as the normal to the plane and a point contained in it:
255
256<informalexample>
257<programlisting>
258    G4ErrorPlaneSurfaceTarget(const G4Normal3D &amp;n,
259                              const G4Point3D &amp;p);
260</programlisting>
261</informalexample>
262
263or as three points contained in it:
264
265<informalexample>
266<programlisting>
267    G4ErrorPlaneSurfaceTarget(const G4Point3D &amp;p1,
268                              const G4Point3D &amp;p2,
269                              const G4Point3D &amp;p3);
270</programlisting>
271</informalexample>
272</sect4>
273
274<!-- ******************* Section (Level#4) ****************** -->
275<sect4 id="sect.G4Error.Target.SurfaceTarget.CylSurfaceTarget">
276<title>
277Cylindrical surface target
278</title>
279
280<literal>G4ErrorCylSurfaceTarget</literal> implements an infinite-length cylindrical surface (a cylinder without end-caps).
281The surface can be given as the radius, the translation and the rotation
282
283<informalexample>
284<programlisting>
285    G4ErrorCylSurfaceTarget( const G4double&amp; radius,
286                             const G4ThreeVector&amp; trans=G4ThreeVector(),
287                             const G4RotationMatrix&amp; rotm=G4RotationMatrix() );
288</programlisting>
289</informalexample>
290
291or as the radius and the affine transformation
292
293<informalexample>
294<programlisting>
295    G4ErrorCylSurfaceTarget( const G4double&amp; radius,
296                             const G4AffineTransform&amp; trans );
297</programlisting>
298</informalexample>
299
300</sect4>
301
302</sect3>
303
304<!-- ******************* Section (Level#3) ****************** -->
305<sect3 id="sect.G4Error.Target.VolumeTarget">
306<title>
307Geometry volume target
308</title>
309
310<para>
311When the user chooses a <literal>G4ErrorGeomVolumeTarget</literal> as target, the track is propagated until the surface of a GEANT4 volume is reached. User can choose if the track will be stopped only when the track enters the volume, only when the track exits the volume or in both cases.
312</para>
313
314The object has to be instantiated giving the name of a logical volume existing in the geometry:
315
316<informalexample>
317<programlisting>
318  G4ErrorGeomVolumeTarget( const G4String&amp; name );
319</programlisting>
320</informalexample>
321
322</sect3>
323
324<!-- ******************* Section (Level#3) ****************** -->
325<sect3 id="sect.G4Error.Target.TrkLenTarget">
326<title>
327Track Length target
328</title>
329
330<para>
331When the user chooses a <literal>G4ErrorTrackLengthTarget</literal> as target, the track is propagated until the given track length is reached.
332</para>
333
334The object has to be instantiated giving the value of the track length:
335<informalexample>
336<programlisting>
337  G4ErrorTrackLengthTarget(const G4double maxTrkLength );
338</programlisting>
339</informalexample>
340
341<para>
342It is implemented as a <literal>G4VDiscreteProcess</literal> and it limits the step in <literal>PostStepGetPhysicalInteractionLength</literal>.
343To ease its use, the process is registered to all particles in the constructor.
344</para>
345
346</sect3>
347
348</sect2>
349
350
351<!-- ******************* Section (Level#2) ****************** -->
352<sect2 id="sect.G4Error.Propagation">
353<title>
354Managing the track propagation
355</title>
356
357<para>
358The user needs to propagate just one track, so there is no need of run and events. neither of <literal>G4VPrimaryGeneratorAction</literal>.
359
360<literal>G4ErrorPropagator</literal> creates a track from the information given in the <literal>G4ErrorTrajState</literal> and manages the step propagation.
361The propagation is done by the standard GEANT4 methods, invoking <literal>G4SteppingManager::Stepping() </literal> to propagate each step.
362</para>
363
364<para>
365After one step is propagated, <literal>G4ErrorPropagator</literal> takes cares of propagating the track errors for this step, what is done by <literal>G4ErrorTrajStateFree::PropagateError()</literal>.
366The equations of error propagation are only implemented in the representation of  <literal>G4ErrorTrajStateFree</literal>.
367Therefore if the user has provided instead a <literal>G4ErrorTrajStateOnSurface</literal> object, it will be transformed into a <literal>G4ErrorTrajStateFree</literal> at the beginning of tracking, and at the end it is converted back into <literal>G4ErrorTrajStateOnSurface</literal> on the target surface (on the normal plane to the surface at the final point).
368</para>
369
370<para>
371The user <literal>G4VUserTrackingAction::PreUserTrackingAction( const G4Track* )</literal> and  <literal>G4VUserTrackingAction::PreUserTrackingAction( const G4Track* )</literal> are also invoked at the beginning and at the end of the track propagation.
372</para>
373
374<para>
375<literal>G4ErrorPropagator</literal> stops the tracking when one of the three conditions is true:
376
377<itemizedlist>
378<listitem> Energy is exhausted </listitem>
379<listitem> World boundary is reached </listitem>
380<listitem> User-defined target is reached  </listitem>
381</itemizedlist>
382
383In case the defined target is not reached, <literal>G4ErrorPropagator::Propagate()</literal> returns a negative value.
384</para>
385
386<para>
387The propagation of a trajectory state until a user defined target can be done by invoking the method of <literal>G4ErrorPropagatorManager</literal>
388<informalexample>
389<programlisting>
390  G4int Propagate( G4ErrorTrajState* currentTS, const G4ErrorTarget* target,
391                   G4ErrorMode mode = G4ErrorMode_PropForwards );
392</programlisting>
393</informalexample>
394</para>
395
396You can get the pointer to the only instance of <literal>G4ErrorPropagatorManager</literal> with
397
398<informalexample>
399<programlisting>
400  G4ErrorPropagatorManager* g4emgr = G4ErrorPropagatorManager::GetErrorPropagatorManager();
401</programlisting>
402</informalexample>
403
404Another possibility is to invoke the propagation step by step, returning control to the user after each step. This can be done with the method
405
406<informalexample>
407<programlisting>
408  G4int PropagateOneStep( G4ErrorTrajState* currentTS,
409                          G4ErrorMode mode = G4ErrorMode_PropForwards );
410</programlisting>
411</informalexample>
412
413In this case you should register the target first with the command
414
415<informalexample>
416<programlisting>
417  G4ErrorPropagatorData::GetG4ErrorPropagatorData()->SetTarget( theG4eTarget );
418</programlisting>
419</informalexample>
420 
421
422<!-- ******************* Section (Level#3) ****************** -->
423<sect3 id="sect.G4Error.Propagation.ErrorPropagation">
424<title>
425Error propagation
426</title>
427
428As in the GEANT3-based GEANE package, the error propagation is based on the equations of the European Muon Collaboration, that take into account:
429
430<itemizedlist>
431<listitem> Error from curved trajectory in magnetic field </listitem>
432<listitem> Error from multiple scattering </listitem>
433<listitem> Error from ionization </listitem>
434</itemizedlist>
435
436<para>
437The formulas assume propagation along an helix.
438This means that it is necessary to make steps small enough to assure magnetic field constantness and not too big energy loss.
439</para>
440
441</sect3>
442</sect2>
443
444<!-- ******************* Section (Level#2) ****************** -->
445<sect2 id="sect.G4Error.StepLimits">
446<title>
447Limiting the step
448</title>
449
450There are three ways to limit the step.
451The first one is by using a fixed length value. This can be set by invoking the user command :
452
453<informalexample>
454<programlisting>
455G4UImanager::GetUIpointer()->ApplyCommand("/geant4e/limits/stepLength MY_VALUE MY_UNIT");
456</programlisting>
457</informalexample>
458
459The second one is by setting the maximum percentage of energy loss in the step (or energy gain is propagation is backwards).
460This can be set by invoking the user command :
461
462<informalexample>
463<programlisting>
464G4UImanager::GetUIpointer()->ApplyCommand("/geant4e/limits/energyLoss MY_VALUE");
465</programlisting>
466</informalexample>
467
468The last one is by setting the maximum difference between the value of the magnetic field at the beginning and at the end of the step. Indeed what is limited is the curvature, or exactly the value of the magnetic field divided by the value of the momentum transversal to the field.
469This can be set by invoking the user command :
470
471<informalexample>
472<programlisting>
473G4UImanager::GetUIpointer()->ApplyCommand("/geant4e/limits/magField MY_VALUE");
474</programlisting>
475</informalexample>
476
477The classes that limit the step are implemented as GEANT4 processes. Therefore, the invocation of the above-mentioned commands should only be done after the initialization (for example after <literal>G4ErrorPropagatorManager::InitGeant4e()</literal>.
478
479</sect2>
480</sect1>
481
Note: See TracBrowser for help on using the repository browser.