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

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

update ti head

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.8 2010/10/14 18:38:21 vnivanch Exp $
28// GEANT4 tag $Name: xrays-V09-03-05 $
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#include "G4EmProcessSubType.hh"
48
49///////////////////////////////////////////////////////////////////////
50//
51//  Constructor
52//
53
54G4SynchrotronRadiation::G4SynchrotronRadiation(const G4String& processName,
55  G4ProcessType type):G4VDiscreteProcess (processName, type),
56  theGamma (G4Gamma::Gamma() ),
57  theElectron ( G4Electron::Electron() ),
58  thePositron ( G4Positron::Positron() )
59{
60  G4TransportationManager* transportMgr = 
61    G4TransportationManager::GetTransportationManager();
62
63  fFieldPropagator = transportMgr->GetPropagatorInField();
64
65  fLambdaConst = std::sqrt(3.0)*electron_mass_c2/
66                           (2.5*fine_structure_const*eplus*c_light) ;
67  fEnergyConst = 1.5*c_light*c_light*eplus*hbar_Planck/electron_mass_c2  ;
68
69  SetProcessSubType(fSynchrotronRadiation);
70  verboseLevel=1;
71}
72
73/////////////////////////////////////////////////////////////////////////
74//
75// Destructor
76//
77
78G4SynchrotronRadiation::~G4SynchrotronRadiation()
79{}
80
81/////////////////////////////// METHODS /////////////////////////////////
82//
83//
84// Production of synchrotron X-ray photon
85// GEANT4 internal units.
86//
87
88
89G4double
90G4SynchrotronRadiation::GetMeanFreePath( const G4Track& trackData,
91                                               G4double,
92                                               G4ForceCondition* condition)
93{
94  // gives the MeanFreePath in GEANT4 internal units
95  G4double MeanFreePath;
96
97  const G4DynamicParticle* aDynamicParticle = trackData.GetDynamicParticle();
98
99  *condition = NotForced ;
100
101  G4double gamma = aDynamicParticle->GetTotalEnergy()/
102                   aDynamicParticle->GetMass();
103
104  G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
105
106  if ( gamma < 1.0e3 )  MeanFreePath = DBL_MAX;
107  else
108  {
109
110    G4ThreeVector  FieldValue;
111    const G4Field*   pField = 0;
112
113    G4FieldManager* fieldMgr=0;
114    G4bool          fieldExertsForce = false;
115
116    if( (particleCharge != 0.0) )
117    {
118      fieldMgr = fFieldPropagator->FindAndSetFieldManager( trackData.GetVolume() );
119
120      if ( fieldMgr != 0 )
121      {
122        // If the field manager has no field, there is no field !
123
124        fieldExertsForce = ( fieldMgr->GetDetectorField() != 0 );
125      }
126    }
127    if ( fieldExertsForce )
128    {
129      pField = fieldMgr->GetDetectorField() ;
130      G4ThreeVector  globPosition = trackData.GetPosition();
131
132      G4double  globPosVec[3], FieldValueVec[3];
133
134      globPosVec[0] = globPosition.x();
135      globPosVec[1] = globPosition.y();
136      globPosVec[2] = globPosition.z();
137
138      pField->GetFieldValue( globPosVec, FieldValueVec );
139
140      FieldValue = G4ThreeVector( FieldValueVec[0],
141                                   FieldValueVec[1],
142                                   FieldValueVec[2]   );
143
144
145
146      G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
147      G4ThreeVector unitMcrossB  = FieldValue.cross(unitMomentum) ;
148      G4double perpB             = unitMcrossB.mag() ;
149
150      if( perpB > 0.0 ) MeanFreePath = fLambdaConst/perpB;
151      else              MeanFreePath = DBL_MAX;
152
153      static G4bool FirstTime=true;
154      if(verboseLevel > 0 && FirstTime)
155      {
156        G4cout << "G4SynchrotronRadiation::GetMeanFreePath :" << '\n' << std::setprecision(4)
157          << "  MeanFreePath = " << G4BestUnit(MeanFreePath, "Length")
158              << G4endl;
159        if(verboseLevel > 1)
160        {
161          G4ThreeVector pvec=aDynamicParticle->GetMomentum();
162          G4double Btot=FieldValue.getR();
163          G4double ptot=pvec.getR();
164          G4double rho= ptot / (MeV * c_light * Btot ); // full bending radius
165          G4double Theta=unitMomentum.theta(FieldValue); // angle between particle and field
166          G4cout
167                << "  B = " << Btot/tesla << " Tesla"
168            << "  perpB = " << perpB/tesla << " Tesla"
169            << "  Theta = " << Theta << " std::sin(Theta)=" << std::sin(Theta) << '\n'
170            << "  ptot  = " << G4BestUnit(ptot,"Energy")
171            << "  rho   = " << G4BestUnit(rho,"Length")
172                << G4endl;
173        }
174            FirstTime=false;
175      }
176    }
177    else  MeanFreePath = DBL_MAX;
178
179
180  }
181
182  return MeanFreePath;
183}
184
185////////////////////////////////////////////////////////////////////////////////
186//
187//
188
189G4VParticleChange* 
190G4SynchrotronRadiation::PostStepDoIt(const G4Track& trackData,
191                                     const G4Step&  stepData      )
192
193{
194  aParticleChange.Initialize(trackData);
195
196  const G4DynamicParticle* aDynamicParticle=trackData.GetDynamicParticle();
197
198  G4double gamma = aDynamicParticle->GetTotalEnergy()/
199                   (aDynamicParticle->GetMass()              );
200
201  if(gamma <= 1.0e3 ) 
202  {
203    return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
204  }
205  G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
206
207  G4ThreeVector  FieldValue;
208  const G4Field*   pField = 0 ;
209
210  G4FieldManager* fieldMgr=0;
211  G4bool          fieldExertsForce = false;
212
213  if( (particleCharge != 0.0) )
214  {
215    fieldMgr = fFieldPropagator->FindAndSetFieldManager( trackData.GetVolume() );
216    if ( fieldMgr != 0 )
217    {
218      // If the field manager has no field, there is no field !
219
220      fieldExertsForce = ( fieldMgr->GetDetectorField() != 0 );
221    }
222  }
223  if ( fieldExertsForce )
224  {
225    pField = fieldMgr->GetDetectorField() ;
226    G4ThreeVector  globPosition = trackData.GetPosition() ;
227    G4double  globPosVec[3], FieldValueVec[3] ;
228    globPosVec[0] = globPosition.x() ;
229    globPosVec[1] = globPosition.y() ;
230    globPosVec[2] = globPosition.z() ;
231
232    pField->GetFieldValue( globPosVec, FieldValueVec ) ;
233    FieldValue = G4ThreeVector( FieldValueVec[0],
234                                   FieldValueVec[1],
235                                   FieldValueVec[2]   );
236
237    G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
238    G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum);
239    G4double perpB = unitMcrossB.mag() ;
240    if(perpB > 0.0)
241    {
242      // M-C of synchrotron photon energy
243
244      G4double energyOfSR = GetRandomEnergySR(gamma,perpB);
245
246      // check against insufficient energy
247
248      if( energyOfSR <= 0.0 )
249      {
250        return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
251      }
252      G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
253      G4ParticleMomentum
254      particleDirection = aDynamicParticle->GetMomentumDirection();
255
256      // M-C of its direction
257
258      G4double Teta = G4UniformRand()/gamma ;    // Very roughly
259
260      G4double Phi  = twopi * G4UniformRand() ;
261
262      G4double dirx = std::sin(Teta)*std::cos(Phi) ,
263               diry = std::sin(Teta)*std::sin(Phi) ,
264               dirz = std::cos(Teta) ;
265
266      G4ThreeVector gammaDirection ( dirx, diry, dirz);
267      gammaDirection.rotateUz(particleDirection);
268
269      // polarization of new gamma
270
271      // G4double sx =  std::cos(Teta)*std::cos(Phi);
272      // G4double sy =  std::cos(Teta)*std::sin(Phi);
273      // G4double sz = -std::sin(Teta);
274
275      G4ThreeVector gammaPolarization = FieldValue.cross(gammaDirection);
276      gammaPolarization = gammaPolarization.unit();
277
278      // (sx, sy, sz);
279      // gammaPolarization.rotateUz(particleDirection);
280
281      // create G4DynamicParticle object for the SR photon
282
283      G4DynamicParticle* aGamma= new G4DynamicParticle ( theGamma,
284                                                         gammaDirection,
285                                                         energyOfSR      );
286      aGamma->SetPolarization( gammaPolarization.x(),
287                               gammaPolarization.y(),
288                               gammaPolarization.z() );
289
290
291      aParticleChange.SetNumberOfSecondaries(1);
292      aParticleChange.AddSecondary(aGamma);
293
294      // Update the incident particle
295
296      G4double newKinEnergy = kineticEnergy - energyOfSR ;
297      aParticleChange.ProposeLocalEnergyDeposit (0.);
298
299      if (newKinEnergy > 0.)
300      {
301        aParticleChange.ProposeMomentumDirection( particleDirection );
302        aParticleChange.ProposeEnergy( newKinEnergy );
303      }
304      else
305      {
306        aParticleChange.ProposeEnergy( 0. );
307      }
308    }
309  }
310  return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
311}
312
313
314/////////////////////////////////////////////////////////////////////////////////
315//
316//
317
318G4double G4SynchrotronRadiation::InvSynFracInt(G4double x)
319// direct generation
320{
321  // from 0 to 0.7
322  const G4double aa1=0  ,aa2=0.7;
323  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  const G4double aa3=0.9132260271183847;
334  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  const G4double aa4=2.4444485538746025480,aa5=9.3830728608909477079;
345  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  const G4double aa6=33.122936966163038145;
355  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=-std::log(1-x);
369        return y*Chebyshev(aa4,aa5,cheb3,ncheb3,y);
370  }
371  else
372  { G4double y=-std::log(1-x);
373        return y*Chebyshev(aa5,aa6,cheb4,ncheb4,y);
374  }
375}
376
377G4double G4SynchrotronRadiation::GetRandomEnergySR(G4double gamma, G4double perpB)
378{
379
380  G4double Ecr=fEnergyConst*gamma*gamma*perpB;
381
382  static G4bool FirstTime=true;
383  if(verboseLevel > 0 && FirstTime)
384  { G4double Emean=8./(15.*std::sqrt(3.))*Ecr; // mean photon energy
385    G4double E_rms=std::sqrt(211./675.)*Ecr; // rms of photon energy distribution
386    G4cout << "G4SynchrotronRadiation::GetRandomEnergySR :" << '\n' << std::setprecision(4)
387        << "  Ecr   = "    << G4BestUnit(Ecr,"Energy") << '\n'
388        << "  Emean = "    << G4BestUnit(Emean,"Energy") << '\n'
389        << "  E_rms = "    << G4BestUnit(E_rms,"Energy") << G4endl;
390    FirstTime=false;
391  }
392
393  G4double energySR=Ecr*InvSynFracInt(G4UniformRand());
394  return energySR;
395}
396
397
398void G4SynchrotronRadiation::BuildPhysicsTable(const G4ParticleDefinition& part)
399{
400  if(0 < verboseLevel && &part==theElectron ) PrintInfoDefinition();
401}
402
403void G4SynchrotronRadiation::PrintInfoDefinition() // not yet called, usually called from BuildPhysicsTable
404{
405  G4String comments ="Incoherent Synchrotron Radiation\n";
406  G4cout << G4endl << GetProcessName() << ":  " << comments
407         << "      good description for long magnets at all energies" << G4endl;
408}
409
410///////////////////// end of G4SynchrotronRadiation.cc
Note: See TracBrowser for help on using the repository browser.