CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HFGflash.cc
Go to the documentation of this file.
2 
3 #include "G4VPhysicalVolume.hh"
4 #include "G4Step.hh"
5 #include "G4Track.hh"
6 #include "G4Navigator.hh"
7 #include "G4NavigationHistory.hh"
8 #include "CLHEP/Units/PhysicalConstants.h"
9 #include "CLHEP/Units/SystemOfUnits.h"
10 
11 #include "Randomize.hh"
12 #include "G4TransportationManager.hh"
13 #include "G4VPhysicalVolume.hh"
14 #include "G4LogicalVolume.hh"
15 #include "G4VSensitiveDetector.hh"
16 #include "G4EventManager.hh"
17 #include "G4SteppingManager.hh"
18 #include "G4FastTrack.hh"
19 #include "G4ParticleTable.hh"
20 
21 #include "CLHEP/GenericFunctions/IncompleteGamma.hh"
22 
26 
29 
30 #include <math.h>
31 
32 //#define DebugLog
33 
35 
36  edm::ParameterSet m_HF = p.getParameter<edm::ParameterSet>("HFGflash");
37  theBField = m_HF.getUntrackedParameter<double>("BField", 3.8);
38  theWatcherOn = m_HF.getUntrackedParameter<bool>("WatcherOn",true);
39  theFillHisto = m_HF.getUntrackedParameter<bool>("FillHisto",true);
40  edm::LogInfo("HFShower") << "HFGFlash:: Set B-Field to " << theBField
41  << ", WatcherOn to " << theWatcherOn
42  << " and FillHisto to " << theFillHisto;
43 
45  theGflashStep = new G4Step();
46  theGflashNavigator = new G4Navigator();
47  // theGflashNavigator = 0;
48  theGflashTouchableHandle = new G4TouchableHistory();
49 
50 #ifdef DebugLog
51  if (theFillHisto) {
53  if ( tfile.isAvailable() ) {
54  TFileDirectory showerDir = tfile->mkdir("GflashEMShowerProfile");
55 
56  em_incE = showerDir.make<TH1F>("em_incE","Incoming energy (GeV)",500,0,500.);
57  em_ssp_rho = showerDir.make<TH1F>("em_ssp_rho","Shower starting position;#rho (cm);Number of Events",100,100.0,200.0);
58  em_ssp_z = showerDir.make<TH1F>("em_ssp_z","Shower starting position;z (cm);Number of Events",2000,0.0,2000.0);
59  em_long = showerDir.make<TH1F>("em_long","Longitudinal Profile;Radiation Length;Number of Spots",800,800.0,1600.0);
60  em_lateral = showerDir.make<TH1F>("em_lateral","Lateral Profile;Radiation Length;Moliere Radius",100,0.0,5.0);
61  em_2d = showerDir.make<TH2F>("em_2d","Lateral Profile vs. Shower Depth;Radiation Length;Moliere Radius",800,800.0,1600.0,100,0.0,5.0);
62  em_long_sd = showerDir.make<TH1F>("em_long_sd","Longitudinal Profile in Sensitive Detector;Radiation Length;Number of Spots",800,800.0,1600.0);
63  em_lateral_sd = showerDir.make<TH1F>("em_lateral_sd","Lateral Profile vs. Shower Depth in Sensitive Detector;Radiation Length;Moliere Radius",100,0.0,5.0);
64  em_2d_sd = showerDir.make<TH2F>("em_2d_sd","Lateral Profile vs. Shower Depth in Sensitive Detector;Radiation Length;Moliere Radius",800,800.0,1600.0,100,0.0,5.0);
65  em_ze = showerDir.make<TH2F>("em_ze","Profile vs. Energy;Radiation Length;Moliere Radius",800,800.0,1600.0,1000,0.0,1.0);
66  em_ratio = showerDir.make<TH2F>("em_ratio","Profile vs. Energy;Radiation Length;Moliere Radius",800,800.0,1600.0,1000,0.0,100.0);
67  em_ratio_selected = showerDir.make<TH2F>("em_ratio_selected","Profile vs. Energy;Radiation Length;Moliere Radius",800,800.0,1600.0,1000,0.0,100.0);
68  em_nSpots_sd = showerDir.make<TH1F>("em_nSpots_sd","Number of Gflash Spots in Sensitive Detector;Number of Spots;Number of Events",100,0.0,100);
69  em_ze_ratio = showerDir.make<TH1F>("em_ze_ratio","Ratio of Energy and Z Position",1000,0.0,0.001);
70  } else {
71  theFillHisto = false;
72  edm::LogInfo("HFShower") << "HFGFlash::No file is available for saving"
73  << " histos so the flag is set to false";
74  }
75  }
76 #endif
78 
79 }
80 
82  if (theHelix) delete theHelix;
83  if (theGflashStep) delete theGflashStep;
85 }
86 
87 
88 std::vector<HFGflash::Hit> HFGflash::gfParameterization(G4Step * aStep,bool & ok,bool onlyLong) {
89  double tempZCalo = 26; // Gflash::Z[jCalorimeter]
90  double hfcriticalEnergy = 0.021; // Gflash::criticalEnergy
91 
92  std::vector<HFGflash::Hit> hit;
93  HFGflash::Hit oneHit;
94 
95  G4StepPoint * preStepPoint = aStep->GetPreStepPoint();
96  //G4StepPoint * postStepPoint = aStep->GetPostStepPoint();
97  G4Track * track = aStep->GetTrack();
98  // Get Z-direction
99  const G4DynamicParticle *aParticle = track->GetDynamicParticle();
100  G4ThreeVector momDir = aParticle->GetMomentumDirection();
101 
102  G4ThreeVector hitPoint = preStepPoint->GetPosition();
103  G4String partType = track->GetDefinition()->GetParticleName();
104  // int parCode = track->GetDefinition()->GetPDGEncoding();
105 
106  // This part of code is copied from the original GFlash Fortran code.
107  // reference : hep-ex/0001020v1
108 
109  const G4double energyCutoff = 1;
110  const G4int maxNumberOfSpots = 10000000;
111 
112  G4ThreeVector showerStartingPosition = track->GetPosition()/cm;
113  G4ThreeVector showerMomentum = preStepPoint->GetMomentum()/GeV;
115 
116  G4double logEinc = std::log((preStepPoint->GetTotalEnergy())/GeV);
117 
118  G4double y = ((preStepPoint->GetTotalEnergy())/GeV) / hfcriticalEnergy; // y = E/Ec, criticalEnergy is in GeV
119  G4double logY = std::log(y);
120 
121 
122  G4double nSpots = 93.0 * std::log(tempZCalo) * ((preStepPoint->GetTotalEnergy())/GeV); // total number of spot due linearization
123  if(preStepPoint->GetTotalEnergy()/GeV < 1.6) nSpots = 140.4 * std::log(tempZCalo) * std::pow(((preStepPoint->GetTotalEnergy())/GeV),0.876);
124 
125 
126  // // implementing magnetic field effects
127  double charge = track->GetStep()->GetPreStepPoint()->GetCharge();
128  theHelix->initializeTrajectory(showerMomentum,showerStartingPosition,charge,theBField);
129 
130  G4double pathLength0 = theHelix->getPathLengthAtZ(showerStartingPosition.getZ());
131  G4double pathLength = pathLength0; // this will grow along the shower development
132 
133  //--- 2.2 Fix intrinsic properties of em. showers.
134 
135  G4double fluctuatedTmax = std::log(logY - 0.7157);
136  G4double fluctuatedAlpha= std::log(0.7996 +(0.4581 + 1.8628/tempZCalo)*logY);
137 
138  G4double sigmaTmax = 1.0/( -1.4 + 1.26 * logY);
139  G4double sigmaAlpha = 1.0/( -0.58 + 0.86 * logY);
140  G4double rho = 0.705 - 0.023 * logY;
141  G4double sqrtPL = std::sqrt((1.0+rho)/2.0);
142  G4double sqrtLE = std::sqrt((1.0-rho)/2.0);
143 
144  G4double norm1 = G4RandGauss::shoot();
145  G4double norm2 = G4RandGauss::shoot();
146  G4double tempTmax = fluctuatedTmax + sigmaTmax*(sqrtPL*norm1 + sqrtLE*norm2);
147  G4double tempAlpha = fluctuatedAlpha + sigmaAlpha*(sqrtPL*norm1 - sqrtLE*norm2);
148 
149  // tmax, alpha, beta : parameters of gamma distribution
150  G4double tmax = std::exp(tempTmax);
151  G4double alpha = std::exp(tempAlpha);
152  G4double beta = (alpha - 1.0)/tmax;
153 
154  if (!alpha) return hit;
155  if (!beta) return hit;
156  if (alpha < 0.00001) return hit;
157  if (beta < 0.00001) return hit;
158 
159  // spot fluctuations are added to tmax, alpha, beta
160  G4double averageTmax = logY-0.858;
161  G4double averageAlpha = 0.21+(0.492+2.38/tempZCalo)*logY;
162  G4double spotTmax = averageTmax * (0.698 + .00212*tempZCalo);
163  G4double spotAlpha= averageAlpha * (0.639 + .00334*tempZCalo);
164  G4double spotBeta = (spotAlpha-1.0)/spotTmax;
165 
166  if (!spotAlpha) return hit;
167  if (!spotBeta) return hit;
168  if (spotAlpha < 0.00001) return hit;
169  if (spotBeta < 0.00001) return hit;
170 
171 #ifdef DebugLog
172  LogDebug("HFShower") << "Incoming energy = " << ((preStepPoint->GetTotalEnergy())/GeV) << " Position (rho,z) = (" << showerStartingPosition.rho() << ", " << showerStartingPosition.z() << ")";
173 
174  if(theFillHisto) {
175  em_incE->Fill(((preStepPoint->GetTotalEnergy())/GeV));
176  em_ssp_rho->Fill(showerStartingPosition.rho());
177  em_ssp_z->Fill(std::abs(showerStartingPosition.z()));
178  }
179 #endif
180  // parameters for lateral distribution and fluctuation
181  G4double z1=0.0251+0.00319*logEinc;
182  G4double z2=0.1162-0.000381*tempZCalo;
183 
184  G4double k1=0.659 - 0.00309 * tempZCalo;
185  G4double k2=0.645;
186  G4double k3=-2.59;
187  G4double k4=0.3585+ 0.0421*logEinc;
188 
189  G4double p1=2.623 -0.00094*tempZCalo;
190  G4double p2=0.401 +0.00187*tempZCalo;
191  G4double p3=1.313 -0.0686*logEinc;
192 
193  // // @@@ dwjang, intial tuning by comparing 20-150GeV TB data
194  // // the width of energy response is not yet tuned.
195  G4double e25Scale = 1.03551;
196  z1 *= 9.76972e-01 - 3.85026e-01 * std::tanh(1.82790e+00*std::log(((preStepPoint->GetTotalEnergy())/GeV)) - 3.66237e+00);
197  p1 *= 0.96;
198 
199  G4double stepLengthLeft = 10000;
200  G4int nSpots_sd = 0; // count total number of spots in SD
201  G4double zInX0 = 0.0; // shower depth in X0 unit
202  G4double deltaZInX0 = 0.0; // segment of depth in X0 unit
203  G4double deltaZ = 0.0; // segment of depth in cm
204  G4double stepLengthLeftInX0 = 0.0; // step length left in X0 unit
205 
206  const G4double divisionStepInX0 = 0.1; //step size in X0 unit
207  G4double energy = ((preStepPoint->GetTotalEnergy())/GeV); // energy left in GeV
208 
209  Genfun::IncompleteGamma gammaDist;
210 
211  G4double energyInGamma = 0.0; // energy in a specific depth(z) according to Gamma distribution
212  G4double preEnergyInGamma = 0.0; // energy calculated in a previous depth
213  G4double sigmaInGamma = 0.; // sigma of energy in a specific depth(z) according to Gamma distribution
214  G4double preSigmaInGamma = 0.0; // sigma of energy in a previous depth
215 
216  //energy segment in Gamma distribution of shower in each step
217  G4double deltaEnergy =0.0 ; // energy in deltaZ
218  G4int spotCounter = 0; // keep track of number of spots generated
219 
220  //step increment along the shower direction
221  G4double deltaStep = 0.0;
222 
223  // Uniqueness of G4Step is important otherwise hits won't be created.
224  G4double timeGlobal = track->GetStep()->GetPreStepPoint()->GetGlobalTime();
225 
226  // this needs to be deleted manually at the end of this loop.
227  // theGflashNavigator = new G4Navigator();
228  theGflashNavigator->SetWorldVolume(G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume());
229 
230  // // loop for longitudinal integration
231 
232 #ifdef DebugLog
233  LogDebug("HFShower") << " Energy = " << energy << " Step Length Left = " << stepLengthLeft;
234 #endif
235  while(energy > 0.0 && stepLengthLeft > 0.0) {
236  stepLengthLeftInX0 = stepLengthLeft / Gflash::radLength[jCalorimeter];
237 
238  if ( stepLengthLeftInX0 < divisionStepInX0 ) {
239  deltaZInX0 = stepLengthLeftInX0;
240  deltaZ = deltaZInX0 * Gflash::radLength[jCalorimeter];
241  stepLengthLeft = 0.0;
242  } else {
243  deltaZInX0 = divisionStepInX0;
244  deltaZ = deltaZInX0 * Gflash::radLength[jCalorimeter];
245  stepLengthLeft -= deltaZ;
246  }
247 
248  zInX0 += deltaZInX0;
249 
250 #ifdef DebugLog
251  LogDebug("HFShower") << " zInX0 = " << zInX0 << " spotBeta*zInX0 = " << spotBeta*zInX0;
252 #endif
253  if ((!zInX0) || (!spotBeta*zInX0) || (zInX0 < 0.01) ||
254  (spotBeta*zInX0 < 0.00001) || (!zInX0*beta) || (zInX0*beta < 0.00001))
255  return hit;
256 
257  G4int nSpotsInStep = 0;
258 
259 #ifdef DebugLog
260  LogDebug("HFShower") << " Energy - Energy Cut off = " << energy - energyCutoff;
261 #endif
262 
263  if ( energy > energyCutoff ) {
264  preEnergyInGamma = energyInGamma;
265  gammaDist.a().setValue(alpha); //alpha
266 
267  energyInGamma = gammaDist(beta*zInX0); //beta
268  G4double energyInDeltaZ = energyInGamma - preEnergyInGamma;
269  deltaEnergy = std::min(energy,((preStepPoint->GetTotalEnergy())/GeV)*energyInDeltaZ);
270 
271  preSigmaInGamma = sigmaInGamma;
272  gammaDist.a().setValue(spotAlpha); //alpha spot
273  sigmaInGamma = gammaDist(spotBeta*zInX0); //beta spot
274  nSpotsInStep = std::max(1,int(nSpots * (sigmaInGamma - preSigmaInGamma)));
275  } else {
276  deltaEnergy = energy;
277  preSigmaInGamma = sigmaInGamma;
278  nSpotsInStep = std::max(1,int(nSpots * (1.0 - preSigmaInGamma)));
279  }
280 
281  if ( deltaEnergy > energy || (energy-deltaEnergy) < energyCutoff ) deltaEnergy = energy;
282 
283  energy -= deltaEnergy;
284 
285  if ( spotCounter+nSpotsInStep > maxNumberOfSpots ) {
286  nSpotsInStep = maxNumberOfSpots - spotCounter;
287  if (nSpotsInStep < 1) nSpotsInStep = 1;
288  }
289 
290 
291  // // It begins with 0.5 of deltaZ and then icreases by 1 deltaZ
292  deltaStep += 0.5*deltaZ;
293  pathLength += deltaStep;
294  deltaStep = 0.5*deltaZ;
295 
296 
297  //lateral shape and fluctuations for homogenous calo.
298  G4double tScale = tmax *alpha/(alpha-1.0) * (std::exp(fluctuatedAlpha)-1.0)/std::exp(fluctuatedAlpha);
299  G4double tau = std::min(10.0,(zInX0 - 0.5*deltaZInX0)/tScale);
300  G4double rCore = z1 + z2 * tau;
301  G4double rTail = k1 *( std::exp(k3*(tau-k2)) + std::exp(k4*(tau-k2))); // @@ check RT3 sign
302  G4double p23 = (p2 - tau)/p3;
303  G4double probabilityWeight = p1 * std::exp( p23 - std::exp(p23) );
304 
305 
306  // Deposition of spots according to lateral distr.
307  // Apply absolute energy scale
308  // Convert into MeV unit
309  G4double emSpotEnergy = deltaEnergy / nSpotsInStep * e25Scale * GeV;
310 
311 
312 #ifdef DebugLog
313  LogDebug("HFShower") << " nSpotsInStep = " << nSpotsInStep;
314 #endif
315 
316 
317 
318  for (G4int ispot = 0 ; ispot < nSpotsInStep ; ispot++) {
319  spotCounter++;
320  G4double u1 = G4UniformRand();
321  G4double u2 = G4UniformRand();
322  G4double rInRM = 0.0;
323 
324  if (u1 < probabilityWeight) {
325  rInRM = rCore * std::sqrt( u2/(1.0-u2) );
326  } else {
327  rInRM = rTail * std::sqrt( u2/(1.0-u2) );
328  }
329 
330  G4double rShower = rInRM * Gflash::rMoliere[jCalorimeter];
331 
332  //Uniform & random rotation of spot along the azimuthal angle
333  G4double azimuthalAngle = twopi*G4UniformRand();
334 
335  //Compute global position of generated spots with taking into account magnetic field
336  //Divide deltaZ into nSpotsInStep and give a spot a global position
337  G4double incrementPath = (deltaZ/nSpotsInStep)*(ispot+0.5 - 0.5*nSpotsInStep);
338 
339  // trajectoryPoint give a spot an imaginary point along the shower development
340  GflashTrajectoryPoint trajectoryPoint;
341  theHelix->getGflashTrajectoryPoint(trajectoryPoint,pathLength+incrementPath);
342 
343 
344  G4ThreeVector SpotPosition0 = trajectoryPoint.getPosition() + rShower*std::cos(azimuthalAngle)*trajectoryPoint.getOrthogonalUnitVector() + rShower*std::sin(azimuthalAngle)*trajectoryPoint.getCrossUnitVector();
345 
346  // Convert into mm unit
347  SpotPosition0 *= cm;
348 
349 
350  //---------------------------------------------------
351  // fill a fake step to send it to hit maker
352  //---------------------------------------------------
353 
354  // to make a different time for each fake step. (0.03 nsec is corresponding to 1cm step size)
355  timeGlobal += 0.0001*nanosecond;
356 
357  //fill equivalent changes to a (fake) step associated with a spot
358 
359  G4double zInX0_spot = std::abs(pathLength+incrementPath - pathLength0)/Gflash::radLength[jCalorimeter];
360 
361 #ifdef DebugLog
362  LogDebug("HFShower") << "zInX0_spot,emSpotEnergy/GeV =" << zInX0_spot << " , " << emSpotEnergy/GeV << "emSpotEnergy/GeV =" << emSpotEnergy/GeV;
363 #endif
364 
365  if ((!zInX0_spot) || (zInX0_spot < 0)) continue;
366  if ((!emSpotEnergy/GeV) || (emSpotEnergy < 0)) continue;
367  if ((!rShower/Gflash::rMoliere[jCalorimeter]) || (rShower/Gflash::rMoliere[jCalorimeter] < 0)) continue;
368 
369  oneHit.depth = 1;
370 
371 #ifdef DebugLog
372  if (theFillHisto) {
373  em_long->Fill(SpotPosition0.z()/cm,emSpotEnergy/GeV);
374  em_lateral->Fill(rShower/Gflash::rMoliere[jCalorimeter],emSpotEnergy/GeV);
375  em_2d->Fill(SpotPosition0.z()/cm,rShower/Gflash::rMoliere[jCalorimeter],emSpotEnergy/GeV);
376  }
377 #endif
378 
379  if(SpotPosition0 == 0) continue;
380 
381  double energyratio = emSpotEnergy/(preStepPoint->GetTotalEnergy()/(nSpots*e25Scale));
382 
383  if (emSpotEnergy/GeV < 0.0001) continue;
384  if (energyratio > 80) continue;
385 
386  double zshift =0;
387  if(SpotPosition0.z() > 0) zshift = 18;
388  if(SpotPosition0.z() < 0) zshift = -18;
389 
390 
391  G4ThreeVector gfshift(0,0,zshift*(pow(100,0.1)/pow(preStepPoint->GetTotalEnergy()/GeV,0.1)));
392 
393  G4ThreeVector SpotPosition = gfshift + SpotPosition0;
394 
395  double LengthWeight = std::fabs(std::pow(SpotPosition0.z()/11370,1));
396  if (G4UniformRand()> 0.0021 * energyratio * LengthWeight) continue;
397 
398  oneHit.position = SpotPosition;
399  oneHit.time = timeGlobal;
400  oneHit.edep = emSpotEnergy/GeV;
401  hit.push_back(oneHit);
402  nSpots_sd++;
403 
404  } // end of for spot iteration
405 
406  } // end of while for longitudinal integration
407 #ifdef DebugLog
408  if (theFillHisto) {
409  em_nSpots_sd->Fill(nSpots_sd);
410  }
411 #endif
412  // delete theGflashNavigator;
413  return hit;
414 }
#define LogDebug(id)
const double beta
std::vector< Hit > gfParameterization(G4Step *aStep, bool &ok, bool onlyLong=false)
Definition: HFGflash.cc:88
virtual ~HFGflash()
Definition: HFGflash.cc:81
T getParameter(std::string const &) const
G4Step * theGflashStep
Definition: HFGflash.h:50
T getUntrackedParameter(std::string const &, T const &) const
TH2F * em_ze
Definition: HFGflash.h:69
G4Navigator * theGflashNavigator
Definition: HFGflash.h:51
float alpha
Definition: AMPTWrapper.h:95
double getPathLengthAtZ(double z) const
Gflash3Vector getCrossUnitVector()
TH1F * em_ssp_z
Definition: HFGflash.h:67
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Definition: DDAxes.h:10
#define abs(x)
Definition: mlp_lapack.h:159
#define min(a, b)
Definition: mlp_lapack.h:161
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
Gflash::CalorimeterNumber jCalorimeter
Definition: HFGflash.h:54
G4double theBField
Definition: HFGflash.h:58
double charge(const std::vector< uint8_t > &Ampls)
bool theFillHisto
Definition: HFGflash.h:57
const double rMoliere[kNumberCalorimeter]
const T & max(const T &a, const T &b)
T sqrt(T t)
Definition: SSEVec.h:28
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
TH2F * em_ratio_selected
Definition: HFGflash.h:69
bool isAvailable() const
Definition: Service.h:47
G4TouchableHandle theGflashTouchableHandle
Definition: HFGflash.h:52
double edep
Definition: HFGflash.h:41
bool theWatcherOn
Definition: HFGflash.h:56
double p2[4]
Definition: TauolaWrapper.h:90
TH2F * em_ratio
Definition: HFGflash.h:69
TH1F * em_ze_ratio
Definition: HFGflash.h:68
HFGflash(edm::ParameterSet const &p)
Definition: HFGflash.cc:34
Gflash3Vector getOrthogonalUnitVector()
void getGflashTrajectoryPoint(GflashTrajectoryPoint &point, double s) const
static const double tmax[3]
TH1F * em_ssp_rho
Definition: HFGflash.h:67
const double radLength[kNumberCalorimeter]
Log< T >::type log(const T &t)
Definition: Log.h:22
double time
Definition: HFGflash.h:40
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
TH1F * em_lateral_sd
Definition: HFGflash.h:68
TH1F * em_nSpots_sd
Definition: HFGflash.h:68
TH1F * em_lateral
Definition: HFGflash.h:67
TH1F * em_long
Definition: HFGflash.h:67
Gflash3Vector & getPosition()
double p1[4]
Definition: TauolaWrapper.h:89
G4ThreeVector position
Definition: HFGflash.h:38
T * make() const
make new ROOT object
TH1F * em_incE
Definition: HFGflash.h:67
GflashTrajectory * theHelix
Definition: HFGflash.h:49
TH1F * em_long_sd
Definition: HFGflash.h:68
TH2F * em_2d
Definition: HFGflash.h:69
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
double p3[4]
Definition: TauolaWrapper.h:91
TH2F * em_2d_sd
Definition: HFGflash.h:69
void initializeTrajectory(const HepGeom::Vector3D< double > &, const HepGeom::Point3D< double > &, double q, double Field)