source: trunk/source/processes/electromagnetic/xrays/src/G4SynchrotronRadiation.cc @ 961

Last change on this file since 961 was 961, checked in by garnier, 15 years ago

update processes

File size: 15.2 KB
Line 
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: G4SynchrotronRadiation.cc,v 1.5 2006/06/29 19:56:15 gunter Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
29//
30// --------------------------------------------------------------
31//      GEANT 4 class implementation file
32//      CERN Geneva Switzerland
33//
34//      History: first implementation,
35//      21-5-98 V.Grichine
36//      28-05-01, V.Ivanchenko minor changes to provide ANSI -wall compilation
37//      04.03.05, V.Grichine: get local field interface     
38//      18-05-06 H. Burkhardt: Energy spectrum from function rather than table
39//                   
40//
41//
42//
43///////////////////////////////////////////////////////////////////////////
44
45#include "G4SynchrotronRadiation.hh"
46// #include "G4Integrator.hh"
47#include "G4UnitsTable.hh"
48
49using namespace std;
50
51///////////////////////////////////////////////////////////////////////
52//
53//  Constructor
54//
55
56G4SynchrotronRadiation::G4SynchrotronRadiation(const G4String& processName,
57  G4ProcessType type):G4VDiscreteProcess (processName, type),
58  theGamma (G4Gamma::Gamma() ),
59  theElectron ( G4Electron::Electron() ),
60  thePositron ( G4Positron::Positron() )
61{
62  G4TransportationManager* transportMgr = 
63    G4TransportationManager::GetTransportationManager();
64
65  fFieldPropagator = transportMgr->GetPropagatorInField();
66
67  fLambdaConst = sqrt(3.0)*electron_mass_c2/
68                           (2.5*fine_structure_const*eplus*c_light) ;
69  fEnergyConst = 1.5*c_light*c_light*eplus*hbar_Planck/electron_mass_c2  ;
70  verboseLevel=1;
71}
72
73/////////////////////////////////////////////////////////////////////////
74//
75// Destructor
76//
77
78G4SynchrotronRadiation::~G4SynchrotronRadiation()
79{
80     ;
81}
82
83/////////////////////////////// METHODS /////////////////////////////////
84//
85//
86// Production of synchrotron X-ray photon
87// GEANT4 internal units.
88//
89
90
91G4double
92G4SynchrotronRadiation::GetMeanFreePath( const G4Track& trackData,
93                                               G4double,
94                                               G4ForceCondition* condition)
95{
96  // gives the MeanFreePath in GEANT4 internal units
97  G4double MeanFreePath;
98
99  const G4DynamicParticle* aDynamicParticle = trackData.GetDynamicParticle();
100
101  *condition = NotForced ;
102
103  G4double gamma = aDynamicParticle->GetTotalEnergy()/
104                   aDynamicParticle->GetMass();
105
106  G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
107
108  if ( gamma < 1.0e3 )  MeanFreePath = DBL_MAX;
109  else
110  {
111
112    G4ThreeVector  FieldValue;
113    const G4Field*   pField = 0;
114
115    G4FieldManager* fieldMgr=0;
116    G4bool          fieldExertsForce = false;
117
118    if( (particleCharge != 0.0) )
119    {
120      fieldMgr = fFieldPropagator->FindAndSetFieldManager( trackData.GetVolume() );
121
122      if ( fieldMgr != 0 )
123      {
124        // If the field manager has no field, there is no field !
125
126        fieldExertsForce = ( fieldMgr->GetDetectorField() != 0 );
127      }
128    }
129    if ( fieldExertsForce )
130    {
131      pField = fieldMgr->GetDetectorField() ;
132      G4ThreeVector  globPosition = trackData.GetPosition();
133
134      G4double  globPosVec[3], FieldValueVec[3];
135
136      globPosVec[0] = globPosition.x();
137      globPosVec[1] = globPosition.y();
138      globPosVec[2] = globPosition.z();
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' << std::setprecision(4)
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 << " sin(Theta)=" << 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
191G4VParticleChange* 
192G4SynchrotronRadiation::PostStepDoIt(const G4Track& trackData,
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[3], FieldValueVec[3] ;
230    globPosVec[0] = globPosition.x() ;
231    globPosVec[1] = globPosition.y() ;
232    globPosVec[2] = globPosition.z() ;
233
234    pField->GetFieldValue( globPosVec, FieldValueVec ) ;
235    FieldValue = G4ThreeVector( FieldValueVec[0],
236                                   FieldValueVec[1],
237                                   FieldValueVec[2]   );
238
239    G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
240    G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum);
241    G4double perpB = unitMcrossB.mag() ;
242    if(perpB > 0.0)
243    {
244      // M-C of synchrotron photon energy
245
246      G4double energyOfSR = GetRandomEnergySR(gamma,perpB);
247
248      // check against insufficient energy
249
250      if( energyOfSR <= 0.0 )
251      {
252        return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
253      }
254      G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
255      G4ParticleMomentum
256      particleDirection = aDynamicParticle->GetMomentumDirection();
257
258      // M-C of its direction
259
260      G4double Teta = G4UniformRand()/gamma ;    // Very roughly
261
262      G4double Phi  = twopi * G4UniformRand() ;
263
264      G4double dirx = sin(Teta)*cos(Phi) ,
265               diry = sin(Teta)*sin(Phi) ,
266               dirz = cos(Teta) ;
267
268      G4ThreeVector gammaDirection ( dirx, diry, dirz);
269      gammaDirection.rotateUz(particleDirection);
270
271      // polarization of new gamma
272
273      // G4double sx =  cos(Teta)*cos(Phi);
274      // G4double sy =  cos(Teta)*sin(Phi);
275      // G4double sz = -sin(Teta);
276
277      G4ThreeVector gammaPolarization = FieldValue.cross(gammaDirection);
278      gammaPolarization = gammaPolarization.unit();
279
280      // (sx, sy, sz);
281      // gammaPolarization.rotateUz(particleDirection);
282
283      // create G4DynamicParticle object for the SR photon
284
285      G4DynamicParticle* aGamma= new G4DynamicParticle ( theGamma,
286                                                         gammaDirection,
287                                                         energyOfSR      );
288      aGamma->SetPolarization( gammaPolarization.x(),
289                               gammaPolarization.y(),
290                               gammaPolarization.z() );
291
292
293      aParticleChange.SetNumberOfSecondaries(1);
294      aParticleChange.AddSecondary(aGamma);
295
296      // Update the incident particle
297
298      G4double newKinEnergy = kineticEnergy - energyOfSR ;
299      aParticleChange.ProposeLocalEnergyDeposit (0.);
300
301      if (newKinEnergy > 0.)
302      {
303        aParticleChange.ProposeMomentumDirection( particleDirection );
304        aParticleChange.ProposeEnergy( newKinEnergy );
305      }
306      else
307      {
308        aParticleChange.ProposeEnergy( 0. );
309      }
310    }
311  }
312  return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
313}
314
315
316/////////////////////////////////////////////////////////////////////////////////
317//
318//
319
320G4double G4SynchrotronRadiation::InvSynFracInt(G4double x)
321// direct generation
322{
323  // from 0 to 0.7
324  const G4double aa1=0  ,aa2=0.7;
325  const G4int ncheb1=27;
326  static const G4double cheb1[] =
327  { 1.22371665676046468821,0.108956475422163837267,0.0383328524358594396134,0.00759138369340257753721,
328   0.00205712048644963340914,0.000497810783280019308661,0.000130743691810302187818,0.0000338168760220395409734,
329   8.97049680900520817728e-6,2.38685472794452241466e-6,6.41923109149104165049e-7,1.73549898982749277843e-7,
330   4.72145949240790029153e-8,1.29039866111999149636e-8,3.5422080787089834182e-9,9.7594757336403784905e-10,
331   2.6979510184976065731e-10,7.480422622550977077e-11,2.079598176402699913e-11,5.79533622220841193e-12,
332   1.61856011449276096e-12,4.529450993473807e-13,1.2698603951096606e-13,3.566117394511206e-14,1.00301587494091e-14,
333   2.82515346447219e-15,7.9680747949792e-16};
334  //   from 0.7 to 0.9132260271183847
335  const G4double aa3=0.9132260271183847;
336  const G4int ncheb2=27;
337  static const G4double cheb2[] =
338  { 1.1139496701107756,0.3523967429328067,0.0713849171926623,0.01475818043595387,0.003381255637322462,
339        0.0008228057599452224,0.00020785506681254216,0.00005390169253706556,0.000014250571923902464,3.823880733161044e-6,
340        1.0381966089136036e-6,2.8457557457837253e-7,7.86223332179956e-8,2.1866609342508474e-8,6.116186259857143e-9,
341        1.7191233618437565e-9,4.852755117740807e-10,1.3749966961763457e-10,3.908961987062447e-11,1.1146253766895824e-11,
342        3.1868887323415814e-12,9.134319791300977e-13,2.6211077371181566e-13,7.588643377757906e-14,2.1528376972619e-14,
343        6.030906040404772e-15,1.9549163926819867e-15};
344   // Chebyshev with exp/log  scale
345   // a = -Log[1 - SynFracInt[1]]; b = -Log[1 - SynFracInt[7]];
346  const G4double aa4=2.4444485538746025480,aa5=9.3830728608909477079;
347  const G4int ncheb3=28;
348  static const G4double cheb3[] =
349  { 1.2292683840435586977,0.160353449247864455879,-0.0353559911947559448721,0.00776901561223573936985,
350   -0.00165886451971685133259,0.000335719118906954279467,-0.0000617184951079161143187,9.23534039743246708256e-6,
351   -6.06747198795168022842e-7,-3.07934045961999778094e-7,1.98818772614682367781e-7,-8.13909971567720135413e-8,
352   2.84298174969641838618e-8,-9.12829766621316063548e-9,2.77713868004820551077e-9,-8.13032767247834023165e-10,
353   2.31128525568385247392e-10,-6.41796873254200220876e-11,1.74815310473323361543e-11,-4.68653536933392363045e-12,
354   1.24016595805520752748e-12,-3.24839432979935522159e-13,8.44601465226513952994e-14,-2.18647276044246803998e-14,
355   5.65407548745690689978e-15,-1.46553625917463067508e-15,3.82059606377570462276e-16,-1.00457896653436912508e-16};
356  const G4double aa6=33.122936966163038145;
357  const G4int ncheb4=27;
358  static const G4double cheb4[] =
359  {1.69342658227676741765,0.0742766400841232319225,-0.019337880608635717358,0.00516065527473364110491,
360   -0.00139342012990307729473,0.000378549864052022522193,-0.000103167085583785340215,0.0000281543441271412178337,
361   -7.68409742018258198651e-6,2.09543221890204537392e-6,-5.70493140367526282946e-7,1.54961164548564906446e-7,
362   -4.19665599629607704794e-8,1.13239680054166507038e-8,-3.04223563379021441863e-9,8.13073745977562957997e-10,
363   -2.15969415476814981374e-10,5.69472105972525594811e-11,-1.48844799572430829499e-11,3.84901514438304484973e-12,
364   -9.82222575944247161834e-13,2.46468329208292208183e-13,-6.04953826265982691612e-14,1.44055805710671611984e-14,
365   -3.28200813577388740722e-15,6.96566359173765367675e-16,-1.294122794852896275e-16};
366
367  if(x<aa2)      return x*x*x*Chebyshev(aa1,aa2,cheb1,ncheb1,x);
368  else if(x<aa3) return       Chebyshev(aa2,aa3,cheb2,ncheb2,x);
369  else if(x<1-0.0000841363)
370  { G4double y=-log(1-x);
371        return y*Chebyshev(aa4,aa5,cheb3,ncheb3,y);
372  }
373  else
374  { G4double y=-log(1-x);
375        return y*Chebyshev(aa5,aa6,cheb4,ncheb4,y);
376  }
377}
378
379G4double G4SynchrotronRadiation::GetRandomEnergySR(G4double gamma, G4double perpB)
380{
381
382  G4double Ecr=fEnergyConst*gamma*gamma*perpB;
383
384  static G4bool FirstTime=true;
385  if(verboseLevel > 0 && FirstTime)
386  { G4double Emean=8./(15.*sqrt(3.))*Ecr; // mean photon energy
387    G4double E_rms=sqrt(211./675.)*Ecr; // rms of photon energy distribution
388    G4cout << "G4SynchrotronRadiation::GetRandomEnergySR :" << '\n' << std::setprecision(4)
389        << "  Ecr   = "    << G4BestUnit(Ecr,"Energy") << '\n'
390        << "  Emean = "    << G4BestUnit(Emean,"Energy") << '\n'
391        << "  E_rms = "    << G4BestUnit(E_rms,"Energy") << G4endl;
392    FirstTime=false;
393  }
394
395  G4double energySR=Ecr*InvSynFracInt(G4UniformRand());
396  return energySR;
397}
398
399
400void G4SynchrotronRadiation::BuildPhysicsTable(const G4ParticleDefinition& part)
401{
402  if(0 < verboseLevel && &part==theElectron ) PrintInfoDefinition();
403}
404
405void G4SynchrotronRadiation::PrintInfoDefinition() // not yet called, usually called from BuildPhysicsTable
406{
407  G4String comments ="Incoherent Synchrotron Radiation\n";
408  G4cout << G4endl << GetProcessName() << ":  " << comments
409         << "      good description for long magnets at all energies" << G4endl;
410}
411
412///////////////////// end of G4SynchrotronRadiation.cc
Note: See TracBrowser for help on using the repository browser.