source: trunk/source/geometry/magneticfield/src/G4MagErrorStepper.cc @ 1274

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

update from CVS

File size: 4.6 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// $Id: G4MagErrorStepper.cc,v 1.13 2006/06/29 18:24:18 gunter Exp $
28// GEANT4 tag $Name:  $
29//
30// --------------------------------------------------------------------
31
32#include "G4MagErrorStepper.hh"
33#include "G4LineSection.hh"
34
35G4MagErrorStepper::~G4MagErrorStepper()
36{
37   delete[] yMiddle;
38   delete[] dydxMid;
39   delete[] yInitial;
40   delete[] yOneStep;
41}
42
43void
44G4MagErrorStepper::Stepper( const G4double yInput[],
45                            const G4double dydx[],
46                                  G4double hstep,
47                                  G4double yOutput[],
48                                  G4double yError []      )
49{ 
50   const G4int nvar = this->GetNumberOfVariables() ;
51   const G4int maxvar= GetNumberOfStateVariables();
52
53   G4int i;
54   // correction for Richardson Extrapolation.
55   G4double  correction = 1. / ( (1 << IntegratorOrder()) -1 );
56   
57   //  Saving yInput because yInput and yOutput can be aliases for same array
58
59   for(i=0;i<nvar;i++) yInitial[i]=yInput[i];
60   yInitial[7]= yInput[7];    // Copy the time in case ... even if not really needed
61   yMiddle[7] = yInput[7];  // Copy the time from initial value
62   yOneStep[7] = yInput[7]; // As it contributes to final value of yOutput ?
63   // yOutput[7] = yInput[7];  // -> dumb stepper does it too for RK4
64   for(i=nvar;i<maxvar;i++) yOutput[i]=yInput[i];
65   // yError[7] = 0.0;         
66
67   G4double halfStep = hstep * 0.5; 
68
69   // Do two half steps
70
71   DumbStepper  (yInitial,  dydx,   halfStep, yMiddle);
72   RightHandSide(yMiddle, dydxMid);   
73   DumbStepper  (yMiddle, dydxMid, halfStep, yOutput); 
74
75   // Store midpoint, chord calculation
76
77   fMidPoint = G4ThreeVector( yMiddle[0],  yMiddle[1],  yMiddle[2]); 
78
79   // Do a full Step
80   DumbStepper(yInitial, dydx, hstep, yOneStep);
81   for(i=0;i<nvar;i++) {
82      yError [i] = yOutput[i] - yOneStep[i] ;
83      yOutput[i] += yError[i]*correction ;  // Provides accuracy increased
84                                            // by 1 order via the
85                                            // Richardson Extrapolation 
86   }
87
88   fInitialPoint = G4ThreeVector( yInitial[0], yInitial[1], yInitial[2]); 
89   fFinalPoint   = G4ThreeVector( yOutput[0],  yOutput[1],  yOutput[2]); 
90
91   return ;
92}
93
94
95
96G4double
97G4MagErrorStepper::DistChord() const 
98{
99  // Estimate the maximum distance from the curve to the chord
100  //
101  //  We estimate this using the distance of the midpoint to
102  //  chord (the line between
103  //
104  //  Method below is good only for angle deviations < 2 pi,
105  //   This restriction should not a problem for the Runge cutta methods,
106  //   which generally cannot integrate accurately for large angle deviations.
107  G4double distLine, distChord; 
108
109  if (fInitialPoint != fFinalPoint) {
110     distLine= G4LineSection::Distline( fMidPoint, fInitialPoint, fFinalPoint );
111     // This is a class method that gives distance of Mid
112     //  from the Chord between the Initial and Final points.
113
114     distChord = distLine;
115  }else{
116     distChord = (fMidPoint-fInitialPoint).mag();
117  }
118
119  return distChord;
120}
121 
Note: See TracBrowser for help on using the repository browser.