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

G4TCachedMagneticField is a templated implementation for a specialisation of G4MagneticField used to cache the Magnetic Field value. More...

#include <G4TCachedMagneticField.hh>

Inheritance diagram for G4TCachedMagneticField< T_Field >:

Public Member Functions

 G4TCachedMagneticField (T_Field *pTField, G4double distance)
virtual ~G4TCachedMagneticField ()=default
 G4TCachedMagneticField (const G4TCachedMagneticField< T_Field > &rightCMF)
G4TCachedMagneticFieldoperator= (const G4TCachedMagneticField &right)
G4TCachedMagneticFieldClone () const
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 ()
Public Member Functions inherited from G4MagneticField
 G4MagneticField ()
 ~G4MagneticField () override=default
 G4MagneticField (const G4MagneticField &r)
G4MagneticFieldoperator= (const G4MagneticField &p)
G4bool DoesFieldChangeEnergy () const override
G4FieldType GetFieldType () const override
Public Member Functions inherited from G4Field
 G4Field (G4bool gravityOn=false)
virtual ~G4Field ()=default
 G4Field (const G4Field &p)=default
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 >

G4TCachedMagneticField is a templated implementation for a specialisation of G4MagneticField used to cache the Magnetic Field value.

Definition at line 44 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 48 of file G4TCachedMagneticField.hh.

50 , fLastLocation(DBL_MAX, DBL_MAX, DBL_MAX)
51 , fLastValue(DBL_MAX, DBL_MAX, DBL_MAX)
52 , fCountCalls(0)
54 {
55 fpMagneticField = pTField;
56 fDistanceConst = distance;
57
58 this->ClearCounts();
59 }
G4TCachedMagneticField is a templated implementation for a specialisation of G4MagneticField used to ...

Referenced by Clone(), G4TCachedMagneticField(), and operator=().

◆ ~G4TCachedMagneticField()

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

◆ G4TCachedMagneticField() [2/2]

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

Definition at line 63 of file G4TCachedMagneticField.hh.

64 {
65 fpMagneticField = rightCMF.fpMagneticField;
66 fDistanceConst = rightCMF.fDistanceConst;
67 fLastLocation = rightCMF.fLastLocation;
68 fLastValue = rightCMF.fLastValue;
69 this->ClearCounts();
70 }

Member Function Documentation

◆ ClearCounts()

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

Definition at line 132 of file G4TCachedMagneticField.hh.

133 {
134 fCountCalls = 0;
136 }

Referenced by G4TCachedMagneticField(), and G4TCachedMagneticField().

◆ Clone()

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

Interface method to implement cloning, needed by multi-threading. Here issuing a fatal exception, as expecting this to be implemented concretely in derived classes.

Reimplemented from G4Field.

Definition at line 88 of file G4TCachedMagneticField.hh.

89 {
90 G4cout << "Clone is called" << G4endl;
91 // Cannot use copy constructor: I need to clone associated magnetic field
92 T_Field* aF = this->fpMagneticField->T_Field::Clone();
93 G4TCachedMagneticField* cloned = new G4TCachedMagneticField(aF, this->fDistanceConst);
94 cloned->fLastLocation = this->fLastLocation;
95 cloned->fLastValue = this->fLastValue;
96 return cloned;
97 }
G4TCachedMagneticField(T_Field *pTField, G4double distance)

◆ GetConstDistance()

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

Definition at line 127 of file G4TCachedMagneticField.hh.

127{ return fDistanceConst; }

◆ GetCountCalls()

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

Definition at line 130 of file G4TCachedMagneticField.hh.

130{ return fCountCalls; }

◆ GetCountEvaluations()

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

Definition at line 131 of file G4TCachedMagneticField.hh.

131{ return fCountEvaluations; }

◆ GetFieldValue()

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

Given the position time vector 'Point', returns the value of the field in the array 'Bfield'.

Parameters
[in]PointThe position time vector.
[out]BfieldThe field array in output.

Implements G4MagneticField.

Definition at line 106 of file G4TCachedMagneticField.hh.

107 {
109
110 G4double distSq = (newLocation - fLastLocation).mag2();
111 fCountCalls++;
112 if(distSq < fDistanceConst * fDistanceConst)
113 {
114 Bfield[0] = fLastValue.x();
115 Bfield[1] = fLastValue.y();
116 Bfield[2] = fLastValue.z();
117 }
118 else
119 {
120 fpMagneticField->T_Field::GetFieldValue(Point, Bfield);
122 fLastLocation = G4ThreeVector(Point[0], Point[1], Point[2]);
123 fLastValue = G4ThreeVector(Bfield[0], Bfield[1], Bfield[2]);
124 }
125 }
CLHEP::Hep3Vector G4ThreeVector

◆ operator=()

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

Definition at line 72 of file G4TCachedMagneticField.hh.

73 {
74 if(&right == this) { return *this; }
75
76 fpMagneticField = right.fpMagneticField;
77
78 fDistanceConst= right.fDistanceConst;
79 fLastLocation = right.fLastLocation;
80 fLastValue = right.fLastValue;
81
82 fCountCalls = 0;
84
85 return *this;
86 }

◆ ReportStatistics()

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

Definition at line 99 of file G4TCachedMagneticField.hh.

100 {
101 G4cout << " Cached field: " << G4endl
102 << " Number of calls: " << fCountCalls << G4endl
103 << " Number of evaluations : " << fCountEvaluations << G4endl;
104 }

◆ SetConstDistance()

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

Definition at line 128 of file G4TCachedMagneticField.hh.

128{ fDistanceConst = dist; }

Member Data Documentation

◆ fCountCalls

template<class T_Field>
G4int G4TCachedMagneticField< T_Field >::fCountCalls
mutableprotected

◆ fCountEvaluations

template<class T_Field>
G4int G4TCachedMagneticField< T_Field >::fCountEvaluations
protected

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