Geant4  10.00.p02
G4HelixExplicitEuler.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: G4HelixExplicitEuler.cc 66356 2012-12-18 09:02:32Z gcosmo $
28 //
29 //
30 // Helix Explicit Euler: x_1 = x_0 + helix(h)
31 // with helix(h) being a helix piece of length h
32 // most simple approach for solving linear differential equations.
33 // Take the current derivative and add it to the current position.
34 //
35 // W.Wander <wwc@mit.edu> 12/09/97
36 // -------------------------------------------------------------------
37 
38 #include "G4HelixExplicitEuler.hh"
39 #include "G4PhysicalConstants.hh"
40 #include "G4ThreeVector.hh"
41 
42 
44  const G4double*,
45  G4double Step,
46  G4double yOut[7],
47  G4double yErr[])
48 
49 {
50 
51  //Estimation of the Stepping Angle
52 
53  G4ThreeVector Bfld;
54  MagFieldEvaluate(yInput, Bfld);
55 
56  const G4int nvar = 6 ;
57  G4int i;
58  G4double yTemp[7], yIn[7] ;
59  G4ThreeVector Bfld_midpoint;
60  // Saving yInput because yInput and yOut can be aliases for same array
61  for(i=0;i<nvar;i++) yIn[i]=yInput[i];
62 
63  G4double h = Step * 0.5;
64 
65  // Do full step and two half steps
66  G4double yTemp2[7];
67  AdvanceHelix(yIn, Bfld, h, yTemp2,yTemp);
68  MagFieldEvaluate(yTemp2, Bfld_midpoint) ;
69  AdvanceHelix(yTemp2, Bfld_midpoint, h, yOut);
70 
71  // Error estimation
72  for(i=0;i<nvar;i++) {
73  yErr[i] = yOut[i] - yTemp[i] ;
74  }
75 
76 }
77 
79 {
80  // Implementation : must check whether h/R > 2 pi !!
81  // If( h/R < pi) use G4LineSection::DistLine
82  // Else DistChord=R_helix
83  //
84  G4double distChord;
85  G4double Ang_curve=GetAngCurve();
86 
87 
88  if(Ang_curve<=pi){
89  distChord=GetRadHelix()*(1-std::cos(0.5*Ang_curve));
90  }
91  else
92  if(Ang_curve<twopi){
93  distChord=GetRadHelix()*(1+std::cos(0.5*(twopi-Ang_curve)));
94  }
95  else{
96  distChord=2.*GetRadHelix();
97  }
98 
99  return distChord;
100 
101 }
102 void
104  G4ThreeVector Bfld,
105  G4double h,
106  G4double yOut[])
107 {
108 
109  AdvanceHelix(yIn, Bfld, h, yOut);
110 
111 }
void AdvanceHelix(const G4double yIn[], G4ThreeVector Bfld, G4double h, G4double yHelix[], G4double yHelix2[]=0)
CLHEP::Hep3Vector G4ThreeVector
G4double GetRadHelix() const
const G4double pi
void MagFieldEvaluate(const G4double y[], G4ThreeVector &Bfield)
int G4int
Definition: G4Types.hh:78
G4double DistChord() const
void DumbStepper(const G4double y[], G4ThreeVector Bfld, G4double h, G4double yout[])
Definition: Step.hh:41
void Stepper(const G4double y[], const G4double *, G4double h, G4double yout[], G4double yerr[])
double G4double
Definition: G4Types.hh:76
G4double GetAngCurve() const