Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VChannelingFastSimCrystalData Class Referenceabstract

#include <G4VChannelingFastSimCrystalData.hh>

+ Inheritance diagram for G4VChannelingFastSimCrystalData:

Public Member Functions

 G4VChannelingFastSimCrystalData ()
 
virtual ~G4VChannelingFastSimCrystalData ()
 
G4double Ex (G4double x, G4double y)
 electric fields produced by crystal lattice
 
G4double Ey (G4double x, G4double y)
 
G4double ElectronDensity (G4double x, G4double y)
 electron density function
 
G4double MinIonizationEnergy (G4double x, G4double y)
 minimum energy of ionization function
 
G4double NuclearDensity (G4double x, G4double y, G4int ielement)
 nuclear density function (normalized to average nuclear density)
 
G4double GetLindhardAngle (G4double etotal, G4double mass, G4double charge)
 Calculate the value of the Lindhard angle (!!! the value for a straight crystal)
 
G4double GetLindhardAngle ()
 Calculate the value of the Lindhard angle (!!! the value for a straight crystal)
 
G4double GetSimulationStep (G4double tx, G4double ty)
 
G4double GetMaxSimulationStep (G4double etotal, G4double mass, G4double charge)
 Calculate maximal simulation step (standard value for channeling particles)
 
G4double GetBeta ()
 get particle velocity/c
 
G4int GetNelements ()
 
G4int GetModel ()
 
G4double GetBendingAngle ()
 
G4double GetMiscutAngle ()
 
G4double GetCurv (G4double z)
 
G4double GetCUx (G4double z)
 get crystalline undulator wave function
 
G4double GetCUtetax (G4double z)
 get crystalline undulator wave 1st derivative function
 
virtual void SetMaterialProperties (const G4Material *crystal, const G4String &lattice, const G4String &filePath)=0
 
void SetGeometryParameters (const G4LogicalVolume *crystallogic)
 set geometry parameters from current logical volume
 
void SetBendingAngle (G4double tetab, const G4LogicalVolume *crystallogic)
 
void SetMiscutAngle (G4double tetam, const G4LogicalVolume *crystallogic)
 
void SetCrystallineUndulatorParameters (G4double amplitude, G4double period, G4double phase, const G4LogicalVolume *crystallogic)
 
void SetCUParameters (const G4ThreeVector &amplitudePeriodPhase, const G4LogicalVolume *crystallogic)
 
void SetParticleProperties (G4double etotal, G4double mp, G4double charge, const G4String &particleName)
 
virtual G4ThreeVector CoordinatesFromBoxToLattice (const G4ThreeVector &pos0)=0
 
virtual G4ThreeVector CoordinatesFromLatticeToBox (const G4ThreeVector &pos)=0
 
virtual G4ThreeVector ChannelChange (G4double &x, G4double &y, G4double &z)=0
 change the channel if necessary, recalculate x o y
 
G4double GetCorrectionZ ()
 
virtual G4double AngleXFromBoxToLattice (G4double tx, G4double z)=0
 
virtual G4double AngleXFromLatticeToBox (G4double tx, G4double z)=0
 
virtual G4double AngleXShift (G4double z)=0
 auxialiary function to transform the horizontal angle
 
G4ThreeVector CoulombAtomicScattering (G4double effectiveStep, G4double step, G4int ielement)
 multiple and single scattering on screened potential
 
G4ThreeVector CoulombElectronScattering (G4double eMinIonization, G4double electronDensity, G4double step)
 multiple and single scattering on electrons
 
G4double IonizationLosses (G4double dz, G4int ielement)
 ionization losses
 
void SetVerbosity (G4int ver)
 

Protected Attributes

G4ChannelingFastSimInterpolationfElectricFieldX {nullptr}
 classes containing interpolation coefficients
 
G4ChannelingFastSimInterpolationfElectricFieldY {nullptr}
 
G4ChannelingFastSimInterpolationfElectronDensity {nullptr}
 
G4ChannelingFastSimInterpolationfMinIonizationEnergy {nullptr}
 
std::vector< G4ChannelingFastSimInterpolation * > fNucleiDensity
 
G4ThreeVector fHalfDimBoundingBox
 values related to the crystal geometry
 
G4int fBent =0
 
G4double fBendingAngle =0.
 
G4double fBendingR = 0.
 
G4double fBending2R =0.
 
G4double fBendingRsquare =0.
 
G4double fCurv =0.
 
G4double fMiscutAngle = 0.
 
G4double fCosMiscutAngle =1.
 
G4double fSinMiscutAngle =0.
 
G4double fCorrectionZ = 1.
 
G4bool fCU = false
 
G4double fCUAmplitude =0.
 
G4double fCUK =0.
 
G4double fCUPhase =0.
 
G4double fCUAmplitudeK =0.
 
G4double fCUK2 =0.
 
G4int fNelements =1
 values related to the crystal lattice
 
G4int iModel =1
 
G4double fVmax =0
 
G4double fVmax2 =0
 
G4double fVMinCrystal =0
 
G4double fChangeStep =0
 
std::vector< G4doublefI0
 
std::vector< G4doublefRF
 
std::vector< G4doublefTeta10
 angles necessary for multiple and single coulomb scattering
 
std::vector< G4doublefTetamax0
 
std::vector< G4doublefTetamax2
 
std::vector< G4doublefTetamax12
 
std::vector< G4doublefTeta12
 
