Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4SynchrotronRadiation.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$
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// 18-05-06 H. Burkhardt: Energy spectrum from function rather than table
38//
39//
40//
41//
42///////////////////////////////////////////////////////////////////////////
43
46#include "G4SystemOfUnits.hh"
47#include "G4UnitsTable.hh"
48#include "G4EmProcessSubType.hh"
49
50///////////////////////////////////////////////////////////////////////
51//
52// Constructor
53//
54
56 G4ProcessType type):G4VDiscreteProcess (processName, type),
57 theGamma (G4Gamma::Gamma() ),
58 theElectron ( G4Electron::Electron() ),
59 thePositron ( G4Positron::Positron() )
60{
61 G4TransportationManager* transportMgr =
63
64 fFieldPropagator = transportMgr->GetPropagatorInField();
65
66 fLambdaConst = std::sqrt(3.0)*electron_mass_c2/
67 (2.5*fine_structure_const*eplus*c_light);
68 fEnergyConst = 1.5*c_light*c_light*eplus*hbar_Planck/electron_mass_c2 ;
69
72}
73
74/////////////////////////////////////////////////////////////////////////
75//
76// Destructor
77//
78
80{}
81
82/////////////////////////////// METHODS /////////////////////////////////
83//
84//
85// Production of synchrotron X-ray photon
86// GEANT4 internal units.
87//
88
89
94{
95 // gives the MeanFreePath in GEANT4 internal units
96 G4double MeanFreePath;
97
98 const G4DynamicParticle* aDynamicParticle = trackData.GetDynamicParticle();
99
101
102 G4double gamma = aDynamicParticle->GetTotalEnergy()/
103 aDynamicParticle->GetMass();
104
105 G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
106
107 if ( gamma < 1.0e3 ) MeanFreePath = DBL_MAX;
108 else
109 {
110
111 G4ThreeVector FieldValue;
112 const G4Field* pField = 0;
113
114 G4FieldManager* fieldMgr=0;
115 G4bool fieldExertsForce = false;
116
117 if( (particleCharge != 0.0) )
118 {
119 fieldMgr = fFieldPropagator->FindAndSetFieldManager( trackData.GetVolume() );
120
121 if ( fieldMgr != 0 )
122 {
123 // If the field manager has no field, there is no field !
124
125 fieldExertsForce = ( fieldMgr->GetDetectorField() != 0 );
126 }
127 }
128 if ( fieldExertsForce )
129 {
130 pField = fieldMgr->GetDetectorField();
131 G4ThreeVector globPosition = trackData.GetPosition();
132
133 G4double globPosVec[4], FieldValueVec[6];
134
135 globPosVec[0] = globPosition.x();
136 globPosVec[1] = globPosition.y();
137 globPosVec[2] = globPosition.z();
138 globPosVec[3] = trackData.GetGlobalTime();
139
140 pField->GetFieldValue( globPosVec, FieldValueVec );
141
142 FieldValue = G4ThreeVector( FieldValueVec[0],
143 FieldValueVec[1],
144 FieldValueVec[2] );
145
146
147
148 G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
149 G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum);
150 G4double perpB = unitMcrossB.mag();
151
152 if( perpB > 0.0 ) MeanFreePath = fLambdaConst/perpB;
153 else MeanFreePath = DBL_MAX;
154
155 static G4bool FirstTime=true;
156 if(verboseLevel > 0 && FirstTime)
157 {
158 G4cout << "G4SynchrotronRadiation::GetMeanFreePath :" << '\n'
159 << " MeanFreePath = " << G4BestUnit(MeanFreePath, "Length")
160 << G4endl;
161 if(verboseLevel > 1)
162 {
163 G4ThreeVector pvec=aDynamicParticle->GetMomentum();
164 G4double Btot=FieldValue.getR();
165 G4double ptot=pvec.getR();
166 G4double rho= ptot / (MeV * c_light * Btot ); // full bending radius
167 G4double Theta=unitMomentum.theta(FieldValue); // angle between particle and field
168 G4cout
169 << " B = " << Btot/tesla << " Tesla"
170 << " perpB = " << perpB/tesla << " Tesla"
171 << " Theta = " << Theta << " std::sin(Theta)=" << std::sin(Theta) << '\n'
172 << " ptot = " << G4BestUnit(ptot,"Energy")
173 << " rho = " << G4BestUnit(rho,"Length")
174 << G4endl;
175 }
176 FirstTime=false;
177 }
178 }
179 else MeanFreePath = DBL_MAX;
180
181
182 }
183
184 return MeanFreePath;
185}
186
187////////////////////////////////////////////////////////////////////////////////
188//
189//
190
193 const G4Step& stepData )
194
195{
196 aParticleChange.Initialize(trackData);
197
198 const G4DynamicParticle* aDynamicParticle=trackData.GetDynamicParticle();
199
200 G4double gamma = aDynamicParticle->GetTotalEnergy()/
201 (aDynamicParticle->GetMass() );
202
203 if(gamma <= 1.0e3 )
204 {
205 return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
206 }
207 G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
208
209 G4ThreeVector FieldValue;
210 const G4Field* pField = 0;
211
212 G4FieldManager* fieldMgr=0;
213 G4bool fieldExertsForce = false;
214
215 if( (particleCharge != 0.0) )
216 {
217 fieldMgr = fFieldPropagator->FindAndSetFieldManager( trackData.GetVolume() );
218 if ( fieldMgr != 0 )
219 {
220 // If the field manager has no field, there is no field !
221
222 fieldExertsForce = ( fieldMgr->GetDetectorField() != 0 );
223 }
224 }
225 if ( fieldExertsForce )
226 {
227 pField = fieldMgr->GetDetectorField();
228 G4ThreeVector globPosition = trackData.GetPosition();
229 G4double globPosVec[4], FieldValueVec[6];
230 globPosVec[0] = globPosition.x();
231 globPosVec[1] = globPosition.y();
232 globPosVec[2] = globPosition.z();
233 globPosVec[3] = trackData.GetGlobalTime();
234
235 pField->GetFieldValue( globPosVec, FieldValueVec );
236 FieldValue = G4ThreeVector( FieldValueVec[0],
237 FieldValueVec[1],
238 FieldValueVec[2] );
239
240 G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
241 G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum);
242 G4double perpB = unitMcrossB.mag();
243 if(perpB > 0.0)
244 {
245 // M-C of synchrotron photon energy
246
247 G4double energyOfSR = GetRandomEnergySR(gamma,perpB);
248
249 // check against insufficient energy
250
251 if( energyOfSR <= 0.0 )
252 {
253 return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
254 }
255 G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
257 particleDirection = aDynamicParticle->GetMomentumDirection();
258
259 // M-C of its direction, simplified dipole boosted approach
260
261 // G4double Teta, fteta; // = G4UniformRand()/gamma; // Very roughly
262
263 G4double cosTheta, sinTheta, fcos, beta;
264
265 do
266 {
267 cosTheta = 1. - 2.*G4UniformRand();
268 fcos = (1 + cosTheta*cosTheta)*0.5;
269 }
270 while( fcos < G4UniformRand() );
271
272 beta = std::sqrt(1. - 1./(gamma*gamma));
273
274 cosTheta = (cosTheta + beta)/(1. + beta*cosTheta);
275
276 if( cosTheta > 1. ) cosTheta = 1.;
277 if( cosTheta < -1. ) cosTheta = -1.;
278
279 sinTheta = std::sqrt(1. - cosTheta*cosTheta );
280
281 G4double Phi = twopi * G4UniformRand();
282
283 G4double dirx = sinTheta*std::cos(Phi) ,
284 diry = sinTheta*std::sin(Phi) ,
285 dirz = cosTheta;
286
287 G4ThreeVector gammaDirection ( dirx, diry, dirz);
288 gammaDirection.rotateUz(particleDirection);
289
290 // polarization of new gamma
291
292 // G4double sx = std::cos(Teta)*std::cos(Phi);
293 // G4double sy = std::cos(Teta)*std::sin(Phi);
294 // G4double sz = -std::sin(Teta);
295
296 G4ThreeVector gammaPolarization = FieldValue.cross(gammaDirection);
297 gammaPolarization = gammaPolarization.unit();
298
299 // (sx, sy, sz);
300 // gammaPolarization.rotateUz(particleDirection);
301
302 // create G4DynamicParticle object for the SR photon
303
304 G4DynamicParticle* aGamma= new G4DynamicParticle ( theGamma,
305 gammaDirection,
306 energyOfSR );
307 aGamma->SetPolarization( gammaPolarization.x(),
308 gammaPolarization.y(),
309 gammaPolarization.z() );
310
311
314
315 // Update the incident particle
316
317 G4double newKinEnergy = kineticEnergy - energyOfSR;
319
320 if (newKinEnergy > 0.)
321 {
322 aParticleChange.ProposeMomentumDirection( particleDirection );
323 aParticleChange.ProposeEnergy( newKinEnergy );
324 }
325 else
326 {
328 }
329 }
330 }
331 return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
332}
333
334
335/////////////////////////////////////////////////////////////////////////////////
336//
337//
338
340// direct generation
341{
342 // from 0 to 0.7
343 const G4double aa1=0 ,aa2=0.7;
344 const G4int ncheb1=27;
345 static const G4double cheb1[] =
346 { 1.22371665676046468821,0.108956475422163837267,0.0383328524358594396134,0.00759138369340257753721,
347 0.00205712048644963340914,0.000497810783280019308661,0.000130743691810302187818,0.0000338168760220395409734,
348 8.97049680900520817728e-6,2.38685472794452241466e-6,6.41923109149104165049e-7,1.73549898982749277843e-7,
349 4.72145949240790029153e-8,1.29039866111999149636e-8,3.5422080787089834182e-9,9.7594757336403784905e-10,
350 2.6979510184976065731e-10,7.480422622550977077e-11,2.079598176402699913e-11,5.79533622220841193e-12,
351 1.61856011449276096e-12,4.529450993473807e-13,1.2698603951096606e-13,3.566117394511206e-14,1.00301587494091e-14,
352 2.82515346447219e-15,7.9680747949792e-16};
353 // from 0.7 to 0.9132260271183847
354 const G4double aa3=0.9132260271183847;
355 const G4int ncheb2=27;
356 static const G4double cheb2[] =
357 { 1.1139496701107756,0.3523967429328067,0.0713849171926623,0.01475818043595387,0.003381255637322462,
358 0.0008228057599452224,0.00020785506681254216,0.00005390169253706556,0.000014250571923902464,3.823880733161044e-6,
359 1.0381966089136036e-6,2.8457557457837253e-7,7.86223332179956e-8,2.1866609342508474e-8,6.116186259857143e-9,
360 1.7191233618437565e-9,4.852755117740807e-10,1.3749966961763457e-10,3.908961987062447e-11,1.1146253766895824e-11,
361 3.1868887323415814e-12,9.134319791300977e-13,2.6211077371181566e-13,7.588643377757906e-14,2.1528376972619e-14,
362 6.030906040404772e-15,1.9549163926819867e-15};
363 // Chebyshev with exp/log scale
364 // a = -Log[1 - SynFracInt[1]]; b = -Log[1 - SynFracInt[7]];
365 const G4double aa4=2.4444485538746025480,aa5=9.3830728608909477079;
366 const G4int ncheb3=28;
367 static const G4double cheb3[] =
368 { 1.2292683840435586977,0.160353449247864455879,-0.0353559911947559448721,0.00776901561223573936985,
369 -0.00165886451971685133259,0.000335719118906954279467,-0.0000617184951079161143187,9.23534039743246708256e-6,
370 -6.06747198795168022842e-7,-3.07934045961999778094e-7,1.98818772614682367781e-7,-8.13909971567720135413e-8,
371 2.84298174969641838618e-8,-9.12829766621316063548e-9,2.77713868004820551077e-9,-8.13032767247834023165e-10,
372 2.31128525568385247392e-10,-6.41796873254200220876e-11,1.74815310473323361543e-11,-4.68653536933392363045e-12,
373 1.24016595805520752748e-12,-3.24839432979935522159e-13,8.44601465226513952994e-14,-2.18647276044246803998e-14,
374 5.65407548745690689978e-15,-1.46553625917463067508e-15,3.82059606377570462276e-16,-1.00457896653436912508e-16};
375 const G4double aa6=33.122936966163038145;
376 const G4int ncheb4=27;
377 static const G4double cheb4[] =
378 {1.69342658227676741765,0.0742766400841232319225,-0.019337880608635717358,0.00516065527473364110491,
379 -0.00139342012990307729473,0.000378549864052022522193,-0.000103167085583785340215,0.0000281543441271412178337,
380 -7.68409742018258198651e-6,2.09543221890204537392e-6,-5.70493140367526282946e-7,1.54961164548564906446e-7,
381 -4.19665599629607704794e-8,1.13239680054166507038e-8,-3.04223563379021441863e-9,8.13073745977562957997e-10,
382 -2.15969415476814981374e-10,5.69472105972525594811e-11,-1.48844799572430829499e-11,3.84901514438304484973e-12,
383 -9.82222575944247161834e-13,2.46468329208292208183e-13,-6.04953826265982691612e-14,1.44055805710671611984e-14,
384 -3.28200813577388740722e-15,6.96566359173765367675e-16,-1.294122794852896275e-16};
385
386 if(x<aa2) return x*x*x*Chebyshev(aa1,aa2,cheb1,ncheb1,x);
387 else if(x<aa3) return Chebyshev(aa2,aa3,cheb2,ncheb2,x);
388 else if(x<1-0.0000841363)
389 { G4double y=-std::log(1-x);
390 return y*Chebyshev(aa4,aa5,cheb3,ncheb3,y);
391 }
392 else
393 { G4double y=-std::log(1-x);
394 return y*Chebyshev(aa5,aa6,cheb4,ncheb4,y);
395 }
396}
397
399{
400
401 G4double Ecr=fEnergyConst*gamma*gamma*perpB;
402
403 static G4bool FirstTime=true;
404 if(verboseLevel > 0 && FirstTime)
405 { G4double Emean=8./(15.*std::sqrt(3.))*Ecr; // mean photon energy
406 G4double E_rms=std::sqrt(211./675.)*Ecr; // rms of photon energy distribution
407 G4int prec = G4cout.precision();
408 G4cout << "G4SynchrotronRadiation::GetRandomEnergySR :" << '\n' << std::setprecision(4)
409 << " Ecr = " << G4BestUnit(Ecr,"Energy") << '\n'
410 << " Emean = " << G4BestUnit(Emean,"Energy") << '\n'
411 << " E_rms = " << G4BestUnit(E_rms,"Energy") << G4endl;
412 FirstTime=false;
413 G4cout.precision(prec);
414 }
415
416 G4double energySR=Ecr*InvSynFracInt(G4UniformRand());
417 return energySR;
418}
419
420
422{
423 if(0 < verboseLevel && &part==theElectron ) PrintInfoDefinition();
424}
425
426void G4SynchrotronRadiation::PrintInfoDefinition() // not yet called, usually called from BuildPhysicsTable
427{
428 G4String comments ="Incoherent Synchrotron Radiation\n";
429 G4cout << G4endl << GetProcessName() << ": " << comments
430 << " good description for long magnets at all energies" << G4endl;
431}
432
433///////////////////// end of G4SynchrotronRadiation.cc
@ fSynchrotronRadiation
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ NotForced
G4ProcessType
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4double fcos(G4double arg)
#define G4UniformRand()
Definition: Randomize.hh:53
double z() const
Hep3Vector unit() const
double theta() const
double x() const
double getR() const
double y() const
Hep3Vector cross(const Hep3Vector &) const
double mag() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
G4double GetMass() const
const G4ThreeVector & GetMomentumDirection() const
void SetPolarization(G4double polX, G4double polY, G4double polZ)
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
G4ThreeVector GetMomentum() const
const G4Field * GetDetectorField() const
virtual void GetFieldValue(const double Point[4], double *fieldArr) const =0
void AddSecondary(G4Track *aSecondary)
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
virtual void Initialize(const G4Track &)
G4double GetPDGCharge() const
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
Definition: G4Step.hh:78
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &Step)
void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetRandomEnergySR(G4double, G4double)
G4SynchrotronRadiation(const G4String &pName="SynRad", G4ProcessType type=fElectromagnetic)
G4double Chebyshev(G4double a, G4double b, const G4double c[], G4int n, G4double x)
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4double InvSynFracInt(G4double x)
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
const G4DynamicParticle * GetDynamicParticle() const
static G4TransportationManager * GetTransportationManager()
G4PropagatorInField * GetPropagatorInField() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
G4int verboseLevel
Definition: G4VProcess.hh:368
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379
#define DBL_MAX
Definition: templates.hh:83