Geant4  10.03
G4FieldManager.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 // $Id: G4FieldManager.cc 93661 2015-10-28 09:47:47Z gcosmo $
28 //
29 // -------------------------------------------------------------------
30 
31 #include "G4FieldManager.hh"
32 #include "G4Field.hh"
33 #include "G4MagneticField.hh"
34 #include "G4ChordFinder.hh"
35 #include "G4FieldManagerStore.hh"
36 
38  G4ChordFinder *pChordFinder,
39  G4bool fieldChangesEnergy
40  )
41  : fDetectorField(detectorField),
42  fChordFinder(pChordFinder),
43  fAllocatedChordFinder(false),
44  fEpsilonMinDefault(5.0e-5),
45  fEpsilonMaxDefault(0.001),
46  fDefault_Delta_One_Step_Value(0.01), // mm
47  fDefault_Delta_Intersection_Val(0.001), // mm
48  fEpsilonMin( fEpsilonMinDefault ),
49  fEpsilonMax( fEpsilonMaxDefault)
50 {
53  if ( detectorField )
55  else
56  fFieldChangesEnergy= fieldChangesEnergy;
57 
58  // Add to store
60 }
61 
63  : fDetectorField(detectorField), fAllocatedChordFinder(true),
64  fEpsilonMinDefault(5.0e-5),
65  fEpsilonMaxDefault(0.001),
66  fFieldChangesEnergy(false),
67  fDefault_Delta_One_Step_Value(0.01), // mm
68  fDefault_Delta_Intersection_Val(0.001), // mm
69  fEpsilonMin( fEpsilonMinDefault ),
70  fEpsilonMax( fEpsilonMaxDefault)
71 {
72  fChordFinder= new G4ChordFinder( detectorField );
75  // Add to store
77 }
78 
80 {
81  G4Field* aField = 0;
82  G4FieldManager* aFM = 0;
83  G4ChordFinder* aCF = 0;
84  try {
85  if ( this->fDetectorField )
86  aField = this->fDetectorField->Clone();
87 
88  //Create a new field manager, note that we do not set any chordfinder now.
89  aFM = new G4FieldManager( aField , 0 , this->fFieldChangesEnergy );
90 
91  //Check if orignally we have the fAllocatedChordFinder variable set, in case, call chord
92  //constructor
93  if ( this->fAllocatedChordFinder )
94  {
95  aFM->CreateChordFinder( dynamic_cast<G4MagneticField*>(aField) );
96  }
97  else
98  {
99  //Chord was specified by user, should we clone?
100  //TODO: For the moment copy pointer, to be understood if cloning of ChordFinder is needed
101  aCF = this->fChordFinder;/*->Clone*/
102  aFM->fChordFinder = aCF;
103  }
104  //Copy values of other variables
105  aFM->fEpsilonMax = this->fEpsilonMax;
106  aFM->fEpsilonMin = this->fEpsilonMin;
111  //TODO: Should we really add to the store the cloned FM? Who will use this?
112  }
113  catch ( ... )
114  {
115  //Failed creating clone: probably user did not implement Clone method
116  //in derived classes?
117  //Perform clean-up after ourselves...
118  delete aField;
119  delete aFM;
120  delete aCF;
121  throw;
122  }
123  return aFM;
124 }
125 
127 {
128  // Default is to do nothing!
129  ;
130 }
131 
133 {
134  if( fAllocatedChordFinder ){
135  delete fChordFinder;
136  }
138 }
139 
140 void
142 {
143  if ( fAllocatedChordFinder )
144  delete fChordFinder;
145  fChordFinder= new G4ChordFinder( detectorMagField );
146  fAllocatedChordFinder= true;
147 }
148 
150 {
151  fDetectorField= pDetectorField;
152 
153  if ( pDetectorField )
154  fFieldChangesEnergy= pDetectorField->DoesFieldChangeEnergy();
155  else
156  fFieldChangesEnergy= false; // No field
157 
158  return false;
159 }
160 
virtual G4FieldManager * Clone() const
G4bool SetDetectorField(G4Field *detectorField)
G4double fDelta_One_Step_Value
G4double fEpsilonMin
G4double fDefault_Delta_Intersection_Val
virtual G4Field * Clone() const
Definition: G4Field.cc:54
virtual G4bool DoesFieldChangeEnergy() const =0
G4double fDefault_Delta_One_Step_Value
virtual void ConfigureForTrack(const G4Track *)
G4FieldManager(G4Field *detectorField=0, G4ChordFinder *pChordFinder=0, G4bool b=true)
bool G4bool
Definition: G4Types.hh:79
G4double fEpsilonMax
G4ChordFinder * fChordFinder
virtual ~G4FieldManager()
G4bool fFieldChangesEnergy
static void Register(G4FieldManager *pVolume)
G4Field * fDetectorField
void CreateChordFinder(G4MagneticField *detectorMagField)
G4double fDelta_Intersection_Val
G4bool fAllocatedChordFinder
static void DeRegister(G4FieldManager *pVolume)