Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// GEANT 4 class implementation file
28//
29// History: first implementation,
30// 21-5-98 V.Grichine
31// 28-05-01, V.Ivanchenko minor changes to provide ANSI -wall compilation
32// 04.03.05, V.Grichine: get local field interface
33// 19-05-06, V.Ivanchenko rename from G4SynchrotronRadiation
34//
35///////////////////////////////////////////////////////////////////////////
36
38
39#include "G4EmProcessSubType.hh"
40#include "G4Field.hh"
41#include "G4FieldManager.hh"
42#include "G4Integrator.hh"
45#include "G4SystemOfUnits.hh"
47
48const G4double G4SynchrotronRadiationInMat::fIntegralProbabilityOfSR[200] = {
49 1.000000e+00, 9.428859e-01, 9.094095e-01, 8.813971e-01, 8.565154e-01,
50 8.337008e-01, 8.124961e-01, 7.925217e-01, 7.735517e-01, 7.554561e-01,
51 7.381233e-01, 7.214521e-01, 7.053634e-01, 6.898006e-01, 6.747219e-01,
52 6.600922e-01, 6.458793e-01, 6.320533e-01, 6.185872e-01, 6.054579e-01,
53 5.926459e-01, 5.801347e-01, 5.679103e-01, 5.559604e-01, 5.442736e-01,
54 5.328395e-01, 5.216482e-01, 5.106904e-01, 4.999575e-01, 4.894415e-01,
55 4.791351e-01, 4.690316e-01, 4.591249e-01, 4.494094e-01, 4.398800e-01,
56 4.305320e-01, 4.213608e-01, 4.123623e-01, 4.035325e-01, 3.948676e-01,
57 3.863639e-01, 3.780179e-01, 3.698262e-01, 3.617858e-01, 3.538933e-01,
58 3.461460e-01, 3.385411e-01, 3.310757e-01, 3.237474e-01, 3.165536e-01,
59 3.094921e-01, 3.025605e-01, 2.957566e-01, 2.890784e-01, 2.825237e-01,
60 2.760907e-01, 2.697773e-01, 2.635817e-01, 2.575020e-01, 2.515365e-01,
61 2.456834e-01, 2.399409e-01, 2.343074e-01, 2.287812e-01, 2.233607e-01,
62 2.180442e-01, 2.128303e-01, 2.077174e-01, 2.027040e-01, 1.977885e-01,
63 1.929696e-01, 1.882457e-01, 1.836155e-01, 1.790775e-01, 1.746305e-01,
64 1.702730e-01, 1.660036e-01, 1.618212e-01, 1.577243e-01, 1.537117e-01,
65 1.497822e-01, 1.459344e-01, 1.421671e-01, 1.384791e-01, 1.348691e-01,
66 1.313360e-01, 1.278785e-01, 1.244956e-01, 1.211859e-01, 1.179483e-01,
67 1.147818e-01, 1.116850e-01, 1.086570e-01, 1.056966e-01, 1.028026e-01,
68 9.997405e-02, 9.720975e-02, 9.450865e-02, 9.186969e-02, 8.929179e-02,
69 8.677391e-02, 8.431501e-02, 8.191406e-02, 7.957003e-02, 7.728192e-02,
70 7.504872e-02, 7.286944e-02, 7.074311e-02, 6.866874e-02, 6.664538e-02,
71 6.467208e-02, 6.274790e-02, 6.087191e-02, 5.904317e-02, 5.726079e-02,
72 5.552387e-02, 5.383150e-02, 5.218282e-02, 5.057695e-02, 4.901302e-02,
73 4.749020e-02, 4.600763e-02, 4.456450e-02, 4.315997e-02, 4.179325e-02,
74 4.046353e-02, 3.917002e-02, 3.791195e-02, 3.668855e-02, 3.549906e-02,
75 3.434274e-02, 3.321884e-02, 3.212665e-02, 3.106544e-02, 3.003452e-02,
76 2.903319e-02, 2.806076e-02, 2.711656e-02, 2.619993e-02, 2.531021e-02,
77 2.444677e-02, 2.360897e-02, 2.279620e-02, 2.200783e-02, 2.124327e-02,
78 2.050194e-02, 1.978324e-02, 1.908662e-02, 1.841151e-02, 1.775735e-02,
79 1.712363e-02, 1.650979e-02, 1.591533e-02, 1.533973e-02, 1.478250e-02,
80 1.424314e-02, 1.372117e-02, 1.321613e-02, 1.272755e-02, 1.225498e-02,
81 1.179798e-02, 1.135611e-02, 1.092896e-02, 1.051609e-02, 1.011712e-02,
82 9.731635e-03, 9.359254e-03, 8.999595e-03, 8.652287e-03, 8.316967e-03,
83 7.993280e-03, 7.680879e-03, 7.379426e-03, 7.088591e-03, 6.808051e-03,
84 6.537491e-03, 6.276605e-03, 6.025092e-03, 5.782661e-03, 5.549027e-03,
85 5.323912e-03, 5.107045e-03, 4.898164e-03, 4.697011e-03, 4.503336e-03,
86 4.316896e-03, 4.137454e-03, 3.964780e-03, 3.798649e-03, 3.638843e-03,
87 3.485150e-03, 3.337364e-03, 3.195284e-03, 3.058715e-03, 2.927469e-03,
88 2.801361e-03, 2.680213e-03, 2.563852e-03, 2.452110e-03, 2.344824e-03
89};
90
91///////////////////////////////////////////////////////////////////////
92// Constructor
94 const G4String& processName, G4ProcessType type)
95 : G4VDiscreteProcess(processName, type)
96 , theGamma(G4Gamma::Gamma())
97 , theElectron(G4Electron::Electron())
98 , thePositron(G4Positron::Positron())
99 , LowestKineticEnergy(10. * keV)
100 , fAlpha(0.0)
101 , fRootNumber(80)
102 , fVerboseLevel(verboseLevel)
103{
104 G4TransportationManager* transportMgr =
106
107 fFieldPropagator = transportMgr->GetPropagatorInField();
108 secID = G4PhysicsModelCatalog::GetModelID("model_SynchrotronRadiation");
110 CutInRange = GammaCutInKineticEnergyNow = ElectronCutInKineticEnergyNow =
111 PositronCutInKineticEnergyNow = ParticleCutInKineticEnergyNow = fKsi =
112 fPsiGamma = fEta = fOrderAngleK = 0.0;
113}
114
115/////////////////////////////////////////////////////////////////////////
116// Destructor
118
120 const G4ParticleDefinition& particle)
121{
122 return ((&particle == (const G4ParticleDefinition*) theElectron) ||
123 (&particle == (const G4ParticleDefinition*) thePositron));
124}
125
127
129
130// Production of synchrotron X-ray photon
131// Geant4 internal units.
133 const G4Track& trackData, G4double, G4ForceCondition* condition)
134{
135 // gives the MeanFreePath in GEANT4 internal units
136 G4double MeanFreePath;
137
138 const G4DynamicParticle* aDynamicParticle = trackData.GetDynamicParticle();
139
141
142 G4double gamma =
143 aDynamicParticle->GetTotalEnergy() / aDynamicParticle->GetMass();
144
145 G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
146
147 G4double KineticEnergy = aDynamicParticle->GetKineticEnergy();
148
149 if(KineticEnergy < LowestKineticEnergy || gamma < 1.0e3)
150 MeanFreePath = DBL_MAX;
151 else
152 {
153 G4ThreeVector FieldValue;
154 const G4Field* pField = nullptr;
155
156 G4FieldManager* fieldMgr = nullptr;
157 G4bool fieldExertsForce = false;
158
159 if((particleCharge != 0.0))
160 {
161 fieldMgr =
162 fFieldPropagator->FindAndSetFieldManager(trackData.GetVolume());
163
164 if(fieldMgr != nullptr)
165 {
166 // If the field manager has no field, there is no field !
167 fieldExertsForce = (fieldMgr->GetDetectorField() != nullptr);
168 }
169 }
170 if(fieldExertsForce)
171 {
172 pField = fieldMgr->GetDetectorField();
173 G4ThreeVector globPosition = trackData.GetPosition();
174
175 G4double globPosVec[4], FieldValueVec[6];
176
177 globPosVec[0] = globPosition.x();
178 globPosVec[1] = globPosition.y();
179 globPosVec[2] = globPosition.z();
180 globPosVec[3] = trackData.GetGlobalTime();
181
182 pField->GetFieldValue(globPosVec, FieldValueVec);
183
184 FieldValue =
185 G4ThreeVector(FieldValueVec[0], FieldValueVec[1], FieldValueVec[2]);
186
187 G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
188 G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum);
189 G4double perpB = unitMcrossB.mag();
190 G4double beta = aDynamicParticle->GetTotalMomentum() /
191 (aDynamicParticle->GetTotalEnergy());
192
193 if(perpB > 0.0)
194 MeanFreePath = fLambdaConst * beta / perpB;
195 else
196 MeanFreePath = DBL_MAX;
197 }
198 else
199 MeanFreePath = DBL_MAX;
200 }
201 if(fVerboseLevel > 0)
202 {
203 G4cout << "G4SynchrotronRadiationInMat::MeanFreePath = " << MeanFreePath / m
204 << " m" << G4endl;
205 }
206 return MeanFreePath;
207}
208
209////////////////////////////////////////////////////////////////////////////////
211 const G4Track& trackData, const G4Step& stepData)
212
213{
214 aParticleChange.Initialize(trackData);
215
216 const G4DynamicParticle* aDynamicParticle = trackData.GetDynamicParticle();
217
218 G4double gamma =
219 aDynamicParticle->GetTotalEnergy() / (aDynamicParticle->GetMass());
220
221 if(gamma <= 1.0e3)
222 {
223 return G4VDiscreteProcess::PostStepDoIt(trackData, stepData);
224 }
225 G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
226
227 G4ThreeVector FieldValue;
228 const G4Field* pField = nullptr;
229
230 G4FieldManager* fieldMgr = nullptr;
231 G4bool fieldExertsForce = false;
232
233 if((particleCharge != 0.0))
234 {
235 fieldMgr = fFieldPropagator->FindAndSetFieldManager(trackData.GetVolume());
236 if(fieldMgr != nullptr)
237 {
238 // If the field manager has no field, there is no field !
239 fieldExertsForce = (fieldMgr->GetDetectorField() != nullptr);
240 }
241 }
242 if(fieldExertsForce)
243 {
244 pField = fieldMgr->GetDetectorField();
245 G4ThreeVector globPosition = trackData.GetPosition();
246 G4double globPosVec[4], FieldValueVec[6];
247 globPosVec[0] = globPosition.x();
248 globPosVec[1] = globPosition.y();
249 globPosVec[2] = globPosition.z();
250 globPosVec[3] = trackData.GetGlobalTime();
251
252 pField->GetFieldValue(globPosVec, FieldValueVec);
253 FieldValue =
254 G4ThreeVector(FieldValueVec[0], FieldValueVec[1], FieldValueVec[2]);
255
256 G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
257 G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum);
258 G4double perpB = unitMcrossB.mag();
259 if(perpB > 0.0)
260 {
261 // M-C of synchrotron photon energy
262 G4double energyOfSR = GetRandomEnergySR(gamma, perpB);
263
264 if(fVerboseLevel > 0)
265 {
266 G4cout << "SR photon energy = " << energyOfSR / keV << " keV" << G4endl;
267 }
268
269 // check against insufficient energy
270 if(energyOfSR <= 0.0)
271 {
272 return G4VDiscreteProcess::PostStepDoIt(trackData, stepData);
273 }
274 G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
275 G4ParticleMomentum particleDirection =
276 aDynamicParticle->GetMomentumDirection();
277
278 // M-C of its direction, simplified dipole busted approach
279 G4double cosTheta, sinTheta, fcos, beta;
280
281 do
282 {
283 cosTheta = 1. - 2. * G4UniformRand();
284 fcos = (1 + cosTheta * cosTheta) * 0.5;
285 }
286 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
287 while(fcos < G4UniformRand());
288
289 beta = std::sqrt(1. - 1. / (gamma * gamma));
290
291 cosTheta = (cosTheta + beta) / (1. + beta * cosTheta);
292
293 if(cosTheta > 1.)
294 cosTheta = 1.;
295 if(cosTheta < -1.)
296 cosTheta = -1.;
297
298 sinTheta = std::sqrt(1. - cosTheta * cosTheta);
299
300 G4double Phi = twopi * G4UniformRand();
301
302 G4double dirx = sinTheta * std::cos(Phi);
303 G4double diry = sinTheta * std::sin(Phi);
304 G4double dirz = cosTheta;
305
306 G4ThreeVector gammaDirection(dirx, diry, dirz);
307 gammaDirection.rotateUz(particleDirection);
308
309 // polarization of new gamma
310 G4ThreeVector gammaPolarization = FieldValue.cross(gammaDirection);
311 gammaPolarization = gammaPolarization.unit();
312
313 // create G4DynamicParticle object for the SR photon
314 auto aGamma =
315 new G4DynamicParticle(G4Gamma::Gamma(), gammaDirection, energyOfSR);
316 aGamma->SetPolarization(gammaPolarization.x(), gammaPolarization.y(),
317 gammaPolarization.z());
318
320
321 // Update the incident particle
322 G4double newKinEnergy = kineticEnergy - energyOfSR;
323
324 if(newKinEnergy > 0.)
325 {
326 aParticleChange.ProposeMomentumDirection(particleDirection);
327 aParticleChange.ProposeEnergy(newKinEnergy);
329 }
330 else
331 {
334 G4double charge = aDynamicParticle->GetDefinition()->GetPDGCharge();
335 if(charge < 0.)
336 {
338 }
339 else
340 {
342 }
343 }
344
345 // Create the G4Track
346 G4Track* aSecondaryTrack = new G4Track(aGamma, trackData.GetGlobalTime(), trackData.GetPosition());
347 aSecondaryTrack->SetTouchableHandle(stepData.GetPostStepPoint()->GetTouchableHandle());
348 aSecondaryTrack->SetParentID(trackData.GetTrackID());
349 aSecondaryTrack->SetCreatorModelID(secID);
350 aParticleChange.AddSecondary(aSecondaryTrack);
351 }
352 else
353 {
354 return G4VDiscreteProcess::PostStepDoIt(trackData, stepData);
355 }
356 }
357 return G4VDiscreteProcess::PostStepDoIt(trackData, stepData);
358}
359
361 const G4Step&)
362{
363 G4int i;
364 G4double energyOfSR = -1.0;
365
366 const G4DynamicParticle* aDynamicParticle = trackData.GetDynamicParticle();
367
368 G4double gamma =
369 aDynamicParticle->GetTotalEnergy() / (aDynamicParticle->GetMass());
370
371 G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
372
373 G4ThreeVector FieldValue;
374 const G4Field* pField = nullptr;
375
376 G4FieldManager* fieldMgr = nullptr;
377 G4bool fieldExertsForce = false;
378
379 if((particleCharge != 0.0))
380 {
381 fieldMgr = fFieldPropagator->FindAndSetFieldManager(trackData.GetVolume());
382 if(fieldMgr != nullptr)
383 {
384 // If the field manager has no field, there is no field !
385 fieldExertsForce = (fieldMgr->GetDetectorField() != nullptr);
386 }
387 }
388 if(fieldExertsForce)
389 {
390 pField = fieldMgr->GetDetectorField();
391 G4ThreeVector globPosition = trackData.GetPosition();
392 G4double globPosVec[3], FieldValueVec[3];
393
394 globPosVec[0] = globPosition.x();
395 globPosVec[1] = globPosition.y();
396 globPosVec[2] = globPosition.z();
397
398 pField->GetFieldValue(globPosVec, FieldValueVec);
399 FieldValue =
400 G4ThreeVector(FieldValueVec[0], FieldValueVec[1], FieldValueVec[2]);
401
402 G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
403 G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum);
404 G4double perpB = unitMcrossB.mag();
405 if(perpB > 0.0)
406 {
407 // M-C of synchrotron photon energy
408 G4double random = G4UniformRand();
409 for(i = 0; i < 200; ++i)
410 {
411 if(random >= fIntegralProbabilityOfSR[i])
412 break;
413 }
414 energyOfSR = 0.0001 * i * i * fEnergyConst * gamma * gamma * perpB;
415
416 // check against insufficient energy
417 if(energyOfSR <= 0.0)
418 {
419 return -1.0;
420 }
421 }
422 else
423 {
424 return -1.0;
425 }
426 }
427 return energyOfSR;
428}
429
430/////////////////////////////////////////////////////////////////////////////////
432 G4double perpB)
433{
434 G4int i;
435 static constexpr G4int iMax = 200;
436 G4double energySR, random, position;
437
438 random = G4UniformRand();
439
440 for(i = 0; i < iMax; ++i)
441 {
442 if(random >= fIntegralProbabilityOfSR[i])
443 break;
444 }
445 if(i <= 0)
447 else if(i >= iMax)
448 position = G4double(iMax);
449 else
450 position = i + G4UniformRand();
451
452 energySR =
453 0.0001 * position * position * fEnergyConst * gamma * gamma * perpB;
454
455 if(energySR < 0.)
456 energySR = 0.;
457
458 return energySR;
459}
460
461/////////////////////////////////////////////////////////////////////////
463{
464 G4double result, hypCos2, hypCos = std::cosh(t);
465
466 hypCos2 = hypCos * hypCos;
467 result = std::cosh(5. * t / 3.) * std::exp(t - fKsi * hypCos); // fKsi > 0. !
468 result /= hypCos2;
469 return result;
470}
471
472///////////////////////////////////////////////////////////////////////////
473// return the probability to emit SR photon with relative energy
474// energy/energy_c >= ksi
475// for ksi <= 0. P = 1., however the method works for ksi > 0 only!
477{
478 if(ksi <= 0.)
479 return 1.0;
480 fKsi = ksi; // should be > 0. !
481 G4int n;
482 G4double result, a;
483
484 a = fAlpha; // always = 0.
485 n = fRootNumber; // around default = 80
486
489 integral;
490
491 result = integral.Laguerre(
493
494 result *= 3. / 5. / pi;
495
496 return result;
497}
498
499/////////////////////////////////////////////////////////////////////////
500// return an auxiliary function for K_5/3 integral representation
502{
503 G4double result, hypCos = std::cosh(t);
504
505 result = std::cosh(5. * t / 3.) * std::exp(t - fKsi * hypCos); // fKsi > 0. !
506 result /= hypCos;
507 return result;
508}
509
510///////////////////////////////////////////////////////////////////////////
511// return the probability to emit SR photon energy with relative energy
512// energy/energy_c >= ksi
513// for ksi <= 0. P = 1., however the method works for ksi > 0 only!
515{
516 if(ksi <= 0.)
517 return 1.0;
518 fKsi = ksi; // should be > 0. !
519 G4int n;
520 G4double result, a;
521
522 a = fAlpha; // always = 0.
523 n = fRootNumber; // around default = 80
524
527 integral;
528
529 result = integral.Laguerre(
531
532 result *= 9. * std::sqrt(3.) * ksi / 8. / pi;
533
534 return result;
535}
536
537/////////////////////////////////////////////////////////////////////////////
539{
540 G4double result, hypCos = std::cosh(t);
541
542 result =
543 std::cosh(fOrderAngleK * t) * std::exp(t - fEta * hypCos); // fEta > 0. !
544 result /= hypCos;
545 return result;
546}
547
548//////////////////////////////////////////////////////////////////////////
549// Return K 1/3 or 2/3 for angular distribution
551{
552 fEta = eta; // should be > 0. !
553 G4int n;
554 G4double result, a;
555
556 a = fAlpha; // always = 0.
557 n = fRootNumber; // around default = 80
558
561 integral;
562
563 result = integral.Laguerre(
565
566 return result;
567}
568
569/////////////////////////////////////////////////////////////////////////
570// Relative angle diff distribution for given fKsi, which is set externally
572{
573 G4double result, funK, funK2, gpsi2 = gpsi * gpsi;
574
575 fPsiGamma = gpsi;
576 fEta = 0.5 * fKsi * (1. + gpsi2) * std::sqrt(1. + gpsi2);
577
578 fOrderAngleK = 1. / 3.;
579 funK = GetAngleK(fEta);
580 funK2 = funK * funK;
581
582 result = gpsi2 * funK2 / (1. + gpsi2);
583
584 fOrderAngleK = 2. / 3.;
585 funK = GetAngleK(fEta);
586 funK2 = funK * funK;
587
588 result += funK2;
589 result *= (1. + gpsi2) * fKsi;
590
591 return result;
592}
@ fSynchrotronRadiation
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ NotForced
const G4double fAlpha
G4ProcessType
CLHEP::Hep3Vector G4ThreeVector
@ fStopAndKill
@ fStopButAlive
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
double z() const
Hep3Vector unit() const
double x() const
double y() const
Hep3Vector cross(const Hep3Vector &) const
double mag() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
G4double GetMass() const
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
G4double GetTotalMomentum() const
const G4Field * GetDetectorField() const
virtual void GetFieldValue(const G4double Point[4], G4double *fieldArr) const =0
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
void AddSecondary(G4Track *aSecondary)
void Initialize(const G4Track &) override
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double GetPDGCharge() const
static G4int GetModelID(const G4int modelIndex)
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
const G4TouchableHandle & GetTouchableHandle() const
Definition: G4Step.hh:62
G4StepPoint * GetPostStepPoint() const
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &Step) override
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4SynchrotronRadiationInMat(const G4String &processName="SynchrotronRadiation", G4ProcessType type=fElectromagnetic)
G4double GetRandomEnergySR(G4double, G4double)
G4double GetPhotonEnergy(const G4Track &trackData, const G4Step &stepData)
G4bool IsApplicable(const G4ParticleDefinition &) override
G4int GetTrackID() const
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
const G4DynamicParticle * GetDynamicParticle() const
void SetCreatorModelID(const G4int id)
void SetParentID(const G4int aValue)
static G4TransportationManager * GetTransportationManager()
G4PropagatorInField * GetPropagatorInField() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:331
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
#define DBL_MAX
Definition: templates.hh:62