source: trunk/Documentation/geant4/UserDocumentation/UsersGuides/ForApplicationDeveloper/html/ch04s03.html@ 901

Last change on this file since 901 was 901, checked in by garnier, 17 years ago

Add Geant4 Documentation at 8.12.2008

File size: 28.5 KB
Line 
1<html><head><meta http-equiv="Content-Type" content="text/html; charset=ISO-8859-1"><title>4.3.  Electromagnetic Field</title><link rel="stylesheet" href="../xml/XSLCustomizationLayer/G4HTMLStylesheet.css" type="text/css"><meta name="generator" content="DocBook XSL Stylesheets V1.71.1"><link rel="start" href="index.html" title="Geant4 User's Guide for Application Developers"><link rel="up" href="ch04.html" title="Chapter 4.  Detector Definition and Response"><link rel="prev" href="ch04s02.html" title="4.2.  Material"><link rel="next" href="ch04s04.html" title="4.4.  Hits"><script language="JavaScript">
2function remote_win(fName)
3{
4 var url = "AllResources/Detector/geometry.src/" + fName;
5 RemoteWin=window.open(url,"","resizable=no,toolbar=0,location=0,directories=0,status=0,menubar=0,scrollbars=0,copyhistory=0,width=520,height=520")
6 RemoteWin.creator=self
7}
8</script></head><body bgcolor="white" text="black" link="#0000FF" vlink="#840084" alink="#0000FF"><div class="navheader"><table width="100%" summary="Navigation header"><tr><th colspan="3" align="center">4.3. 
9Electromagnetic Field
10</th></tr><tr><td width="20%" align="left"><a accesskey="p" href="ch04s02.html"><img src="AllResources/IconsGIF/prev.gif" alt="Prev"></a> </td><th width="60%" align="center">Chapter 4. 
11Detector Definition and Response
12</th><td width="20%" align="right"> <a accesskey="n" href="ch04s04.html"><img src="AllResources/IconsGIF/next.gif" alt="Next"></a></td></tr></table><hr></div><div class="sect1" lang="en"><div class="titlepage"><div><div><h2 class="title" style="clear: both"><a name="sect.EMField"></a>4.3. 
13Electromagnetic Field
14</h2></div></div></div><div class="sect2" lang="en"><div class="titlepage"><div><div><h3 class="title"><a name="sect.EMField.Overview"></a>4.3.1. 
15An Overview of Propagation in a Field
16</h3></div></div></div><p>
17Geant4 is capable of describing and propagating in a variety of
18fields. Magnetic fields, electric fields and electromagnetic,
19uniform or non-uniform, can specified for a Geant4 setup. The
20propagation of tracks inside them can be performed to a
21user-defined accuracy.
22</p><p>
23In order to propagate a track inside a field, the equation of
24motion of the particle in the field is integrated. In general, this
25is done using a Runge-Kutta method for the integration of ordinary
26differential equations. However, for specific cases where an
27analytical solution is known, it is possible to utilize this
28instead. Several Runge-Kutta methods are available, suitable for
29different conditions. In specific cases (such as a uniform field
30where the analytical solution is known) different solvers can also
31be used. In addition, when an approximate analytical solution is
32known, it is possible to utilize it in an iterative manner in order
33to converge to the solution to the precision required. This latter
34method is currently implemented and can be used particularly well
35for magnetic fields that are almost uniform.
36</p><p>
37Once a method is chosen that calculates the track's propagation
38in a specific field, the curved path is broken up into linear chord
39segments. These chord segments are determined so that they closely
40approximate the curved path. The chords are then used to
41interrogate the Navigator as to whether the track has crossed a
42volume boundary. Several parameters are available to adjust the
43accuracy of the integration and the subsequent interrogation of the
44model geometry.
45</p><p>
46How closely the set of chords approximates a curved trajectory
47is governed by a parameter called the <span class="emphasis"><em>miss distance</em></span> (also
48called the <span class="emphasis"><em>chord distance</em></span> ). This is an upper bound for the
49value of the sagitta - the distance between the 'real' curved
50trajectory and the approximate linear trajectory of the chord. By
51setting this parameter, the user can control the precision of the
52volume interrogation. Every attempt has been made to ensure that
53all volume interrogations will be made to an accuracy within this
54<span class="emphasis"><em>miss distance</em></span>.
55
56</p><div class="figure"><a name="fig.EMField_1"></a><div class="figure-contents"><div class="mediaobject" align="center"><img src="./AllResources/Detector/electroMagneticField.src/MissDistance.jpg" align="middle" alt="Miss Distance"></div></div><p class="title"><b>Figure 4.6. 
57 The curved trajectory will be approximated by chords,
58 so that the maximum estimated distance between curve and chord is
59 less than the the <span class="emphasis"><em>miss distance</em></span>.
60</b></p></div><p><br class="figure-break">
61</p><p>
62In addition to the <span class="emphasis"><em>miss distance</em></span> there are two more
63parameters which the user can set in order to adjust the accuracy
64(and performance) of tracking in a field. In particular these
65parameters govern the accuracy of the intersection with a volume
66boundary and the accuracy of the integration of other steps. As
67such they play an important role for tracking.
68</p><p>
69The <span class="emphasis"><em>delta intersection</em></span> parameter is the accuracy to which
70an intersection with a volume boundary is calculated. If a
71candidate boundary intersection is estimated to have a precision
72better than this, it is accepted. This parameter is especially
73important because it is used to limit a bias that our algorithm
74(for boundary crossing in a field) exhibits. This algorithm
75calculates the intersection with a volume boundary using a chord
76between two points on the curved particle trajectory. As such, the
77intersection point is always on the 'inside' of the curve. By
78setting a value for this parameter that is much smaller than some
79acceptable error, the user can limit the effect of this bias on,
80for example, the future estimation of the reconstructed particle
81momentum.
82
83</p><div class="figure"><a name="fig.EMField_2"></a><div class="figure-contents"><div class="mediaobject" align="center"><img src="./AllResources/Detector/electroMagneticField.src/IntersectionError.jpg" align="middle" alt="Figure: The distance between the calculated chord intersection point and a computed curve point."></div></div><p class="title"><b>Figure 4.7. 
84 The distance between the calculated chord intersection
85 point C and a computed curve point D is used to determine whether C
86 is an accurate representation of the intersection of the curved
87 path ADB with a volume boundary. Here CD is likely too large, and a
88 new intersection on the chord AD will be calculated.
89</b></p></div><p><br class="figure-break">
90</p><p>
91The <span class="emphasis"><em>delta one step</em></span> parameter is the accuracy for the
92endpoint of 'ordinary' integration steps, those which do not
93intersect a volume boundary. This parameter is a limit on the
94estimated error of the endpoint of each physics step. It can be
95seen as akin to a statistical uncertainty and is not expected to
96contribute any systematic behavior to physical quantities. In
97contrast, the bias addressed by <span class="emphasis"><em>delta intersection</em></span> is
98clearly correlated with potential systematic errors in the momentum
99of reconstructed tracks. Thus very strict limits on the
100intersection parameter should be used in tracking detectors or
101wherever the intersections are used to reconstruct a track's
102momentum.
103</p><p>
104<span class="emphasis"><em>Delta intersection</em></span> and <span class="emphasis"><em>delta one step</em></span> are
105parameters of the Field Manager; the user can set them according to
106the demands of his application. Because it is possible to use more
107than one field manager, different values can be set for different
108detector regions.
109</p><p>
110Note that reasonable values for the two parameters are strongly
111coupled: it does not make sense to request an accuracy of 1 nm for
112<span class="emphasis"><em>delta intersection</em></span> and accept 100 &amp;#956m for the
113<span class="emphasis"><em>delta one step</em></span> error value. Nevertheless <span class="emphasis"><em>delta
114intersection</em></span> is the more important of the two. It is
115recommended that these parameters should not differ significantly -
116certainly not by more than an order of magnitude.
117</p></div><div class="sect2" lang="en"><div class="titlepage"><div><div><h3 class="title"><a name="sect.EMField.Pract"></a>4.3.2. 
118Practical Aspects
119</h3></div></div></div><div class="sect3" lang="en"><div class="titlepage"><div><div><h4 class="title"><a name="sect.EMField.Pract.MagField"></a>4.3.2.1. 
120Creating a Magnetic Field for a Detector
121</h4></div></div></div><p>
122The simplest way to define a field for a detector involves the
123following steps:
124
125</p><div class="orderedlist"><ol type="1" compact><li><p>
126 create a field:
127 </p><div class="informalexample"><pre class="programlisting">
128 G4UniformMagField* magField
129 = new G4UniformMagField(G4ThreeVector(0.,0.,fieldValue));
130 </pre></div><p>
131 </p></li><li><p>
132 set it as the default field:
133 </p><div class="informalexample"><pre class="programlisting">
134 G4FieldManager* fieldMgr
135 = G4TransportationManager::GetTransportationManager()
136 -&gt;GetFieldManager();
137 fieldMgr-&gt;SetDetectorField(magField);
138 </pre></div><p>
139 </p></li><li><p>
140 create the objects which calculate the trajectory:
141 </p><div class="informalexample"><pre class="programlisting">
142 fieldMgr-&gt;CreateChordFinder(magField);
143 </pre></div><p>
144</p></li></ol></div><p>
145</p><p>
146To change the accuracy of volume intersection use the
147<code class="literal">SetDeltaChord</code> method:
148
149</p><div class="informalexample"><pre class="programlisting">
150 fieldMgr-&gt;GetChordFinder()-&gt;SetDeltaChord( G4double newValue);
151</pre></div><p>
152</p></div><div class="sect3" lang="en"><div class="titlepage"><div><div><h4 class="title"><a name="sect.EMField.Pract.PartVol"></a>4.3.2.2. 
153Creating a Field for a Part of the Volume Hierarchy
154</h4></div></div></div><p>
155It is possible to create a field for a part of the detector. In
156particular it can describe the field (with pointer fEmField, for
157example) inside a logical volume and all its daughters. This can be
158done by simply creating a <code class="literal">G4FieldManager</code> and attaching it
159to a logical volume (with pointer, logicVolumeWithField, for
160example) or set of logical volumes.
161
162</p><div class="informalexample"><pre class="programlisting">
163 G4bool allLocal = true;
164 logicVolumeWithField-&gt;SetFieldManager(localFieldManager, allLocal);
165</pre></div><p>
166</p><p>
167Using the second parameter to <code class="literal">SetFieldManager</code> you choose
168whether daughter volumes of this logical volume will also be given this new
169field. If it has the value <code class="literal">true</code>, the field will be
170assigned also to its daughters, and all their sub-volumes.
171Else, if it is <code class="literal">false</code>, it will be copied only to those
172daughter volumes which do not have a field manager already.
173</p></div><div class="sect3" lang="en"><div class="titlepage"><div><div><h4 class="title"><a name="sect.EMField.Pract.Electric"></a>4.3.2.3. 
174Creating an Electric or Electromagnetic Field
175</h4></div></div></div><p>
176The design and implementation of the <span class="emphasis"><em>Field</em></span> category
177allows and enables the use of an electric or combined
178electromagnetic field. These fields can also vary with time, as can
179magnetic fields.
180</p><p>
181Source listing <a href="ch04s03.html#programlist_EMField_1" title="Example 4.12. 
182How to define a uniform electric field for the whole of a detector,
183extracted from example in examples/extended/field/field02 .
184">Example 4.12</a> shows how to
185define a uniform electric field for the whole of a detector.
186
187</p><div class="example"><a name="programlist_EMField_1"></a><p class="title"><b>Example 4.12. 
188How to define a uniform electric field for the whole of a detector,
189extracted from example in examples/extended/field/field02 .
190</b></p><div class="example-contents"><pre class="programlisting">
191 // in the header file (or first)
192 #include "G4EqMagElectricField.hh"
193 #include "G4UniformElectricField.hh"
194
195 ...
196 G4ElectricField* fEMfield;
197 G4EqMagElectricField* fEquation;
198 G4MagIntegratorStepper* fStepper;
199 G4FieldManager* fFieldMgr;
200 G4double fMinStep ;
201 G4ChordFinder* fChordFinder ;
202
203 // in the source file
204
205 {
206 fEMfield = new G4UniformElectricField(
207 G4ThreeVector(0.0,100000.0*kilovolt/cm,0.0));
208
209 // Create an equation of motion for this field
210 fEquation = new G4EqMagElectricField(fEMfield);
211
212 G4int nvar = 8;
213 fStepper = new G4ClassicalRK4( fEquation, nvar );
214
215 // Get the global field manager
216 fFieldManager= G4TransportationManager::GetTransportationManager()-&gt;
217 GetFieldManager();
218 // Set this field to the global field manager
219 fFieldManager-&gt;SetDetectorField(fEMfield );
220
221 fMinStep = 0.010*mm ; // minimal step of 10 microns
222
223 fIntgrDriver = new G4MagInt_Driver(fMinStep,
224 fStepper,
225 fStepper-&gt;GetNumberOfVariables() );
226
227 fChordFinder = new G4ChordFinder(fIntgrDriver);
228 fFieldManager-&gt;SetChordFinder( fChordFinder );
229
230 }
231</pre></div></div><p><br class="example-break">
232</p><p>
233An example with an electric field is
234examples/extended/field/field02, where the class
235F02ElectricFieldSetup demonstrates how to set these and other
236parameters, and how to choose different Integration Steppers.
237</p><p>
238The user can also create their own type of field, inheriting from
239<code class="literal">G4VField</code>, and an associated Equation of Motion class
240(inheriting from <code class="literal">G4EqRhs</code>) to simulate other types of
241fields.
242</p></div><div class="sect3" lang="en"><div class="titlepage"><div><div><h4 class="title"><a name="sect.EMField.Pract.Stepper"></a>4.3.2.4. 
243Choosing a Stepper
244</h4></div></div></div><p>
245Runge-Kutta integration is used to compute the motion of a
246charged track in a general field. There are many general steppers
247from which to choose, of low and high order, and specialized
248steppers for pure magnetic fields. By default, Geant4 uses the
249classical fourth-order Runge-Kutta stepper, which is general
250purpose and robust. If the field is known to have specific
251properties, lower or higher order steppers can be used to obtain
252the same quality results using fewer computing cycles.
253</p><p>
254In particular, if the field is calculated from a field map, a
255lower order stepper is recommended. The less smooth the field is,
256the lower the order of the stepper that should be used. The choice
257of lower order steppers includes the third order stepper
258<code class="literal">G4SimpleHeum</code>, the second order <code class="literal">G4ImplicitEuler</code>
259and <code class="literal">G4SimpleRunge</code>, and the first order
260<code class="literal">G4ExplicitEuler</code>. A first order stepper would be useful
261only for very rough fields. For somewhat smooth fields
262(intermediate), the choice between second and third order steppers
263should be made by trial and error. Trying a few different types of
264steppers for a particular field or application is suggested if
265maximum performance is a goal.
266</p><p>
267The choice of stepper depends on the type of field: magnetic or
268general. A general field can be an electric or electromagnetic
269field, it can be a magnetic field or a user-defined field (which
270requires a user-defined equation of motion.) For a general field
271several steppers are available as alternatives to the default
272(<code class="literal">G4ClassicalRK4</code>):
273
274</p><div class="informalexample"><pre class="programlisting">
275 G4int nvar = 8; // To integrate time &amp; energy
276 // in addition to position, momentum
277 G4EqMagElectricField* fEquation= new G4EqMagElectricField(fEMfield);
278
279 fStepper = new G4SimpleHeum( fEquation, nvar );
280 // 3rd order, a good alternative to ClassicalRK
281 fStepper = new G4SimpleRunge( fEquation, nvar );
282 // 2nd order, for less smooth fields
283 fStepper = new G4CashKarpRKF45( fEquation );
284 // 4/5th order for very smooth fields
285</pre></div><p>
286</p><p>
287Specialized steppers for pure magnetic fields are also
288available. They take into account the fact that a local trajectory
289in a slowly varying field will not vary significantly from a helix.
290Combining this in with a variation the Runge-Kutta method can
291provide higher accuracy at lower computational cost when large
292steps are possible.
293
294</p><div class="informalexample"><pre class="programlisting">
295 G4Mag_UsualEqRhs*
296 fEquation = new G4Mag_UsualEqRhs(fMagneticField);
297 fStepper = new G4HelixImplicitEuler( fEquation );
298 // Note that for magnetic field that do not vary with time,
299 // the default number of variables suffices.
300
301 // or ..
302 fStepper = new G4HelixExplicitEuler( fEquation );
303 fStepper = new G4HelixSimpleRunge( fEquation );
304</pre></div><p>
305</p><p>
306You can choose an alternative stepper either when the field
307manager is constructed or later. At the construction of the
308ChordFinder it is an optional argument:
309
310</p><div class="informalexample"><pre class="programlisting">
311 G4ChordFinder( G4MagneticField* itsMagField,
312 G4double stepMinimum = 1.0e-2 * mm,
313 G4MagIntegratorStepper* pItsStepper = 0 );
314</pre></div><p>
315</p><p>
316To change the stepper at a later time use
317
318</p><div class="informalexample"><pre class="programlisting">
319 pChordFinder-&gt;GetIntegrationDriver()
320 -&gt;RenewStepperAndAdjust( newStepper );
321</pre></div><p>
322</p></div><div class="sect3" lang="en"><div class="titlepage"><div><div><h4 class="title"><a name="sect.EMField.Pract.Adjust"></a>4.3.2.5. 
323How to Adjust the Accuracy of Propagation
324</h4></div></div></div><p>
325In order to obtain a particular accuracy in tracking particles
326through an electromagnetic field, it is necessary to adjust the
327parameters of the field propagation module. In the following
328section, some of these additional parameters are discussed.
329</p><p>
330When integration is used to calculate the trajectory, it is
331necessary to determine an acceptable level of numerical imprecision
332in order to get performant simulation with acceptable errors. The
333parameters in Geant4 tell the field module what level of
334integration inaccuracy is acceptable.
335</p><p>
336In all quantities which are integrated (position, momentum,
337energy) there will be errors. Here, however, we focus on the error
338in two key quantities: the position and the momentum. (The error in
339the energy will come from the momentum integration).
340</p><p>
341Three parameters exist which are relevant to the integration
342accuracy. DeltaOneStep is a distance and is roughly the position
343error which is acceptable in an integration step. Since many
344integration steps may be required for a single physics step,
345DeltaOneStep should be a fraction of the average physics step size.
346The next two parameters impose a further limit on the relative
347error of the position/momentum inaccuracy. EpsilonMin and
348EpsilonMax impose a minimum and maximum on this relative error -
349and take precedence over DeltaOneStep. (Note: if you set
350EpsilonMin=EpsilonMax=your-value, then all steps will be made to
351this relative precision.
352
353</p><div class="example"><a name="programlist_EMField_2"></a><p class="title"><b>Example 4.13. 
354How to set accuracy parameters for the 'global' field of the setup.
355</b></p><div class="example-contents"><pre class="programlisting">
356 G4FieldManager *globalFieldManager;
357
358 G4TransportationManager *transportMgr=
359 G4TransportationManager::GetTransportationManager();
360
361 globalFieldManager = transportMgr-&gt;GetFieldManager();
362 // Relative accuracy values:
363 G4double minEps= 1.0e-5; // Minimum &amp; value for smallest steps
364 G4double maxEps= 1.0e-4; // Maximum &amp; value for largest steps
365
366 globalFieldManager-&gt;SetMinimumEpsilonStep( minEps );
367 globalFieldManager-&gt;SetMaximumEpsilonStep( maxEps );
368 globalFieldManager-&gt;SetDeltaOneStep( 0.5e-3 * mm ); // 0.5 micrometer
369
370 G4cout &lt;&lt; "EpsilonStep: set min= " &lt;&lt; minEps &lt;&lt; " max= " &lt;&lt; maxEps &lt;&lt; G4endl;
371</pre></div></div><p><br class="example-break">
372</p><p>
373We note that the relevant parameters above limit the inaccuracy
374in each step. The final inaccuracy due to the full trajectory will
375accumulate!
376</p><p>
377The exact point at which a track crosses a boundary is also
378calculated with finite accuracy. To limit this inaccuracy, a
379parameter called DeltaIntersection is used. This is a maximum for
380the inaccuracy of a single boundary crossing. Thus the accuracy of
381the position of the track after a number of boundary crossings is
382directly proportional to the number of boundaries.
383</p></div><div class="sect3" lang="en"><div class="titlepage"><div><div><h4 class="title"><a name="sect.EMField.Pract.DiffAcc"></a>4.3.2.6. 
384Choosing different accuracies for the same volume
385</h4></div></div></div><p>
386It is possible to create a FieldManager which has different
387properties for particles of different momenta (or depending on
388other parameters of a track). This is useful, for example, in
389obtaining high accuracy for 'important' tracks (e.g. muons) and
390accept less accuracy in tracking others (e.g. electrons). To use
391this, you must create your own field manager which uses the
392method
393
394</p><div class="informalexample"><pre class="programlisting">
395 void ConfigureForTrack( const G4Track * );
396</pre></div><p>
397
398to configure itself using the parameters of the current track.
399An example of this will be available in
400examples/extended/field05.
401</p></div><div class="sect3" lang="en"><div class="titlepage"><div><div><h4 class="title"><a name="sect.EMField.Pract.ParaScal"></a>4.3.2.7. 
402Parameters that must scale with problem size
403</h4></div></div></div><p>
404The default settings of this module are for problems with the
405physical size of a typical high energy physics setup, that is,
406distances smaller than about one kilometer. A few parameters are
407necessary to carry this information to the magnetic field module,
408and must typically be rescaled for problems of vastly different
409sizes in order to get reasonable performance and robustness. Two of
410these parameters are the maximum acceptable step and the minimum
411step size.
412</p><p>
413The <span class="bold"><strong>maximum acceptable step</strong></span> should be set
414to a distance larger than the biggest reasonable step. If the apparatus
415in a setup has a diameter of two meters, a likely maximum acceptable
416steplength would be 10 meters. A particle could then take large
417spiral steps, but would not attempt to take, for example, a
4181000-meter-long step in the case of a very low-density material.
419Similarly, for problems of a planetary scale, such as the earth
420with its radius of roughly 6400 km, a maximum acceptabe steplength
421of a few times this value would be reasonable.
422</p><p>
423An upper limit for the size of a step is a parameter of
424<code class="literal">G4PropagatorInField</code>, and can be set by calling its
425<code class="literal">SetLargestAcceptableStep</code> method.
426</p><p>
427The <span class="bold"><strong>minimum step size</strong></span> is used during
428integration to limit the amount of work in difficult cases. It is
429possible that strong fields or integration problems can force the
430integrator to try very small steps; this parameter stops them from
431becoming unnecessarily small.
432</p><p>
433Trial steps smaller than this parameter will be treated with
434less accuracy, and may even be ignored, depending on the
435situation.
436</p><p>
437The minimum step size is a parameter of the MagInt_Driver, but
438can be set in the contstructor of G4ChordFinder, as in the source
439listing above.
440</p></div><div class="sect3" lang="en"><div class="titlepage"><div><div><h4 class="title"><a name="sect.EMField.Pract.Known"></a>4.3.2.8. 
441Known Issues
442</h4></div></div></div><p>
443Currently it is computationally expensive to change the <span class="emphasis"><em>miss
444distance</em></span> to very small values, as it causes tracks to be
445limited to curved sections whose 'bend' is smaller than this value.
446(The bend is the distance of the mid-point from the chord between
447endpoints.) For tracks with small curvature (typically low momentum
448particles in strong fields) this can cause a large number of
449steps
450
451</p><div class="itemizedlist"><ul type="disc" compact><li><p>
452 even in areas where there are no volumes to intersect
453 (something that is expected to be addressed in future development,
454 in which the safety will be utilized to partially alleviate this
455 limitation)
456 </p></li><li><p>
457 especially in a region near a volume boundary (in which case it
458 is necessary in order to discover whether a track might intersect a
459 volume for only a short distance.)
460 </p></li></ul></div><p>
461</p><p>
462Requiring such precision at the intersection is clearly expensive,
463and new development would be necessary to minimize the expense.
464</p><p>
465By contrast, changing the intersection parameter is less
466computationally expensive. It causes further calculation for only a
467fraction of the steps, in particular those that intersect a volume
468boundary.
469</p></div></div><div class="sect2" lang="en"><div class="titlepage"><div><div><h3 class="title"><a name="sect.EMField.Spin"></a>4.3.3. 
470Spin Tracking
471</h3></div></div></div><p>
472The effects of a particle's motion on the precession of its spin
473angular momentum in slowly varying external fields are simulated.
474The relativistic equation of motion for spin is known as the BMT
475equation. The equation demonstrates a remarkable property; in a
476purely magnetic field, in vacuum, and neglecting small anomalous
477magnetic moments, the particle's spin precesses in such a manner
478that the longitudinal polarization remains a constant, whatever the
479motion of the particle. But when the particle interacts with
480electric fields of the medium and multiple scatters, the spin,
481which is related to the particle's magnetic moment, does not
482participate, and the need thus arises to propagate it independent
483of the momentum vector. In the case of a polarized muon beam, for
484example, it is important to predict the muon's spin direction at
485decay-time in order to simulate the decay electron (Michel)
486distribution correctly.
487</p><p>
488In order to track the spin of a particle in a magnetic field,
489you need to code the following:
490
491</p><div class="orderedlist"><ol type="1" compact><li><p>
492 in your DetectorConstruction
493
494 </p><div class="informalexample"><pre class="programlisting">
495#include "G4Mag_SpinEqRhs.hh"
496
497 G4Mag_EqRhs* fEquation = new G4Mag_SpinEqRhs(magField);
498
499 G4MagIntegratorStepper* pStepper = new G4ClassicalRK4(fEquation,12);
500 <span class="bold"><strong>notice the 12 </strong></span>
501 </pre></div><p>
502 </p></li><li><p>
503 in your PrimaryGenerator
504
505 </p><div class="informalexample"><pre class="programlisting">
506 particleGun-&gt;SetParticlePolarization(G4ThreeVector p)
507 </pre></div><p>
508
509 for example:
510
511 </p><div class="informalexample"><pre class="programlisting">
512 particleGun-&gt;
513 SetParticlePolarization(-(particleGun-&gt;GetParticleMomentumDirection()));
514
515 // or
516 particleGun-&gt;
517 SetParticlePolarization(particleGun-&gt;GetParticleMomentumDirection()
518 .cross(G4ThreeVector(0.,1.,0.)));
519 </pre></div><p>
520
521 where you set the initial spin direction.
522 </p></li></ol></div><p>
523</p><p>
524While the G4Mag_SpinEqRhs class constructor
525
526</p><div class="informalexample"><pre class="programlisting">
527 G4Mag_SpinEqRhs::G4Mag_SpinEqRhs( G4MagneticField* MagField )
528 : G4Mag_EqRhs( MagField )
529 {
530 anomaly = 1.165923e-3;
531 }
532</pre></div><p>
533
534sets the muon anomaly by default, the class also comes with the
535public method:
536
537</p><div class="informalexample"><pre class="programlisting">
538 inline void SetAnomaly(G4double a) { anomaly = a; }
539</pre></div><p>
540
541with which you can set the magnetic anomaly to any value you
542require.
543</p><p>
544For the moment, the code is written such that field tracking of
545the spin is done only for particles with non-zero charge. Please,
546see the Forum posting:
547
548<code class="literal">
549http://geant4-hn.slac.stanford.edu:5090/HyperNews/public/get/emfields/88/3/1.html
550</code>
551
552for modifications the user is required to make to facilitate
553neutron spin tracking.
554</p></div></div><div class="navfooter"><hr><table width="100%" summary="Navigation footer"><tr><td width="40%" align="left"><a accesskey="p" href="ch04s02.html"><img src="AllResources/IconsGIF/prev.gif" alt="Prev"></a> </td><td width="20%" align="center"><a accesskey="u" href="ch04.html"><img src="AllResources/IconsGIF/up.gif" alt="Up"></a></td><td width="40%" align="right"> <a accesskey="n" href="ch04s04.html"><img src="AllResources/IconsGIF/next.gif" alt="Next"></a></td></tr><tr><td width="40%" align="left" valign="top">4.2. 
555Material
556 </td><td width="20%" align="center"><a accesskey="h" href="index.html"><img src="AllResources/IconsGIF/home.gif" alt="Home"></a></td><td width="40%" align="right" valign="top"> 4.4. 
557Hits
558</td></tr></table></div></body></html>
Note: See TracBrowser for help on using the repository browser.