std::vector< G4doublefK20
 coefficients necessary for multiple and single coulomb scattering
 
std::vector< G4doublefK2
 
std::vector< G4doublefK40
 
G4double fK30 =0
 
G4double fK3 =0
 
std::vector< G4doublefKD
 
std::vector< G4doublefLogPlasmaEdI0
 
std::vector< G4doublefPu11
 coefficients for multiple scattering suppression
 
std::vector< G4doublefPzu11
 
std::vector< G4doublefBB
 
std::vector< G4doublefE1XBbb
 
std::vector< G4doublefBBDEXP
 
G4int fVerbosity = 1
 

Detailed Description

Definition at line 60 of file G4VChannelingFastSimCrystalData.hh.

Constructor & Destructor Documentation

◆ G4VChannelingFastSimCrystalData()

G4VChannelingFastSimCrystalData::G4VChannelingFastSimCrystalData ( )

Definition at line 36 of file G4VChannelingFastSimCrystalData.cc.

37{
38
39}

◆ ~G4VChannelingFastSimCrystalData()

G4VChannelingFastSimCrystalData::~G4VChannelingFastSimCrystalData ( )
virtual

Definition at line 43 of file G4VChannelingFastSimCrystalData.cc.

43{;}

Member Function Documentation

◆ AngleXFromBoxToLattice()

virtual G4double G4VChannelingFastSimCrystalData::AngleXFromBoxToLattice ( G4double tx,
G4double z )
pure virtual

calculate the horizontal angle in the co-rotating reference system within a channel (periodic cell) (connected with crystal planes/axes either bent or straight)

Implemented in G4ChannelingFastSimCrystalData.

◆ AngleXFromLatticeToBox()

virtual G4double G4VChannelingFastSimCrystalData::AngleXFromLatticeToBox ( G4double tx,
G4double z )
pure virtual

calculate the horizontal angle in the Box reference system (connected with the bounding box of the volume)

Implemented in G4ChannelingFastSimCrystalData.

◆ AngleXShift()

virtual G4double G4VChannelingFastSimCrystalData::AngleXShift ( G4double z)
pure virtual

auxialiary function to transform the horizontal angle

Implemented in G4ChannelingFastSimCrystalData.

◆ ChannelChange()

virtual G4ThreeVector G4VChannelingFastSimCrystalData::ChannelChange ( G4double & x,
G4double & y,
G4double & z )
pure virtual

change the channel if necessary, recalculate x o y

Implemented in G4ChannelingFastSimCrystalData.

◆ CoordinatesFromBoxToLattice()

virtual G4ThreeVector G4VChannelingFastSimCrystalData::CoordinatesFromBoxToLattice ( const G4ThreeVector & pos0)
pure virtual

calculate the coordinates in the co-rotating reference system within a channel (periodic cell) (connected with crystal planes/axes either bent or straight)

Implemented in G4ChannelingFastSimCrystalData.

◆ CoordinatesFromLatticeToBox()

virtual G4ThreeVector G4VChannelingFastSimCrystalData::CoordinatesFromLatticeToBox ( const G4ThreeVector & pos)
pure virtual

calculate the coordinates in the Box reference system (connected with the bounding box of the volume)

Implemented in G4ChannelingFastSimCrystalData.

◆ CoulombAtomicScattering()

G4ThreeVector G4VChannelingFastSimCrystalData::CoulombAtomicScattering ( G4double effectiveStep,
G4double step,
G4int ielement )

multiple and single scattering on screened potential

Definition at line 340 of file G4VChannelingFastSimCrystalData.cc.

