Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4SynchrotronRadiation Class Reference

#include <G4SynchrotronRadiation.hh>

+ Inheritance 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 &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
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 ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
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 G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
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
 
virtual void ProcessDescription (std::ostream &outfile) 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)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 64 of file G4SynchrotronRadiation.hh.

Constructor & Destructor Documentation

◆ G4SynchrotronRadiation()

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 theManager = G4LossTableManager::Instance();
73 theManager->Register(this);
74}
@ fSynchrotronRadiation
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4LossTableManager * Instance()
void Register(G4VEnergyLossProcess *p)
void SetAngularGenerator(G4VEmAngularDistribution *p)
static G4TransportationManager * GetTransportationManager()
G4PropagatorInField * GetPropagatorInField() const
G4int verboseLevel
Definition: G4VProcess.hh:356
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406

◆ ~G4SynchrotronRadiation()

G4SynchrotronRadiation::~G4SynchrotronRadiation ( )
virtual

Definition at line 81 of file G4SynchrotronRadiation.cc.

82{
83 delete genAngle;
84 theManager->DeRegister(this);
85}
void DeRegister(G4VEnergyLossProcess *p)

Member Function Documentation

◆ BuildPhysicsTable()

void G4SynchrotronRadiation::BuildPhysicsTable ( const G4ParticleDefinition part)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 408 of file G4SynchrotronRadiation.cc.

409{
411 // same for all particles, print only for one (electron)
412}
static G4Electron * Electron()
Definition: G4Electron.cc:93

◆ Chebyshev()

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

Definition at line 114 of file G4SynchrotronRadiation.hh.

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

Referenced by InvSynFracInt().

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 112 of file G4SynchrotronRadiation.cc.

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

◆ GetPhotonEnergy()

G4double G4SynchrotronRadiation::GetPhotonEnergy ( const G4Track trackData,
const G4Step stepData 
)

◆ GetRandomEnergySR()

G4double G4SynchrotronRadiation::GetRandomEnergySR ( G4double  gamma,
G4double  perpB,
G4double  mass_c2 
)

Definition at line 377 of file G4SynchrotronRadiation.cc.

379{
380
381 static const G4double fEnergyConst = 1.5*c_light*c_light*eplus*hbar_Planck;
382 G4double Ecr=fEnergyConst*gamma*gamma*perpB/mass_c2;
383
384 if(verboseLevel > 0 && FirstTime1)
385 {
386 // mean and rms of photon energy
387 G4double Emean=8./(15.*std::sqrt(3.))*Ecr;
388 G4double E_rms=std::sqrt(211./675.)*Ecr;
389 G4int prec = G4cout.precision();
390 G4cout << "G4SynchrotronRadiation::GetRandomEnergySR :" << '\n'
391 << std::setprecision(4)
392 << " Ecr = " << G4BestUnit(Ecr,"Energy") << '\n'
393 << " Emean = " << G4BestUnit(Emean,"Energy") << '\n'
394 << " E_rms = " << G4BestUnit(E_rms,"Energy") << G4endl;
395 FirstTime1=false;
396 G4cout.precision(prec);
397 }
398
399 G4double energySR=Ecr*InvSynFracInt(G4UniformRand());
400 return energySR;
401}
#define G4UniformRand()
Definition: Randomize.hh:52
G4double InvSynFracInt(G4double x)

Referenced by PostStepDoIt().

◆ InvSynFracInt()

G4double G4SynchrotronRadiation::InvSynFracInt ( G4double  x)

Definition at line 318 of file G4SynchrotronRadiation.cc.

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

Referenced by GetRandomEnergySR().

◆ IsApplicable()

G4bool G4SynchrotronRadiation::IsApplicable ( const G4ParticleDefinition particle)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 100 of file G4SynchrotronRadiation.cc.

101{
102 return (particle.GetPDGCharge() != 0.0 && !particle.IsShortLived());
103}

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 215 of file G4SynchrotronRadiation.cc.

