Geant4_10
G4SynchrotronRadiationInMat.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: G4SynchrotronRadiationInMat.cc 68037 2013-03-13 14:15:08Z gcosmo $
28 //
29 // --------------------------------------------------------------
30 // GEANT 4 class implementation file
31 // CERN Geneva Switzerland
32 //
33 // History: first implementation,
34 // 21-5-98 V.Grichine
35 // 28-05-01, V.Ivanchenko minor changes to provide ANSI -wall compilation
36 // 04.03.05, V.Grichine: get local field interface
37 // 19-05-06, V.Ivanchenko rename from G4SynchrotronRadiation
38 //
39 //
41 
43 #include "G4PhysicalConstants.hh"
44 #include "G4SystemOfUnits.hh"
45 #include "G4Integrator.hh"
46 #include "G4EmProcessSubType.hh"
47 
49 //
50 // Constant for calculation of mean free path
51 //
52 
53 const G4double
54 G4SynchrotronRadiationInMat::fLambdaConst = std::sqrt(3.0)*electron_mass_c2/
56 
58 //
59 // Constant for calculation of characterictic energy
60 //
61 
62 const G4double
63 G4SynchrotronRadiationInMat::fEnergyConst = 1.5*c_light*c_light*eplus*hbar_Planck/
65 
67 //
68 // Array of integral probability of synchrotron photons:
69 //
70 // the corresponding energy = 0.0001*i*i*(characteristic energy)
71 //
72 
73 const G4double
74 G4SynchrotronRadiationInMat::fIntegralProbabilityOfSR[200] =
75 {
76  1.000000e+00, 9.428859e-01, 9.094095e-01, 8.813971e-01, 8.565154e-01,
77  8.337008e-01, 8.124961e-01, 7.925217e-01, 7.735517e-01, 7.554561e-01,
78  7.381233e-01, 7.214521e-01, 7.053634e-01, 6.898006e-01, 6.747219e-01,
79  6.600922e-01, 6.458793e-01, 6.320533e-01, 6.185872e-01, 6.054579e-01,
80  5.926459e-01, 5.801347e-01, 5.679103e-01, 5.559604e-01, 5.442736e-01,
81  5.328395e-01, 5.216482e-01, 5.106904e-01, 4.999575e-01, 4.894415e-01,
82  4.791351e-01, 4.690316e-01, 4.591249e-01, 4.494094e-01, 4.398800e-01,
83  4.305320e-01, 4.213608e-01, 4.123623e-01, 4.035325e-01, 3.948676e-01,
84  3.863639e-01, 3.780179e-01, 3.698262e-01, 3.617858e-01, 3.538933e-01,
85  3.461460e-01, 3.385411e-01, 3.310757e-01, 3.237474e-01, 3.165536e-01,
86  3.094921e-01, 3.025605e-01, 2.957566e-01, 2.890784e-01, 2.825237e-01,
87  2.760907e-01, 2.697773e-01, 2.635817e-01, 2.575020e-01, 2.515365e-01,
88  2.456834e-01, 2.399409e-01, 2.343074e-01, 2.287812e-01, 2.233607e-01,
89  2.180442e-01, 2.128303e-01, 2.077174e-01, 2.027040e-01, 1.977885e-01,
90  1.929696e-01, 1.882457e-01, 1.836155e-01, 1.790775e-01, 1.746305e-01,
91  1.702730e-01, 1.660036e-01, 1.618212e-01, 1.577243e-01, 1.537117e-01,
92  1.497822e-01, 1.459344e-01, 1.421671e-01, 1.384791e-01, 1.348691e-01,
93  1.313360e-01, 1.278785e-01, 1.244956e-01, 1.211859e-01, 1.179483e-01,
94  1.147818e-01, 1.116850e-01, 1.086570e-01, 1.056966e-01, 1.028026e-01,
95  9.997405e-02, 9.720975e-02, 9.450865e-02, 9.186969e-02, 8.929179e-02,
96  8.677391e-02, 8.431501e-02, 8.191406e-02, 7.957003e-02, 7.728192e-02,
97  7.504872e-02, 7.286944e-02, 7.074311e-02, 6.866874e-02, 6.664538e-02,
98  6.467208e-02, 6.274790e-02, 6.087191e-02, 5.904317e-02, 5.726079e-02,
99  5.552387e-02, 5.383150e-02, 5.218282e-02, 5.057695e-02, 4.901302e-02,
100  4.749020e-02, 4.600763e-02, 4.456450e-02, 4.315997e-02, 4.179325e-02,
101  4.046353e-02, 3.917002e-02, 3.791195e-02, 3.668855e-02, 3.549906e-02,
102  3.434274e-02, 3.321884e-02, 3.212665e-02, 3.106544e-02, 3.003452e-02,
103  2.903319e-02, 2.806076e-02, 2.711656e-02, 2.619993e-02, 2.531021e-02,
104  2.444677e-02, 2.360897e-02, 2.279620e-02, 2.200783e-02, 2.124327e-02,
105  2.050194e-02, 1.978324e-02, 1.908662e-02, 1.841151e-02, 1.775735e-02,
106  1.712363e-02, 1.650979e-02, 1.591533e-02, 1.533973e-02, 1.478250e-02,
107  1.424314e-02, 1.372117e-02, 1.321613e-02, 1.272755e-02, 1.225498e-02,
108  1.179798e-02, 1.135611e-02, 1.092896e-02, 1.051609e-02, 1.011712e-02,
109  9.731635e-03, 9.359254e-03, 8.999595e-03, 8.652287e-03, 8.316967e-03,
110  7.993280e-03, 7.680879e-03, 7.379426e-03, 7.088591e-03, 6.808051e-03,
111  6.537491e-03, 6.276605e-03, 6.025092e-03, 5.782661e-03, 5.549027e-03,
112  5.323912e-03, 5.107045e-03, 4.898164e-03, 4.697011e-03, 4.503336e-03,
113  4.316896e-03, 4.137454e-03, 3.964780e-03, 3.798649e-03, 3.638843e-03,
114  3.485150e-03, 3.337364e-03, 3.195284e-03, 3.058715e-03, 2.927469e-03,
115  2.801361e-03, 2.680213e-03, 2.563852e-03, 2.452110e-03, 2.344824e-03
116 };
117 
119 //
120 // Constructor
121 //
122 
124  G4ProcessType type):G4VDiscreteProcess (processName, type),
125  LowestKineticEnergy (10.*keV),
126  HighestKineticEnergy (100.*TeV),
127  TotBin(200),
128  theGamma (G4Gamma::Gamma() ),
129  theElectron ( G4Electron::Electron() ),
130  thePositron ( G4Positron::Positron() ),
131  GammaCutInKineticEnergy(0),
132  ElectronCutInKineticEnergy(0),
133  PositronCutInKineticEnergy(0),
134  ParticleCutInKineticEnergy(0),
135  fAlpha(0.0), fRootNumber(80),
136  fVerboseLevel( verboseLevel )
137 {
139 
140  fFieldPropagator = transportMgr->GetPropagatorInField();
142  CutInRange = GammaCutInKineticEnergyNow = ElectronCutInKineticEnergyNow =
143  PositronCutInKineticEnergyNow = ParticleCutInKineticEnergyNow = fKsi =
144  fPsiGamma = fEta = fOrderAngleK = 0.0;
145 }
146 
148 //
149 // Destructor
150 //
151 
153 {}
154 
155 
156 G4bool
158 {
159 
160  return ( ( &particle == (const G4ParticleDefinition *)theElectron ) ||
161  ( &particle == (const G4ParticleDefinition *)thePositron ) );
162 
163 }
164 
166 //
167 //
168 // Production of synchrotron X-ray photon
169 // GEANT4 internal units.
170 //
171 
172 
173 G4double
175  G4double,
177 {
178  // gives the MeanFreePath in GEANT4 internal units
179  G4double MeanFreePath;
180 
181  const G4DynamicParticle* aDynamicParticle = trackData.GetDynamicParticle();
182  // G4Material* aMaterial = trackData.GetMaterial();
183 
184  //G4bool isOutRange ;
185 
186  *condition = NotForced ;
187 
188  G4double gamma = aDynamicParticle->GetTotalEnergy()/
189  aDynamicParticle->GetMass();
190 
191  G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
192 
193  G4double KineticEnergy = aDynamicParticle->GetKineticEnergy();
194 
195  if ( KineticEnergy < LowestKineticEnergy || gamma < 1.0e3 ) MeanFreePath = DBL_MAX;
196  else
197  {
198 
199  G4ThreeVector FieldValue;
200  const G4Field* pField = 0;
201 
202  G4FieldManager* fieldMgr=0;
203  G4bool fieldExertsForce = false;
204 
205  if( (particleCharge != 0.0) )
206  {
207  fieldMgr = fFieldPropagator->FindAndSetFieldManager( trackData.GetVolume() );
208 
209  if ( fieldMgr != 0 )
210  {
211  // If the field manager has no field, there is no field !
212 
213  fieldExertsForce = ( fieldMgr->GetDetectorField() != 0 );
214  }
215  }
216  if ( fieldExertsForce )
217  {
218  pField = fieldMgr->GetDetectorField() ;
219  G4ThreeVector globPosition = trackData.GetPosition();
220 
221  G4double globPosVec[4], FieldValueVec[6];
222 
223  globPosVec[0] = globPosition.x();
224  globPosVec[1] = globPosition.y();
225  globPosVec[2] = globPosition.z();
226  globPosVec[3] = trackData.GetGlobalTime();
227 
228  pField->GetFieldValue( globPosVec, FieldValueVec );
229 
230  FieldValue = G4ThreeVector( FieldValueVec[0],
231  FieldValueVec[1],
232  FieldValueVec[2] );
233 
234 
235 
236  G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
237  G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum) ;
238  G4double perpB = unitMcrossB.mag() ;
239  G4double beta = aDynamicParticle->GetTotalMomentum()/
240  (aDynamicParticle->GetTotalEnergy() );
241 
242  if( perpB > 0.0 ) MeanFreePath = fLambdaConst*beta/perpB;
243  else MeanFreePath = DBL_MAX;
244  }
245  else MeanFreePath = DBL_MAX;
246  }
247  if(fVerboseLevel > 0)
248  {
249  G4cout<<"G4SynchrotronRadiationInMat::MeanFreePath = "<<MeanFreePath/m<<" m"<<G4endl;
250  }
251  return MeanFreePath;
252 }
253 
255 //
256 //
257 
260  const G4Step& stepData )
261 
262 {
263  aParticleChange.Initialize(trackData);
264 
265  const G4DynamicParticle* aDynamicParticle=trackData.GetDynamicParticle();
266 
267  G4double gamma = aDynamicParticle->GetTotalEnergy()/
268  (aDynamicParticle->GetMass() );
269 
270  if(gamma <= 1.0e3 )
271  {
272  return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
273  }
274  G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
275 
276  G4ThreeVector FieldValue;
277  const G4Field* pField = 0 ;
278 
279  G4FieldManager* fieldMgr=0;
280  G4bool fieldExertsForce = false;
281 
282  if( (particleCharge != 0.0) )
283  {
284  fieldMgr = fFieldPropagator->FindAndSetFieldManager( trackData.GetVolume() );
285  if ( fieldMgr != 0 )
286  {
287  // If the field manager has no field, there is no field !
288 
289  fieldExertsForce = ( fieldMgr->GetDetectorField() != 0 );
290  }
291  }
292  if ( fieldExertsForce )
293  {
294  pField = fieldMgr->GetDetectorField() ;
295  G4ThreeVector globPosition = trackData.GetPosition() ;
296  G4double globPosVec[4], FieldValueVec[6] ;
297  globPosVec[0] = globPosition.x() ;
298  globPosVec[1] = globPosition.y() ;
299  globPosVec[2] = globPosition.z() ;
300  globPosVec[3] = trackData.GetGlobalTime();
301 
302  pField->GetFieldValue( globPosVec, FieldValueVec ) ;
303  FieldValue = G4ThreeVector( FieldValueVec[0],
304  FieldValueVec[1],
305  FieldValueVec[2] );
306 
307  G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
308  G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum);
309  G4double perpB = unitMcrossB.mag() ;
310  if(perpB > 0.0)
311  {
312  // M-C of synchrotron photon energy
313 
314  G4double energyOfSR = GetRandomEnergySR(gamma,perpB);
315 
316  if(fVerboseLevel > 0)
317  {
318  G4cout<<"SR photon energy = "<<energyOfSR/keV<<" keV"<<G4endl;
319  }
320  // check against insufficient energy
321 
322  if( energyOfSR <= 0.0 )
323  {
324  return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
325  }
326  G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
328  particleDirection = aDynamicParticle->GetMomentumDirection();
329 
330  // M-C of its direction, simplified dipole busted approach
331 
332  // G4double Teta = G4UniformRand()/gamma ; // Very roughly
333 
334  G4double cosTheta, sinTheta, fcos, beta;
335 
336  do
337  {
338  cosTheta = 1. - 2.*G4UniformRand();
339  fcos = (1 + cosTheta*cosTheta)*0.5;
340  }
341  while( fcos < G4UniformRand() );
342 
343  beta = std::sqrt(1. - 1./(gamma*gamma));
344 
345  cosTheta = (cosTheta + beta)/(1. + beta*cosTheta);
346 
347  if( cosTheta > 1. ) cosTheta = 1.;
348  if( cosTheta < -1. ) cosTheta = -1.;
349 
350  sinTheta = std::sqrt(1. - cosTheta*cosTheta );
351 
352  G4double Phi = twopi * G4UniformRand() ;
353 
354  G4double dirx = sinTheta*std::cos(Phi) ,
355  diry = sinTheta*std::sin(Phi) ,
356  dirz = cosTheta;
357 
358  G4ThreeVector gammaDirection ( dirx, diry, dirz);
359  gammaDirection.rotateUz(particleDirection);
360 
361  // polarization of new gamma
362 
363  // G4double sx = std::cos(Teta)*std::cos(Phi);
364  // G4double sy = std::cos(Teta)*std::sin(Phi);
365  // G4double sz = -std::sin(Teta);
366 
367  G4ThreeVector gammaPolarization = FieldValue.cross(gammaDirection);
368  gammaPolarization = gammaPolarization.unit();
369 
370  // (sx, sy, sz);
371  // gammaPolarization.rotateUz(particleDirection);
372 
373  // create G4DynamicParticle object for the SR photon
374 
376  gammaDirection,
377  energyOfSR );
378  aGamma->SetPolarization( gammaPolarization.x(),
379  gammaPolarization.y(),
380  gammaPolarization.z() );
381 
382 
384  aParticleChange.AddSecondary(aGamma);
385 
386  // Update the incident particle
387 
388  G4double newKinEnergy = kineticEnergy - energyOfSR ;
389 
390  if (newKinEnergy > 0.)
391  {
392  aParticleChange.ProposeMomentumDirection( particleDirection );
393  aParticleChange.ProposeEnergy( newKinEnergy );
395  }
396  else
397  {
400  G4double charge = aDynamicParticle->GetDefinition()->GetPDGCharge();
401  if (charge<0.)
402  {
404  }
405  else
406  {
408  }
409  }
410  }
411  else
412  {
413  return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
414  }
415  }
416  return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
417 }
418 
419 
420 G4double
422  const G4Step& )
423 
424 {
425  G4int i ;
426  G4double energyOfSR = -1.0 ;
427  //G4Material* aMaterial=trackData.GetMaterial() ;
428 
429  const G4DynamicParticle* aDynamicParticle=trackData.GetDynamicParticle();
430 
431  G4double gamma = aDynamicParticle->GetTotalEnergy()/
432  (aDynamicParticle->GetMass() ) ;
433 
434  G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
435 
436  G4ThreeVector FieldValue;
437  const G4Field* pField = 0 ;
438 
439  G4FieldManager* fieldMgr=0;
440  G4bool fieldExertsForce = false;
441 
442  if( (particleCharge != 0.0) )
443  {
444  fieldMgr = fFieldPropagator->FindAndSetFieldManager( trackData.GetVolume() );
445  if ( fieldMgr != 0 )
446  {
447  // If the field manager has no field, there is no field !
448 
449  fieldExertsForce = ( fieldMgr->GetDetectorField() != 0 );
450  }
451  }
452  if ( fieldExertsForce )
453  {
454  pField = fieldMgr->GetDetectorField();
455  G4ThreeVector globPosition = trackData.GetPosition();
456  G4double globPosVec[3], FieldValueVec[3];
457 
458  globPosVec[0] = globPosition.x();
459  globPosVec[1] = globPosition.y();
460  globPosVec[2] = globPosition.z();
461 
462  pField->GetFieldValue( globPosVec, FieldValueVec );
463  FieldValue = G4ThreeVector( FieldValueVec[0],
464  FieldValueVec[1],
465  FieldValueVec[2] );
466 
467  G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
468  G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum) ;
469  G4double perpB = unitMcrossB.mag();
470  if( perpB > 0.0 )
471  {
472  // M-C of synchrotron photon energy
473 
474  G4double random = G4UniformRand() ;
475  for(i=0;i<200;i++)
476  {
477  if(random >= fIntegralProbabilityOfSR[i]) break ;
478  }
479  energyOfSR = 0.0001*i*i*fEnergyConst*gamma*gamma*perpB ;
480 
481  // check against insufficient energy
482 
483  if(energyOfSR <= 0.0)
484  {
485  return -1.0 ;
486  }
487  //G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
488  //G4ParticleMomentum
489  //particleDirection = aDynamicParticle->GetMomentumDirection();
490 
491  // Gamma production cut in this material
492  //G4double
493  //gammaEnergyCut = (G4Gamma::GetCutsInEnergy())[aMaterial->GetIndex()];
494 
495  // SR photon has energy more than the current material cut
496  // M-C of its direction
497 
498  //G4double Teta = G4UniformRand()/gamma ; // Very roughly
499 
500  //G4double Phi = twopi * G4UniformRand() ;
501  }
502  else
503  {
504  return -1.0 ;
505  }
506  }
507  return energyOfSR ;
508 }
509 
511 //
512 //
513 
515 {
516  G4int i, iMax;
517  G4double energySR, random, position;
518 
519  iMax = 200;
520  random = G4UniformRand();
521 
522  for( i = 0; i < iMax; i++ )
523  {
524  if( random >= fIntegralProbabilityOfSR[i] ) break;
525  }
526  if(i <= 0 ) position = G4UniformRand(); // 0.
527  else if( i>= iMax) position = G4double(iMax);
528  else position = i + G4UniformRand(); // -1
529  //
530  // it was in initial implementation:
531  // energyOfSR = 0.0001*i*i*fEnergyConst*gamma*gamma*perpB ;
532 
533  energySR = 0.0001*position*position*fEnergyConst*gamma*gamma*perpB;
534 
535  if( energySR < 0. ) energySR = 0.;
536 
537  return energySR;
538 }
539 
541 //
542 // return
543 
545 {
546  G4double result, hypCos2, hypCos=std::cosh(t);
547 
548  hypCos2 = hypCos*hypCos;
549  result = std::cosh(5.*t/3.)*std::exp(t-fKsi*hypCos); // fKsi > 0. !
550  result /= hypCos2;
551  return result;
552 }
553 
555 //
556 // return the probability to emit SR photon with relative energy
557 // energy/energy_c >= ksi
558 // for ksi <= 0. P = 1., however the method works for ksi > 0 only!
559 
561 {
562  if (ksi <= 0.) return 1.0;
563  fKsi = ksi; // should be > 0. !
564  G4int n;
565  G4double result, a;
566 
567  a = fAlpha; // always = 0.
568  n = fRootNumber; // around default = 80
569 
571 
572  result = integral.Laguerre(this,
574 
575  result *= 3./5./pi;
576 
577  return result;
578 }
579 
581 //
582 // return an auxiliary function for K_5/3 integral representation
583 
585 {
586  G4double result, hypCos=std::cosh(t);
587 
588  result = std::cosh(5.*t/3.)*std::exp(t - fKsi*hypCos); // fKsi > 0. !
589  result /= hypCos;
590  return result;
591 }
592 
594 //
595 // return the probability to emit SR photon energy with relative energy
596 // energy/energy_c >= ksi
597 // for ksi <= 0. P = 1., however the method works for ksi > 0 only!
598 
600 {
601  if (ksi <= 0.) return 1.0;
602  fKsi = ksi; // should be > 0. !
603  G4int n;
604  G4double result, a;
605 
606  a = fAlpha; // always = 0.
607  n = fRootNumber; // around default = 80
608 
610 
611  result = integral.Laguerre(this,
613 
614  result *= 9.*std::sqrt(3.)*ksi/8./pi;
615 
616  return result;
617 }
618 
620 //
621 //
622 
624 {
625  G4double result, hypCos=std::cosh(t);
626 
627  result = std::cosh(fOrderAngleK*t)*std::exp(t - fEta*hypCos); // fEta > 0. !
628  result /= hypCos;
629  return result;
630 }
631 
633 //
634 // Return K 1/3 or 2/3 for angular distribution
635 
637 {
638  fEta = eta; // should be > 0. !
639  G4int n;
640  G4double result, a;
641 
642  a = fAlpha; // always = 0.
643  n = fRootNumber; // around default = 80
644 
646 
647  result = integral.Laguerre(this,
649 
650  return result;
651 }
652 
654 //
655 // Relative angle diff distribution for given fKsi, which is set externally
656 
658 {
659  G4double result, funK, funK2, gpsi2 = gpsi*gpsi;
660 
661  fPsiGamma = gpsi;
662  fEta = 0.5*fKsi*(1 + gpsi2)*std::sqrt(1 + gpsi2);
663 
664  fOrderAngleK = 1./3.;
665  funK = GetAngleK(fEta);
666  funK2 = funK*funK;
667 
668  result = gpsi2*funK2/(1 + gpsi2);
669 
670  fOrderAngleK = 2./3.;
671  funK = GetAngleK(fEta);
672  funK2 = funK*funK;
673 
674  result += funK2;
675  result *= (1 + gpsi2)*fKsi;
676 
677  return result;
678 }
679 
680 
682 
G4double condition(const G4ErrorSymMatrix &m)
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
tuple a
Definition: test.py:11
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4double GetTotalEnergy() const
virtual void GetFieldValue(const double Point[4], double *fieldArr) const =0
double x() const
const G4DynamicParticle * GetDynamicParticle() const
const G4ThreeVector & GetPosition() const
G4double G4NeutronHPJENDLHEData::G4double result
G4ParticleDefinition * GetDefinition() const
G4double fcos(G4double arg)
int G4int
Definition: G4Types.hh:78
double z() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double GetTotalMomentum() const
Char_t n[5]
#define G4UniformRand()
Definition: Randomize.hh:87
G4SynchrotronRadiationInMat(const G4String &processName="SynchrotronRadiation", G4ProcessType type=fElectromagnetic)
G4GLOB_DLL std::ostream G4cout
G4double GetPhotonEnergy(const G4Track &trackData, const G4Step &stepData)
G4double GetMass() const
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
float electron_mass_c2
Definition: hepunit.py:274
Definition: G4Step.hh:76
G4double GetGlobalTime() const
void SetPolarization(G4double polX, G4double polY, G4double polZ)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
static G4TransportationManager * GetTransportationManager()
virtual void Initialize(const G4Track &)
void SetNumberOfSecondaries(G4int totSecondaries)
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &Step)
Hep3Vector unit() const
int position
Definition: filter.cc:7
double y() const
void ProposeEnergy(G4double finalEnergy)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
G4VPhysicalVolume * GetVolume() const
void AddSecondary(G4Track *aSecondary)
#define G4endl
Definition: G4ios.hh:61
Hep3Vector cross(const Hep3Vector &) const
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4bool IsApplicable(const G4ParticleDefinition &)
const G4Field * GetDetectorField() const
G4ForceCondition
G4double GetPDGCharge() const
G4double GetRandomEnergySR(G4double, G4double)
G4PropagatorInField * GetPropagatorInField() const
double mag() const
#define DBL_MAX
Definition: templates.hh:83
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
float c_light
Definition: hepunit.py:257
G4ProcessType
G4double Laguerre(T &typeT, F f, G4double alpha, G4int n)