344{
345 G4double tx = 0.;//horizontal scattering angle
346 G4double ty = 0.;//vertical scattering angle
347
348 G4double ksi=0.1;
349
350// calculation of the teta2-minimal possible angle of a single scattering
351 G4double e1=fK2[ielement]*effectiveStep; //for high speed of a program
352// (real formula is (4*pi*fN0*wpl(x)*dz*fZ1*zz2*alpha*hdc/fPV)**2)
353 G4double teta122=fTetamax12[ielement]/(ksi*fTetamax12[ielement]/e1+1.);
354 // teta122=fTeta12+teta22=teta1^2+teta2^2
355
356 G4double teta22;
357 G4double t;
358// if the angle of a single scattering is less teta1 - minimal possible
359// angle of coulomb scattering defining by the electron shielding than
360// multiple scattering by both nuclei and electrons and electrons will not
361// occur => minimal possible angle of a single scattering is equal to teta1
362 if (teta122<=fTeta12[ielement]*1.000125)
363 {
364 teta22=0.;
365 teta122=fTeta12[ielement];
366 }
367 else
368 {
369 teta22=teta122-fTeta12[ielement];
370 G4double aa=teta22/fTeta12[ielement];
371 G4double aa1=1.+aa;
372
373// crystal, with scattering suppression
374 G4double tetamsi=e1*(G4Log(aa1)+
375 (1.-std::exp(-aa*fBB[ielement]))/aa1+
376 fBBDEXP[ielement]*
377 (expint(fBB[ielement]*aa1)-fE1XBbb[ielement]));
378
379// sumilation of multiple coulomb scattering by nuclei and electrons
380// for high speed of a program, real formula is
381// 4*pi*fN0*wpl(x)*dz*(fZ1*zz2*alpha*hdc/fPV)**2*
382// *(ln(1+a)+(1-exp(-a*b))/(1+a)+(1+b)*exp(b)*(E1XB(b*(1+a))-E1XB(b)))
383
384 ksi=G4UniformRand();
385 t=std::sqrt(-tetamsi*G4Log(ksi));
386
387 ksi=G4UniformRand();
388
389 tx+=t*std::cos(CLHEP::twopi*ksi);
390 ty+=t*std::sin(CLHEP::twopi*ksi);
391
392 }
393// simulation of single coulomb scattering by nuclei (with screened potential)
394 G4double zss=0.;
395 G4double dzss=step;
396
397// (calculation of a distance, at which another single scattering can happen)
398 ksi=G4UniformRand();
399
400 zss=-G4Log(ksi)*step/(e1*(1./teta122-1./fTetamax12[ielement]));
401 G4double tt;
402
403// At some step several single scattering can occur.
404// So,if the distance of the next scattering is less than the step,
405// another scattering can occur. If the distance of the next scattering
406// is less than the difference between the step and the distance of
407// the previous scattering, another scattering can occur. And so on, and so on.
408// In the cycle we simulate each of them. The cycle is finished, when
409// the remaining part of step is less than a distance of the next single scattering.
410//********************************************
411// if at a step a single scattering occurs
412 while (zss<dzss)
413 {
414
415// simulation by Monte-Carlo of angles of single scattering
416 ksi=G4UniformRand();
417
418 tt=fTetamax12[ielement]/(1.+ksi*(fTetamax2[ielement]-teta22)/teta122)-
419 fTeta12[ielement];
420
421 ksi=G4UniformRand();
422
423// suppression of incoherent scattering by the atomic correlations in crystals
424 t=fPzu11[ielement]*tt;
425 t=std::exp(-t);
426
427 if (t<ksi) //if scattering takes place
428 {
429 //scattering angle
430 t=std::sqrt(tt);
431 ksi=G4UniformRand();
432
433 tx+=t*std::cos(CLHEP::twopi*ksi);
434 ty+=t*std::sin(CLHEP::twopi*ksi);
435 }
436
437 dzss-=zss;
438// (calculation of a distance, at which another single scattering can happen)
439 ksi=G4UniformRand();
440
441 zss=-G4Log(ksi)*step/(e1*(1./teta122-1./fTetamax12[ielement]));
442 }
443//********************************************
444 return G4ThreeVector(tx,ty,0.);
445}
G4double G4Log(G4double x)
Definition G4Log.hh:227
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
#define G4UniformRand()
Definition Randomize.hh:52

◆ CoulombElectronScattering()

G4ThreeVector G4VChannelingFastSimCrystalData::CoulombElectronScattering ( G4double eMinIonization,
G4double electronDensity,
G4double step )

multiple and single scattering on electrons

Definition at line 449 of file G4VChannelingFastSimCrystalData.cc.

453{
454
455 G4double zss=0.;
456 G4double dzss=step;
457 G4double ksi = 0.;
458
459 G4double tx = 0.;//horizontal scattering angle
460 G4double ty = 0.;//vertical scattering angle
461 G4double eloss = 0.;//energy loss
462
463 // eMinIonization - minimal energy transfered to electron
464 // a cut to reduce the number of calls of electron scattering
465 // is needed only at low density regions, in many cases does not do anything at all
466 if (eMinIonization<0.5*eV){eMinIonization=0.5*eV;}
467
468 // single scattering on electrons routine
469 if ((eMinIonization<fTmax)&&(electronDensity>DBL_EPSILON))
470 {
471
472// (calculation of a distance, at which another single scattering can happen)
473// simulation of scattering length (by the same way single scattering by nucleus
474 ksi=G4UniformRand();
475
476 zss=-1.0*G4Log(ksi)/(fK3*electronDensity)/(1./eMinIonization-1./fTmax);
477
478//********************************************
479// if at a step a single scattering occur
480 while (zss<dzss)
481 {
482// simulation by Monte-Carlo of angles of single scattering
483 ksi=G4UniformRand();
484
485// energy transfered to electron
486 G4double e1=eMinIonization/(1.-ksi*(1.-eMinIonization/fTmax));
487
488// scattering angle
489 G4double t=0;
490 if(fTmax-e1>DBL_EPSILON) //to be sure e1<fTmax
491 {
492 t=std::sqrt(2.*CLHEP::electron_mass_c2*e1*(1-e1/fTmax))/fPz;
493 }
494
495 // energy losses
496 eloss=e1;
497 ksi=G4UniformRand();
498
499 tx+=t*std::cos(CLHEP::twopi*ksi);
500 ty+=t*std::sin(CLHEP::twopi*ksi);
501
502 dzss-=zss;
503// (calculation of a distance, at which another single scattering can happen)
504// simulation of scattering length
505// (by the same way single scattering by nucleus
506 ksi=G4UniformRand();
507
508 zss=-1.0*G4Log(ksi)/(fK3*electronDensity)/(1./eMinIonization-1./fTmax);
509 }
510//********************************************
511 }
512 return G4ThreeVector(tx,ty,eloss);
513}
#define DBL_EPSILON
Definition templates.hh:66

◆ ElectronDensity()

G4double G4VChannelingFastSimCrystalData::ElectronDensity ( G4double x,
G4double y )
inline

electron density function

Definition at line 71 of file G4VChannelingFastSimCrystalData.hh.

72 {
73 G4double nel0=fElectronDensity->GetIF(x,y);
74 if(nel0<0.) {nel0=0.;}//exception, errors of interpolation functions
75 return nel0;
76 }
G4ChannelingFastSimInterpolation * fElectronDensity

◆ Ex()

G4double G4VChannelingFastSimCrystalData::Ex ( G4double x,
G4double y )
inline

electric fields produced by crystal lattice

Definition at line 67 of file G4VChannelingFastSimCrystalData.hh.

67{return (fElectricFieldX->GetIF(x,y))*(-fZ2/fPV);}
G4ChannelingFastSimInterpolation * fElectricFieldX
classes containing interpolation coefficients

◆ Ey()

G4double G4VChannelingFastSimCrystalData::Ey ( G4double x,
G4double y )
inline

Definition at line 68 of file G4VChannelingFastSimCrystalData.hh.

68{return (fElectricFieldY->GetIF(x,y))*(-fZ2/fPV);}
G4ChannelingFastSimInterpolation * fElectricFieldY

◆ GetBendingAngle()

G4double G4VChannelingFastSimCrystalData::GetBendingAngle ( )
inline

get bending angle of the crystal planes/axes (default BendingAngle=0 => straight crystal);

Definition at line 106 of file G4VChannelingFastSimCrystalData.hh.

◆ GetBeta()

G4double G4VChannelingFastSimCrystalData::GetBeta ( )
inline

get particle velocity/c

Definition at line 99 of file G4VChannelingFastSimCrystalData.hh.

99{return fBeta;}

◆ GetCorrectionZ()

G4double G4VChannelingFastSimCrystalData::GetCorrectionZ ( )
inline

return correction of the longitudinal coordinate (along current plane/axis vs "central plane/axis")

Definition at line 176 of file G4VChannelingFastSimCrystalData.hh.

◆ GetCurv()

G4double G4VChannelingFastSimCrystalData::GetCurv ( G4double z)
inline

get crystal curvature for crystalline undulator the curvature is a function, otherwise it's a constant

Definition at line 114 of file G4VChannelingFastSimCrystalData.hh.

◆ GetCUtetax()

G4double G4VChannelingFastSimCrystalData::GetCUtetax ( G4double z)
inline

◆ GetCUx()

G4double G4VChannelingFastSimCrystalData::GetCUx ( G4double z)
inline

◆ GetLindhardAngle() [1/2]

G4double G4VChannelingFastSimCrystalData::GetLindhardAngle ( )

Calculate the value of the Lindhard angle (!!! the value for a straight crystal)

Definition at line 293 of file G4VChannelingFastSimCrystalData.cc.

294{
295 return fTetaL; //return the Lindhard angle value calculated in SetParticleProperties
296}

Referenced by GetMaxSimulationStep().

◆ GetLindhardAngle() [2/2]

G4double G4VChannelingFastSimCrystalData::GetLindhardAngle ( G4double etotal,
G4double mass,
G4double charge )

Calculate the value of the Lindhard angle (!!! the value for a straight crystal)

Definition at line 282 of file G4VChannelingFastSimCrystalData.cc.

285{
286 G4double pv0 = etotal-mass*mass/etotal;
287 return std::sqrt(2*std::abs(charge)*fVmax/pv0); //Calculate the value of the Lindhard angle
288 //(!!! the value for a straight crystal)
289}

◆ GetMaxSimulationStep()

G4double G4VChannelingFastSimCrystalData::GetMaxSimulationStep ( G4double etotal,
G4double mass,
G4double charge )

Calculate maximal simulation step (standard value for channeling particles)

Definition at line 330 of file G4VChannelingFastSimCrystalData.cc.

333{
334 //standard value of step for channeling particles which is the maximal possible step
335 return fChangeStep/GetLindhardAngle(etotal, mass, charge);
336}
G4double GetLindhardAngle()
Calculate the value of the Lindhard angle (!!! the value for a straight crystal)

◆ GetMiscutAngle()

G4double G4VChannelingFastSimCrystalData::GetMiscutAngle ( )
inline

fBendingAngle MAY BE NOT THE SAME AS THE BENDING ANGLE OF THE CRYSTAL VOLUME: THE VOLUME OF A BENT CRYSTAL MAY BE G4Box, while the planes/axes inside may be bent

Definition at line 110 of file G4VChannelingFastSimCrystalData.hh.

◆ GetModel()

G4int G4VChannelingFastSimCrystalData::GetModel ( )
inline

Definition at line 102 of file G4VChannelingFastSimCrystalData.hh.

102{return iModel;}//=1 for planes, =2 for axes

◆ GetNelements()

G4int G4VChannelingFastSimCrystalData::GetNelements ( )
inline

Definition at line 101 of file G4VChannelingFastSimCrystalData.hh.

101{return fNelements;}
G4int fNelements
values related to the crystal lattice

◆ GetSimulationStep()

G4double G4VChannelingFastSimCrystalData::GetSimulationStep ( G4double tx,
G4double ty )

Calculate simulation step (standard value for channeling particles and reduced value for overbarrier particles)

Definition at line 300 of file G4VChannelingFastSimCrystalData.cc.

301{
302 G4double simulationstep;
303 //find angle of particle w.r.t. the plane or axis
304 G4double angle=0.;
305 if (iModel==1)//1D model
306 {
307 angle = std::abs(tx);
308 }
309 else if (iModel==2)//2D model
310 {
311 angle = std::sqrt(tx*tx+ty*ty);
312 }
313
314 //compare this angle with the Lindhard angle
315 if (angle<fTetaL)
316 {
317 simulationstep = fChannelingStep;
318 }
319 else
320 {
321 simulationstep = fChangeStep;
322 if (angle > 0.0) { simulationstep /= angle; }
323 }
324
325 return simulationstep;
326}

◆ IonizationLosses()

G4double G4VChannelingFastSimCrystalData::IonizationLosses ( G4double dz,
G4int ielement )

ionization losses

Definition at line 517 of file G4VChannelingFastSimCrystalData.cc.

519{
520 //amorphous part of ionization losses
521
522 G4double elosses = 0.;
523 // 1/2 already taken into account in fKD
524
525 G4double loge = G4Log(fMe2Gamma*fGamma*fV2/fI0[ielement]);
526 G4double delta= 2*(G4Log(fBeta*fGamma)+fLogPlasmaEdI0[ielement]-0.5);
527 if(delta<0){delta=0;}
528 loge-=delta;
529 if(fParticleName=="e-")
530 {
531 loge+=(-G4Log(2.) + 1
532 -(2*fGamma - 1)/fGamma/fGamma*G4Log(2.) +
533 1/8*((fGamma - 1)/fGamma)*((fGamma - 1)/fGamma));
534 }
535 else if(fParticleName=="e+")
536 {
537 loge+=(-fV2/12*(11 + 14/(fGamma + 1) + 10/(fGamma + 1)/(fGamma + 1) +
538 4/(fGamma + 1)/(fGamma + 1)/(fGamma + 1)));
539 }
540 else
541 {
542 loge-=fV2;
543 }
544 elosses=fZ2*fZ2*fKD[ielement]/fV2*loge*dz;
545
546 return elosses;}

◆ MinIonizationEnergy()

G4double G4VChannelingFastSimCrystalData::MinIonizationEnergy ( G4double x,
G4double y )
inline

minimum energy of ionization function

Definition at line 78 of file G4VChannelingFastSimCrystalData.hh.

79 {return fMinIonizationEnergy->GetIF(x,y);}
G4ChannelingFastSimInterpolation * fMinIonizationEnergy

◆ NuclearDensity()

G4double G4VChannelingFastSimCrystalData::NuclearDensity ( G4double x,
G4double y,
G4int ielement )
inline

nuclear density function (normalized to average nuclear density)

Definition at line 81 of file G4VChannelingFastSimCrystalData.hh.

82 {return std::abs(fNucleiDensity[ielement]->GetIF(x,y));}
std::vector< G4ChannelingFastSimInterpolation * > fNucleiDensity

◆ SetBendingAngle()

void G4VChannelingFastSimCrystalData::SetBendingAngle ( G4double tetab,
const G4LogicalVolume * crystallogic )

set bending angle of the crystal planes/axes (default fBendingAngle=0 => straight crystal); only non-negative values! crystal is bent in the positive direction of x

Definition at line 71 of file G4VChannelingFastSimCrystalData.cc.

73{
74 G4int crystalID = crystallogic->GetInstanceID();
75
76 //set the bending angle for this logical volume
77 fMapBendingAngle[crystalID]=tetab;
78
79 G4ThreeVector limboxmin;//minimal limits of the box bounding the logical volume
80 G4ThreeVector limboxmax;//maximal limits of the box bounding the logical volume
81 //save the limits of the box bounding the logical volume
82 crystallogic->GetSolid()->BoundingLimits(limboxmin,limboxmax);
83
84 //bounding box half dimensions
85 fHalfDimBoundingBox = (limboxmax-limboxmin)/2.;
86
87 G4double lcr = limboxmax.getZ()-limboxmin.getZ();//crystal thickness
88
89 fBendingAngle=std::abs(tetab);
90 if (fBendingAngle<0.000001)//no bending less then 1 urad
91 {
93 {
94 G4cout << "Channeling model: volume " << crystallogic->GetName() << G4endl;
95 G4cout << "Warning: bending angle is lower than 1 urad => set to 0" << G4endl;
96 }
97
98 fBent=0;
100 fBendingR=0.;//just for convenience (infinity in reality)
101 fBending2R=0.;
103 fCurv=0.;
104
105 fCorrectionZ = 1.;
106 }
107 else
108 {
109 fBent=1;
113 fCurv=1./fBendingR;
114
115 if (tetab<0.)
116 {
117 G4cout << "Channeling model: volume " << crystallogic->GetName() << G4endl;
118 G4cout << "Warning: bending angle is negative => set to be positive" << G4endl;
119 }
120 }
121
122}
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
double getZ() const
G4VSolid * GetSolid() const
G4int GetInstanceID() const
const G4String & GetName() const
G4ThreeVector fHalfDimBoundingBox
values related to the crystal geometry
virtual void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition G4VSolid.cc:680

Referenced by SetCUParameters(), and SetGeometryParameters().

◆ SetCrystallineUndulatorParameters()

void G4VChannelingFastSimCrystalData::SetCrystallineUndulatorParameters ( G4double amplitude,
G4double period,
G4double phase,
const G4LogicalVolume * crystallogic )

set crystalline undulator parameters: amplitude, period and phase (default: all 3 value = 0) function to use in Detector Construction

Definition at line 148 of file G4VChannelingFastSimCrystalData.cc.

