source: trunk/source/geometry/magneticfield/src/G4HelixMixedStepper.cc @ 1202

Last change on this file since 1202 was 831, checked in by garnier, 16 years ago

import all except CVS

File size: 7.0 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// class G4HelixMixedStepper
27//
28// Class description:
29//
30// G4HelixMixedStepper split the Method used for Integration in two:
31//
32// If Stepping Angle ( h / R_curve) < pi/3 
33//        use Stepper for small step(ClassicalRK4 by default)
34// Else use  HelixExplicitEuler Stepper
35//
36// History:
37// Derived from ExactHelicalStepper 18/05/07
38//
39// -------------------------------------------------------------------------
40
41#include "G4HelixMixedStepper.hh"
42#include "G4ClassicalRK4.hh"
43#include "G4CashKarpRKF45.hh"
44#include "G4SimpleRunge.hh"
45#include "G4HelixImplicitEuler.hh"
46#include "G4HelixExplicitEuler.hh"
47#include "G4HelixSimpleRunge.hh"
48#include "G4ExactHelixStepper.hh"
49#include "G4ExplicitEuler.hh"
50#include "G4ImplicitEuler.hh"
51#include "G4SimpleHeum.hh"
52#include "G4RKG3_Stepper.hh"
53
54#include "G4ThreeVector.hh"
55#include "G4LineSection.hh"
56G4HelixMixedStepper::G4HelixMixedStepper(G4Mag_EqRhs *EqRhs,G4int fStepperNumber)
57  : G4MagHelicalStepper(EqRhs)
58   
59{
60   SetVerbose(1); fNumCallsRK4=0; fNumCallsHelix=0;
61   if(!fStepperNumber) fStepperNumber=4;
62   fRK4Stepper =  SetupStepper(EqRhs, fStepperNumber);
63}
64
65
66G4HelixMixedStepper::~G4HelixMixedStepper() {
67     
68     delete(fRK4Stepper);
69     if (fVerbose>0){ PrintCalls();};
70} 
71void G4HelixMixedStepper::Stepper(  const G4double  yInput[7],
72                               const G4double dydx[7],
73                                     G4double Step,
74                                     G4double yOut[7],
75                                     G4double yErr[])
76
77{
78
79 //Estimation of the Stepping Angle
80
81  G4ThreeVector Bfld;
82  MagFieldEvaluate(yInput, Bfld); 
83
84  G4double Bmag = Bfld.mag();
85  const G4double *pIn = yInput+3;
86  G4ThreeVector initVelocity= G4ThreeVector( pIn[0], pIn[1], pIn[2]);
87  G4double      velocityVal = initVelocity.mag();
88  G4double R_1; 
89  G4double Ang_curve;
90
91   R_1=std::abs(GetInverseCurve(velocityVal,Bmag));
92   Ang_curve=R_1*Step;
93   SetAngCurve(Ang_curve);
94   SetCurve(std::abs(1/R_1));
95   
96
97   if(Ang_curve<0.33*pi){
98     fNumCallsRK4++;   
99     fRK4Stepper->Stepper(yInput,dydx,Step,yOut,yErr);
100   
101
102   }
103    else{
104      fNumCallsHelix++;
105      const G4int nvar = 6 ;
106      G4int i;
107      G4double      yTemp[7], yIn[7] ;
108      G4double yTemp2[7];
109      G4ThreeVector  Bfld_midpoint;
110    //  Saving yInput because yInput and yOut can be aliases for same array
111        for(i=0;i<nvar;i++) yIn[i]=yInput[i];
112
113      G4double h = Step * 0.5;
114     // Do two half steps and full step
115          AdvanceHelix(yIn,   Bfld,  h, yTemp,yTemp2);
116          MagFieldEvaluate(yTemp, Bfld_midpoint) ;     
117          AdvanceHelix(yTemp, Bfld_midpoint, h, yOut);
118     // Error estimation
119          for(i=0;i<nvar;i++) {
120          yErr[i] = yOut[i] - yTemp2[i] ;
121         
122          }
123    }
124
125
126
127
128}
129
130void
131G4HelixMixedStepper::DumbStepper( const G4double  yIn[],
132                                   G4ThreeVector   Bfld,
133                                   G4double        h,
134                                   G4double        yOut[])
135{
136 
137   
138       AdvanceHelix(yIn, Bfld, h, yOut);
139
140   
141               
142} 
143
144G4double G4HelixMixedStepper::DistChord()   const 
145{
146  // Implementation : must check whether h/R > 2 pi  !!
147  //   If( h/R <  pi) use G4LineSection::DistLine
148  //   Else           DistChord=R_helix
149  //
150  G4double distChord;
151  G4double Ang_curve=GetAngCurve();
152
153     
154         if(Ang_curve<=pi){
155           distChord=GetRadHelix()*(1-std::cos(0.5*Ang_curve));
156         }
157         else 
158         if(Ang_curve<twopi){
159           distChord=GetRadHelix()*(1+std::cos(0.5*(twopi-Ang_curve)));
160         }
161         else{
162          distChord=2.*GetRadHelix(); 
163         }
164
165   
166
167  return distChord;
168 
169}
170// ---------------------------------------------------------------------------
171void G4HelixMixedStepper::PrintCalls()
172{
173  G4cout<<"In HelixMixedStepper::Number of calls to smallStepStepper = "<<fNumCallsRK4
174        <<"  and Number of calls to Helix = "<<fNumCallsHelix<<G4endl;
175}
176
177
178
179G4MagIntegratorStepper* G4HelixMixedStepper:: SetupStepper(G4Mag_EqRhs* pE, G4int StepperNumber)
180{
181  G4MagIntegratorStepper* pStepper;
182  if (fVerbose>0)G4cout<<"In G4HelixMixedStepper Stepper for small steps is "; 
183  switch ( StepperNumber )
184    {
185     case 0: pStepper = new G4ExplicitEuler( pE ); if (fVerbose>0)G4cout<<"G4ExplicitEuler"<<G4endl; break;
186     case 1: pStepper = new G4ImplicitEuler( pE ); if (fVerbose>0)G4cout<<"G4ImplicitEuler"<<G4endl; break;
187     case 2: pStepper = new G4SimpleRunge( pE ); if (fVerbose>0)G4cout<<"G4SimpleRunge"<<G4endl; break;
188     case 3: pStepper = new G4SimpleHeum( pE );  if (fVerbose>0)G4cout<<"G4SimpleHeum"<<G4endl;break;
189     case 4: pStepper = new G4ClassicalRK4( pE ); if (fVerbose>0)G4cout<<"G4ClassicalRK4"<<G4endl; break;
190     case 5: pStepper = new G4HelixExplicitEuler( pE ); if (fVerbose>0)G4cout<<"G4HelixExplicitEuler"<<G4endl; break;
191     case 6: pStepper = new G4HelixImplicitEuler( pE ); if (fVerbose>0)G4cout<<"G4HelixImplicitEuler"<<G4endl; break;
192     case 7: pStepper = new G4HelixSimpleRunge( pE ); if (fVerbose>0)G4cout<<"G4HelixSimpleRunge"<<G4endl; break;
193     case 8: pStepper = new G4CashKarpRKF45( pE );    if (fVerbose>0)G4cout<<"G4CashKarpRKF45"<<G4endl; break;
194     case 9: pStepper = new G4ExactHelixStepper( pE );  if (fVerbose>0)G4cout<<"G4ExactHelixStepper"<<G4endl;   break;
195     case 10: pStepper = new G4RKG3_Stepper( pE );  if (fVerbose>0)G4cout<<"G4RKG3_Stepper"<<G4endl;   break;
196     
197      default: pStepper = new G4ClassicalRK4( pE );G4cout<<"Default G4ClassicalRK4"<<G4endl; break;
198     
199    }
200  return pStepper;
201}
Note: See TracBrowser for help on using the repository browser.