G4HelixExplicitEuler.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 //
00027 // $Id: G4HelixExplicitEuler.cc 69786 2013-05-15 09:38:51Z gcosmo $
00028 //
00029 //
00030 //  Helix Explicit Euler: x_1 = x_0 + helix(h)
00031 //  with helix(h) being a helix piece of length h
00032 //  most simple approach for solving linear differential equations.
00033 //  Take the current derivative and add it to the current position.
00034 //
00035 //  W.Wander <wwc@mit.edu> 12/09/97 
00036 // -------------------------------------------------------------------
00037 
00038 #include "G4HelixExplicitEuler.hh"
00039 #include "G4PhysicalConstants.hh"
00040 #include "G4ThreeVector.hh"
00041 
00042 
00043 void G4HelixExplicitEuler::Stepper(  const G4double  yInput[7],
00044                                const G4double*,
00045                                      G4double Step,
00046                                      G4double yOut[7],
00047                                      G4double yErr[])
00048 
00049 {
00050 
00051  //Estimation of the Stepping Angle
00052 
00053   G4ThreeVector Bfld;
00054   MagFieldEvaluate(yInput, Bfld); 
00055   
00056   const G4int nvar = 6 ;
00057   G4int i;
00058   G4double      yTemp[7], yIn[7] ;
00059   G4ThreeVector  Bfld_midpoint;
00060   //  Saving yInput because yInput and yOut can be aliases for same array
00061         for(i=0;i<nvar;i++) yIn[i]=yInput[i];
00062      
00063         G4double h = Step * 0.5;
00064  
00065      // Do full step and two half steps
00066         G4double yTemp2[7];
00067         AdvanceHelix(yIn,   Bfld,  h, yTemp2,yTemp);
00068         MagFieldEvaluate(yTemp2, Bfld_midpoint) ;     
00069         AdvanceHelix(yTemp2, Bfld_midpoint, h, yOut);
00070     
00071      // Error estimation
00072         for(i=0;i<nvar;i++) {
00073          yErr[i] = yOut[i] - yTemp[i] ;
00074        }
00075     
00076 }
00077 
00078 G4double G4HelixExplicitEuler::DistChord()   const 
00079 {
00080   // Implementation : must check whether h/R > 2 pi  !!
00081   //   If( h/R <  pi) use G4LineSection::DistLine
00082   //   Else           DistChord=R_helix
00083   //
00084   G4double distChord;
00085   G4double Ang_curve=GetAngCurve();
00086 
00087       
00088          if(Ang_curve<=pi){
00089            distChord=GetRadHelix()*(1-std::cos(0.5*Ang_curve));
00090          }
00091          else 
00092          if(Ang_curve<twopi){
00093            distChord=GetRadHelix()*(1+std::cos(0.5*(twopi-Ang_curve)));
00094          }
00095          else{
00096           distChord=2.*GetRadHelix();  
00097          }
00098 
00099   return distChord;
00100   
00101 }
00102 void
00103 G4HelixExplicitEuler::DumbStepper( const G4double  yIn[],
00104                                    G4ThreeVector   Bfld,
00105                                    G4double        h,
00106                                    G4double        yOut[])
00107 {
00108     
00109        AdvanceHelix(yIn, Bfld, h, yOut);
00110                
00111 }  

Generated on Mon May 27 17:48:29 2013 for Geant4 by  doxygen 1.4.7