Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TCachedMagneticField< T_Field > Class Template Reference

#include <G4TCachedMagneticField.hh>

+ Inheritance diagram for G4TCachedMagneticField< T_Field >:

Public Member Functions

 G4TCachedMagneticField (T_Field *pTField, G4double distance)
 
 G4TCachedMagneticField (const G4TCachedMagneticField< T_Field > &rightCMF)
 
G4TCachedMagneticFieldClone () const
 
virtual ~G4TCachedMagneticField ()
 
void ReportStatistics ()
 
virtual void GetFieldValue (const G4double Point[4], G4double *Bfield) const
 
G4double GetConstDistance () const
 
void SetConstDistance (G4double dist)
 
G4int GetCountCalls () const
 
G4int GetCountEvaluations () const
 
void ClearCounts ()
 
G4TCachedMagneticFieldoperator= (const G4TCachedMagneticField &right)
 
- Public Member Functions inherited from G4MagneticField
 G4MagneticField ()
 
 ~G4MagneticField () override
 
 G4MagneticField (const G4MagneticField &r)
 
G4MagneticFieldoperator= (const G4MagneticField &p)
 
G4bool DoesFieldChangeEnergy () const override
 
- Public Member Functions inherited from G4Field
 G4Field (G4bool gravityOn=false)
 
 G4Field (const G4Field &)
 
virtual ~G4Field ()
 
G4Fieldoperator= (const G4Field &p)
 
G4bool IsGravityActive () const
 
void SetGravityActive (G4bool OnOffFlag)
 

Protected Attributes

G4int fCountCalls
 
G4int fCountEvaluations
 

Additional Inherited Members

- Static Public Attributes inherited from G4Field
static constexpr G4int MAX_NUMBER_OF_COMPONENTS = 24
 

Detailed Description

template<class T_Field>
class G4TCachedMagneticField< T_Field >

Definition at line 35 of file G4TCachedMagneticField.hh.

Constructor & Destructor Documentation

◆ G4TCachedMagneticField() [1/2]

template<class T_Field >
G4TCachedMagneticField< T_Field >::G4TCachedMagneticField ( T_Field * pTField,
G4double distance )
inline

Definition at line 38 of file G4TCachedMagneticField.hh.

40 , fLastLocation(DBL_MAX, DBL_MAX, DBL_MAX)
41 , fLastValue(DBL_MAX, DBL_MAX, DBL_MAX)
42 , fCountCalls(0)
44 {
45 fpMagneticField = pTField;
46 fDistanceConst = distance;
47
48 // G4cout << " Cached-B-Field constructor> Distance = " << distance <<
49 // G4endl;
50 this->ClearCounts();
51 }
#define DBL_MAX
Definition templates.hh:62

Referenced by G4TCachedMagneticField< T_Field >::Clone().

◆ G4TCachedMagneticField() [2/2]

template<class T_Field >
G4TCachedMagneticField< T_Field >::G4TCachedMagneticField ( const G4TCachedMagneticField< T_Field > & rightCMF)
inline

Definition at line 53 of file G4TCachedMagneticField.hh.

54 {
55 fpMagneticField = rightCMF.fpMagneticField;
56 fDistanceConst = rightCMF.fDistanceConst;
57 fLastLocation = rightCMF.fLastLocation;
58 fLastValue = rightCMF.fLastValue;
59 this->ClearCounts();
60 }

◆ ~G4TCachedMagneticField()

template<class T_Field >
virtual G4TCachedMagneticField< T_Field >::~G4TCachedMagneticField ( )
inlinevirtual

Definition at line 75 of file G4TCachedMagneticField.hh.

75{ ; }

Member Function Documentation

◆ ClearCounts()

template<class T_Field >
void G4TCachedMagneticField< T_Field >::ClearCounts ( )
inline

◆ Clone()

template<class T_Field >
G4TCachedMagneticField * G4TCachedMagneticField< T_Field >::Clone ( ) const
inlinevirtual

Reimplemented from G4Field.

Definition at line 62 of file G4TCachedMagneticField.hh.

63 {
64 G4cout << "Clone is called" << G4endl;
65 // Cannot use copy constructor: I need to clone the associated magnetif
66 // field
67 T_Field* aF = this->fpMagneticField->T_Field::Clone();
69 new G4TCachedMagneticField(aF, this->fDistanceConst);
70 cloned->fLastLocation = this->fLastLocation;
71 cloned->fLastValue = this->fLastValue;
72 return cloned;
73 }
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4TCachedMagneticField(T_Field *pTField, G4double distance)

◆ GetConstDistance()

template<class T_Field >
G4double G4TCachedMagneticField< T_Field >::GetConstDistance ( ) const
inline

Definition at line 111 of file G4TCachedMagneticField.hh.

111{ return fDistanceConst; }

◆ GetCountCalls()

template<class T_Field >
G4int G4TCachedMagneticField< T_Field >::GetCountCalls ( ) const
inline

Definition at line 114 of file G4TCachedMagneticField.hh.

114{ return fCountCalls; }

◆ GetCountEvaluations()

template<class T_Field >
G4int G4TCachedMagneticField< T_Field >::GetCountEvaluations ( ) const
inline

Definition at line 115 of file G4TCachedMagneticField.hh.

115{ return fCountEvaluations; }

◆ GetFieldValue()

template<class T_Field >
virtual void G4TCachedMagneticField< T_Field >::GetFieldValue ( const G4double Point[4],
G4double * Bfield ) const
inlinevirtual

Implements G4MagneticField.

Definition at line 85 of file G4TCachedMagneticField.hh.

86 {
87 G4ThreeVector newLocation(Point[0], Point[1], Point[2]);
88
89 G4double distSq = (newLocation - fLastLocation).mag2();
91 if(distSq < fDistanceConst * fDistanceConst)
92 {
93 Bfield[0] = fLastValue.x();
94 Bfield[1] = fLastValue.y();
95 Bfield[2] = fLastValue.z();
96 }
97 else
98 {
99 // G4CachedMagneticField* thisNonC=
100 // const_cast<G4CachedMagneticField*>(this);
101 fpMagneticField->T_Field::GetFieldValue(Point, Bfield);
102 // G4cout << " Evaluating. " << G4endl;
104 // thisNonC->
105 fLastLocation = G4ThreeVector(Point[0], Point[1], Point[2]);
106 // thisNonC->
107 fLastValue = G4ThreeVector(Bfield[0], Bfield[1], Bfield[2]);
108 }
109 }
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
double z() const
double x() const
double y() const

◆ operator=()

template<class T_Field >
G4TCachedMagneticField & G4TCachedMagneticField< T_Field >::operator= ( const G4TCachedMagneticField< T_Field > & right)
inline

Definition at line 122 of file G4TCachedMagneticField.hh.

123 {
124 if(&right == this)
125 return *this;
126
127 fpMagneticField = right.fpMagneticField;
128
129 fDistanceConst= right.fDistanceConst;
130 fLastLocation = right.fLastLocation;
131 fLastValue = right.fLastValue;
132
133 fCountCalls = 0;
135
136 return *this;
137 }

◆ ReportStatistics()

template<class T_Field >
void G4TCachedMagneticField< T_Field >::ReportStatistics ( )
inline

Definition at line 78 of file G4TCachedMagneticField.hh.

79 {
80 G4cout << " Cached field: " << G4endl
81 << " Number of calls: " << fCountCalls << G4endl
82 << " Number of evaluations : " << fCountEvaluations << G4endl;
83 }

◆ SetConstDistance()

template<class T_Field >
void G4TCachedMagneticField< T_Field >::SetConstDistance ( G4double dist)
inline

Definition at line 112 of file G4TCachedMagneticField.hh.

112{ fDistanceConst = dist; }

Member Data Documentation

◆ fCountCalls

◆ fCountEvaluations


The documentation for this class was generated from the following file: