G4ErrorMagFieldLimitProcess Class Reference

#include <G4ErrorMagFieldLimitProcess.hh>

Inheritance diagram for G4ErrorMagFieldLimitProcess:

G4VErrorLimitProcess G4VDiscreteProcess G4VProcess

Public Member Functions

 G4ErrorMagFieldLimitProcess (const G4String &processName="G4ErrorMagFieldLimit")
 ~G4ErrorMagFieldLimitProcess ()
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)

Detailed Description

Definition at line 54 of file G4ErrorMagFieldLimitProcess.hh.


Constructor & Destructor Documentation

G4ErrorMagFieldLimitProcess::G4ErrorMagFieldLimitProcess ( const G4String processName = "G4ErrorMagFieldLimit"  ) 

Definition at line 46 of file G4ErrorMagFieldLimitProcess.cc.

References G4VErrorLimitProcess::theStepLimit.

00047   : G4VErrorLimitProcess(processName) 
00048 {
00049   theStepLimit = kInfinity;
00050 }

G4ErrorMagFieldLimitProcess::~G4ErrorMagFieldLimitProcess (  ) 

Definition at line 54 of file G4ErrorMagFieldLimitProcess.cc.

00055 { }


Member Function Documentation

G4double G4ErrorMagFieldLimitProcess::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
) [virtual]

Implements G4VErrorLimitProcess.

Definition at line 60 of file G4ErrorMagFieldLimitProcess.cc.

References G4cout, G4endl, G4FieldManager::GetDetectorField(), G4TransportationManager::GetFieldManager(), G4Field::GetFieldValue(), G4Track::GetMomentum(), G4Track::GetPosition(), G4TransportationManager::GetTransportationManager(), NotForced, G4VErrorLimitProcess::theStepLength, G4VErrorLimitProcess::theStepLimit, and G4ErrorPropagatorData::verbose().

00062 {
00063   *condition = NotForced;
00064   const G4Field* field =
00065     G4TransportationManager::GetTransportationManager()->GetFieldManager()
00066     ->GetDetectorField();
00067 
00068   theStepLength = kInfinity;
00069   if( field != 0 ) {
00070     G4ThreeVector trkPosi = aTrack.GetPosition();
00071     G4double pos1[3];
00072        pos1[0] = trkPosi.x(); pos1[1] = trkPosi.y(); pos1[2] = trkPosi.z();
00073 
00074     G4double h1[3];
00075     field->GetFieldValue( pos1, h1 );
00076 
00077     G4ThreeVector BVec(h1[0],h1[1],h1[2]);
00078     G4double pmag = aTrack.GetMomentum().mag();
00079     G4double BPerpMom = BVec.cross( G4ThreeVector(pmag,0.,0.) ).mag() / pmag;
00080 
00081     theStepLength = theStepLimit * pmag / BPerpMom; 
00082 #ifdef G4VERBOSE
00083   if(G4ErrorPropagatorData::verbose() >= 3 ) { 
00084     G4cout <<  "G4ErrorMagFieldLimitProcess:: stepLength "
00085            << theStepLength << " B " << BPerpMom << " BVec " << BVec
00086            << " pmag " << pmag << G4endl;
00087   }
00088 #endif
00089   }
00090 
00091   return theStepLength;
00092 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:51:55 2013 for Geant4 by  doxygen 1.4.7