1 | <!-- ******************************************************** --> |
---|
2 | <!-- --> |
---|
3 | <!-- [History] --> |
---|
4 | <!-- Changed by: Katsuya Amako, 30-Nov-1998 --> |
---|
5 | <!-- Converted with LaTeX2HTML 96.1-g (July 19, 1996) by --> |
---|
6 | <!-- Nikos Drakos (nikos@cbl.leeds.ac.uk), --> |
---|
7 | <!-- CBLU, University of Leeds --> |
---|
8 | <!-- Changed by: John APOSTOLAKIS, 7-Dec-1998 --> |
---|
9 | <!-- Proof read by: Joe Chuma, 29-Jun-1999 --> |
---|
10 | <!-- Changed by: Dennis Wright, 12-Dec-2002 --> |
---|
11 | <!-- Corrections for field changes: John Apostolakis, --> |
---|
12 | <!-- 17-Jun-2005 --> |
---|
13 | <!-- Corrections, changes by: John APOSTOLAKIS, 7-Jul-2005 --> |
---|
14 | <!-- Converted to DocBook: Katsuya Amako, Aug-2006 --> |
---|
15 | <!-- --> |
---|
16 | <!-- ******************************************************** --> |
---|
17 | |
---|
18 | |
---|
19 | <!-- ******************* Section (Level#1) ****************** --> |
---|
20 | <sect1 id="sect.EMField"> |
---|
21 | <title> |
---|
22 | Electromagnetic Field |
---|
23 | </title> |
---|
24 | |
---|
25 | |
---|
26 | <!-- ******************* Section (Level#2) ****************** --> |
---|
27 | <sect2 id="sect.EMField.Overview"> |
---|
28 | <title> |
---|
29 | An Overview of Propagation in a Field |
---|
30 | </title> |
---|
31 | |
---|
32 | <para> |
---|
33 | Geant4 is capable of describing and propagating in a variety of |
---|
34 | fields. Magnetic fields, electric fields and electromagnetic, |
---|
35 | uniform or non-uniform, can specified for a Geant4 setup. The |
---|
36 | propagation of tracks inside them can be performed to a |
---|
37 | user-defined accuracy. |
---|
38 | </para> |
---|
39 | |
---|
40 | <para> |
---|
41 | In order to propagate a track inside a field, the equation of |
---|
42 | motion of the particle in the field is integrated. In general, this |
---|
43 | is done using a Runge-Kutta method for the integration of ordinary |
---|
44 | differential equations. However, for specific cases where an |
---|
45 | analytical solution is known, it is possible to utilize this |
---|
46 | instead. Several Runge-Kutta methods are available, suitable for |
---|
47 | different conditions. In specific cases (such as a uniform field |
---|
48 | where the analytical solution is known) different solvers can also |
---|
49 | be used. In addition, when an approximate analytical solution is |
---|
50 | known, it is possible to utilize it in an iterative manner in order |
---|
51 | to converge to the solution to the precision required. This latter |
---|
52 | method is currently implemented and can be used particularly well |
---|
53 | for magnetic fields that are almost uniform. |
---|
54 | </para> |
---|
55 | |
---|
56 | <para> |
---|
57 | Once a method is chosen that calculates the track's propagation |
---|
58 | in a specific field, the curved path is broken up into linear chord |
---|
59 | segments. These chord segments are determined so that they closely |
---|
60 | approximate the curved path. The chords are then used to |
---|
61 | interrogate the Navigator as to whether the track has crossed a |
---|
62 | volume boundary. Several parameters are available to adjust the |
---|
63 | accuracy of the integration and the subsequent interrogation of the |
---|
64 | model geometry. |
---|
65 | </para> |
---|
66 | |
---|
67 | <para> |
---|
68 | How closely the set of chords approximates a curved trajectory |
---|
69 | is governed by a parameter called the <emphasis>miss distance</emphasis> (also |
---|
70 | called the <emphasis>chord distance</emphasis> ). This is an upper bound for the |
---|
71 | value of the sagitta - the distance between the 'real' curved |
---|
72 | trajectory and the approximate linear trajectory of the chord. By |
---|
73 | setting this parameter, the user can control the precision of the |
---|
74 | volume interrogation. Every attempt has been made to ensure that |
---|
75 | all volume interrogations will be made to an accuracy within this |
---|
76 | <emphasis>miss distance</emphasis>. |
---|
77 | |
---|
78 | <figure id="fig.EMField_1"> |
---|
79 | <title> |
---|
80 | The curved trajectory will be approximated by chords, |
---|
81 | so that the maximum estimated distance between curve and chord is |
---|
82 | less than the the <emphasis>miss distance</emphasis>. |
---|
83 | </title> |
---|
84 | <mediaobject> |
---|
85 | <imageobject role="fo"> |
---|
86 | <imagedata fileref="./AllResources/Detector/electroMagneticField.src/MissDistance.jpg" |
---|
87 | format="JPG" contentwidth="12.0cm" align="center" /> |
---|
88 | </imageobject> |
---|
89 | <imageobject role="html"> |
---|
90 | <imagedata fileref="./AllResources/Detector/electroMagneticField.src/MissDistance.jpg" |
---|
91 | format="JPG" align="center" /> |
---|
92 | </imageobject> |
---|
93 | <textobject> |
---|
94 | <phrase>Miss Distance</phrase> |
---|
95 | </textobject> |
---|
96 | </mediaobject> |
---|
97 | </figure> |
---|
98 | </para> |
---|
99 | |
---|
100 | <para> |
---|
101 | In addition to the <emphasis>miss distance</emphasis> there are two more |
---|
102 | parameters which the user can set in order to adjust the accuracy |
---|
103 | (and performance) of tracking in a field. In particular these |
---|
104 | parameters govern the accuracy of the intersection with a volume |
---|
105 | boundary and the accuracy of the integration of other steps. As |
---|
106 | such they play an important role for tracking. |
---|
107 | </para> |
---|
108 | |
---|
109 | <para> |
---|
110 | The <emphasis>delta intersection</emphasis> parameter is the accuracy to which |
---|
111 | an intersection with a volume boundary is calculated. If a |
---|
112 | candidate boundary intersection is estimated to have a precision |
---|
113 | better than this, it is accepted. This parameter is especially |
---|
114 | important because it is used to limit a bias that our algorithm |
---|
115 | (for boundary crossing in a field) exhibits. This algorithm |
---|
116 | calculates the intersection with a volume boundary using a chord |
---|
117 | between two points on the curved particle trajectory. As such, the |
---|
118 | intersection point is always on the 'inside' of the curve. By |
---|
119 | setting a value for this parameter that is much smaller than some |
---|
120 | acceptable error, the user can limit the effect of this bias on, |
---|
121 | for example, the future estimation of the reconstructed particle |
---|
122 | momentum. |
---|
123 | |
---|
124 | <figure id="fig.EMField_2"> |
---|
125 | <title> |
---|
126 | The distance between the calculated chord intersection |
---|
127 | point C and a computed curve point D is used to determine whether C |
---|
128 | is an accurate representation of the intersection of the curved |
---|
129 | path ADB with a volume boundary. Here CD is likely too large, and a |
---|
130 | new intersection on the chord AD will be calculated. |
---|
131 | </title> |
---|
132 | <mediaobject> |
---|
133 | <imageobject role="fo"> |
---|
134 | <imagedata fileref="./AllResources/Detector/electroMagneticField.src/IntersectionError.jpg" |
---|
135 | format="JPG" contentdepth="5.0cm" align="center" /> |
---|
136 | </imageobject> |
---|
137 | <imageobject role="html"> |
---|
138 | <imagedata fileref="./AllResources/Detector/electroMagneticField.src/IntersectionError.jpg" |
---|
139 | format="JPG" align="center" /> |
---|
140 | </imageobject> |
---|
141 | <textobject> |
---|
142 | <phrase>Figure: The distance between the calculated chord intersection |
---|
143 | point and a computed curve point.</phrase> |
---|
144 | </textobject> |
---|
145 | </mediaobject> |
---|
146 | </figure> |
---|
147 | </para> |
---|
148 | |
---|
149 | <para> |
---|
150 | The <emphasis>delta one step</emphasis> parameter is the accuracy for the |
---|
151 | endpoint of 'ordinary' integration steps, those which do not |
---|
152 | intersect a volume boundary. This parameter is a limit on the |
---|
153 | estimated error of the endpoint of each physics step. It can be |
---|
154 | seen as akin to a statistical uncertainty and is not expected to |
---|
155 | contribute any systematic behavior to physical quantities. In |
---|
156 | contrast, the bias addressed by <emphasis>delta intersection</emphasis> is |
---|
157 | clearly correlated with potential systematic errors in the momentum |
---|
158 | of reconstructed tracks. Thus very strict limits on the |
---|
159 | intersection parameter should be used in tracking detectors or |
---|
160 | wherever the intersections are used to reconstruct a track's |
---|
161 | momentum. |
---|
162 | </para> |
---|
163 | |
---|
164 | <para> |
---|
165 | <emphasis>Delta intersection</emphasis> and <emphasis>delta one step</emphasis> are |
---|
166 | parameters of the Field Manager; the user can set them according to |
---|
167 | the demands of his application. Because it is possible to use more |
---|
168 | than one field manager, different values can be set for different |
---|
169 | detector regions. |
---|
170 | </para> |
---|
171 | |
---|
172 | <para> |
---|
173 | Note that reasonable values for the two parameters are strongly |
---|
174 | coupled: it does not make sense to request an accuracy of 1 nm for |
---|
175 | <emphasis>delta intersection</emphasis> and accept 100 &#956m for the |
---|
176 | <emphasis>delta one step</emphasis> error value. Nevertheless <emphasis>delta |
---|
177 | intersection</emphasis> is the more important of the two. It is |
---|
178 | recommended that these parameters should not differ significantly - |
---|
179 | certainly not by more than an order of magnitude. |
---|
180 | </para> |
---|
181 | |
---|
182 | </sect2> |
---|
183 | |
---|
184 | |
---|
185 | <!-- ******************* Section (Level#2) ****************** --> |
---|
186 | <sect2 id="sect.EMField.Pract"> |
---|
187 | <title> |
---|
188 | Practical Aspects |
---|
189 | </title> |
---|
190 | |
---|
191 | <!-- ******************* Section (Level#3) ****************** --> |
---|
192 | <sect3 id="sect.EMField.Pract.MagField"> |
---|
193 | <title> |
---|
194 | Creating a Magnetic Field for a Detector |
---|
195 | </title> |
---|
196 | |
---|
197 | <para> |
---|
198 | The simplest way to define a field for a detector involves the |
---|
199 | following steps: |
---|
200 | |
---|
201 | <orderedlist spacing="compact"> |
---|
202 | <listitem><para> |
---|
203 | create a field: |
---|
204 | <informalexample> |
---|
205 | <programlisting> |
---|
206 | G4UniformMagField* magField |
---|
207 | = new G4UniformMagField(G4ThreeVector(0.,0.,fieldValue)); |
---|
208 | </programlisting> |
---|
209 | </informalexample> |
---|
210 | </para></listitem> |
---|
211 | <listitem><para> |
---|
212 | set it as the default field: |
---|
213 | <informalexample> |
---|
214 | <programlisting> |
---|
215 | G4FieldManager* fieldMgr |
---|
216 | = G4TransportationManager::GetTransportationManager() |
---|
217 | ->GetFieldManager(); |
---|
218 | fieldMgr->SetDetectorField(magField); |
---|
219 | </programlisting> |
---|
220 | </informalexample> |
---|
221 | </para></listitem> |
---|
222 | <listitem><para> |
---|
223 | create the objects which calculate the trajectory: |
---|
224 | <informalexample> |
---|
225 | <programlisting> |
---|
226 | fieldMgr->CreateChordFinder(magField); |
---|
227 | </programlisting> |
---|
228 | </informalexample> |
---|
229 | </para></listitem> |
---|
230 | </orderedlist> |
---|
231 | </para> |
---|
232 | |
---|
233 | <para> |
---|
234 | To change the accuracy of volume intersection use the |
---|
235 | <literal>SetDeltaChord</literal> method: |
---|
236 | |
---|
237 | <informalexample> |
---|
238 | <programlisting> |
---|
239 | fieldMgr->GetChordFinder()->SetDeltaChord( G4double newValue); |
---|
240 | </programlisting> |
---|
241 | </informalexample> |
---|
242 | </para> |
---|
243 | |
---|
244 | </sect3> |
---|
245 | |
---|
246 | <!-- ******************* Section (Level#3) ****************** --> |
---|
247 | <sect3 id="sect.EMField.Pract.PartVol"> |
---|
248 | <title> |
---|
249 | Creating a Field for a Part of the Volume Hierarchy |
---|
250 | </title> |
---|
251 | |
---|
252 | <para> |
---|
253 | It is possible to create a field for a part of the detector. In |
---|
254 | particular it can describe the field (with pointer fEmField, for |
---|
255 | example) inside a logical volume and all its daughters. This can be |
---|
256 | done by simply creating a <literal>G4FieldManager</literal> and attaching it |
---|
257 | to a logical volume (with pointer, logicVolumeWithField, for |
---|
258 | example) or set of logical volumes. |
---|
259 | |
---|
260 | <informalexample> |
---|
261 | <programlisting> |
---|
262 | G4bool allLocal = true; |
---|
263 | logicVolumeWithField->SetFieldManager(localFieldManager, allLocal); |
---|
264 | </programlisting> |
---|
265 | </informalexample> |
---|
266 | </para> |
---|
267 | |
---|
268 | <para> |
---|
269 | Using the second parameter to <literal>SetFieldManager</literal> you choose |
---|
270 | whether daughter volumes of this logical volume will also be given this new |
---|
271 | field. If it has the value <literal>true</literal>, the field will be |
---|
272 | assigned also to its daughters, and all their sub-volumes. |
---|
273 | Else, if it is <literal>false</literal>, it will be copied only to those |
---|
274 | daughter volumes which do not have a field manager already. |
---|
275 | </para> |
---|
276 | |
---|
277 | </sect3> |
---|
278 | |
---|
279 | <!-- ******************* Section (Level#3) ****************** --> |
---|
280 | <sect3 id="sect.EMField.Pract.Electric"> |
---|
281 | <title> |
---|
282 | Creating an Electric or Electromagnetic Field |
---|
283 | </title> |
---|
284 | |
---|
285 | <para> |
---|
286 | The design and implementation of the <emphasis>Field</emphasis> category |
---|
287 | allows and enables the use of an electric or combined |
---|
288 | electromagnetic field. These fields can also vary with time, as can |
---|
289 | magnetic fields. |
---|
290 | </para> |
---|
291 | |
---|
292 | <para> |
---|
293 | Source listing <xref linkend="programlist_EMField_1" /> shows how to |
---|
294 | define a uniform electric field for the whole of a detector. |
---|
295 | |
---|
296 | <example id="programlist_EMField_1"> |
---|
297 | <title> |
---|
298 | How to define a uniform electric field for the whole of a detector, |
---|
299 | extracted from example in examples/extended/field/field02 . |
---|
300 | </title> |
---|
301 | |
---|
302 | <programlisting> |
---|
303 | // in the header file (or first) |
---|
304 | #include "G4EqMagElectricField.hh" |
---|
305 | #include "G4UniformElectricField.hh" |
---|
306 | |
---|
307 | ... |
---|
308 | G4ElectricField* fEMfield; |
---|
309 | G4EqMagElectricField* fEquation; |
---|
310 | G4MagIntegratorStepper* fStepper; |
---|
311 | G4FieldManager* fFieldMgr; |
---|
312 | G4double fMinStep ; |
---|
313 | G4ChordFinder* fChordFinder ; |
---|
314 | |
---|
315 | // in the source file |
---|
316 | |
---|
317 | { |
---|
318 | fEMfield = new G4UniformElectricField( |
---|
319 | G4ThreeVector(0.0,100000.0*kilovolt/cm,0.0)); |
---|
320 | |
---|
321 | // Create an equation of motion for this field |
---|
322 | fEquation = new G4EqMagElectricField(fEMfield); |
---|
323 | |
---|
324 | G4int nvar = 8; |
---|
325 | fStepper = new G4ClassicalRK4( fEquation, nvar ); |
---|
326 | |
---|
327 | // Get the global field manager |
---|
328 | fFieldManager= G4TransportationManager::GetTransportationManager()-> |
---|
329 | GetFieldManager(); |
---|
330 | // Set this field to the global field manager |
---|
331 | fFieldManager->SetDetectorField(fEMfield ); |
---|
332 | |
---|
333 | fMinStep = 0.010*mm ; // minimal step of 10 microns |
---|
334 | |
---|
335 | fIntgrDriver = new G4MagInt_Driver(fMinStep, |
---|
336 | fStepper, |
---|
337 | fStepper->GetNumberOfVariables() ); |
---|
338 | |
---|
339 | fChordFinder = new G4ChordFinder(fIntgrDriver); |
---|
340 | fFieldManager->SetChordFinder( fChordFinder ); |
---|
341 | |
---|
342 | } |
---|
343 | </programlisting> |
---|
344 | </example> |
---|
345 | </para> |
---|
346 | |
---|
347 | <para> |
---|
348 | An example with an electric field is |
---|
349 | examples/extended/field/field02, where the class |
---|
350 | F02ElectricFieldSetup demonstrates how to set these and other |
---|
351 | parameters, and how to choose different Integration Steppers. |
---|
352 | </para> |
---|
353 | |
---|
354 | <para> |
---|
355 | The user can also create their own type of field, inheriting from |
---|
356 | <literal>G4VField</literal>, and an associated Equation of Motion class |
---|
357 | (inheriting from <literal>G4EqRhs</literal>) to simulate other types of |
---|
358 | fields. |
---|
359 | </para> |
---|
360 | |
---|
361 | </sect3> |
---|
362 | |
---|
363 | <!-- ******************* Section (Level#3) ****************** --> |
---|
364 | <sect3 id="sect.EMField.Pract.Stepper"> |
---|
365 | <title> |
---|
366 | Choosing a Stepper |
---|
367 | </title> |
---|
368 | |
---|
369 | <para> |
---|
370 | Runge-Kutta integration is used to compute the motion of a |
---|
371 | charged track in a general field. There are many general steppers |
---|
372 | from which to choose, of low and high order, and specialized |
---|
373 | steppers for pure magnetic fields. By default, Geant4 uses the |
---|
374 | classical fourth-order Runge-Kutta stepper, which is general |
---|
375 | purpose and robust. If the field is known to have specific |
---|
376 | properties, lower or higher order steppers can be used to obtain |
---|
377 | the same quality results using fewer computing cycles. |
---|
378 | </para> |
---|
379 | |
---|
380 | <para> |
---|
381 | In particular, if the field is calculated from a field map, a |
---|
382 | lower order stepper is recommended. The less smooth the field is, |
---|
383 | the lower the order of the stepper that should be used. The choice |
---|
384 | of lower order steppers includes the third order stepper |
---|
385 | <literal>G4SimpleHeum</literal>, the second order <literal>G4ImplicitEuler</literal> |
---|
386 | and <literal>G4SimpleRunge</literal>, and the first order |
---|
387 | <literal>G4ExplicitEuler</literal>. A first order stepper would be useful |
---|
388 | only for very rough fields. For somewhat smooth fields |
---|
389 | (intermediate), the choice between second and third order steppers |
---|
390 | should be made by trial and error. Trying a few different types of |
---|
391 | steppers for a particular field or application is suggested if |
---|
392 | maximum performance is a goal. |
---|
393 | </para> |
---|
394 | |
---|
395 | <para> |
---|
396 | The choice of stepper depends on the type of field: magnetic or |
---|
397 | general. A general field can be an electric or electromagnetic |
---|
398 | field, it can be a magnetic field or a user-defined field (which |
---|
399 | requires a user-defined equation of motion.) For a general field |
---|
400 | several steppers are available as alternatives to the default |
---|
401 | (<literal>G4ClassicalRK4</literal>): |
---|
402 | |
---|
403 | <informalexample> |
---|
404 | <programlisting> |
---|
405 | G4int nvar = 8; // To integrate time & energy |
---|
406 | // in addition to position, momentum |
---|
407 | G4EqMagElectricField* fEquation= new G4EqMagElectricField(fEMfield); |
---|
408 | |
---|
409 | fStepper = new G4SimpleHeum( fEquation, nvar ); |
---|
410 | // 3rd order, a good alternative to ClassicalRK |
---|
411 | fStepper = new G4SimpleRunge( fEquation, nvar ); |
---|
412 | // 2nd order, for less smooth fields |
---|
413 | fStepper = new G4CashKarpRKF45( fEquation ); |
---|
414 | // 4/5th order for very smooth fields |
---|
415 | </programlisting> |
---|
416 | </informalexample> |
---|
417 | </para> |
---|
418 | |
---|
419 | <para> |
---|
420 | Specialized steppers for pure magnetic fields are also |
---|
421 | available. They take into account the fact that a local trajectory |
---|
422 | in a slowly varying field will not vary significantly from a helix. |
---|
423 | Combining this in with a variation the Runge-Kutta method can |
---|
424 | provide higher accuracy at lower computational cost when large |
---|
425 | steps are possible. |
---|
426 | |
---|
427 | <informalexample> |
---|
428 | <programlisting> |
---|
429 | G4Mag_UsualEqRhs* |
---|
430 | fEquation = new G4Mag_UsualEqRhs(fMagneticField); |
---|
431 | fStepper = new G4HelixImplicitEuler( fEquation ); |
---|
432 | // Note that for magnetic field that do not vary with time, |
---|
433 | // the default number of variables suffices. |
---|
434 | |
---|
435 | // or .. |
---|
436 | fStepper = new G4HelixExplicitEuler( fEquation ); |
---|
437 | fStepper = new G4HelixSimpleRunge( fEquation ); |
---|
438 | </programlisting> |
---|
439 | </informalexample> |
---|
440 | </para> |
---|
441 | |
---|
442 | <para> |
---|
443 | You can choose an alternative stepper either when the field |
---|
444 | manager is constructed or later. At the construction of the |
---|
445 | ChordFinder it is an optional argument: |
---|
446 | |
---|
447 | <informalexample> |
---|
448 | <programlisting> |
---|
449 | G4ChordFinder( G4MagneticField* itsMagField, |
---|
450 | G4double stepMinimum = 1.0e-2 * mm, |
---|
451 | G4MagIntegratorStepper* pItsStepper = 0 ); |
---|
452 | </programlisting> |
---|
453 | </informalexample> |
---|
454 | </para> |
---|
455 | |
---|
456 | <para> |
---|
457 | To change the stepper at a later time use |
---|
458 | |
---|
459 | <informalexample> |
---|
460 | <programlisting> |
---|
461 | pChordFinder->GetIntegrationDriver() |
---|
462 | ->RenewStepperAndAdjust( newStepper ); |
---|
463 | </programlisting> |
---|
464 | </informalexample> |
---|
465 | </para> |
---|
466 | |
---|
467 | </sect3> |
---|
468 | |
---|
469 | <!-- ******************* Section (Level#3) ****************** --> |
---|
470 | <sect3 id="sect.EMField.Pract.Adjust"> |
---|
471 | <title> |
---|
472 | How to Adjust the Accuracy of Propagation |
---|
473 | </title> |
---|
474 | |
---|
475 | <para> |
---|
476 | In order to obtain a particular accuracy in tracking particles |
---|
477 | through an electromagnetic field, it is necessary to adjust the |
---|
478 | parameters of the field propagation module. In the following |
---|
479 | section, some of these additional parameters are discussed. |
---|
480 | </para> |
---|
481 | |
---|
482 | <para> |
---|
483 | When integration is used to calculate the trajectory, it is |
---|
484 | necessary to determine an acceptable level of numerical imprecision |
---|
485 | in order to get performant simulation with acceptable errors. The |
---|
486 | parameters in Geant4 tell the field module what level of |
---|
487 | integration inaccuracy is acceptable. |
---|
488 | </para> |
---|
489 | |
---|
490 | <para> |
---|
491 | In all quantities which are integrated (position, momentum, |
---|
492 | energy) there will be errors. Here, however, we focus on the error |
---|
493 | in two key quantities: the position and the momentum. (The error in |
---|
494 | the energy will come from the momentum integration). |
---|
495 | </para> |
---|
496 | |
---|
497 | <para> |
---|
498 | Three parameters exist which are relevant to the integration |
---|
499 | accuracy. DeltaOneStep is a distance and is roughly the position |
---|
500 | error which is acceptable in an integration step. Since many |
---|
501 | integration steps may be required for a single physics step, |
---|
502 | DeltaOneStep should be a fraction of the average physics step size. |
---|
503 | The next two parameters impose a further limit on the relative |
---|
504 | error of the position/momentum inaccuracy. EpsilonMin and |
---|
505 | EpsilonMax impose a minimum and maximum on this relative error - |
---|
506 | and take precedence over DeltaOneStep. (Note: if you set |
---|
507 | EpsilonMin=EpsilonMax=your-value, then all steps will be made to |
---|
508 | this relative precision. |
---|
509 | |
---|
510 | <example id="programlist_EMField_2"> |
---|
511 | <title> |
---|
512 | How to set accuracy parameters for the 'global' field of the setup. |
---|
513 | </title> |
---|
514 | |
---|
515 | <programlisting> |
---|
516 | G4FieldManager *globalFieldManager; |
---|
517 | |
---|
518 | G4TransportationManager *transportMgr= |
---|
519 | G4TransportationManager::GetTransportationManager(); |
---|
520 | |
---|
521 | globalFieldManager = transportMgr->GetFieldManager(); |
---|
522 | // Relative accuracy values: |
---|
523 | G4double minEps= 1.0e-5; // Minimum & value for smallest steps |
---|
524 | G4double maxEps= 1.0e-4; // Maximum & value for largest steps |
---|
525 | |
---|
526 | globalFieldManager->SetMinimumEpsilonStep( minEps ); |
---|
527 | globalFieldManager->SetMaximumEpsilonStep( maxEps ); |
---|
528 | globalFieldManager->SetDeltaOneStep( 0.5e-3 * mm ); // 0.5 micrometer |
---|
529 | |
---|
530 | G4cout << "EpsilonStep: set min= " << minEps << " max= " << maxEps << G4endl; |
---|
531 | </programlisting> |
---|
532 | </example> |
---|
533 | </para> |
---|
534 | |
---|
535 | <para> |
---|
536 | We note that the relevant parameters above limit the inaccuracy |
---|
537 | in each step. The final inaccuracy due to the full trajectory will |
---|
538 | accumulate! |
---|
539 | </para> |
---|
540 | |
---|
541 | <para> |
---|
542 | The exact point at which a track crosses a boundary is also |
---|
543 | calculated with finite accuracy. To limit this inaccuracy, a |
---|
544 | parameter called DeltaIntersection is used. This is a maximum for |
---|
545 | the inaccuracy of a single boundary crossing. Thus the accuracy of |
---|
546 | the position of the track after a number of boundary crossings is |
---|
547 | directly proportional to the number of boundaries. |
---|
548 | </para> |
---|
549 | |
---|
550 | </sect3> |
---|
551 | |
---|
552 | <!-- ******************* Section (Level#3) ****************** --> |
---|
553 | <sect3 id="sect.EMField.Pract.DiffAcc"> |
---|
554 | <title> |
---|
555 | Choosing different accuracies for the same volume |
---|
556 | </title> |
---|
557 | |
---|
558 | <para> |
---|
559 | It is possible to create a FieldManager which has different |
---|
560 | properties for particles of different momenta (or depending on |
---|
561 | other parameters of a track). This is useful, for example, in |
---|
562 | obtaining high accuracy for 'important' tracks (e.g. muons) and |
---|
563 | accept less accuracy in tracking others (e.g. electrons). To use |
---|
564 | this, you must create your own field manager which uses the |
---|
565 | method |
---|
566 | |
---|
567 | <informalexample> |
---|
568 | <programlisting> |
---|
569 | void ConfigureForTrack( const G4Track * ); |
---|
570 | </programlisting> |
---|
571 | </informalexample> |
---|
572 | |
---|
573 | to configure itself using the parameters of the current track. |
---|
574 | An example of this will be available in |
---|
575 | examples/extended/field05. |
---|
576 | </para> |
---|
577 | |
---|
578 | </sect3> |
---|
579 | |
---|
580 | <!-- ******************* Section (Level#3) ****************** --> |
---|
581 | <sect3 id="sect.EMField.Pract.ParaScal"> |
---|
582 | <title> |
---|
583 | Parameters that must scale with problem size |
---|
584 | </title> |
---|
585 | |
---|
586 | <para> |
---|
587 | The default settings of this module are for problems with the |
---|
588 | physical size of a typical high energy physics setup, that is, |
---|
589 | distances smaller than about one kilometer. A few parameters are |
---|
590 | necessary to carry this information to the magnetic field module, |
---|
591 | and must typically be rescaled for problems of vastly different |
---|
592 | sizes in order to get reasonable performance and robustness. Two of |
---|
593 | these parameters are the maximum acceptable step and the minimum |
---|
594 | step size. |
---|
595 | </para> |
---|
596 | |
---|
597 | |
---|
598 | <para> |
---|
599 | The <emphasis role="bold">maximum acceptable step</emphasis> should be set |
---|
600 | to a distance larger than the biggest reasonable step. If the apparatus |
---|
601 | in a setup has a diameter of two meters, a likely maximum acceptable |
---|
602 | steplength would be 10 meters. A particle could then take large |
---|
603 | spiral steps, but would not attempt to take, for example, a |
---|
604 | 1000-meter-long step in the case of a very low-density material. |
---|
605 | Similarly, for problems of a planetary scale, such as the earth |
---|
606 | with its radius of roughly 6400 km, a maximum acceptabe steplength |
---|
607 | of a few times this value would be reasonable. |
---|
608 | </para> |
---|
609 | |
---|
610 | <para> |
---|
611 | An upper limit for the size of a step is a parameter of |
---|
612 | <literal>G4PropagatorInField</literal>, and can be set by calling its |
---|
613 | <literal>SetLargestAcceptableStep</literal> method. |
---|
614 | </para> |
---|
615 | |
---|
616 | <para> |
---|
617 | The <emphasis role="bold">minimum step size</emphasis> is used during |
---|
618 | integration to limit the amount of work in difficult cases. It is |
---|
619 | possible that strong fields or integration problems can force the |
---|
620 | integrator to try very small steps; this parameter stops them from |
---|
621 | becoming unnecessarily small. |
---|
622 | </para> |
---|
623 | |
---|
624 | <para> |
---|
625 | Trial steps smaller than this parameter will be treated with |
---|
626 | less accuracy, and may even be ignored, depending on the |
---|
627 | situation. |
---|
628 | </para> |
---|
629 | |
---|
630 | <para> |
---|
631 | The minimum step size is a parameter of the MagInt_Driver, but |
---|
632 | can be set in the contstructor of G4ChordFinder, as in the source |
---|
633 | listing above. |
---|
634 | </para> |
---|
635 | |
---|
636 | </sect3> |
---|
637 | |
---|
638 | <!-- ******************* Section (Level#3) ****************** --> |
---|
639 | <sect3 id="sect.EMField.Pract.Known"> |
---|
640 | <title> |
---|
641 | Known Issues |
---|
642 | </title> |
---|
643 | |
---|
644 | <para> |
---|
645 | Currently it is computationally expensive to change the <emphasis>miss |
---|
646 | distance</emphasis> to very small values, as it causes tracks to be |
---|
647 | limited to curved sections whose 'bend' is smaller than this value. |
---|
648 | (The bend is the distance of the mid-point from the chord between |
---|
649 | endpoints.) For tracks with small curvature (typically low momentum |
---|
650 | particles in strong fields) this can cause a large number of |
---|
651 | steps |
---|
652 | |
---|
653 | <itemizedlist spacing="compact"> |
---|
654 | <listitem><para> |
---|
655 | even in areas where there are no volumes to intersect |
---|
656 | (something that is expected to be addressed in future development, |
---|
657 | in which the safety will be utilized to partially alleviate this |
---|
658 | limitation) |
---|
659 | </para></listitem> |
---|
660 | <listitem><para> |
---|
661 | especially in a region near a volume boundary (in which case it |
---|
662 | is necessary in order to discover whether a track might intersect a |
---|
663 | volume for only a short distance.) |
---|
664 | </para></listitem> |
---|
665 | </itemizedlist> |
---|
666 | </para> |
---|
667 | |
---|
668 | <para> |
---|
669 | Requiring such precision at the intersection is clearly expensive, |
---|
670 | and new development would be necessary to minimize the expense. |
---|
671 | </para> |
---|
672 | |
---|
673 | <para> |
---|
674 | By contrast, changing the intersection parameter is less |
---|
675 | computationally expensive. It causes further calculation for only a |
---|
676 | fraction of the steps, in particular those that intersect a volume |
---|
677 | boundary. |
---|
678 | </para> |
---|
679 | |
---|
680 | </sect3> |
---|
681 | </sect2> |
---|
682 | |
---|
683 | |
---|
684 | <!-- ******************* Section (Level#2) ****************** --> |
---|
685 | <sect2 id="sect.EMField.Spin"> |
---|
686 | <title> |
---|
687 | Spin Tracking |
---|
688 | </title> |
---|
689 | |
---|
690 | <para> |
---|
691 | The effects of a particle's motion on the precession of its spin |
---|
692 | angular momentum in slowly varying external fields are simulated. |
---|
693 | The relativistic equation of motion for spin is known as the BMT |
---|
694 | equation. The equation demonstrates a remarkable property; in a |
---|
695 | purely magnetic field, in vacuum, and neglecting small anomalous |
---|
696 | magnetic moments, the particle's spin precesses in such a manner |
---|
697 | that the longitudinal polarization remains a constant, whatever the |
---|
698 | motion of the particle. But when the particle interacts with |
---|
699 | electric fields of the medium and multiple scatters, the spin, |
---|
700 | which is related to the particle's magnetic moment, does not |
---|
701 | participate, and the need thus arises to propagate it independent |
---|
702 | of the momentum vector. In the case of a polarized muon beam, for |
---|
703 | example, it is important to predict the muon's spin direction at |
---|
704 | decay-time in order to simulate the decay electron (Michel) |
---|
705 | distribution correctly. |
---|
706 | </para> |
---|
707 | |
---|
708 | <para> |
---|
709 | In order to track the spin of a particle in a magnetic field, |
---|
710 | you need to code the following: |
---|
711 | |
---|
712 | <orderedlist spacing="compact"> |
---|
713 | <listitem><para> |
---|
714 | in your DetectorConstruction |
---|
715 | |
---|
716 | <informalexample> |
---|
717 | <programlisting> |
---|
718 | #include "G4Mag_SpinEqRhs.hh" |
---|
719 | |
---|
720 | G4Mag_EqRhs* fEquation = new G4Mag_SpinEqRhs(magField); |
---|
721 | |
---|
722 | G4MagIntegratorStepper* pStepper = new G4ClassicalRK4(fEquation,12); |
---|
723 | <emphasis role="bold">notice the 12 </emphasis> |
---|
724 | </programlisting> |
---|
725 | </informalexample> |
---|
726 | </para></listitem> |
---|
727 | <listitem><para> |
---|
728 | in your PrimaryGenerator |
---|
729 | |
---|
730 | <informalexample> |
---|
731 | <programlisting> |
---|
732 | particleGun->SetParticlePolarization(G4ThreeVector p) |
---|
733 | </programlisting> |
---|
734 | </informalexample> |
---|
735 | |
---|
736 | for example: |
---|
737 | |
---|
738 | <informalexample> |
---|
739 | <programlisting> |
---|
740 | particleGun-> |
---|
741 | SetParticlePolarization(-(particleGun->GetParticleMomentumDirection())); |
---|
742 | |
---|
743 | // or |
---|
744 | particleGun-> |
---|
745 | SetParticlePolarization(particleGun->GetParticleMomentumDirection() |
---|
746 | .cross(G4ThreeVector(0.,1.,0.))); |
---|
747 | </programlisting> |
---|
748 | </informalexample> |
---|
749 | |
---|
750 | where you set the initial spin direction. |
---|
751 | </para></listitem> |
---|
752 | </orderedlist> |
---|
753 | </para> |
---|
754 | |
---|
755 | <para> |
---|
756 | While the G4Mag_SpinEqRhs class constructor |
---|
757 | |
---|
758 | <informalexample> |
---|
759 | <programlisting> |
---|
760 | G4Mag_SpinEqRhs::G4Mag_SpinEqRhs( G4MagneticField* MagField ) |
---|
761 | : G4Mag_EqRhs( MagField ) |
---|
762 | { |
---|
763 | anomaly = 1.165923e-3; |
---|
764 | } |
---|
765 | </programlisting> |
---|
766 | </informalexample> |
---|
767 | |
---|
768 | sets the muon anomaly by default, the class also comes with the |
---|
769 | public method: |
---|
770 | |
---|
771 | <informalexample> |
---|
772 | <programlisting> |
---|
773 | inline void SetAnomaly(G4double a) { anomaly = a; } |
---|
774 | </programlisting> |
---|
775 | </informalexample> |
---|
776 | |
---|
777 | with which you can set the magnetic anomaly to any value you |
---|
778 | require. |
---|
779 | </para> |
---|
780 | |
---|
781 | |
---|
782 | <para> |
---|
783 | For the moment, the code is written such that field tracking of |
---|
784 | the spin is done only for particles with non-zero charge. Please, |
---|
785 | see the Forum posting: |
---|
786 | |
---|
787 | <literal> |
---|
788 | http://geant4-hn.slac.stanford.edu:5090/HyperNews/public/get/emfields/88/3/1.html |
---|
789 | </literal> |
---|
790 | |
---|
791 | for modifications the user is required to make to facilitate |
---|
792 | neutron spin tracking. |
---|
793 | </para> |
---|
794 | |
---|
795 | </sect2> |
---|
796 | </sect1> |
---|