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

Last change on this file since 1337 was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

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