153{
154 if (amplitude<DBL_EPSILON||period<DBL_EPSILON)
155 {
156 amplitude = 0.;
157 period=0.;
158 phase=0.;
159 G4cout << "Channeling model: volume " << crystallogic->GetName() << G4endl;
160 G4cout << "Warning: The crystalline undulator parameters are out of range "
161 "=> the crystalline undulator mode switched off" << G4endl;
162 }
163
164 SetCUParameters(G4ThreeVector(amplitude,period,phase),crystallogic);
165}
void SetCUParameters(const G4ThreeVector &amplitudePeriodPhase, const G4LogicalVolume *crystallogic)

◆ SetCUParameters()

void G4VChannelingFastSimCrystalData::SetCUParameters ( const G4ThreeVector & amplitudePeriodPhase,
const G4LogicalVolume * crystallogic )

set crystalline undulator parameters (internal function of the model) for convenience we put amplitude, period and phase in a G4ThreeVector

Definition at line 169 of file G4VChannelingFastSimCrystalData.cc.

172{
173 G4int crystalID = crystallogic->GetInstanceID();
174
175 //set the crystalline undulator parameters for this logical volume
176 fMapCUAmplitudePeriodPhase[crystalID]=amplitudePeriodPhase;
177 fCUAmplitude=amplitudePeriodPhase.x();
178 G4double period = amplitudePeriodPhase.y();
179 fCUPhase = amplitudePeriodPhase.z();
180
181 //if the amplidude of the crystalline undulator is 0 => no undulator
183 {
184 //crystalline undulator flag
185 fCU = true;
186
187 fCUK = CLHEP::twopi/period;
188
190 {
191 //bent and periodically bent crystal are not compatible
192 SetBendingAngle(0,crystallogic);
193
194 G4cout << "Channeling model: volume " << crystallogic->GetName() << G4endl;
195 G4cout << "Warning: crystalline undulator is not compatible with "
196 "a bent crystal mode => setting bending angle to 0." << G4endl;
197 }
198 }
199 else
200 {
201 fCU = false;
202 fCUAmplitude = 0.;
203 fCUK = 0.;
204 fCUPhase = 0.;
205 fMapCUAmplitudePeriodPhase[crystalID] = G4ThreeVector(0.,0.,0.);
206 }
207
208 fCUK2 = fCUK*fCUK;
210}
double z() const
double x() const
double y() const
void SetBendingAngle(G4double tetab, const G4LogicalVolume *crystallogic)

Referenced by SetCrystallineUndulatorParameters(), and SetGeometryParameters().

◆ SetGeometryParameters()

void G4VChannelingFastSimCrystalData::SetGeometryParameters ( const G4LogicalVolume * crystallogic)

set geometry parameters from current logical volume

Definition at line 47 of file G4VChannelingFastSimCrystalData.cc.

49{
50 G4int crystalID = crystallogic->GetInstanceID();
51
52 //set bending angle if the it exists in the list, otherwise default = 0
53 (fMapBendingAngle.count(crystalID) > 0)
54 ? SetBendingAngle(fMapBendingAngle[crystalID],crystallogic)
55 : SetBendingAngle(0.,crystallogic);
56
57 //set miscut angle if the it exists in the list, otherwise default = 0
58 (fMapMiscutAngle.count(crystalID) > 0)
59 ? SetMiscutAngle(fMapMiscutAngle[crystalID],crystallogic)
60 : SetMiscutAngle(0.,crystallogic);
61
62 //set crystalline undulator parameters if they exist in the list,
63 //otherwise default = G4ThreeVector(0,0,0).
64 (fMapCUAmplitudePeriodPhase.count(crystalID) > 0)
65 ? SetCUParameters(fMapCUAmplitudePeriodPhase[crystalID],crystallogic)
66 : SetCUParameters(G4ThreeVector(0.,0.,0.),crystallogic);
67}
void SetMiscutAngle(G4double tetam, const G4LogicalVolume *crystallogic)

◆ SetMaterialProperties()

virtual void G4VChannelingFastSimCrystalData::SetMaterialProperties ( const G4Material * crystal,
const G4String & lattice,
const G4String & filePath )
pure virtual

find and upload crystal lattice input files, calculate all the basic values (to do only once)

Implemented in G4ChannelingFastSimCrystalData.

◆ SetMiscutAngle()

void G4VChannelingFastSimCrystalData::SetMiscutAngle ( G4double tetam,
const G4LogicalVolume * crystallogic )

fBendingAngle MAY BE NOT THE SAME AS THE BENDING ANGLE OF THE CRYSTAL VOLUME THE VOLUME OF A BENT CRYSTAL MAY BE G4Box, while the planes/axes inside may be bent set miscut angle (default fMiscutAngle=0), acceptable range +-1 mrad, otherwise geometry routines may be unstable

Definition at line 126 of file G4VChannelingFastSimCrystalData.cc.

128{
129 G4int crystalID = crystallogic->GetInstanceID();
130
131 //set the bending angle for this logical volume
132 fMapMiscutAngle[crystalID]=tetam;
133
134 // fMiscutAngle>0: rotation of xz coordinate planes clockwise in the xz plane
135 fMiscutAngle=tetam;
136 if (std::abs(tetam)>1.*mrad)
137 {
138 G4cout << "Channeling model: volume " << crystallogic->GetName() << G4endl;
139 G4cout << "Warning: miscut angle is higher than 1 mrad => " << G4endl;
140 G4cout << "coordinate transformation routines may be unstable" << G4endl;
141 }
144}

Referenced by SetGeometryParameters().

◆ SetParticleProperties()

