Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ExactHelixStepper.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 //
27 // $Id: G4ExactHelixStepper.cc 66356 2012-12-18 09:02:32Z gcosmo $
28 //
29 // Helix a-la-Explicity Euler: x_1 = x_0 + helix(h)
30 // with helix(h) being a helix piece of length h
31 // simplest approach for solving linear differential equations.
32 // Take the current derivative and add it to the current position.
33 //
34 // As the field is assumed constant, an error is not calculated.
35 //
36 // Author: J. Apostolakis, 28 Jan 2005
37 // Implementation adapted from ExplicitEuler of W.Wander
38 // -------------------------------------------------------------------
39 
40 #include "G4ExactHelixStepper.hh"
41 #include "G4PhysicalConstants.hh"
42 #include "G4ThreeVector.hh"
43 #include "G4LineSection.hh"
44 
46  : G4MagHelicalStepper(EqRhs),
47  fBfieldValue(DBL_MAX, DBL_MAX, DBL_MAX),
48  fPtrMagEqOfMot(EqRhs)
49 {
50  ;
51 }
52 
54 
55 void
57  const G4double*,
58  G4double hstep,
59  G4double yOut[],
60  G4double yErr[] )
61 {
62  const G4int nvar = 6;
63 
64  G4int i;
65  G4ThreeVector Bfld_value;
66 
67  MagFieldEvaluate(yInput, Bfld_value);
68  AdvanceHelix(yInput, Bfld_value, hstep, yOut);
69 
70  // We are assuming a constant field: helix is exact
71  //
72  for(i=0;i<nvar;i++)
73  {
74  yErr[i] = 0.0 ;
75  }
76 
77  fBfieldValue=Bfld_value;
78 }
79 
80 void
82  G4ThreeVector Bfld,
83  G4double h,
84  G4double yOut[])
85 {
86  // Assuming a constant field: solution is a helix
87 
88  AdvanceHelix(yIn, Bfld, h, yOut);
89 
90  G4Exception("G4ExactHelixStepper::DumbStepper",
91  "GeomField0002", FatalException,
92  "Should not be called. Stepper must do all the work." );
93 }
94 
95 
96 // ---------------------------------------------------------------------------
97 
100 {
101  // Implementation : must check whether h/R > pi !!
102  // If( h/R < pi) DistChord=h/2*std::tan(Ang_curve/4)
103  // Else DistChord=R_helix
104 
105  G4double distChord;
106  G4double Ang_curve=GetAngCurve();
107 
108  if (Ang_curve<=pi)
109  {
110  distChord=GetRadHelix()*(1-std::cos(0.5*Ang_curve));
111  }
112  else if(Ang_curve<twopi)
113  {
114  distChord=GetRadHelix()*(1+std::cos(0.5*(twopi-Ang_curve)));
115  }
116  else
117  {
118  distChord=2.*GetRadHelix();
119  }
120 
121  return distChord;
122 }
123 
124 G4int
126 {
127  return 1;
128 }
void AdvanceHelix(const G4double yIn[], G4ThreeVector Bfld, G4double h, G4double yHelix[], G4double yHelix2[]=0)
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])
G4double GetRadHelix() const
void MagFieldEvaluate(const G4double y[], G4ThreeVector &Bfield)
G4double DistChord() const
int G4int
Definition: G4Types.hh:78
virtual G4int IntegratorOrder() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void DumbStepper(const G4double y[], G4ThreeVector Bfld, G4double h, G4double yout[])
double G4double
Definition: G4Types.hh:76
G4double GetAngCurve() const
#define DBL_MAX
Definition: templates.hh:83
G4ExactHelixStepper(G4Mag_EqRhs *EqRhs)