Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VChannelingFastSimCrystalData.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// Author: Alexei Sytov
27// Co-author: Gianfranco PaternĂ² (modifications & testing)
28// On the base of the CRYSTALRAD realization of scattering model:
29// A. I. Sytov, V. V. Tikhomirov, and L. Bandiera PRAB 22, 064601 (2019)
30
32#include "G4Log.hh"
33
35
36//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
37
39
40//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
41
43 (const G4LogicalVolume *crystallogic)
44{
45 G4int crystalID = crystallogic->GetInstanceID();
46
47 //set bending angle if the it exists in the list, otherwise default = 0
48 (fMapBendingAngle.count(crystalID) > 0)
49 ? SetBendingAngle(fMapBendingAngle[crystalID],crystallogic)
50 : SetBendingAngle(0.,crystallogic);
51
52 //set miscut angle if the it exists in the list, otherwise default = 0
53 (fMapMiscutAngle.count(crystalID) > 0)
54 ? SetMiscutAngle(fMapMiscutAngle[crystalID],crystallogic)
55 : SetMiscutAngle(0.,crystallogic);
56
57 //set crystalline undulator parameters if they exist in the list,
58 //otherwise default = G4ThreeVector(0,0,0).
59 (fMapCUAmplitudePeriodPhase.count(crystalID) > 0)
60 ? SetCUParameters(fMapCUAmplitudePeriodPhase[crystalID],crystallogic)
61 : SetCUParameters(G4ThreeVector(0.,0.,0.),crystallogic);
62}
63
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
65
67 const G4LogicalVolume* crystallogic)
68{
69 G4int crystalID = crystallogic->GetInstanceID();
70
71 //set the bending angle for this logical volume
72 fMapBendingAngle[crystalID]=tetab;
73
74 G4ThreeVector limboxmin;//minimal limits of the box bounding the logical volume
75 G4ThreeVector limboxmax;//maximal limits of the box bounding the logical volume
76 //save the limits of the box bounding the logical volume
77 crystallogic->GetSolid()->BoundingLimits(limboxmin,limboxmax);
78
79 //bounding box half dimensions
80 fHalfDimBoundingBox = (limboxmax-limboxmin)/2.;
81
82 G4double lcr = limboxmax.getZ()-limboxmin.getZ();//crystal thickness
83
84 fBendingAngle=std::abs(tetab);
85 if (fBendingAngle<0.000001)//no bending less then 1 urad
86 {
88 {
89 G4cout << "Channeling model: volume " << crystallogic->GetName() << G4endl;
90 G4cout << "Warning: bending angle is lower than 1 urad => set to 0" << G4endl;
91 }
92
93 fBent=0;
95 fBendingR=0.;//just for convenience (infinity in reality)
96 fBending2R=0.;
98 fCurv=0.;
99
100 fCorrectionZ = 1.;
101 }
102 else
103 {
104 fBent=1;
108 fCurv=1./fBendingR;
109
110 if (tetab<0.)
111 {
112 G4cout << "Channeling model: volume " << crystallogic->GetName() << G4endl;
113 G4cout << "Warning: bending angle is negative => set to be positive" << G4endl;
114 }
115 }
116
117}
118
119//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
120
122 const G4LogicalVolume *crystallogic)
123{
124 G4int crystalID = crystallogic->GetInstanceID();
125
126 //set the bending angle for this logical volume
127 fMapMiscutAngle[crystalID]=tetam;
128
129 // fMiscutAngle>0: rotation of xz coordinate planes clockwise in the xz plane
130 fMiscutAngle=tetam;
131 if (std::abs(tetam)>1.*CLHEP::mrad)
132 {
133 G4cout << "Channeling model: volume " << crystallogic->GetName() << G4endl;
134 G4cout << "Warning: miscut angle is higher than 1 mrad => " << G4endl;
135 G4cout << "coordinate transformation routines may be unstable" << G4endl;
136 }
139}
140
141//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
142
144 G4double amplitude,
145 G4double period,
146 G4double phase,
147 const G4LogicalVolume *crystallogic)
148{
149 if (amplitude<DBL_EPSILON||period<DBL_EPSILON)
150 {
151 amplitude = 0.;
152 period=0.;
153 phase=0.;
154 G4cout << "Channeling model: volume " << crystallogic->GetName() << G4endl;
155 G4cout << "Warning: The crystalline undulator parameters are out of range "
156 "=> the crystalline undulator mode switched off" << G4endl;
157 }
158
159 SetCUParameters(G4ThreeVector(amplitude,period,phase),crystallogic);
160}
161
162//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
163
165 const G4LogicalVolume *crystallogic,
166 const G4String &filename)
167{
168 G4int crystalID = crystallogic->GetInstanceID();
169 fMapImportCrystalGeometry[crystalID] = true;
170
171 G4PhysicsLinearVector vec(true); //CU geometry data x=f(z) - transverse coordinate vs longitudinal
172 G4PhysicsLinearVector vecD(true); // = first derivative of vec
173 G4PhysicsLinearVector vecDD(true);// = second derivative of vec
174
175 //open input file for CU data
176 std::ifstream vfilein;
177 vfilein.open(filename);
178
179 if (!vfilein) {
180 G4String message = "Input file " +
181 filename +
182 " is not found!";
183 G4Exception("SetCrystallineUndulatorParameters",// Origin of the exception
184 "001", // Unique error code
185 FatalException, // Terminate the program
186 message);
187 }
188
189 G4cout << "Importing the Crystal geometry from the file: "
190 << filename << G4endl;
191
192 //reading from file
193 //CAUTION: all the values in the file in mm!!!
194 vec.Retrieve(vfilein,true);
195
196 vfilein.close();
197
198 //check the start coordinate (10.*DBL_EPSILON to add more flexibility, if not equal, no problem)
199 if(std::abs(vec.GetMinEnergy()) > 10.*DBL_EPSILON)
200 {
201 G4Exception("SetCrystallineUndulatorParameters",// Origin of the exception
202 "001", // Unique error code
203 FatalException, // Terminate the program
204 "The interpolation of the crystalline undulator does not start from 0.! "
205 "The program will terminate.");
206
207 }
208
209 //check the end coordinate (1.e8*DBL_EPSILON - single precision check)
210 if(std::abs(vec.GetMaxEnergy() - 2.*fHalfDimBoundingBox.z()) > 100000000.*DBL_EPSILON)
211 {
212 G4Exception("SetCrystallineUndulatorParameters",// Origin of the exception
213 "001", // Unique error code
214 FatalException, // Terminate the program
215 "The interpolation of the crystalline undulator does not end with at"
216 "the boundary of the volume => check the crystal thickness in z! "
217 "The program will terminate.");
218 }
219
220 // Get the total length of vectors
221 std::size_t inodes = vec.GetVectorLength();
222
223 //check if the interpolation nodes are equidistant (otherwise )
224 for(std::size_t i=1; i<inodes-1; i++)
225 {
226 if(std::abs((vec.Energy(i+1)-vec.Energy(i))-
227 (vec.Energy(i) -vec.Energy(i-1)))>100000000.*DBL_EPSILON)//single precision
228 {
229 G4String message = "The interpolation nodes in the file " + filename +
230 " are not equidistant! The program will terminate.";
231 G4Exception("SetCrystallineUndulatorParameters",// Origin of the exception
232 "001", // Unique error code
233 FatalException, // Terminate the program
234 message);
235 }
236 }
237
238 //necessary to use spline interpolation
240
241 //setting up the vectors (to have the same arguments)
242 vecD = vec;
243 vecDD = vec;
244
245 //calculating the derivatives
246 for(std::size_t i=1; i<inodes-1; i++)
247 {
248 vecD.PutValue(i,(vec[i+1]-vec[i-1])/(vec.Energy(i+1)-vec.Energy(i-1)));
249 vecDD.PutValue(i,((vec[i+1]-vec[i])/(vec.Energy(i+1)-vec.Energy(i)) -
250 (vec[i]-vec[i-1])/(vec.Energy(i)-vec.Energy(i-1)))*
251 2./(vec.Energy(i+1)-vec.Energy(i-1)));
252 }
253
254 //end points have the same values as neighbouring
255 vecD.PutValue(0,vecD[1]);
256 vecD.PutValue(inodes-1,vecD[inodes-2]);
257 vecDD.PutValue(0,vecDD[1]);
258 vecDD.PutValue(inodes-1,vecDD[inodes-2]);
259
260 //necessary to use spline interpolation
263
264 fVecCUx.push_back(vec);
265 fVecCUtetax.push_back(vecD);
266 fVecCUCurv.push_back(vecDD);
267
268 //check if the fVecCUx was already set up for this logical volume
269 if(fMapCUID.count(crystalID) > 0)
270 {
271 G4Exception("SetCrystallineUndulatorParameters",// Origin of the exception
272 "001", // Unique error code
273 FatalException, // Terminate the program
274 "It is not allowed to set up the crystalline undulator geometry "
275 "the second time for the same logical volume! The program will terminate.");
276 }
277
278 //save the number of current element in fVecCUx, fVecCUtetax, fVecCUCurv
279 //to call them later using crystalID
280 fMapCUID[crystalID] = (G4int)fVecCUx.size()-1;
281
282 SetCUParameters(G4ThreeVector(1.,1.,0.),crystallogic);// input G4ThreeVector(1.,1.,0.)
283 // is needed just to setup
284 // the undulator in a normal way =>
285 // first 2 arguments > 0.
286 // values do not matter
287}
288
289//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
290
292 const G4ThreeVector &amplitudePeriodPhase,
293 const G4LogicalVolume *crystallogic)
294{
295 G4int crystalID = crystallogic->GetInstanceID();
296
297 //set the crystalline undulator parameters for this logical volume
298 fMapCUAmplitudePeriodPhase[crystalID]=amplitudePeriodPhase;
299 fCUAmplitude=amplitudePeriodPhase.x();
300 G4double period = amplitudePeriodPhase.y();
301 fCUPhase = amplitudePeriodPhase.z();
302
303 // flag of custom crystal geometry uploaded from a file
304 fImportCrystalGeometry = (fMapImportCrystalGeometry.count(crystalID) > 0) ?
305 fMapImportCrystalGeometry[crystalID] :
306 false;
307
308 //if the amplidude of the crystalline undulator is 0 => no undulator
310 {
311 //crystalline undulator flag
312 fCU = true;
313
314 fCUK = CLHEP::twopi/period;
315
317 {
318 //bent and periodically bent crystal are not compatible
319 SetBendingAngle(0,crystallogic);
320
321 G4cout << "Channeling model: volume " << crystallogic->GetName() << G4endl;
322 G4cout << "Warning: crystalline undulator is not compatible with "
323 "a bent crystal mode => setting bending angle to 0." << G4endl;
324 }
325
326 //setup crystal internal geometry data to access the correct elements of
327 //fVecCUx, fVecCUtetax, fVecCUCurv
328 if(fImportCrystalGeometry)
329 {
330 fCUID = fMapCUID[crystalID];
331 }
332 }
333 else
334 {
335 fCU = false;
336 fCUAmplitude = 0.;
337 fCUK = 0.;
338 fCUPhase = 0.;
339 fMapCUAmplitudePeriodPhase[crystalID] = G4ThreeVector(0.,0.,0.);
340 }
341
342 fCUK2 = fCUK*fCUK;
344}
345
346//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
347
349 G4double mass,
350 G4double charge,
351 const G4String& particleName)
352{
353 G4double teta1;
354 fZ2=charge;
355 G4double zz22=fZ2*fZ2;
356 fParticleName=particleName;
357
358// particle momentum and energy
359 G4double t=etotal*etotal-mass*mass; // economy of operations
360 fPz=std::sqrt(t); // momentum of particle
361 fPV=t/etotal; // pv
362 fBeta=fPz/etotal; // velocity/c
363 fTetaL = std::sqrt(std::abs(fZ2)*fVmax2/fPV); //Lindhard angle
364 fChannelingStep = fChangeStep/fTetaL; //standard simulation step
365
366// Energy losses
367 fV2 = fBeta*fBeta; // particle (velocity/c)^2
368 fGamma = etotal/mass; // Lorentz factor
369 fMe2Gamma = 2*CLHEP::electron_mass_c2*fGamma;
370// max ionization losses
371 fTmax = fMe2Gamma*fGamma*fV2/
372 (CLHEP::electron_mass_c2/mass*CLHEP::electron_mass_c2/mass +
373 1. + fMe2Gamma/mass);
374// max ionization losses for electrons
375 if(fParticleName=="e-"){fTmax/=2;}
376
377 for(G4int i=0; i<fNelements; i++)
378 {
379
380// minimal scattering angle by coulomb scattering on nuclei
381// defining by shielding by electrons
382// teta1=hdc/(fPz*fRF)*DSQRT(1.13D0+3.76D0*(alpha*fZ1*fZ2/fBeta)**2){ev*cm/(eV*cm)}
383 teta1=fTeta10[i]*std::sqrt(1.13+fK40[i]*zz22/fV2); // /fPz later to speed up
384 // the calculations
385
386// the coefficient for multiple scattering
387 fBB[i]=teta1*teta1*fPu11[i];
388 fE1XBbb[i]=expint(fBB[i]);
389 fBBDEXP[i]=(1.+fBB[i])*std::exp(fBB[i]);
390// necessary for suppression of incoherent scattering
391// by the atomic correlations in crystals for single scattering on nucleus
392// (screened atomic potential): EXP(-(fPz*teta*fU1)**2)=EXP(-fPzu11*teta**2).GE.ksi
393// =>no scattering
394 fPzu11[i]=fPu11[i]*fPz*fPz;
395
396 teta1=teta1/fPz; //
397 fTeta12[i]=teta1*teta1;
398// maximal scattering angle by coulomb scattering on nuclei
399// defining by nucleus radius
400// tetamax=hc/(fPz*1.D-6*fR0*fAN**(1.D0/3.D0))// {Mev*fermi/(MeV*fermi)}
401 G4double tetamax=fTetamax0[i]/fPz;
402 fTetamax2[i]=tetamax*tetamax;
404
405// a coefficient in a formula for scattering (for high speed of simulation)
406// fK2=(fZ2)**2*alphahbarc2*4.*pi*fN0*(fZ1/fPV)**2
407 fK2[i]=fK20[i]*zz22/fPV/fPV;
408 }
409
410// fK3=(fZ2)**2*alphahbarc2*pi/electron_mass_c2/(fV2)**2
411 fK3=fK30*zz22/fV2;
412}
413
414//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
415
417 G4double mass,
418 G4double charge)
419{
420 G4double pv0 = etotal-mass*mass/etotal;
421 return std::sqrt(2*std::abs(charge)*fVmax/pv0); //Calculate the value of the Lindhard angle
422 //(!!! the value for a straight crystal)
423}
424
425//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
426
428{
429 return fTetaL; //return the Lindhard angle value calculated in SetParticleProperties
430}
431
432//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
433
435{
436 G4double simulationstep;
437 //find angle of particle w.r.t. the plane or axis
438 G4double angle=0.;
439 if (iModel==1)//1D model
440 {
441 angle = std::abs(tx);
442 }
443 else if (iModel==2)//2D model
444 {
445 angle = std::sqrt(tx*tx+ty*ty);
446 }
447
448 //compare this angle with the Lindhard angle
449 if (angle<fTetaL)
450 {
451 simulationstep = fChannelingStep;
452 }
453 else
454 {
455 simulationstep = fChangeStep;
456 if (angle > 0.0) { simulationstep /= angle; }
457 }
458
459 return simulationstep;
460}
461
462//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
463
465 G4double mass,
466 G4double charge)
467{
468 //standard value of step for channeling particles which is the maximal possible step
469 return fChangeStep/GetLindhardAngle(etotal, mass, charge);
470}
471
472//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
473
475 G4double effectiveStep,
476 G4double step,
477 G4int ielement)
478{
479 G4double tx = 0.;//horizontal scattering angle
480 G4double ty = 0.;//vertical scattering angle
481
482 G4double ksi=0.1;
483
484// calculation of the teta2-minimal possible angle of a single scattering
485 G4double e1=fK2[ielement]*effectiveStep; //for high speed of a program
486// (real formula is (4*pi*fN0*wpl(x)*dz*fZ1*zz2*alpha*hdc/fPV)**2)
487 G4double teta122=fTetamax12[ielement]/(ksi*fTetamax12[ielement]/e1+1.);
488 // teta122=fTeta12+teta22=teta1^2+teta2^2
489
490 G4double teta22;
491 G4double t;
492// if the angle of a single scattering is less teta1 - minimal possible
493// angle of coulomb scattering defining by the electron shielding than
494// multiple scattering by both nuclei and electrons and electrons will not
495// occur => minimal possible angle of a single scattering is equal to teta1
496 if (teta122<=fTeta12[ielement]*1.000125)
497 {
498 teta22=0.;
499 teta122=fTeta12[ielement];
500 }
501 else
502 {
503 teta22=teta122-fTeta12[ielement];
504 G4double aa=teta22/fTeta12[ielement];
505 G4double aa1=1.+aa;
506
507// crystal, with scattering suppression
508 G4double tetamsi=e1*(G4Log(aa1)+
509 (1.-std::exp(-aa*fBB[ielement]))/aa1+
510 fBBDEXP[ielement]*
511 (expint(fBB[ielement]*aa1)-fE1XBbb[ielement]));
512
513// sumilation of multiple coulomb scattering by nuclei and electrons
514// for high speed of a program, real formula is
515// 4*pi*fN0*wpl(x)*dz*(fZ1*zz2*alpha*hdc/fPV)**2*
516// *(ln(1+a)+(1-exp(-a*b))/(1+a)+(1+b)*exp(b)*(E1XB(b*(1+a))-E1XB(b)))
517
518 ksi=G4UniformRand();
519 t=std::sqrt(-tetamsi*G4Log(ksi));
520
521 ksi=G4UniformRand();
522
523 tx+=t*std::cos(CLHEP::twopi*ksi);
524 ty+=t*std::sin(CLHEP::twopi*ksi);
525
526 }
527// simulation of single coulomb scattering by nuclei (with screened potential)
528 G4double zss=0.;
529 G4double dzss=step;
530
531// (calculation of a distance, at which another single scattering can happen)
532 ksi=G4UniformRand();
533
534 zss=-G4Log(ksi)*step/(e1*(1./teta122-1./fTetamax12[ielement]));
535 G4double tt;
536
537// At some step several single scattering can occur.
538// So,if the distance of the next scattering is less than the step,
539// another scattering can occur. If the distance of the next scattering
540// is less than the difference between the step and the distance of
541// the previous scattering, another scattering can occur. And so on, and so on.
542// In the cycle we simulate each of them. The cycle is finished, when
543// the remaining part of step is less than a distance of the next single scattering.
544//********************************************
545// if at a step a single scattering occurs
546 while (zss<dzss)
547 {
548
549// simulation by Monte-Carlo of angles of single scattering
550 ksi=G4UniformRand();
551
552 tt=fTetamax12[ielement]/(1.+ksi*(fTetamax2[ielement]-teta22)/teta122)-
553 fTeta12[ielement];
554
555 ksi=G4UniformRand();
556
557// suppression of incoherent scattering by the atomic correlations in crystals
558 t=fPzu11[ielement]*tt;
559 t=std::exp(-t);
560
561 if (t<ksi) //if scattering takes place
562 {
563 //scattering angle
564 t=std::sqrt(tt);
565 ksi=G4UniformRand();
566
567 tx+=t*std::cos(CLHEP::twopi*ksi);
568 ty+=t*std::sin(CLHEP::twopi*ksi);
569 }
570
571 dzss-=zss;
572// (calculation of a distance, at which another single scattering can happen)
573 ksi=G4UniformRand();
574
575 zss=-G4Log(ksi)*step/(e1*(1./teta122-1./fTetamax12[ielement]));
576 }
577//********************************************
578 return G4ThreeVector(tx,ty,0.);
579}
580
581//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
582
584 G4double eMinIonization,
585 G4double electronDensity,
586 G4double step)
587{
588
589 G4double zss=0.;
590 G4double dzss=step;
591 G4double ksi = 0.;
592
593 G4double tx = 0.;//horizontal scattering angle
594 G4double ty = 0.;//vertical scattering angle
595 G4double eloss = 0.;//energy loss
596
597 // eMinIonization - minimal energy transfered to electron
598 // a cut to reduce the number of calls of electron scattering
599 // is needed only at low density regions, in many cases does not do anything at all
600 if (eMinIonization<0.5*CLHEP::eV){eMinIonization=0.5*CLHEP::eV;}
601
602 // single scattering on electrons routine
603 if ((eMinIonization<fTmax)&&(electronDensity>DBL_EPSILON))
604 {
605
606// (calculation of a distance, at which another single scattering can happen)
607// simulation of scattering length (by the same way single scattering by nucleus
608 ksi=G4UniformRand();
609
610 zss=-1.0*G4Log(ksi)/(fK3*electronDensity)/(1./eMinIonization-1./fTmax);
611
612//********************************************
613// if at a step a single scattering occur
614 while (zss<dzss)
615 {
616// simulation by Monte-Carlo of angles of single scattering
617 ksi=G4UniformRand();
618
619// energy transfered to electron
620 G4double e1=eMinIonization/(1.-ksi*(1.-eMinIonization/fTmax));
621
622// scattering angle
623 G4double t=0;
624 if(fTmax-e1>DBL_EPSILON) //to be sure e1<fTmax
625 {
626 t=std::sqrt(2.*CLHEP::electron_mass_c2*e1*(1-e1/fTmax))/fPz;
627 }
628
629 // energy losses
630 eloss=e1;
631 ksi=G4UniformRand();
632
633 tx+=t*std::cos(CLHEP::twopi*ksi);
634 ty+=t*std::sin(CLHEP::twopi*ksi);
635
636 dzss-=zss;
637// (calculation of a distance, at which another single scattering can happen)
638// simulation of scattering length
639// (by the same way single scattering by nucleus
640 ksi=G4UniformRand();
641
642 zss=-1.0*G4Log(ksi)/(fK3*electronDensity)/(1./eMinIonization-1./fTmax);
643 }
644//********************************************
645 }
646 return G4ThreeVector(tx,ty,eloss);
647}
648
649//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
650
652 G4int ielement)
653{
654 //amorphous part of ionization losses
655
656 G4double elosses = 0.;
657 // 1/2 already taken into account in fKD
658
659 G4double loge = G4Log(fMe2Gamma*fGamma*fV2/fI0[ielement]);
660 G4double delta= 2.*(G4Log(fBeta*fGamma)+fLogPlasmaEdI0[ielement]-0.5);
661 if(delta<0.){delta=0.;}
662 loge-=delta;
663 if(fParticleName=="e-")
664 {
665 loge+=(-G4Log(2.) + 1.
666 -(2.*fGamma - 1.)/fGamma/fGamma*G4Log(2.) +
667 1./8.*((fGamma - 1.)/fGamma)*((fGamma - 1.)/fGamma));
668 }
669 else if(fParticleName=="e+")
670 {
671 loge+=(-fV2/12.*(11. + 14./(fGamma + 1.) + 10./(fGamma + 1.)/(fGamma + 1.) +
672 4./(fGamma + 1.)/(fGamma + 1.)/(fGamma + 1.)));
673 }
674 else
675 {
676 loge-=fV2;
677 }
678 elosses=fZ2*fZ2*fKD[ielement]/fV2*loge*dz;
679
680 return elosses;}
681
682//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
683
684G4double G4VChannelingFastSimCrystalData::expint(G4double X)
685{
686// ============================================
687// Purpose: Compute exponential integral E1(x)
688// Input : x --- Argument of E1(x)
689// Output: E1 --- E1(x)
690// ============================================
691
692G4double E1, R, T, T0;
693G4int M;
694
695if (X==0)
696{
697 E1=1.e300;
698}
699else if (X<=1.)
700{
701 E1=1.;
702 R=1.;
703
704
705 for(int K=1; K<=25; K++)
706 {
707 R=-R*K*X/std::pow(K+1.,2.);
708 E1=E1+R;
709 if (std::abs(R)<=std::abs(E1)*1.0e-15) {break;}
710 }
711
712 E1=-0.5772156649015328-G4Log(X)+X*E1;
713}
714else
715{
716 M=20+std::trunc(80.0/X);
717 T0=0.;
718
719 for(int K=M; K>=1; K--)
720 {
721 T0=K/(1.0+K/(X+T0));
722 }
723
724 T=1.0/(X+T0);
725 E1=std::exp(-X)*T;
726}
727
728return E1;
729}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4double G4Log(G4double x)
Definition G4Log.hh:169
#define M(row, col)
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
Definition of the G4VChannelingFastSimCrystalData class The class contains the data and properties re...
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
double z() const
double getZ() const
double x() const
double y() const
G4LogicalVolume represents a leaf node or unpositioned subtree in the geometry hierarchy....
G4VSolid * GetSolid() const
G4int GetInstanceID() const
const G4String & GetName() const
void PutValue(const std::size_t index, const G4double value)
G4double GetMaxEnergy() const
G4double Energy(const std::size_t index) const
G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
G4double GetMinEnergy() const
std::size_t GetVectorLength() const
void FillSecondDerivatives(const G4SplineType=G4SplineType::Base, const G4double dir1=0.0, const G4double dir2=0.0)
G4double IonizationLosses(G4double dz, G4int ielement)
ionization losses
std::vector< G4double > fK20
coefficients necessary for multiple and single coulomb scattering
G4double GetLindhardAngle()
Calculate the value of the Lindhard angle (!!! the value for a straight crystal).
void SetGeometryParameters(const G4LogicalVolume *crystallogic)
set geometry parameters from current logical volume
void SetCUParameters(const G4ThreeVector &amplitudePeriodPhase, const G4LogicalVolume *crystallogic)
G4double GetMaxSimulationStep(G4double etotal, G4double mass, G4double charge)
Calculate maximal simulation step (standard value for channeling particles).
void SetCrystallineUndulatorParameters(G4double amplitude, G4double period, G4double phase, const G4LogicalVolume *crystallogic)
G4ThreeVector fHalfDimBoundingBox
values related to the crystal geometry
G4double GetSimulationStep(G4double tx, G4double ty)
void SetBendingAngle(G4double tetab, const G4LogicalVolume *crystallogic)
G4int fNelements
values related to the crystal lattice
G4ThreeVector CoulombElectronScattering(G4double eMinIonization, G4double electronDensity, G4double step)
multiple and single scattering on electrons
std::vector< G4double > fTeta10
angles necessary for multiple and single coulomb scattering
std::vector< G4double > fPu11
coefficients for multiple scattering suppression
G4ThreeVector CoulombAtomicScattering(G4double effectiveStep, G4double step, G4int ielement)
multiple and single scattering on screened potential
void SetMiscutAngle(G4double tetam, const G4LogicalVolume *crystallogic)
void SetParticleProperties(G4double etotal, G4double mp, G4double charge, const G4String &particleName)
virtual void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition G4VSolid.cc:691
#define DBL_EPSILON
Definition templates.hh:66