Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4SynchrotronRadiation Class Reference

#include <G4SynchrotronRadiation.hh>

Inheritance diagram for G4SynchrotronRadiation:
Collaboration diagram for G4SynchrotronRadiation:

Public Member Functions

 G4SynchrotronRadiation (const G4String &pName="SynRad", G4ProcessType type=fElectromagnetic)
 
virtual ~G4SynchrotronRadiation ()
 
virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &Step) override
 
G4double GetPhotonEnergy (const G4Track &trackData, const G4Step &stepData)
 
G4double GetRandomEnergySR (G4double, G4double, G4double)
 
G4double InvSynFracInt (G4double x)
 
G4double Chebyshev (G4double a, G4double b, const G4double c[], G4int n, G4double x)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &) override
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &) override
 
virtual void PrintInfoDefinition ()
 
void SetAngularGenerator (G4VEmAngularDistribution *p)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 64 of file G4SynchrotronRadiation.hh.

Constructor & Destructor Documentation

G4SynchrotronRadiation::G4SynchrotronRadiation ( const G4String pName = "SynRad",
G4ProcessType  type = fElectromagnetic 
)
explicit

Definition at line 57 of file G4SynchrotronRadiation.cc.

58  :G4VDiscreteProcess (processName, type),
59  theGamma (G4Gamma::Gamma() )
60 {
61  G4TransportationManager* transportMgr =
63 
64  fFieldPropagator = transportMgr->GetPropagatorInField();
65 
67  verboseLevel = 1;
68  FirstTime = true;
69  FirstTime1 = true;
70  genAngle = nullptr;
72 }
G4int verboseLevel
Definition: G4VProcess.hh:368
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4TransportationManager * GetTransportationManager()
G4PropagatorInField * GetPropagatorInField() const
void SetAngularGenerator(G4VEmAngularDistribution *p)

Here is the call graph for this function:

G4SynchrotronRadiation::~G4SynchrotronRadiation ( )
virtual

Definition at line 79 of file G4SynchrotronRadiation.cc.

80 {
81  delete genAngle;
82 }

Member Function Documentation

void G4SynchrotronRadiation::BuildPhysicsTable ( const G4ParticleDefinition part)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 405 of file G4SynchrotronRadiation.cc.

406 {
408  // same for all particles, print only for one (electron)
409 }
G4int verboseLevel
Definition: G4VProcess.hh:368
static G4Electron * Electron()
Definition: G4Electron.cc:94

Here is the call graph for this function:

G4double G4SynchrotronRadiation::Chebyshev ( G4double  a,
G4double  b,
const G4double  c[],
G4int  n,
G4double  x 
)
inline

Definition at line 113 of file G4SynchrotronRadiation.hh.

115 {
116  G4double y;
117  G4double y2=2.0*(y=(2.0*x-a-b)/(b-a)); // Change of variable.
118  G4double d=0.,dd=0.;
119  for (G4int j=n-1;j>=1;--j) // Clenshaw's recurrence.
120  { G4double sv=d;
121  d=y2*d-dd+c[j];
122  dd=sv;
123  }
124  return y*d-dd+0.5*c[0];
125 }
int G4int
Definition: G4Types.hh:78
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

G4double G4SynchrotronRadiation::GetMeanFreePath ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
overridevirtual

Implements G4VDiscreteProcess.

Definition at line 109 of file G4SynchrotronRadiation.cc.