void G4VChannelingFastSimCrystalData::SetParticleProperties ( G4double etotal,
G4double mp,
G4double charge,
const G4String & particleName )

recalculate all the important values (to do both at the trajectory start and after energy loss)

Definition at line 214 of file G4VChannelingFastSimCrystalData.cc.

218{
219 G4double teta1;
220 fZ2=charge;
221 G4double zz22=fZ2*fZ2;
222 fParticleName=particleName;
223
224// particle momentum and energy
225 G4double t=etotal*etotal-mass*mass; // economy of operations
226 fPz=std::sqrt(t); // momentum of particle
227 fPV=t/etotal; // pv
228 fBeta=fPz/etotal; // velocity/c
229 fTetaL = std::sqrt(std::abs(fZ2)*fVmax2/fPV); //Lindhard angle
230 fChannelingStep = fChangeStep/fTetaL; //standard simulation step
231
232// Energy losses
233 fV2 = fBeta*fBeta; // particle (velocity/c)^2
234 fGamma = etotal/mass; // Lorentz factor
235 fMe2Gamma = 2*CLHEP::electron_mass_c2*fGamma;
236// max ionization losses
237 fTmax = fMe2Gamma*fGamma*fV2/
238 (CLHEP::electron_mass_c2/mass*CLHEP::electron_mass_c2/mass +
239 1. + fMe2Gamma/mass);
240// max ionization losses for electrons
241 if(fParticleName=="e-"){fTmax/=2;}
242
243 for(G4int i=0; i<fNelements; i++)
244 {
245
246// minimal scattering angle by coulomb scattering on nuclei
247// defining by shielding by electrons
248// teta1=hdc/(fPz*fRF)*DSQRT(1.13D0+3.76D0*(alpha*fZ1*fZ2/fBeta)**2){ev*cm/(eV*cm)}
249 teta1=fTeta10[i]*std::sqrt(1.13+fK40[i]*zz22/fV2); // /fPz later to speed up
250 // the calculations
251
252// the coefficient for multiple scattering
253 fBB[i]=teta1*teta1*fPu11[i];
254 fE1XBbb[i]=expint(fBB[i]);
255 fBBDEXP[i]=(1.+fBB[i])*std::exp(fBB[i]);
256// necessary for suppression of incoherent scattering
257// by the atomic correlations in crystals for single scattering on nucleus
258// (screened atomic potential): EXP(-(fPz*teta*fU1)**2)=EXP(-fPzu11*teta**2).GE.ksi
259// =>no scattering
260 fPzu11[i]=fPu11[i]*fPz*fPz;
261
262 teta1=teta1/fPz; //
263 fTeta12[i]=teta1*teta1;
264// maximal scattering angle by coulomb scattering on nuclei
265// defining by nucleus radius
266// tetamax=hc/(fPz*1.D-6*fR0*fAN**(1.D0/3.D0))// {Mev*fermi/(MeV*fermi)}
267 G4double tetamax=fTetamax0[i]/fPz;
268 fTetamax2[i]=tetamax*tetamax;
270
271// a coefficient in a formula for scattering (for high speed of simulation)
272// fK2=(fZ2)**2*alphahbarc2*4.*pi*fN0*(fZ1/fPV)**2
273 fK2[i]=fK20[i]*zz22/fPV/fPV;
274 }
275
276// fK3=(fZ2)**2*alphahbarc2*pi/electron_mass_c2/(fV2)**2
277 fK3=fK30*zz22/fV2;
278}
std::vector< G4double > fK20
coefficients necessary for multiple and single coulomb scattering
std::vector< G4double > fTeta10
angles necessary for multiple and single coulomb scattering
std::vector< G4double > fPu11
coefficients for multiple scattering suppression

◆ SetVerbosity()

void G4VChannelingFastSimCrystalData::SetVerbosity ( G4int ver)
inline

Member Data Documentation

◆ fBB

std::vector<G4double> G4VChannelingFastSimCrystalData::fBB
protected

◆ fBBDEXP

std::vector<G4double> G4VChannelingFastSimCrystalData::fBBDEXP
protected

◆ fBending2R

◆ fBendingAngle

G4double G4VChannelingFastSimCrystalData::fBendingAngle =0.
protected

◆ fBendingR

◆ fBendingRsquare

G4double G4VChannelingFastSimCrystalData::fBendingRsquare =0.
protected

◆ fBent

◆ fChangeStep

G4double G4VChannelingFastSimCrystalData::fChangeStep =0
protected

◆ fCorrectionZ

◆ fCosMiscutAngle

◆ fCU

◆ fCUAmplitude

G4double G4VChannelingFastSimCrystalData::fCUAmplitude =0.
protected

Definition at line 243 of file G4VChannelingFastSimCrystalData.hh.

Referenced by GetCUx(), and SetCUParameters().

◆ fCUAmplitudeK

G4double G4VChannelingFastSimCrystalData::fCUAmplitudeK =0.
protected

Definition at line 246 of file G4VChannelingFastSimCrystalData.hh.

Referenced by GetCUtetax(), and SetCUParameters().

◆ fCUK

G4double G4VChannelingFastSimCrystalData::fCUK =0.
protected

Definition at line 244 of file G4VChannelingFastSimCrystalData.hh.

Referenced by GetCUtetax(), GetCUx(), and SetCUParameters().

◆ fCUK2

G4double G4VChannelingFastSimCrystalData::fCUK2 =0.
protected

Definition at line 247 of file G4VChannelingFastSimCrystalData.hh.

Referenced by GetCurv(), and SetCUParameters().

◆ fCUPhase

G4double G4VChannelingFastSimCrystalData::fCUPhase =0.
protected

Definition at line 245 of file G4VChannelingFastSimCrystalData.hh.

Referenced by GetCUtetax(), GetCUx(), and SetCUParameters().

◆ fCurv

G4double G4VChannelingFastSimCrystalData::fCurv =0.
protected

◆ fE1XBbb

std::vector<G4double> G4VChannelingFastSimCrystalData::fE1XBbb
protected

◆ fElectricFieldX

G4ChannelingFastSimInterpolation* G4VChannelingFastSimCrystalData::fElectricFieldX {nullptr}
protected

classes containing interpolation coefficients

Definition at line 207 of file G4VChannelingFastSimCrystalData.hh.

207{nullptr};

Referenced by Ex(), and G4ChannelingFastSimCrystalData::SetMaterialProperties().

◆ fElectricFieldY

G4ChannelingFastSimInterpolation* G4VChannelingFastSimCrystalData::fElectricFieldY {nullptr}
protected

Definition at line 209 of file G4VChannelingFastSimCrystalData.hh.

209{nullptr};

Referenced by Ey(), and G4ChannelingFastSimCrystalData::SetMaterialProperties().

◆ fElectronDensity

G4ChannelingFastSimInterpolation* G4VChannelingFastSimCrystalData::fElectronDensity {nullptr}
protected

◆ fHalfDimBoundingBox

G4ThreeVector G4VChannelingFastSimCrystalData::fHalfDimBoundingBox
protected

◆ fI0

std::vector<G4double> G4VChannelingFastSimCrystalData::fI0
protected

◆ fK2

std::vector<G4double> G4VChannelingFastSimCrystalData::fK2
protected

◆ fK20

std::vector<G4double> G4VChannelingFastSimCrystalData::fK20
protected

coefficients necessary for multiple and single coulomb scattering

Definition at line 279 of file G4VChannelingFastSimCrystalData.hh.

Referenced by G4ChannelingFastSimCrystalData::SetMaterialProperties(), and SetParticleProperties().

◆ fK3

G4double G4VChannelingFastSimCrystalData::fK3 =0
protected

◆ fK30

G4double G4VChannelingFastSimCrystalData::fK30 =0
protected

◆ fK40

std::vector<G4double> G4VChannelingFastSimCrystalData::fK40
protected

◆ fKD

std::vector<G4double> G4VChannelingFastSimCrystalData::fKD
protected

◆ fLogPlasmaEdI0

std::vector<G4double> G4VChannelingFastSimCrystalData::fLogPlasmaEdI0
protected

◆ fMinIonizationEnergy

G4ChannelingFastSimInterpolation* G4VChannelingFastSimCrystalData::fMinIonizationEnergy {nullptr}
protected

◆ fMiscutAngle

G4double G4VChannelingFastSimCrystalData::fMiscutAngle = 0.
protected

◆ fNelements

G4int G4VChannelingFastSimCrystalData::fNelements =1
protected

values related to the crystal lattice

Definition at line 250 of file G4VChannelingFastSimCrystalData.hh.

Referenced by GetNelements(), G4ChannelingFastSimCrystalData::SetMaterialProperties(), and SetParticleProperties().

◆ fNucleiDensity

std::vector<G4ChannelingFastSimInterpolation*> G4VChannelingFastSimCrystalData::fNucleiDensity
protected

◆ fPu11

std::vector<G4double> G4VChannelingFastSimCrystalData::fPu11
protected

coefficients for multiple scattering suppression

Definition at line 290 of file G4VChannelingFastSimCrystalData.hh.

Referenced by G4ChannelingFastSimCrystalData::SetMaterialProperties(), and SetParticleProperties().

◆ fPzu11

std::vector<G4double> G4VChannelingFastSimCrystalData::fPzu11
protected

◆ fRF

std::vector<G4double> G4VChannelingFastSimCrystalData::fRF
protected

◆ fSinMiscutAngle

◆ fTeta10

std::vector<G4double> G4VChannelingFastSimCrystalData::fTeta10
protected

angles necessary for multiple and single coulomb scattering

Definition at line 270 of file G4VChannelingFastSimCrystalData.hh.

Referenced by G4ChannelingFastSimCrystalData::SetMaterialProperties(), and SetParticleProperties().

◆ fTeta12

std::vector<G4double> G4VChannelingFastSimCrystalData::fTeta12
protected

◆ fTetamax0

std::vector<G4double> G4VChannelingFastSimCrystalData::fTetamax0
protected

◆ fTetamax12

std::vector<G4double> G4VChannelingFastSimCrystalData::fTetamax12
protected

◆ fTetamax2

std::vector<G4double> G4VChannelingFastSimCrystalData::fTetamax2
protected

◆ fVerbosity

G4int G4VChannelingFastSimCrystalData::fVerbosity = 1
protected

◆ fVmax

G4double G4VChannelingFastSimCrystalData::fVmax =0
protected

◆ fVmax2

G4double G4VChannelingFastSimCrystalData::fVmax2 =0
protected

◆ fVMinCrystal

G4double G4VChannelingFastSimCrystalData::fVMinCrystal =0
protected

◆ iModel

G4int G4VChannelingFastSimCrystalData::iModel =1
protected

The documentation for this class was generated from the following files: