Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4SynchrotronRadiation.cc
이 파일의 문서화 페이지로 가기
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 97385 2016-06-02 09:59:53Z gcosmo $
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 //
43 
45 #include "G4PhysicalConstants.hh"
46 #include "G4SystemOfUnits.hh"
47 #include "G4UnitsTable.hh"
48 #include "G4EmProcessSubType.hh"
49 #include "G4DipBustGenerator.hh"
50 #include "G4Log.hh"
51 
53 //
54 // Constructor
55 //
56 
58  G4ProcessType type):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 }
73 
75 //
76 // Destructor
77 //
78 
80 {
81  delete genAngle;
82 }
83 
85 //
86 
87 void
89 {
90  if(p != genAngle) {
91  delete genAngle;
92  genAngle = p;
93  }
94 }
95 
96 G4bool
98 {
99  return (particle.GetPDGCharge() != 0.0 && !particle.IsShortLived());
100 }
101 
103 //
104 // Production of synchrotron X-ray photon
105 // GEANT4 internal units.
106 //
107 
108 G4double
110  G4double,
112 {
113  // gives the MeanFreePath in GEANT4 internal units
114  G4double MeanFreePath = DBL_MAX;
115 
116  const G4DynamicParticle* aDynamicParticle = trackData.GetDynamicParticle();
117 
118  *condition = NotForced;
119 
120  G4double gamma = aDynamicParticle->GetTotalEnergy()/
121  aDynamicParticle->GetMass();
122 
123  G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
124 
125  if ( gamma < 1.0e3 || 0.0 == particleCharge) { MeanFreePath = DBL_MAX; }
126  else
127  {
128 
129  G4ThreeVector FieldValue;
130  const G4Field* pField = nullptr;
131  G4bool fieldExertsForce = false;
132 
133 
134  G4FieldManager* fieldMgr =
136 
137  if ( fieldMgr != nullptr )
138  {
139  // If the field manager has no field, there is no field !
140 
141  fieldExertsForce = ( fieldMgr->GetDetectorField() != nullptr );
142  }
143 
144  if ( fieldExertsForce )
145  {
146  pField = fieldMgr->GetDetectorField();
147  G4ThreeVector globPosition = trackData.GetPosition();
148 
149  G4double globPosVec[4], FieldValueVec[6];
150 
151  globPosVec[0] = globPosition.x();
152  globPosVec[1] = globPosition.y();
153  globPosVec[2] = globPosition.z();
154  globPosVec[3] = trackData.GetGlobalTime();
155 
156  pField->GetFieldValue( globPosVec, FieldValueVec );
157 
158  FieldValue = G4ThreeVector( FieldValueVec[0],
159  FieldValueVec[1],
160  FieldValueVec[2] );
161 
162  G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
163  G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum);
164  G4double perpB = unitMcrossB.mag();
165 
166  static const G4double fLambdaConst = std::sqrt(3.0)*eplus/
168 
169  if( perpB > 0.0 )
170  {
171  MeanFreePath =
172  fLambdaConst*aDynamicParticle->GetDefinition()->GetPDGMass()
173  /(perpB*particleCharge*particleCharge);
174  }
175  if(verboseLevel > 0 && FirstTime)
176  {
177  G4cout << "G4SynchrotronRadiation::GetMeanFreePath "
178  << " for particle "
179  << aDynamicParticle->GetDefinition()->GetParticleName()
180  << ":" << '\n' //hbunew
181  << " MeanFreePath = " << G4BestUnit(MeanFreePath, "Length")
182  << G4endl;
183  if(verboseLevel > 1)
184  {
185  G4ThreeVector pvec = aDynamicParticle->GetMomentum();
186  G4double Btot = FieldValue.getR();
187  G4double ptot = pvec.getR();
188  G4double rho = ptot / (MeV * c_light * Btot );
189  // full bending radius
190  G4double Theta=unitMomentum.theta(FieldValue);
191  // angle between particle and field
192  G4cout << " B = " << Btot/tesla << " Tesla"
193  << " perpB = " << perpB/tesla << " Tesla"
194  << " Theta = " << Theta << " std::sin(Theta)="
195  << std::sin(Theta) << '\n'
196  << " ptot = " << G4BestUnit(ptot,"Energy")
197  << " rho = " << G4BestUnit(rho,"Length")
198  << G4endl;
199  }
200  FirstTime=false;
201  }
202  }
203  }
204  return MeanFreePath;
205 }
206 
208 //
209 //
210 
213  const G4Step& stepData )
214 
215 {
216  aParticleChange.Initialize(trackData);
217 
218  const G4DynamicParticle* aDynamicParticle = trackData.GetDynamicParticle();
219 
220  G4double gamma = aDynamicParticle->GetTotalEnergy()/
221  (aDynamicParticle->GetDefinition()->GetPDGMass());
222 
223  G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
224  if(gamma <= 1.0e3 || 0.0 == particleCharge)
225  {
226  return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
227  }
228 
229  G4ThreeVector FieldValue;
230  const G4Field* pField = nullptr;
231 
232  G4bool fieldExertsForce = false;
233  G4FieldManager* fieldMgr =
235 
236  if ( fieldMgr != nullptr )
237  {
238  // If the field manager has no field, there is no field !
239  fieldExertsForce = ( fieldMgr->GetDetectorField() != nullptr );
240  }
241 
242  if ( fieldExertsForce )
243  {
244  pField = fieldMgr->GetDetectorField();
245  G4ThreeVector globPosition = trackData.GetPosition();
246  G4double globPosVec[4], FieldValueVec[6];
247  globPosVec[0] = globPosition.x();
248  globPosVec[1] = globPosition.y();
249  globPosVec[2] = globPosition.z();
250  globPosVec[3] = trackData.GetGlobalTime();
251 
252  pField->GetFieldValue( globPosVec, FieldValueVec );
253  FieldValue = G4ThreeVector( FieldValueVec[0],
254  FieldValueVec[1],
255  FieldValueVec[2] );
256 
257  G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
258  G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum);
259  G4double perpB = unitMcrossB.mag();
260  if(perpB > 0.0)
261  {
262  // M-C of synchrotron photon energy
263 
264  G4double energyOfSR =
265  GetRandomEnergySR(gamma,perpB,
266  aDynamicParticle->GetDefinition()->GetPDGMass());
267 
268  // check against insufficient energy
269 
270  if( energyOfSR <= 0.0 )
271  {
272  return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
273  }
274  G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
275  G4ThreeVector gammaDirection =
276  genAngle->SampleDirection(aDynamicParticle,
277  energyOfSR, 1, 0);
278 
279  G4ThreeVector gammaPolarization = FieldValue.cross(gammaDirection);
280  gammaPolarization = gammaPolarization.unit();
281 
282  // create G4DynamicParticle object for the SR photon
283 
285  gammaDirection,
286  energyOfSR );
287  aGamma->SetPolarization( gammaPolarization.x(),
288  gammaPolarization.y(),
289  gammaPolarization.z() );
290 
293 
294  // Update the incident particle
295 
296  G4double newKinEnergy = kineticEnergy - energyOfSR;
297 
298  if (newKinEnergy > 0.)
299  {
300  aParticleChange.ProposeEnergy( newKinEnergy );
301  }
302  else
303  {
305  }
306  }
307  }
308  return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
309 }
310 
312 //
313 //
314 
316 // direct generation
317 {
318  // from 0 to 0.7
319  static const G4double aa1=0 ,aa2=0.7;
320  static const G4int ncheb1=27;
321  static const G4double cheb1[] =
322  { 1.22371665676046468821,0.108956475422163837267,0.0383328524358594396134,0.00759138369340257753721,
323  0.00205712048644963340914,0.000497810783280019308661,0.000130743691810302187818,0.0000338168760220395409734,
324  8.97049680900520817728e-6,2.38685472794452241466e-6,6.41923109149104165049e-7,1.73549898982749277843e-7,
325  4.72145949240790029153e-8,1.29039866111999149636e-8,3.5422080787089834182e-9,9.7594757336403784905e-10,
326  2.6979510184976065731e-10,7.480422622550977077e-11,2.079598176402699913e-11,5.79533622220841193e-12,
327  1.61856011449276096e-12,4.529450993473807e-13,1.2698603951096606e-13,3.566117394511206e-14,1.00301587494091e-14,
328  2.82515346447219e-15,7.9680747949792e-16};
329  // from 0.7 to 0.9132260271183847
330  static const G4double aa3=0.9132260271183847;
331  static const G4int ncheb2=27;
332  static const G4double cheb2[] =
333  { 1.1139496701107756,0.3523967429328067,0.0713849171926623,0.01475818043595387,0.003381255637322462,
334  0.0008228057599452224,0.00020785506681254216,0.00005390169253706556,0.000014250571923902464,3.823880733161044e-6,
335  1.0381966089136036e-6,2.8457557457837253e-7,7.86223332179956e-8,2.1866609342508474e-8,6.116186259857143e-9,
336  1.7191233618437565e-9,4.852755117740807e-10,1.3749966961763457e-10,3.908961987062447e-11,1.1146253766895824e-11,
337  3.1868887323415814e-12,9.134319791300977e-13,2.6211077371181566e-13,7.588643377757906e-14,2.1528376972619e-14,
338  6.030906040404772e-15,1.9549163926819867e-15};
339  // Chebyshev with exp/log scale
340  // a = -Log[1 - SynFracInt[1]]; b = -Log[1 - SynFracInt[7]];
341  static const G4double aa4=2.4444485538746025480,aa5=9.3830728608909477079;
342  static const G4int ncheb3=28;
343  static const G4double cheb3[] =
344  { 1.2292683840435586977,0.160353449247864455879,-0.0353559911947559448721,0.00776901561223573936985,
345  -0.00165886451971685133259,0.000335719118906954279467,-0.0000617184951079161143187,9.23534039743246708256e-6,
346  -6.06747198795168022842e-7,-3.07934045961999778094e-7,1.98818772614682367781e-7,-8.13909971567720135413e-8,
347  2.84298174969641838618e-8,-9.12829766621316063548e-9,2.77713868004820551077e-9,-8.13032767247834023165e-10,
348  2.31128525568385247392e-10,-6.41796873254200220876e-11,1.74815310473323361543e-11,-4.68653536933392363045e-12,
349  1.24016595805520752748e-12,-3.24839432979935522159e-13,8.44601465226513952994e-14,-2.18647276044246803998e-14,
350  5.65407548745690689978e-15,-1.46553625917463067508e-15,3.82059606377570462276e-16,-1.00457896653436912508e-16};
351  static const G4double aa6=33.122936966163038145;
352  static const G4int ncheb4=27;
353  static const G4double cheb4[] =
354  {1.69342658227676741765,0.0742766400841232319225,-0.019337880608635717358,0.00516065527473364110491,
355  -0.00139342012990307729473,0.000378549864052022522193,-0.000103167085583785340215,0.0000281543441271412178337,
356  -7.68409742018258198651e-6,2.09543221890204537392e-6,-5.70493140367526282946e-7,1.54961164548564906446e-7,
357  -4.19665599629607704794e-8,1.13239680054166507038e-8,-3.04223563379021441863e-9,8.13073745977562957997e-10,
358  -2.15969415476814981374e-10,5.69472105972525594811e-11,-1.48844799572430829499e-11,3.84901514438304484973e-12,
359  -9.82222575944247161834e-13,2.46468329208292208183e-13,-6.04953826265982691612e-14,1.44055805710671611984e-14,
360  -3.28200813577388740722e-15,6.96566359173765367675e-16,-1.294122794852896275e-16};
361 
362  if(x<aa2) return x*x*x*Chebyshev(aa1,aa2,cheb1,ncheb1,x);
363  else if(x<aa3) return Chebyshev(aa2,aa3,cheb2,ncheb2,x);
364  else if(x<1-0.0000841363)
365  { G4double y=-G4Log(1-x);
366  return y*Chebyshev(aa4,aa5,cheb3,ncheb3,y);
367  }
368  else
369  { G4double y=-G4Log(1-x);
370  return y*Chebyshev(aa5,aa6,cheb4,ncheb4,y);
371  }
372 }
373 
375  G4double gamma, G4double perpB, G4double mass_c2)
376 {
377 
378  static const G4double fEnergyConst = 1.5*c_light*c_light*eplus*hbar_Planck;
379  G4double Ecr=fEnergyConst*gamma*gamma*perpB/mass_c2;
380 
381  if(verboseLevel > 0 && FirstTime1)
382  {
383  // mean and rms of photon energy
384  G4double Emean=8./(15.*std::sqrt(3.))*Ecr;
385  G4double E_rms=std::sqrt(211./675.)*Ecr;
386  G4int prec = G4cout.precision();
387  G4cout << "G4SynchrotronRadiation::GetRandomEnergySR :" << '\n'
388  << std::setprecision(4)
389  << " Ecr = " << G4BestUnit(Ecr,"Energy") << '\n'
390  << " Emean = " << G4BestUnit(Emean,"Energy") << '\n'
391  << " E_rms = " << G4BestUnit(E_rms,"Energy") << G4endl;
392  FirstTime1=false;
393  G4cout.precision(prec);
394  }
395 
396  G4double energySR=Ecr*InvSynFracInt(G4UniformRand());
397  return energySR;
398 }
399 
401 //
402 //
403 
404 void
406 {
408  // same for all particles, print only for one (electron)
409 }
410 
412 //
413 //
414 
416 // not yet called, usually called from BuildPhysicsTable
417 {
418  G4String comments ="Incoherent Synchrotron Radiation\n";
419  G4cout << G4endl << GetProcessName() << ": " << comments
420  << " good description for long magnets at all energies"
421  << G4endl;
422 }
423 
Float_t x
Definition: compare.C:6
void SetAngularGenerator(G4VEmAngularDistribution *p)
static const double prec
Definition: RanecuEngine.cc:58
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
CLHEP::Hep3Vector G4ThreeVector
void SetPolarization(G4double polX, G4double polY, G4double polZ)
static constexpr double hbar_Planck
static constexpr double MeV
Definition: G4SIunits.hh:214
G4ParticleDefinition * theGamma
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &Step) override
#define G4endl
Definition: G4ios.hh:61
Float_t y
Definition: compare.C:6
const G4ThreeVector & GetMomentumDirection() const
const char * p
Definition: xmltok.h:285
const G4String & GetParticleName() const
G4double GetPDGCharge() const
double z() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
virtual void GetFieldValue(const double Point[4], double *fieldArr) const =0
G4double condition(const G4ErrorSymMatrix &m)
G4double GetPDGMass() const
G4PropagatorInField * fFieldPropagator
G4double G4Log(G4double x)
Definition: G4Log.hh:230
double theta() const
void AddSecondary(G4Track *aSecondary)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
G4VEmAngularDistribution * genAngle
G4SynchrotronRadiation(const G4String &pName="SynRad", G4ProcessType type=fElectromagnetic)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
TString part[npart]
Definition: Style.C:32
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
G4ParticleDefinition * GetDefinition() const
virtual void Initialize(const G4Track &)
G4ProcessType
G4double Chebyshev(G4double a, G4double b, const G4double c[], G4int n, G4double x)
G4double InvSynFracInt(G4double x)
virtual G4bool IsApplicable(const G4ParticleDefinition &) override
#define G4UniformRand()
Definition: Randomize.hh:53
G4double GetGlobalTime() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:411
const G4ThreeVector & GetPosition() const
const G4Field * GetDetectorField() const
Definition: G4Step.hh:76
void ProposeEnergy(G4double finalEnergy)
static G4Electron * Electron()
Definition: G4Electron.cc:94
Hep3Vector unit() const
Hep3Vector cross(const Hep3Vector &) const
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
static G4TransportationManager * GetTransportationManager()
G4ThreeVector GetMomentum() const
static constexpr double eplus
Definition: G4SIunits.hh:199
double mag() const
double getR() const
int G4int
Definition: G4Types.hh:78
static constexpr double c_light
G4int verboseLevel
Definition: G4VProcess.hh:371
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:435
G4ForceCondition
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
double x() const
G4PropagatorInField * GetPropagatorInField() const
virtual void BuildPhysicsTable(const G4ParticleDefinition &) override
G4double GetMass() const
G4double GetTotalEnergy() const
double y() const
static constexpr double fine_structure_const
#define DBL_MAX
Definition: templates.hh:83
virtual G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
const G4DynamicParticle * GetDynamicParticle() const
static constexpr double tesla
Definition: G4SIunits.hh:268
G4VPhysicalVolume * GetVolume() const
G4double GetRandomEnergySR(G4double, G4double, G4double)
void SetNumberOfSecondaries(G4int totSecondaries)