source: trunk/source/geometry/magneticfield/include/G4RKG3_Stepper.hh @ 1340

Last change on this file since 1340 was 1340, checked in by garnier, 14 years ago

update ti head

File size: 4.1 KB
Line 
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//
28// $Id: G4RKG3_Stepper.hh,v 1.14 2010/07/23 14:13:09 tnikitin Exp $
29// GEANT4 tag $Name: field-V09-03-03 $
30//
31//
32//
33// class G4RKG3_Stepper
34//
35// Class description:
36//
37// Integrator Runga-Kutta Stepper from Geant3.
38//
39// History:
40// - Created. J.Apostolakis, V.Grichine - 30.01.97
41// -------------------------------------------------------------------
42
43#ifndef G4RKG3_Stepper_hh
44#define G4RKG3_Stepper_hh
45
46#include "G4Types.hh"
47#include "G4MagIntegratorStepper.hh"
48#include "G4ThreeVector.hh"
49
50class G4Mag_EqRhs;
51
52class G4RKG3_Stepper : public G4MagIntegratorStepper
53{
54  public:  // with description
55
56    G4RKG3_Stepper(G4Mag_EqRhs *EqRhs);
57      // Integrate over 6 variables only:  position & velocity.
58      // Not implemented yet !
59
60    ~G4RKG3_Stepper();
61
62    void Stepper( const G4double yIn[],
63                  const G4double dydx[],
64                        G4double h,
65                        G4double yOut[],
66                        G4double yErr[]  );
67      // The method which must be provided, even if less efficient.
68
69    G4double  DistChord() const ;
70 
71    void StepNoErr( const G4double tIn[8],
72                    const G4double dydx[6],
73                          G4double Step,
74                          G4double tOut[8],
75                          G4double B[3] );
76      // Integrator RK Stepper from G3 with only two field evaluation per
77      // Step. It is used in propagation initial Step by small substeps
78      // after solution error and delta geometry considerations.
79      // B[3] is magnetic field which is passed from substep to substep.
80
81    void StepWithEst( const G4double  tIn[8],
82                      const G4double dydx[6],
83                            G4double Step,
84                            G4double tOut[8],
85                            G4double& alpha2,
86                            G4double& beta2,
87                      const G4double B1[3],
88                            G4double B2[3] );
89      // Integrator for RK from G3 with evaluation of error in solution and delta
90      // geometry based on naive similarity with the case of uniform magnetic field.
91      // B1[3] is input  and is the first magnetic field values
92      // B2[3] is output and is the final magnetic field values.
93
94  public:  // without description
95
96    G4int IntegratorOrder() const { return 4; }
97
98  private:
99
100    G4ThreeVector fyInitial,
101                  fyMidPoint,
102                  fyFinal;
103   G4ThreeVector  fpInitial;
104   G4ThreeVector  BfldIn;
105   G4double       hStep;
106 
107   G4Mag_EqRhs*  fPtrMagEqOfMot;
108};
109
110#endif
Note: See TracBrowser for help on using the repository browser.