218{
219 aParticleChange.Initialize(trackData);
220
221 const G4DynamicParticle* aDynamicParticle = trackData.GetDynamicParticle();
222
223 G4double gamma = aDynamicParticle->GetTotalEnergy()/
224 (aDynamicParticle->GetDefinition()->GetPDGMass());
225
226 G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
227 if(gamma <= 1.0e3 || 0.0 == particleCharge)
228 {
229 return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
230 }
231
232 G4ThreeVector FieldValue;
233 const G4Field* pField = nullptr;
234
235 G4bool fieldExertsForce = false;
236 G4FieldManager* fieldMgr =
237 fFieldPropagator->FindAndSetFieldManager(trackData.GetVolume());
238
239 if ( fieldMgr != nullptr )
240 {
241 // If the field manager has no field, there is no field !
242 fieldExertsForce = ( fieldMgr->GetDetectorField() != nullptr );
243 }
244
245 if ( fieldExertsForce )
246 {
247 pField = fieldMgr->GetDetectorField();
248 G4ThreeVector globPosition = trackData.GetPosition();
249 G4double globPosVec[4], FieldValueVec[6];
250 globPosVec[0] = globPosition.x();
251 globPosVec[1] = globPosition.y();
252 globPosVec[2] = globPosition.z();
253 globPosVec[3] = trackData.GetGlobalTime();
254
255 pField->GetFieldValue( globPosVec, FieldValueVec );
256 FieldValue = G4ThreeVector( FieldValueVec[0],
257 FieldValueVec[1],
258 FieldValueVec[2] );
259
260 G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
261 G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum);
262 G4double perpB = unitMcrossB.mag();
263 if(perpB > 0.0)
264 {
265 // M-C of synchrotron photon energy
266
267 G4double energyOfSR =
268 GetRandomEnergySR(gamma,perpB,
269 aDynamicParticle->GetDefinition()->GetPDGMass());
270
271 // check against insufficient energy
272
273 if( energyOfSR <= 0.0 )
274 {
275 return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
276 }
277 G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
278 G4ThreeVector gammaDirection =
279 genAngle->SampleDirection(aDynamicParticle,
280 energyOfSR, 1, 0);
281
282 G4ThreeVector gammaPolarization = FieldValue.cross(gammaDirection);
283 gammaPolarization = gammaPolarization.unit();
284
285 // create G4DynamicParticle object for the SR photon
286
287 G4DynamicParticle* aGamma= new G4DynamicParticle ( theGamma,
288 gammaDirection,
289 energyOfSR );
290 aGamma->SetPolarization( gammaPolarization.x(),
291 gammaPolarization.y(),
292 gammaPolarization.z() );
293
296
297 // Update the incident particle
298
299 G4double newKinEnergy = kineticEnergy - energyOfSR;
300
301 if (newKinEnergy > 0.)
302 {
303 aParticleChange.ProposeEnergy( newKinEnergy );
304 }
305 else
306 {
308 }
309 }
310 }
311 return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
312}
Hep3Vector unit() const
void SetPolarization(const G4ThreeVector &)
G4double GetKineticEnergy() const
void AddSecondary(G4Track *aSecondary)
void ProposeEnergy(G4double finalEnergy)
virtual void Initialize(const G4Track &)
G4double GetRandomEnergySR(G4double, G4double, G4double)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:327

◆ PrintInfoDefinition()

void G4SynchrotronRadiation::PrintInfoDefinition ( )
virtual

Definition at line 418 of file G4SynchrotronRadiation.cc.

420{
421 G4String comments ="Incoherent Synchrotron Radiation\n";
422 G4cout << G4endl << GetProcessName() << ": " << comments
423 << " good description for long magnets at all energies"
424 << G4endl;
425}
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382

Referenced by BuildPhysicsTable().

◆ SetAngularGenerator()

void G4SynchrotronRadiation::SetAngularGenerator ( G4VEmAngularDistribution p)

Definition at line 91 of file G4SynchrotronRadiation.cc.

92{
93 if(p != genAngle) {
94 delete genAngle;
95 genAngle = p;
96 }
97}

Referenced by G4SynchrotronRadiation().


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