112 {
113  // gives the MeanFreePath in GEANT4 internal units
114  G4double MeanFreePath = DBL_MAX;
115 
116  const G4DynamicParticle* aDynamicParticle = trackData.GetDynamicParticle();
117 
118  *condition = NotForced;
119 
120  G4double gamma = aDynamicParticle->GetTotalEnergy()/
121  aDynamicParticle->GetMass();
122 
123  G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
124 
125  if ( gamma < 1.0e3 || 0.0 == particleCharge) { MeanFreePath = DBL_MAX; }
126  else
127  {
128 
129  G4ThreeVector FieldValue;
130  const G4Field* pField = nullptr;
131  G4bool fieldExertsForce = false;
132 
133 
134  G4FieldManager* fieldMgr =
135  fFieldPropagator->FindAndSetFieldManager(trackData.GetVolume());
136 
137  if ( fieldMgr != nullptr )
138  {
139  // If the field manager has no field, there is no field !
140 
141  fieldExertsForce = ( fieldMgr->GetDetectorField() != nullptr );
142  }
143 
144  if ( fieldExertsForce )
145  {
146  pField = fieldMgr->GetDetectorField();
147  G4ThreeVector globPosition = trackData.GetPosition();
148 
149  G4double globPosVec[4], FieldValueVec[6];
150 
151  globPosVec[0] = globPosition.x();
152  globPosVec[1] = globPosition.y();
153  globPosVec[2] = globPosition.z();
154  globPosVec[3] = trackData.GetGlobalTime();
155 
156  pField->GetFieldValue( globPosVec, FieldValueVec );
157 
158  FieldValue = G4ThreeVector( FieldValueVec[0],
159  FieldValueVec[1],
160  FieldValueVec[2] );
161 
162  G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
163  G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum);
164  G4double perpB = unitMcrossB.mag();
165 
166  static const G4double fLambdaConst = std::sqrt(3.0)*eplus/
168 
169  if( perpB > 0.0 )
170  {
171  MeanFreePath =
172  fLambdaConst*aDynamicParticle->GetDefinition()->GetPDGMass()
173  /(perpB*particleCharge*particleCharge);
174  }
175  if(verboseLevel > 0 && FirstTime)
176  {
177  G4cout << "G4SynchrotronRadiation::GetMeanFreePath "
178  << " for particle "
179  << aDynamicParticle->GetDefinition()->GetParticleName()
180  << ":" << '\n' //hbunew
181  << " MeanFreePath = " << G4BestUnit(MeanFreePath, "Length")
182  << G4endl;
183  if(verboseLevel > 1)
184  {
185  G4ThreeVector pvec = aDynamicParticle->GetMomentum();
186  G4double Btot = FieldValue.getR();
187  G4double ptot = pvec.getR();
188  G4double rho = ptot / (MeV * c_light * Btot );
189  // full bending radius
190  G4double Theta=unitMomentum.theta(FieldValue);
191  // angle between particle and field
192  G4cout << " B = " << Btot/tesla << " Tesla"
193  << " perpB = " << perpB/tesla << " Tesla"
194  << " Theta = " << Theta << " std::sin(Theta)="
195  << std::sin(Theta) << '\n'
196  << " ptot = " << G4BestUnit(ptot,"Energy")
197  << " rho = " << G4BestUnit(rho,"Length")
198  << G4endl;
199  }
200  FirstTime=false;
201  }
202  }
203  }
204  return MeanFreePath;
205 }
static constexpr double tesla
Definition: G4SIunits.hh:268
G4double condition(const G4ErrorSymMatrix &m)
double getR() const
G4int verboseLevel
Definition: G4VProcess.hh:368
CLHEP::Hep3Vector G4ThreeVector
G4double GetTotalEnergy() const
virtual void GetFieldValue(const double Point[4], double *fieldArr) const =0
double x() const
G4ParticleDefinition * GetDefinition() const
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
const G4String & GetParticleName() const
double z() const
G4GLOB_DLL std::ostream G4cout
G4double GetMass() const
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
static constexpr double eplus
Definition: G4SIunits.hh:199
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
double theta() const
G4double GetPDGMass() const
static constexpr double c_light
double y() const
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
Hep3Vector cross(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76
const G4Field * GetDetectorField() const
static constexpr double fine_structure_const
G4double GetPDGCharge() const
double mag() const
#define DBL_MAX
Definition: templates.hh:83
G4ThreeVector GetMomentum() const

Here is the call graph for this function:

G4double G4SynchrotronRadiation::GetPhotonEnergy ( const G4Track trackData,
const G4Step stepData 
)
G4double G4SynchrotronRadiation::GetRandomEnergySR ( G4double  gamma,
G4double  perpB,
G4double  mass_c2 
)

Definition at line 374 of file G4SynchrotronRadiation.cc.

376 {
377 
378  static const G4double fEnergyConst = 1.5*c_light*c_light*eplus*hbar_Planck;
379  G4double Ecr=fEnergyConst*gamma*gamma*perpB/mass_c2;
380 
381  if(verboseLevel > 0 && FirstTime1)
382  {
383  // mean and rms of photon energy
384  G4double Emean=8./(15.*std::sqrt(3.))*Ecr;
385  G4double E_rms=std::sqrt(211./675.)*Ecr;
386  G4int prec = G4cout.precision();
387  G4cout << "G4SynchrotronRadiation::GetRandomEnergySR :" << '\n'
388  << std::setprecision(4)
389  << " Ecr = " << G4BestUnit(Ecr,"Energy") << '\n'
390  << " Emean = " << G4BestUnit(Emean,"Energy") << '\n'
391  << " E_rms = " << G4BestUnit(E_rms,"Energy") << G4endl;
392  FirstTime1=false;
393  G4cout.precision(prec);
394  }
395 
396  G4double energySR=Ecr*InvSynFracInt(G4UniformRand());
397  return energySR;
398 }
G4int verboseLevel
Definition: G4VProcess.hh:368
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
static const double prec
Definition: RanecuEngine.cc:58
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
static constexpr double eplus
Definition: G4SIunits.hh:199
static constexpr double c_light
#define G4endl
Definition: G4ios.hh:61
G4double InvSynFracInt(G4double x)
double G4double
Definition: G4Types.hh:76
static constexpr double hbar_Planck

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4SynchrotronRadiation::InvSynFracInt ( G4double  x)

Definition at line 315 of file G4SynchrotronRadiation.cc.

317 {
318  // from 0 to 0.7
319  static const G4double aa1=0 ,aa2=0.7;
320  static const G4int ncheb1=27;
321  static const G4double cheb1[] =
322  { 1.22371665676046468821,0.108956475422163837267,0.0383328524358594396134,0.00759138369340257753721,
323  0.00205712048644963340914,0.000497810783280019308661,0.000130743691810302187818,0.0000338168760220395409734,
324  8.97049680900520817728e-6,2.38685472794452241466e-6,6.41923109149104165049e-7,1.73549898982749277843e-7,
325  4.72145949240790029153e-8,1.29039866111999149636e-8,3.5422080787089834182e-9,9.7594757336403784905e-10,
326  2.6979510184976065731e-10,7.480422622550977077e-11,2.079598176402699913e-11,5.79533622220841193e-12,
327  1.61856011449276096e-12,4.529450993473807e-13,1.2698603951096606e-13,3.566117394511206e-14,1.00301587494091e-14,
328  2.82515346447219e-15,7.9680747949792e-16};
329  // from 0.7 to 0.9132260271183847
330  static const G4double aa3=0.9132260271183847;
331  static const G4int ncheb2=27;
332  static const G4double cheb2[] =
333  { 1.1139496701107756,0.3523967429328067,0.0713849171926623,0.01475818043595387,0.003381255637322462,
334  0.0008228057599452224,0.00020785506681254216,0.00005390169253706556,0.000014250571923902464,3.823880733161044e-6,
335  1.0381966089136036e-6,2.8457557457837253e-7,7.86223332179956e-8,2.1866609342508474e-8,6.116186259857143e-9,
336  1.7191233618437565e-9,4.852755117740807e-10,1.3749966961763457e-10,3.908961987062447e-11,1.1146253766895824e-11,
337  3.1868887323415814e-12,9.134319791300977e-13,2.6211077371181566e-13,7.588643377757906e-14,2.1528376972619e-14,
338  6.030906040404772e-15,1.9549163926819867e-15};
339  // Chebyshev with exp/log scale
340  // a = -Log[1 - SynFracInt[1]]; b = -Log[1 - SynFracInt[7]];
341  static const G4double aa4=2.4444485538746025480,aa5=9.3830728608909477079;
342  static const G4int ncheb3=28;
343  static const G4double cheb3[] =
344  { 1.2292683840435586977,0.160353449247864455879,-0.0353559911947559448721,0.00776901561223573936985,
345  -0.00165886451971685133259,0.000335719118906954279467,-0.0000617184951079161143187,9.23534039743246708256e-6,
346  -6.06747198795168022842e-7,-3.07934045961999778094e-7,1.98818772614682367781e-7,-8.13909971567720135413e-8,
347  2.84298174969641838618e-8,-9.12829766621316063548e-9,2.77713868004820551077e-9,-8.13032767247834023165e-10,
348  2.31128525568385247392e-10,-6.41796873254200220876e-11,1.74815310473323361543e-11,-4.68653536933392363045e-12,
349  1.24016595805520752748e-12,-3.24839432979935522159e-13,8.44601465226513952994e-14,-2.18647276044246803998e-14,
350  5.65407548745690689978e-15,-1.46553625917463067508e-15,3.82059606377570462276e-16,-1.00457896653436912508e-16};
351  static const G4double aa6=33.122936966163038145;
352  static const G4int ncheb4=27;
353  static const G4double cheb4[] =
354  {1.69342658227676741765,0.0742766400841232319225,-0.019337880608635717358,0.00516065527473364110491,
355  -0.00139342012990307729473,0.000378549864052022522193,-0.000103167085583785340215,0.0000281543441271412178337,
356  -7.68409742018258198651e-6,2.09543221890204537392e-6,-5.70493140367526282946e-7,1.54961164548564906446e-7,
357  -4.19665599629607704794e-8,1.13239680054166507038e-8,-3.04223563379021441863e-9,8.13073745977562957997e-10,
358  -2.15969415476814981374e-10,5.69472105972525594811e-11,-1.48844799572430829499e-11,3.84901514438304484973e-12,
359  -9.82222575944247161834e-13,2.46468329208292208183e-13,-6.04953826265982691612e-14,1.44055805710671611984e-14,
360  -3.28200813577388740722e-15,6.96566359173765367675e-16,-1.294122794852896275e-16};
361 
362  if(x<aa2) return x*x*x*Chebyshev(aa1,aa2,cheb1,ncheb1,x);
363  else if(x<aa3) return Chebyshev(aa2,aa3,cheb2,ncheb2,x);
364  else if(x<1-0.0000841363)
365  { G4double y=-G4Log(1-x);
366  return y*Chebyshev(aa4,aa5,cheb3,ncheb3,y);
367  }
368  else
369  { G4double y=-G4Log(1-x);
370  return y*Chebyshev(aa5,aa6,cheb4,ncheb4,y);
371  }
372 }
G4double Chebyshev(G4double a, G4double b, const G4double c[], G4int n, G4double x)
int G4int
Definition: G4Types.hh:78
G4double G4Log(G4double x)
Definition: G4Log.hh:230
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

G4bool G4SynchrotronRadiation::IsApplicable ( const G4ParticleDefinition particle)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 97 of file G4SynchrotronRadiation.cc.

98 {
99  return (particle.GetPDGCharge() != 0.0 && !particle.IsShortLived());
100 }
G4double GetPDGCharge() const

Here is the call graph for this function:

G4VParticleChange * G4SynchrotronRadiation::PostStepDoIt ( const G4Track track,
const G4Step Step 
)
overridevirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 212 of file G4SynchrotronRadiation.cc.

215 {
216  aParticleChange.Initialize(trackData);
217 
218  const G4DynamicParticle* aDynamicParticle = trackData.GetDynamicParticle();
219 
220  G4double gamma = aDynamicParticle->GetTotalEnergy()/
221  (aDynamicParticle->GetDefinition()->GetPDGMass());
222 
223  G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
224  if(gamma <= 1.0e3 || 0.0 == particleCharge)
225  {
226  return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
227  }
228 
229  G4ThreeVector FieldValue;
230  const G4Field* pField = nullptr;
231 
232  G4bool fieldExertsForce = false;
233  G4FieldManager* fieldMgr =
234  fFieldPropagator->FindAndSetFieldManager(trackData.GetVolume());
235 
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 = G4ThreeVector( FieldValueVec[0],
254  FieldValueVec[1],
255  FieldValueVec[2] );
256 
257  G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
258  G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum);
259  G4double perpB = unitMcrossB.mag();
260  if(perpB > 0.0)
261  {
262  // M-C of synchrotron photon energy
263 
264  G4double energyOfSR =
265  GetRandomEnergySR(gamma,perpB,
266  aDynamicParticle->GetDefinition()->GetPDGMass());
267 
268  // check against insufficient energy
269 
270  if( energyOfSR <= 0.0 )
271  {
272  return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
273  }
274  G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
275  G4ThreeVector gammaDirection =
276  genAngle->SampleDirection(aDynamicParticle,
277  energyOfSR, 1, 0);
278 
279  G4ThreeVector gammaPolarization = FieldValue.cross(gammaDirection);
280  gammaPolarization = gammaPolarization.unit();
281 
282  // create G4DynamicParticle object for the SR photon
283 
284  G4DynamicParticle* aGamma= new G4DynamicParticle ( theGamma,
285  gammaDirection,
286  energyOfSR );
287  aGamma->SetPolarization( gammaPolarization.x(),
288  gammaPolarization.y(),
289  gammaPolarization.z() );
290 
293 
294  // Update the incident particle
295 
296  G4double newKinEnergy = kineticEnergy - energyOfSR;
297 
298  if (newKinEnergy > 0.)
299  {
300  aParticleChange.ProposeEnergy( newKinEnergy );
301  }
302  else
303  {
305  }
306  }
307  }
308  return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
309 }
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4double GetTotalEnergy() const
virtual void GetFieldValue(const double Point[4], double *fieldArr) const =0
double x() const
G4ParticleDefinition * GetDefinition() const
double z() const
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
void SetPolarization(G4double polX, G4double polY, G4double polZ)
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
virtual void Initialize(const G4Track &)
G4double GetPDGMass() const
void SetNumberOfSecondaries(G4int totSecondaries)
Hep3Vector unit() const
double y() const
void ProposeEnergy(G4double finalEnergy)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
void AddSecondary(G4Track *aSecondary)
G4double GetRandomEnergySR(G4double, G4double, G4double)
Hep3Vector cross(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76
const G4Field * GetDetectorField() const
G4double GetPDGCharge() const
double mag() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)

Here is the call graph for this function:

void G4SynchrotronRadiation::PrintInfoDefinition ( )
virtual

Definition at line 415 of file G4SynchrotronRadiation.cc.

417 {
418  G4String comments ="Incoherent Synchrotron Radiation\n";
419  G4cout << G4endl << GetProcessName() << ": " << comments
420  << " good description for long magnets at all energies"
421  << G4endl;
422 }
G4GLOB_DLL std::ostream G4cout
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

Here is the caller graph for this function:

void G4SynchrotronRadiation::SetAngularGenerator ( G4VEmAngularDistribution p)

Definition at line 88 of file G4SynchrotronRadiation.cc.

89 {
90  if(p != genAngle) {
91  delete genAngle;
92  genAngle = p;
93  }
94 }
const char * p
Definition: xmltok.h:285

Here is the caller graph for this function:


The documentation for this class was generated from the following files: