CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
FieldStepper.cc
Go to the documentation of this file.
3 
4 #include "G4BogackiShampine45.hh"
5 #include "G4CashKarpRKF45.hh"
6 #include "G4TCashKarpRKF45.hh"
7 #include "G4ClassicalRK4.hh"
8 #include "G4TClassicalRK4.hh"
9 #include "G4DormandPrince745.hh"
10 #include "G4TDormandPrince45.hh"
11 #include "G4HelixExplicitEuler.hh"
12 #include "G4HelixHeum.hh"
13 #include "G4HelixImplicitEuler.hh"
14 #include "G4HelixSimpleRunge.hh"
15 #include "G4ImplicitEuler.hh"
16 #include "G4Mag_UsualEqRhs.hh"
17 #include "G4TMagFieldEquation.hh"
18 #include "G4NystromRK4.hh"
19 #include "G4SimpleHeum.hh"
20 #include "G4SimpleRunge.hh"
21 #include "G4TsitourasRK45.hh"
22 
23 FieldStepper::FieldStepper(G4Mag_UsualEqRhs *eq, double del, const std::string &nam)
24  : G4MagIntegratorStepper(eq, 6), theEquation(eq), theDelta(del) {
25  selectStepper(nam);
26 }
27 
29 
30 void FieldStepper::Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[]) {
31  theStepper->Stepper(y, dydx, h, yout, yerr);
32 }
33 
34 G4double FieldStepper::DistChord() const { return theStepper->DistChord(); }
35 
36 G4int FieldStepper::IntegratorOrder() const { return theStepper->IntegratorOrder(); }
37 
39  if (ss == "G4ClassicalRK4")
40  theStepper = new G4ClassicalRK4(theEquation);
41  else if (ss == "G4TClassicalRK4")
42  theStepper = new G4TClassicalRK4<G4Mag_UsualEqRhs, 8>(theEquation);
43  else if (ss == "G4NystromRK4")
44  theStepper = new G4NystromRK4(theEquation, theDelta);
45  else if (ss == "G4SimpleRunge")
46  theStepper = new G4SimpleRunge(theEquation);
47  else if (ss == "G4SimpleHeum")
48  theStepper = new G4SimpleHeum(theEquation);
49  else if (ss == "G4CashKarpRKF45")
50  theStepper = new G4CashKarpRKF45(theEquation);
51  else if (ss == "G4TCashKarpRKF45")
52  theStepper = new G4TCashKarpRKF45<G4Mag_UsualEqRhs>(theEquation);
53  else if (ss == "G4DormandPrince745")
54  theStepper = new G4DormandPrince745(theEquation);
55  else if (ss == "G4TDormandPrince45")
56  theStepper = new G4TDormandPrince45<G4Mag_UsualEqRhs>(theEquation);
57  else if (ss == "G4BogackiShampine45")
58  theStepper = new G4BogackiShampine45(theEquation);
59  else if (ss == "G4TsitourasRK45")
60  theStepper = new G4TsitourasRK45(theEquation);
61  else if (ss == "G4ImplicitEuler")
62  theStepper = new G4ImplicitEuler(theEquation);
63  else if (ss == "G4HelixExplicitEuler")
64  theStepper = new G4HelixExplicitEuler(theEquation);
65  else if (ss == "G4HelixImplicitEuler")
66  theStepper = new G4HelixImplicitEuler(theEquation);
67  else if (ss == "G4HelixSimpleRunge")
68  theStepper = new G4HelixSimpleRunge(theEquation);
69  else if (ss == "G4HelixHeum")
70  theStepper = new G4HelixHeum(theEquation);
71  else {
72  edm::LogWarning("SimG4CoreMagneticField")
73  << " FieldStepper <" << ss << "> is not known, defaulting to G4ClassicalRK4 ";
74  theStepper = new G4ClassicalRK4(theEquation);
75  }
76  edm::LogVerbatim("SimG4CoreMagneticField") << "### FieldStepper: <" << ss << ">";
77 }
Log< level::Info, true > LogVerbatim
G4int IntegratorOrder() const override
Definition: FieldStepper.cc:36
G4double DistChord() const override
Definition: FieldStepper.cc:34
~FieldStepper() override
Definition: FieldStepper.cc:28
FieldStepper(G4Mag_UsualEqRhs *eq, double del, const std::string &name)
Definition: FieldStepper.cc:23
G4MagIntegratorStepper * theStepper
Definition: FieldStepper.h:21
double theDelta
Definition: FieldStepper.h:23
void selectStepper(const std::string &)
Definition: FieldStepper.cc:38
G4Mag_UsualEqRhs * theEquation
Definition: FieldStepper.h:22
Log< level::Warning, false > LogWarning
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[]) override
Definition: FieldStepper.cc:30