Geant4-11
G4TCachedMagneticField.hh
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#ifndef TCACHED_MAGNETIC_FIELD_DEF
28#define TCACHED_MAGNETIC_FIELD_DEF
29
30#include "G4Types.hh"
31#include "G4ThreeVector.hh"
32#include "G4MagneticField.hh"
33
34template <class T_Field>
36{
37 public: // with description
38 G4TCachedMagneticField(T_Field* pTField, G4double distance)
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 }
52
54 {
58 fLastValue = rightCMF.fLastValue;
59 this->ClearCounts();
60 }
61
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();
70 cloned->fLastLocation = this->fLastLocation;
71 cloned->fLastValue = this->fLastValue;
72 return cloned;
73 }
74
76 // Constructor and destructor. No actions.
77
79 {
80 G4cout << " Cached field: " << G4endl
81 << " Number of calls: " << fCountCalls << G4endl
82 << " Number of evaluations : " << fCountEvaluations << G4endl;
83 }
84
85 virtual void GetFieldValue(const G4double Point[4], G4double* Bfield) const
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 }
110
113
114 G4int GetCountCalls() const { return fCountCalls; }
117 {
118 fCountCalls = 0;
120 }
121
123 {
124 if(&right == this)
125 return *this;
126
128
131 fLastValue = right.fLastValue;
132
133 fCountCalls = 0;
135
136 return *this;
137 }
138
139 private:
141 // When the field is evaluated within this distance it will not change
143 // Caching state
146
147 protected:
149};
150
151#endif /* TCACHED_MAGNETIC_FIELD_DEF */
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double y() const
G4TCachedMagneticField(T_Field *pTField, G4double distance)
G4TCachedMagneticField(const G4TCachedMagneticField< T_Field > &rightCMF)
G4TCachedMagneticField * Clone() const
G4TCachedMagneticField & operator=(const G4TCachedMagneticField &right)
void SetConstDistance(G4double dist)
virtual void GetFieldValue(const G4double Point[4], G4double *Bfield) const
#define DBL_MAX
Definition: templates.hh:62