Loading...
Searching...
No Matches
F01FieldSetup.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26/// \file field/field01/src/F01FieldSetup.cc
27/// \brief Implementation of the F01FieldSetup class
28//
29//
30//
31// User Field setup class implementation.
32//
33//
34//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
35//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
36
37#include "F01FieldSetup.hh"
38#include "F01FieldMessenger.hh"
39
40#include "G4MagneticField.hh"
41#include "G4UniformMagField.hh"
42#include "G4FieldManager.hh"
43#include "G4TransportationManager.hh"
44#include "G4Mag_UsualEqRhs.hh"
45#include "G4MagIntegratorStepper.hh"
46#include "G4ChordFinder.hh"
47
48#include "G4ExplicitEuler.hh"
49#include "G4ImplicitEuler.hh"
50#include "G4SimpleRunge.hh"
51#include "G4SimpleHeum.hh"
52#include "G4ClassicalRK4.hh"
53#include "G4HelixExplicitEuler.hh"
54#include "G4HelixImplicitEuler.hh"
55#include "G4HelixSimpleRunge.hh"
56#include "G4CashKarpRKF45.hh"
57#include "G4RKG3_Stepper.hh"
58#include "G4ConstRK4.hh"
59#include "G4NystromRK4.hh"
60#include "G4HelixMixedStepper.hh"
61#include "G4ExactHelixStepper.hh"
62
63#include "G4BogackiShampine23.hh"
64#include "G4BogackiShampine45.hh"
65#include "G4DormandPrince745.hh"
66#include "G4DormandPrinceRK56.hh"
67#include "G4DormandPrinceRK78.hh"
68#include "G4TsitourasRK45.hh"
69
70#include "G4PhysicalConstants.hh"
71#include "G4SystemOfUnits.hh"
72
73//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
74
80
81// Constructors:
82
83F01FieldSetup::F01FieldSetup(G4ThreeVector fieldVector,
84 G4int stepperNum,
85 G4bool useFSALstepper )
86 : fMagneticField(new G4UniformMagField(fieldVector)),
87 fUseFSALstepper(useFSALstepper),
88 fStepperType(0)
89{
90 G4cout << " F01FieldSetup: magnetic field set to Uniform( "
91 << fieldVector << " ) " << G4endl;
92
93 if( stepperNum == -1000 )
94 {
95 fUseFSALstepper = useFSALstepper;
96 if( !useFSALstepper )
97 fStepperType= 17; // Use Dormand Prince (7) 4/5 as default stepper
98 else
99 fStepperType = 101;
100 }
101 else
102 {
103 fUseFSALstepper = ( stepperNum > 0 );
104 if( stepperNum > 0 )
105 fStepperType = stepperNum;
106 else
107 fStepperType = - stepperNum;
108 }
109
111}
112
113//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
114
116 : fMagneticField(new G4UniformMagField(G4ThreeVector())),
117 fUseFSALstepper(false),
118 fStepperType(17) // Use Dormand Prince (7) 4/5 as default stepper
119{
120 G4cout << " F01FieldSetup: magnetic field set to Uniform( 0.0, 0, 0 ) "
121 << G4endl;
123}
124
125//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
126
128{
130
132
133 fMinStep = 3.0e-3*mm;
134 // minimal step of 1 um is default ==> accept any error for smallersteps!
135 fDeltaOneStep = 1.0e-5*mm;
136 // Errors of this size in an integration sub-step are acceptable
137 // except limited by the relative integration error limits (epsilon_min/max)
138 // Notes: - their initial values are set in the header.
139 // - both this and the eps min/max can be changed using Set methods.
140 fFieldManager = G4TransportationManager::GetTransportationManager()
141 ->GetFieldManager();
142
143 if( fUseFSALstepper ) {
145 }
146 else
147 {
149 // To try the symplectic method (Boris Scheme/Driver) replace the line above
150 // with the one below:
151 // CreateAndSetupBorisDriver();
152 }
153
154 G4cout << " 4/5. Updating eps_min and eps_max in Field Manager." << G4endl;
155 fFieldManager->SetChordFinder( fChordFinder );
156 fFieldManager->SetDetectorField(fMagneticField );
157
158 // For controling the accurancy
159 fFieldManager -> SetMinimumEpsilonStep( fDesiredEpsilonMin ) ;
160 //
161 // const G4double increaseFactor = 3.0 ; // typical rangle 1.0 - 10.0
162 // maxEpsilon must not exceed a ceiling, ideally 0.001 -- above this integration is unreliable
163 // const G4double maxEpsilon = increaseFactor * fDesiredEpsilonMin;
164 // fFieldManager -> SetMaximumEpsilonStep( maxEpsilon );
165
166 // For now, though, let's demonstrate how to extend the ceiling, if needed:
167 if( fDesiredEpsilonMax < G4FieldManager::GetMaxAcceptedEpsilon() )
168 {
169 G4cout << "F01FieldSetup: requesting Max(imum) Epsilon = " << fDesiredEpsilonMax << G4endl;
170 fFieldManager -> SetMaximumEpsilonStep( fDesiredEpsilonMax );
171 }
172 else
173 {
174 G4bool softFail= true;
175 G4bool good= G4FieldManager::SetMaxAcceptedEpsilon(fDesiredEpsilonMax, softFail);
176
177 if( good ) {
178 G4cout << "F01FieldSetup: requested increase Max(imum) Accepted Epsilon to "
179 << fDesiredEpsilonMax << " succeeded." << G4endl;
180 fFieldManager -> SetMaximumEpsilonStep( fDesiredEpsilonMax );
181 } else {
182 G4double maxEpsAllowed = G4FieldManager::GetMaxAcceptedEpsilon();
183 fFieldManager -> SetMaximumEpsilonStep( maxEpsAllowed );
184 G4cout << "F01FieldSetup: requested increase Max(imum) Accepted Epsilon to "
185 << fDesiredEpsilonMax << " FAILED." << G4endl
186 << " Setting it to max allowed = " << maxEpsAllowed << G4endl;
187 }
188 }
189 // To demonstrate that it is now possible to change the maximum accepted epsilon
190
191 // Note: The values of both epsilon parameters must be between
192 // fMaxAcceptedEpsilon = 0.001
193 // to ensure robustness of integration (adequate accuracy of intermediate results)
194 // and (much bigger than)
195 // fMinAcceptedEpsilon ~= 2.2e-13 ( 1000.0 * std::numeric_limits<G4double>::epsilon() )
196 // which even the best integration methods would struggle greatly to achieve.
197
198 G4cout << " Changed FieldManager epsilon values to epsilon_min= "
199 << fFieldManager -> GetMinimumEpsilonStep()
200 << " and epsilon_max= "
201 << fFieldManager -> GetMaximumEpsilonStep() << G4endl;
202}
203
204//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
205
207{
208 delete fMagneticField;
209 delete fChordFinder;
210 delete fStepper;
211 delete fFieldMessenger;
212}
213
214//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
215
217{
218 delete fChordFinder;
219 fChordFinder= nullptr;
220
221 // Update field
222 G4cout << " F01FieldSetup::CreateStepperAndChordFinder() called. " << G4endl
223 << " 1. Creating Stepper." << G4endl;
224
225 SetStepper();
226 G4cout<<"The minimal step is equal to "<<fMinStep/mm<<" mm"<<G4endl;
227
228 G4cout << " 2. Creating ChordFinder." << G4endl;
230
231 G4cout << " 3. Updating Field Manager (chord finder, field-ptr)." << G4endl;
232 fFieldManager->SetChordFinder( fChordFinder );
233 fFieldManager->SetDetectorField(fMagneticField );
234}
235
236//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
237
239{
240// Set stepper according to the stepper type
241
242 delete fStepper;
243
244 switch ( fStepperType )
245 {
246 // The new default in G4 and here ( since G4 10.4 Dec 2017 )
247 case 17:
248 case 457:
249 case 745:
251 G4cout<<"G4DormandPrince745 Stepper is chosen"<<G4endl;
252 break;
253
254 case 0:
256 G4cout<<"G4ExplicitEuler is chosen."<<G4endl;
257 break;
258 case 1:
260 G4cout<<"G4ImplicitEuler is chosen"<<G4endl;
261 break;
262 case 2:
264 G4cout<<"G4SimpleRunge is chosen"<<G4endl;
265 break;
266 case 3:
268 G4cout<<"G4SimpleHeum is chosen"<<G4endl;
269 break;
270 case 4:
272 G4cout<<"G4ClassicalRK4 (default) is chosen"<<G4endl;
273 break;
274 case 5:
276 G4cout<<"G4HelixExplicitEuler is chosen"<<G4endl;
277 break;
278 case 6:
280 G4cout<<"G4HelixImplicitEuler is chosen"<<G4endl;
281 break;
282 case 7:
284 G4cout<<"G4HelixSimpleRunge is chosen"<<G4endl;
285 break;
286 case 8:
288 G4cout<<"G4CashKarpRKF45 is chosen"<<G4endl;
289 break;
290 case 9:
292 G4cout<<"G4RKG3_Stepper is chosen"<<G4endl;
293 break;
294 case 10:
296 G4cout<<"G4ExactHelixStepper is chosen"<<G4endl;
297 break;
298 case 11:
300 G4cout<<"G4HelixMixedStepper is chosen"<<G4endl;
301 break;
302 case 12:
304 G4cout<<"G4ConstRK4 Stepper is chosen"<<G4endl;
305 break;
306 case 13:
307 case 40:
309 G4cout<<" G4NystromRK4 Stepper is chosen"<<G4endl;
310 break;
311 case 14:
312 case 23:
314 G4cout<<"G4BogackiShampine23 Stepper is chosen"<<G4endl;
315 break;
316
317 // Other optimised 4/5th order embedded steppers
318 case 15:
319 case 45:
321 G4cout<<"G4BogackiShampine45 Stepper is chosen"<<G4endl;
322 break;
323
324 // case 145:
325 case kTsitouras45:
327 G4cout<<"G4TsitourasRK45 Stepper is chosen"<<G4endl;
328 break;
329
330 // Higher order embedded steppers - for very smooth fields
331 case 56:
333 G4cout<<"G4DormandPrinceRK56 Stepper is chosen"<<G4endl;
334 break;
335 case 78:
337 G4cout<<"G4DormandPrinceRK78 Stepper is chosen"<<G4endl;
338 break;
339
340 default:
341 // fStepper = new G4ClassicalRK4( fEquation );
342 // G4cout<<"G4ClassicalRK4 Stepper (default) is chosen"<<G4endl;
344 G4cout<<"G4DormandPrince745 (default) Stepper is chosen"<<G4endl;
345 break;
346 }
347}
348
349//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
350
351#include "G4VIntegrationDriver.hh"
352#include "G4FSALIntegrationDriver.hh"
353#include "G4RK547FEq1.hh"
354#include "G4RK547FEq2.hh"
355#include "G4RK547FEq3.hh"
356
359{
360 // using FsalStepperType = G4RK547FEq1;
361 const char *methodName= "F01FieldSetup::CreateFSALStepperAndDriver()";
362 delete fStepper;
363 fStepper = nullptr;
364
365 G4cout << " F01FieldSetup::CreateFSALStepperAndDriver() called. " << G4endl;
366 G4cout << " 1. Creating Stepper." << G4endl;
367 // auto fsalStepper = new FsalStepperType( fEquation );
368 G4RK547FEq1* stepper1 = nullptr;
369 G4RK547FEq2* stepper2 = nullptr;
370 G4RK547FEq3* stepper3 = nullptr;
371
372 G4cout << " 2. Creating FSAL Driver." << G4endl;
373 G4VIntegrationDriver* fsalDriver = nullptr;
374 switch ( fStepperType )
375 {
376 case 1:
377 case 101:
378 stepper1 = new G4RK547FEq1( fEquation );
379 fsalDriver = new G4FSALIntegrationDriver<G4RK547FEq1>( fMinStep, stepper1 );
380 G4cout << " Stepper type '1' is G4RK547FEq1 stepper (in FSAL mode) with FSAL driver. "
381 << G4endl;
382 fStepper = stepper1;
383 stepper1 = nullptr;
384 break;
385
386 case 2:
387 case 102:
388 stepper2= new G4RK547FEq2( fEquation );
389 fsalDriver = new G4FSALIntegrationDriver<G4RK547FEq2>( fMinStep, stepper2 );
390 G4cout << " Stepper type '2' is G4RK547FEq2 stepper (in FSAL mode) with FSAL driver. "
391 << G4endl;
392 fStepper = stepper2;
393 stepper2 = nullptr;
394 break;
395
396 case 3:
397 case 103:
398 stepper3 = new G4RK547FEq3( fEquation );
399 fsalDriver = new G4FSALIntegrationDriver<G4RK547FEq3>( fMinStep, stepper3 );
400 G4cout << " Stepper type '3' is G4RK547FEq3 stepper (in FSAL mode) with FSAL driver. "
401 << G4endl;
402 fStepper = stepper3;
403 stepper3 = nullptr;
404 break;
405
406 default:
407 G4cout << " Warning from " << methodName << " : stepperType (= "
408 << fStepperType << " ) is unknown. " << G4endl
409 << " Using value '1' instead - i.e. G4RK547FEq1 stepper. "
410 << G4endl;
411 stepper1 = new G4RK547FEq1( fEquation );
412 fsalDriver = new G4FSALIntegrationDriver<G4RK547FEq1>( fMinStep, stepper1 );
413 fStepper = stepper1;
414 stepper1 = nullptr;
415 break;
416 }
417
418 delete stepper1; stepper1 = nullptr;
419 delete stepper2; stepper2 = nullptr;
420 delete stepper3; stepper3 = nullptr;
421
422 if( fsalDriver )
423 fStepper = fsalDriver->GetStepper();
424
425 return fsalDriver;
426}
427
429{
430 // using FsalStepperType = G4DormandPrince745; // eventually ?
431 delete fChordFinder;
432 fChordFinder= nullptr;
433
434 G4cout << " F01FieldSetup::CreateFSALStepperAndChordFinder() called. " << G4endl;
435
436 auto FSALdriver= CreateFSALStepperAndDriver();
437 fDriver = FSALdriver;
438 G4cout<<"The minimal step is equal to "<<fMinStep/mm<<" mm"<<G4endl;
439
440 G4cout << " 3. Creating ChordFinder." << G4endl;
441 fChordFinder = new G4ChordFinder( FSALdriver ); // ( fMagneticField, fMinStep, fStepper );
442}
443
444
445//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
446
447void F01FieldSetup::SetFieldZValue(G4double fieldStrength)
448{
449 // Set the value of the Global Field to fieldValue along Z
450
451 SetFieldValue(G4ThreeVector(0, 0, fieldStrength));
452}
453
454//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
455
456void F01FieldSetup::SetFieldValue(G4ThreeVector fieldVector)
457{
458 // Set the value of the Global Field
459
460 delete fMagneticField;
461
462#ifdef G4VERBOSE
463 G4cout << "Setting Field strength to "
464 << fieldVector / gauss << " Gauss." << G4endl;
465#endif
466
467 if (fieldVector != G4ThreeVector(0.,0.,0.))
468 {
469 fMagneticField = new G4UniformMagField(fieldVector);
470 }
471 else
472 {
473#ifdef G4VERBOSE
474 G4cout << " Magnetic field pointer is null." << G4endl;
475#endif
476 // If the new field's value is Zero, signal it as below
477 // so that it is not used for propagation.
478 fMagneticField = nullptr;
479 }
480
481 // Set this as the field of the global Field Manager
482 GetGlobalFieldManager()->SetDetectorField(fMagneticField);
483
484 // Now notify equation of new field
485 fEquation->SetFieldObj( fMagneticField );
486
487}
488
489//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
490
492{
493 // Utility method
494
495 return G4TransportationManager::GetTransportationManager()
496 ->GetFieldManager();
497}
498
499//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
500
501#include "G4BorisScheme.hh"
502#include"G4BorisDriver.hh"
503
504void
506{
507
508 G4cout << " F01FieldSetup::CreateAndSetupBorisDriver() called. " << G4endl;
509 G4cout << " 1. Creating Scheme (Stepper)." << G4endl;
510 auto borisStepr = new G4BorisScheme(fEquation);
511 G4cout << " 2. Creating Driver." << G4endl;
512 auto driver = new G4BorisDriver(fMinStep, borisStepr);
513
514 GetGlobalFieldManager()->SetFieldChangesEnergy(true); // To test energy conservation!!
515
516 G4cout << " 3. Creating ChordFinder." << G4endl;
517 fChordFinder = new G4ChordFinder( driver );
518
519 G4cout << " 4. Updating Field Manager (with ChordFinder, field)." << G4endl;
520 fFieldManager->SetChordFinder( fChordFinder );
521 fFieldManager->SetDetectorField(fMagneticField );
522}
523
524//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
525
F01FieldMessenger allows interactive user control of.
EStepperNumber
@ kDormandPrince56
@ kNystromRK4
@ kTsitouras45
@ kDormandPrince45
@ kDormandPrince78
@ kCashKarp
@ kBogackiShampine45
@ kBogackiShampine23
@ kClassicalRK4
Definition of the F01FieldSetup class.
G4double fDesiredEpsilonMin
G4MagneticField * fMagneticField
void CreateStepperAndChordFinder()
G4FieldManager * fFieldManager
G4double fDeltaOneStep
G4VIntegrationDriver * fDriver
void CreateFSALStepperAndChordFinder()
virtual ~F01FieldSetup()
G4FieldManager * GetGlobalFieldManager()
void SetFieldZValue(G4double fieldValue)
G4Mag_UsualEqRhs * fEquation
G4MagIntegratorStepper * fStepper
G4VIntegrationDriver * CreateFSALStepperAndDriver()
void CreateAndSetupBorisDriver()
G4ChordFinder * fChordFinder
void SetFieldValue(G4ThreeVector fieldVector)
G4double fDesiredEpsilonMax
F01FieldMessenger * fFieldMessenger

Applications | User Support | Publications | Collaboration