Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ErrorMagFieldLimitProcess.cc
Go to the documentation of this file.
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// GEANT 4 class implementation file
29// ------------------------------------------------------------
30//
31
33#include "G4ErrorMessenger.hh"
35#include "G4FieldManager.hh"
36#include "G4Field.hh"
37#include "G4Track.hh"
38
39#ifdef G4VERBOSE
41#endif
42
43//------------------------------------------------------------------------
45G4ErrorMagFieldLimitProcess(const G4String& processName)
46 : G4VErrorLimitProcess(processName)
47{
48 theStepLimit = kInfinity;
49}
50
51
52//------------------------------------------------------------------------
54{ }
55
56
57//------------------------------------------------------------------------
61{
63 const G4Field* field =
66
67 theStepLength = kInfinity;
68 if( field != 0 ) {
69 G4ThreeVector trkPosi = aTrack.GetPosition();
70 G4double pos1[3];
71 pos1[0] = trkPosi.x(); pos1[1] = trkPosi.y(); pos1[2] = trkPosi.z();
72
73 G4double h1[3];
74 field->GetFieldValue( pos1, h1 );
75
76 G4ThreeVector BVec(h1[0],h1[1],h1[2]);
77 G4double pmag = aTrack.GetMomentum().mag();
78 G4double BPerpMom = BVec.cross( G4ThreeVector(pmag,0.,0.) ).mag() / pmag;
79
80 theStepLength = theStepLimit * pmag / BPerpMom;
81#ifdef G4VERBOSE
83 G4cout << "G4ErrorMagFieldLimitProcess:: stepLength "
84 << theStepLength << " B " << BPerpMom << " BVec " << BVec
85 << " pmag " << pmag << G4endl;
86 }
87#endif
88 }
89
90 return theStepLength;
91}
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ NotForced
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double y() const
Hep3Vector cross(const Hep3Vector &) const
double mag() const
G4ErrorMagFieldLimitProcess(const G4String &processName="G4ErrorMagFieldLimit")
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
const G4Field * GetDetectorField() const
virtual void GetFieldValue(const G4double Point[4], G4double *fieldArr) const =0
const G4ThreeVector & GetPosition() const
G4ThreeVector GetMomentum() const
static G4TransportationManager * GetTransportationManager()
G4FieldManager * GetFieldManager() const