ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4VXTRenergyLoss.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4VXTRenergyLoss.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 //
28 // History:
29 // 2001-2002 R&D by V.Grichine
30 // 19.06.03 V. Grichine, modifications in BuildTable for the integration
31 // in respect of angle: range is increased, accuracy is
32 // improved
33 // 28.07.05, P.Gumplinger add G4ProcessType to constructor
34 // 28.09.07, V.Ivanchenko general cleanup without change of algorithms
35 //
36 
37 #include "G4VXTRenergyLoss.hh"
38 
39 #include "G4Timer.hh"
40 #include "G4PhysicalConstants.hh"
41 #include "G4SystemOfUnits.hh"
42 
43 #include "G4Poisson.hh"
44 #include "G4MaterialTable.hh"
45 #include "G4VDiscreteProcess.hh"
46 #include "G4VParticleChange.hh"
47 #include "G4VSolid.hh"
48 
49 #include "G4RotationMatrix.hh"
50 #include "G4ThreeVector.hh"
51 #include "G4AffineTransform.hh"
52 #include "G4SandiaTable.hh"
53 
54 #include "G4PhysicsVector.hh"
55 #include "G4PhysicsFreeVector.hh"
56 #include "G4PhysicsLinearVector.hh"
57 #include "G4EmProcessSubType.hh"
58 
60 //
61 // Constructor, destructor
62 
64  G4Material* foilMat,G4Material* gasMat,
66  G4int n,const G4String& processName,
67  G4ProcessType type) :
68  G4VDiscreteProcess(processName, type),
69  fGammaCutInKineticEnergy(nullptr),
70  fGammaTkinCut(0.0),
71  fAngleDistrTable(nullptr),
72  fEnergyDistrTable(nullptr),
73  fPlatePhotoAbsCof(nullptr),
74  fGasPhotoAbsCof(nullptr),
75  fAngleForEnergyTable(nullptr)
76 {
77  verboseLevel = 1;
79 
80  fPtrGamma = nullptr;
83 
84  // Initialization of local constants
85  fTheMinEnergyTR = 1.0*keV;
86  fTheMaxEnergyTR = 100.0*keV;
87 
88  fTheMaxAngle = 1.0e-2; // 100 mrad
89 
90  fTheMinAngle = 2.5e-5;
91 
92  fBinTR = 200; // 100; // 50;
93 
94  fMinProtonTkin = 100.0*GeV;
95  fMaxProtonTkin = 100.0*TeV;
96  fTotBin = 50; // 100; //
97 
98  // Proton energy vector initialization
101  fTotBin );
102 
105  fBinTR );
106 
108 
110 
111  fEnvelope = anEnvelope;
112 
113  fPlateNumber = n;
114  if(verboseLevel > 0)
115  G4cout<<"### G4VXTRenergyLoss: the number of TR radiator plates = "
116  <<fPlateNumber<<G4endl;
117  if(fPlateNumber == 0)
118  {
119  G4Exception("G4VXTRenergyLoss::G4VXTRenergyLoss()","VXTRELoss01",
120  FatalException,"No plates in X-ray TR radiator");
121  }
122  // default is XTR dEdx, not flux after radiator
123 
124  fExitFlux = false;
125  fAngleRadDistr = true; // false;
126  fCompton = false;
127 
128  fLambda = DBL_MAX;
129  // Mean thicknesses of plates and gas gaps
130 
131  fPlateThick = a;
132  fGasThick = b;
134  if(verboseLevel > 0)
135  G4cout<<"total radiator thickness = "<<fTotalDist/cm<<" cm"<<G4endl;
136 
137  // index of plate material
138  fMatIndex1 = foilMat->GetIndex();
139  if(verboseLevel > 0)
140  G4cout<<"plate material = "<<foilMat->GetName()<<G4endl;
141 
142  // index of gas material
143  fMatIndex2 = gasMat->GetIndex();
144  if(verboseLevel > 0)
145  G4cout<<"gas material = "<<gasMat->GetName()<<G4endl;
146 
147  // plasma energy squared for plate material
148 
150  // fSigma1 = (20.9*eV)*(20.9*eV);
151  if(verboseLevel > 0)
152  G4cout<<"plate plasma energy = "<<std::sqrt(fSigma1)/eV<<" eV"<<G4endl;
153 
154  // plasma energy squared for gas material
155 
157  if(verboseLevel > 0)
158  G4cout<<"gas plasma energy = "<<std::sqrt(fSigma2)/eV<<" eV"<<G4endl;
159 
160  // Compute cofs for preparation of linear photo absorption
161 
164 
166 }
167 
169 
171 {
172  if(fEnvelope) delete fEnvelope;
173  delete fProtonEnergyVector;
174  delete fXTREnergyVector;
175  if(fEnergyDistrTable) {
177  delete fEnergyDistrTable;
178  }
179  if(fAngleRadDistr) {
181  delete fAngleDistrTable;
182  }
185  delete fAngleForEnergyTable;
186  }
187 }
188 
190 //
191 // Returns condition for application of the model depending on particle type
192 
193 
195 {
196  return ( particle.GetPDGCharge() != 0.0 );
197 }
198 
200 //
201 // Calculate step size for XTR process inside raaditor
202 
204  G4double, // previousStepSize,
206 {
207  G4int iTkin, iPlace;
208  G4double lambda, sigma, kinEnergy, mass, gamma;
209  G4double charge, chargeSq, massRatio, TkinScaled;
210  G4double E1,E2,W,W1,W2;
211 
212  *condition = NotForced;
213 
214  if( aTrack.GetVolume()->GetLogicalVolume() != fEnvelope ) lambda = DBL_MAX;
215  else
216  {
217  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
218  kinEnergy = aParticle->GetKineticEnergy();
219  mass = aParticle->GetDefinition()->GetPDGMass();
220  gamma = 1.0 + kinEnergy/mass;
221  if(verboseLevel > 1)
222  {
223  G4cout<<" gamma = "<<gamma<<"; fGamma = "<<fGamma<<G4endl;
224  }
225 
226  if ( std::fabs( gamma - fGamma ) < 0.05*gamma ) lambda = fLambda;
227  else
228  {
229  charge = aParticle->GetDefinition()->GetPDGCharge();
230  chargeSq = charge*charge;
231  massRatio = proton_mass_c2/mass;
232  TkinScaled = kinEnergy*massRatio;
233 
234  for(iTkin = 0; iTkin < fTotBin; iTkin++)
235  {
236  if( TkinScaled < fProtonEnergyVector->GetLowEdgeEnergy(iTkin)) break;
237  }
238  iPlace = iTkin - 1;
239 
240  if(iTkin == 0) lambda = DBL_MAX; // Tkin is too small, neglect of TR photon generation
241  else // general case: Tkin between two vectors of the material
242  {
243  if(iTkin == fTotBin)
244  {
245  sigma = (*(*fEnergyDistrTable)(iPlace))(0)*chargeSq;
246  }
247  else
248  {
249  E1 = fProtonEnergyVector->GetLowEdgeEnergy(iTkin - 1);
251  W = 1.0/(E2 - E1);
252  W1 = (E2 - TkinScaled)*W;
253  W2 = (TkinScaled - E1)*W;
254  sigma = ( (*(*fEnergyDistrTable)(iPlace ))(0)*W1 +
255  (*(*fEnergyDistrTable)(iPlace+1))(0)*W2 )*chargeSq;
256 
257  }
258  if (sigma < DBL_MIN) lambda = DBL_MAX;
259  else lambda = 1./sigma;
260  fLambda = lambda;
261  fGamma = gamma;
262  if(verboseLevel > 1)
263  {
264  G4cout<<" lambda = "<<lambda/mm<<" mm"<<G4endl;
265  }
266  }
267  }
268  }
269  return lambda;
270 }
271 
273 //
274 // Interface for build table from physics list
275 
277 {
278  if(pd.GetPDGCharge() == 0.)
279  {
280  G4Exception("G4VXTRenergyLoss::BuildPhysicsTable", "Notification", JustWarning,
281  "XTR initialisation for neutral particle ?!" );
282  }
284 
285  if ( fAngleRadDistr )
286  {
287  if(verboseLevel > 0)
288  {
289  G4cout<<"Build angle for energy distribution according the current radiator"
290  <<G4endl;
291  }
293  }
294 }
295 
296 
298 //
299 // Build integral energy distribution of XTR photons
300 
302 {
303  G4int iTkin, iTR, iPlace;
304  G4double radiatorCof = 1.0; // for tuning of XTR yield
305  G4double energySum = 0.0;
306 
309 
310  fGammaTkinCut = 0.0;
311 
312  // setting of min/max TR energies
313 
316 
319 
321 
322  G4cout.precision(4);
323  G4Timer timer;
324  timer.Start();
325 
326  if(verboseLevel > 0)
327  {
328  G4cout<<G4endl;
329  G4cout<<"Lorentz Factor"<<"\t"<<"XTR photon number"<<G4endl;
330  G4cout<<G4endl;
331  }
332  for( iTkin = 0; iTkin < fTotBin; iTkin++ ) // Lorentz factor loop
333  {
335  fMaxEnergyTR,
336  fBinTR );
337 
338  fGamma = 1.0 + (fProtonEnergyVector->
339  GetLowEdgeEnergy(iTkin)/proton_mass_c2);
340 
341  fMaxThetaTR = 25.*2500.0/(fGamma*fGamma) ; // theta^2
342 
343  fTheMinAngle = 1.0e-3; // was 5.e-6, e-6 !!!, e-5, e-4
344 
347 
348  energySum = 0.0;
349 
350  energyVector->PutValue(fBinTR-1,energySum);
351 
352  for( iTR = fBinTR - 2; iTR >= 0; iTR-- )
353  {
354  // Legendre96 or Legendre10
355 
356  energySum += radiatorCof*fCofTR*integral.Legendre10(
358  energyVector->GetLowEdgeEnergy(iTR),
359  energyVector->GetLowEdgeEnergy(iTR+1) );
360 
361  energyVector->PutValue(iTR,energySum/fTotalDist);
362  }
363  iPlace = iTkin;
364  fEnergyDistrTable->insertAt(iPlace,energyVector);
365 
366  if(verboseLevel > 0)
367  {
368  G4cout
369  // <<iTkin<<"\t"
370  // <<"fGamma = "
371  <<fGamma<<"\t" // <<" fMaxThetaTR = "<<fMaxThetaTR
372  // <<"sumN = "
373  <<energySum // <<"; sumA = "<<angleSum
374  <<G4endl;
375  }
376  }
377  timer.Stop();
378  G4cout.precision(6);
379  if(verboseLevel > 0)
380  {
381  G4cout<<G4endl;
382  G4cout<<"total time for build X-ray TR energy loss tables = "
383  <<timer.GetUserElapsed()<<" s"<<G4endl;
384  }
385  fGamma = 0.;
386  return;
387 }
388 
390 //
391 // Bank of angle distributions for given energies (slow!)
392 
394 {
395  if( this->GetProcessName() == "TranspRegXTRadiator" ||
396  this->GetProcessName() == "TranspRegXTRmodel" ||
397  this->GetProcessName() == "RegularXTRadiator" ||
398  this->GetProcessName() == "RegularXTRmodel" )
399  {
400  BuildAngleTable();
401  return;
402  }
403  G4int i, iTkin, iTR;
404  G4double angleSum = 0.0;
405 
406 
407  fGammaTkinCut = 0.0;
408 
409  // setting of min/max TR energies
410 
413 
416 
418  fMaxEnergyTR,
419  fBinTR );
420 
422 
423  G4cout.precision(4);
424  G4Timer timer;
425  timer.Start();
426 
427  for( iTkin = 0; iTkin < fTotBin; iTkin++ ) // Lorentz factor loop
428  {
429 
430  fGamma = 1.0 + (fProtonEnergyVector->
431  GetLowEdgeEnergy(iTkin)/proton_mass_c2);
432 
433  fMaxThetaTR = 25*2500.0/(fGamma*fGamma) ; // theta^2
434 
435  fTheMinAngle = 1.0e-3; // was 5.e-6, e-6 !!!, e-5, e-4
436 
439 
441 
442  for( iTR = 0; iTR < fBinTR; iTR++ )
443  {
444  angleSum = 0.0;
445  fEnergy = energyVector->GetLowEdgeEnergy(iTR);
446  G4PhysicsLinearVector* angleVector = new G4PhysicsLinearVector(0.0,
447  fMaxThetaTR,
448  fBinTR );
449 
450  angleVector ->PutValue(fBinTR - 1, angleSum);
451 
452  for( i = fBinTR - 2; i >= 0; i-- )
453  {
454  // Legendre96 or Legendre10
455 
456  angleSum += integral.Legendre10(
458  angleVector->GetLowEdgeEnergy(i),
459  angleVector->GetLowEdgeEnergy(i+1) );
460 
461  angleVector ->PutValue(i, angleSum);
462  }
463  fAngleForEnergyTable->insertAt(iTR, angleVector);
464  }
465  fAngleBank.push_back(fAngleForEnergyTable);
466  }
467  timer.Stop();
468  G4cout.precision(6);
469  if(verboseLevel > 0)
470  {
471  G4cout<<G4endl;
472  G4cout<<"total time for build X-ray TR angle for energy loss tables = "
473  <<timer.GetUserElapsed()<<" s"<<G4endl;
474  }
475  fGamma = 0.;
476  delete energyVector;
477 }
478 
480 //
481 // Build XTR angular distribution at given energy based on the model
482 // of transparent regular radiator
483 
485 {
486  G4int iTkin, iTR;
488 
489  fGammaTkinCut = 0.0;
490 
491  // setting of min/max TR energies
492 
495 
498 
499  G4cout.precision(4);
500  G4Timer timer;
501  timer.Start();
502  if(verboseLevel > 0)
503  {
504  G4cout<<G4endl;
505  G4cout<<"Lorentz Factor"<<"\t"<<"XTR photon number"<<G4endl;
506  G4cout<<G4endl;
507  }
508  for( iTkin = 0; iTkin < fTotBin; iTkin++ ) // Lorentz factor loop
509  {
510 
511  fGamma = 1.0 + (fProtonEnergyVector->
512  GetLowEdgeEnergy(iTkin)/proton_mass_c2);
513 
514  fMaxThetaTR = 25*2500.0/(fGamma*fGamma); // theta^2
515 
516  fTheMinAngle = 1.0e-3; // was 5.e-6, e-6 !!!, e-5, e-4
517 
519  else
520  {
522  }
523 
525 
526  for( iTR = 0; iTR < fBinTR; iTR++ )
527  {
528  // energy = fMinEnergyTR*(iTR+1);
529 
530  energy = fXTREnergyVector->GetLowEdgeEnergy(iTR);
531 
532  G4PhysicsFreeVector* angleVector = GetAngleVector(energy,fBinTR);
533  // G4cout<<G4endl;
534 
535  fAngleForEnergyTable->insertAt(iTR,angleVector);
536  }
537 
538  fAngleBank.push_back(fAngleForEnergyTable);
539  }
540  timer.Stop();
541  G4cout.precision(6);
542  if(verboseLevel > 0)
543  {
544  G4cout<<G4endl;
545  G4cout<<"total time for build XTR angle for given energy tables = "
546  <<timer.GetUserElapsed()<<" s"<<G4endl;
547  }
548  fGamma = 0.;
549 
550  return;
551 }
552 
554 //
555 // Vector of angles and angle integral distributions
556 
558 {
559  G4double theta=0., result, tmp=0., cof1, cof2, cofMin, cofPHC, angleSum = 0.;
560  G4int iTheta, k, /*kMax,*/ kMin;
561 
562  G4PhysicsFreeVector* angleVector = new G4PhysicsFreeVector(n);
563 
564  cofPHC = 4.*pi*hbarc;
565  tmp = (fSigma1 - fSigma2)/cofPHC/energy;
566  cof1 = fPlateThick*tmp;
567  cof2 = fGasThick*tmp;
568 
569  cofMin = energy*(fPlateThick + fGasThick)/fGamma/fGamma;
570  cofMin += (fPlateThick*fSigma1 + fGasThick*fSigma2)/energy;
571  cofMin /= cofPHC;
572 
573  kMin = G4int(cofMin);
574  if (cofMin > kMin) kMin++;
575 
576  //kMax = kMin + fBinTR -1;
577 
578  if(verboseLevel > 2)
579  {
580  G4cout<<"n-1 = "<<n-1<<"; theta = "
581  <<std::sqrt(fMaxThetaTR)*fGamma<<"; tmp = "
582  <<0.
583  <<"; angleSum = "<<angleSum<<G4endl;
584  }
585  // angleVector->PutValue(n-1,fMaxThetaTR, angleSum);
586 
587  for( iTheta = n - 1; iTheta >= 1; iTheta-- )
588  {
589 
590  k = iTheta - 1 + kMin;
591 
592  tmp = pi*fPlateThick*(k + cof2)/(fPlateThick + fGasThick);
593 
594  result = (k - cof1)*(k - cof1)*(k + cof2)*(k + cof2);
595 
596  tmp = std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
597 
598  if( k == kMin && kMin == G4int(cofMin) )
599  {
600  angleSum += 0.5*tmp; // 0.5*std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
601  }
602  else if(iTheta == n-1);
603  else
604  {
605  angleSum += tmp; // std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
606  }
607  theta = std::abs(k-cofMin)*cofPHC/energy/(fPlateThick + fGasThick);
608 
609  if(verboseLevel > 2)
610  {
611  G4cout<<"iTheta = "<<iTheta<<"; k = "<<k<<"; theta = "
612  <<std::sqrt(theta)*fGamma<<"; tmp = "
613  <<tmp // std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result
614  <<"; angleSum = "<<angleSum<<G4endl;
615  }
616  angleVector->PutValue( iTheta, theta, angleSum );
617  }
618  if (theta > 0.)
619  {
620  angleSum += 0.5*tmp;
621  theta = 0.;
622  }
623  if(verboseLevel > 2)
624  {
625  G4cout<<"iTheta = "<<iTheta<<"; theta = "
626  <<std::sqrt(theta)*fGamma<<"; tmp = "
627  <<tmp
628  <<"; angleSum = "<<angleSum<<G4endl;
629  }
630  angleVector->PutValue( iTheta, theta, angleSum );
631 
632  return angleVector;
633 }
634 
636 //
637 // Build XTR angular distribution based on the model of transparent regular radiator
638 
640 {
641  G4int iTkin, iTR, iPlace;
642  G4double radiatorCof = 1.0; // for tuning of XTR yield
643  G4double angleSum;
645 
646  fGammaTkinCut = 0.0;
647 
648  // setting of min/max TR energies
649 
652 
655 
656  G4cout.precision(4);
657  G4Timer timer;
658  timer.Start();
659  if(verboseLevel > 0) {
660  G4cout<<G4endl;
661  G4cout<<"Lorentz Factor"<<"\t"<<"XTR photon number"<<G4endl;
662  G4cout<<G4endl;
663  }
664  for( iTkin = 0; iTkin < fTotBin; iTkin++ ) // Lorentz factor loop
665  {
666 
667  fGamma = 1.0 + (fProtonEnergyVector->
668  GetLowEdgeEnergy(iTkin)/proton_mass_c2);
669 
670  fMaxThetaTR = 25.0/(fGamma*fGamma); // theta^2
671 
672  fTheMinAngle = 1.0e-3; // was 5.e-6, e-6 !!!, e-5, e-4
673 
675  else
676  {
678  }
679  G4PhysicsLinearVector* angleVector = new G4PhysicsLinearVector(0.0,
680  fMaxThetaTR,
681  fBinTR );
682 
683  angleSum = 0.0;
684 
686 
687 
688  angleVector->PutValue(fBinTR-1,angleSum);
689 
690  for( iTR = fBinTR - 2; iTR >= 0; iTR-- )
691  {
692 
693  angleSum += radiatorCof*fCofTR*integral.Legendre96(
695  angleVector->GetLowEdgeEnergy(iTR),
696  angleVector->GetLowEdgeEnergy(iTR+1) );
697 
698  angleVector ->PutValue(iTR,angleSum);
699  }
700  if(verboseLevel > 1) {
701  G4cout
702  // <<iTkin<<"\t"
703  // <<"fGamma = "
704  <<fGamma<<"\t" // <<" fMaxThetaTR = "<<fMaxThetaTR
705  // <<"sumN = "<<energySum // <<"; sumA = "
706  <<angleSum
707  <<G4endl;
708  }
709  iPlace = iTkin;
710  fAngleDistrTable->insertAt(iPlace,angleVector);
711  }
712  timer.Stop();
713  G4cout.precision(6);
714  if(verboseLevel > 0) {
715  G4cout<<G4endl;
716  G4cout<<"total time for build X-ray TR angle tables = "
717  <<timer.GetUserElapsed()<<" s"<<G4endl;
718  }
719  fGamma = 0.;
720 
721  return;
722 }
723 
724 
726 //
727 // The main function which is responsible for the treatment of a particle passage
728 // trough G4Envelope with discrete generation of G4Gamma
729 
731  const G4Step& aStep )
732 {
733  G4int iTkin /*, iPlace*/;
734  G4double energyTR, theta,theta2, phi, dirX, dirY, dirZ;
735 
736 
737  fParticleChange.Initialize(aTrack);
738 
739  if(verboseLevel > 1)
740  {
741  G4cout<<"Start of G4VXTRenergyLoss::PostStepDoIt "<<G4endl;
742  G4cout<<"name of current material = "
743  <<aTrack.GetVolume()->GetLogicalVolume()->GetMaterial()->GetName()<<G4endl;
744  }
745  if( aTrack.GetVolume()->GetLogicalVolume() != fEnvelope )
746  {
747  if(verboseLevel > 0)
748  {
749  G4cout<<"Go out from G4VXTRenergyLoss::PostStepDoIt: wrong volume "<<G4endl;
750  }
751  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
752  }
753  else
754  {
755  G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
756  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
757 
758  // Now we are ready to Generate one TR photon
759 
760  G4double kinEnergy = aParticle->GetKineticEnergy();
761  G4double mass = aParticle->GetDefinition()->GetPDGMass();
762  G4double gamma = 1.0 + kinEnergy/mass;
763 
764  if(verboseLevel > 1 )
765  {
766  G4cout<<"gamma = "<<gamma<<G4endl;
767  }
768  G4double massRatio = proton_mass_c2/mass;
769  G4double TkinScaled = kinEnergy*massRatio;
770  G4ThreeVector position = pPostStepPoint->GetPosition();
771  G4ParticleMomentum direction = aParticle->GetMomentumDirection();
772  G4double startTime = pPostStepPoint->GetGlobalTime();
773 
774  for( iTkin = 0; iTkin < fTotBin; iTkin++ )
775  {
776  if(TkinScaled < fProtonEnergyVector->GetLowEdgeEnergy(iTkin)) break;
777  }
778  //iPlace = iTkin - 1;
779 
780  if(iTkin == 0) // Tkin is too small, neglect of TR photon generation
781  {
782  if( verboseLevel > 0)
783  {
784  G4cout<<"Go out from G4VXTRenergyLoss::PostStepDoIt:iTkin = "<<iTkin<<G4endl;
785  }
786  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
787  }
788  else // general case: Tkin between two vectors of the material
789  {
791 
792  energyTR = GetXTRrandomEnergy(TkinScaled,iTkin);
793 
794  if( verboseLevel > 1)
795  {
796  G4cout<<"energyTR = "<<energyTR/keV<<" keV"<<G4endl;
797  }
798  if ( fAngleRadDistr )
799  {
800  // theta = std::fabs(G4RandGauss::shoot(0.0,pi/gamma));
801 
802  theta2 = GetRandomAngle(energyTR,iTkin);
803 
804  if( theta2 > 0.) theta = std::sqrt(theta2);
805  else theta = 0.; // theta2;
806  }
807  else theta = std::fabs(G4RandGauss::shoot(0.0,pi/gamma));
808 
809  // theta = 0.; // check no spread
810 
811  if( theta >= 0.1 ) theta = 0.1;
812 
813  // G4cout<<" : theta = "<<theta<<endl;
814 
815  phi = twopi*G4UniformRand();
816 
817  dirX = std::sin(theta)*std::cos(phi);
818  dirY = std::sin(theta)*std::sin(phi);
819  dirZ = std::cos(theta);
820 
821  G4ThreeVector directionTR(dirX,dirY,dirZ);
822  directionTR.rotateUz(direction);
823  directionTR.unit();
824 
826  directionTR, energyTR);
827 
828  // A XTR photon is set on the particle track inside the radiator
829  // and is moved to the G4Envelope surface for standard X-ray TR models
830  // only. The case of fExitFlux=true
831 
832  if( fExitFlux )
833  {
834  const G4RotationMatrix* rotM = pPostStepPoint->GetTouchable()->GetRotation();
835  G4ThreeVector transl = pPostStepPoint->GetTouchable()->GetTranslation();
837  transform.Invert();
838  G4ThreeVector localP = transform.TransformPoint(position);
839  G4ThreeVector localV = transform.TransformAxis(directionTR);
840 
841  G4double distance = fEnvelope->GetSolid()->DistanceToOut(localP, localV);
842  if(verboseLevel > 1)
843  {
844  G4cout<<"distance to exit = "<<distance/mm<<" mm"<<G4endl;
845  }
846  position += distance*directionTR;
847  startTime += distance/c_light;
848  }
849  G4Track* aSecondaryTrack = new G4Track( aPhotonTR,
850  startTime, position );
851  aSecondaryTrack->SetTouchableHandle(
853  aSecondaryTrack->SetParentID( aTrack.GetTrackID() );
854 
855  fParticleChange.AddSecondary(aSecondaryTrack);
856  fParticleChange.ProposeEnergy(kinEnergy);
857  }
858  }
859  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
860 }
861 
863 //
864 // This function returns the spectral and angle density of TR quanta
865 // in X-ray energy region generated forward when a relativistic
866 // charged particle crosses interface between two materials.
867 // The high energy small theta approximation is applied.
868 // (matter1 -> matter2, or 2->1)
869 // varAngle =2* (1 - std::cos(theta)) or approximately = theta*theta
870 //
871 
873  G4double gamma,
874  G4double varAngle )
875 {
876  G4complex Z1 = GetPlateComplexFZ(energy,gamma,varAngle);
877  G4complex Z2 = GetGasComplexFZ(energy,gamma,varAngle);
878 
879  G4complex zOut = (Z1 - Z2)*(Z1 - Z2)
880  * (varAngle*energy/hbarc/hbarc);
881  return zOut;
882 
883 }
884 
885 
887 //
888 // For photon energy distribution tables. Integrate first over angle
889 //
890 
892 {
893  G4double result = GetStackFactor(fEnergy,fGamma,varAngle);
894  if(result < 0.0) result = 0.0;
895  return result;
896 }
897 
899 //
900 // For second integration over energy
901 
903 {
904  G4int i, iMax = 8;
905  G4double angleSum = 0.0;
906 
907  G4double lim[8] = { 0.0, 0.01, 0.02, 0.05, 0.1, 0.2, 0.5, 1.0 };
908 
909  for( i = 0; i < iMax; i++ ) lim[i] *= fMaxThetaTR;
910 
912 
913  fEnergy = energy;
914  /*
915  if( fAngleRadDistr && ( fEnergy == fEnergyForAngle ) )
916  {
917  fAngleVector ->PutValue(fBinTR - 1, angleSum);
918 
919  for( i = fBinTR - 2; i >= 0; i-- )
920  {
921 
922  angleSum += integral.Legendre10(
923  this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
924  fAngleVector->GetLowEdgeEnergy(i),
925  fAngleVector->GetLowEdgeEnergy(i+1) );
926 
927  fAngleVector ->PutValue(i, angleSum);
928  }
929  }
930  else
931  */
932  {
933  for( i = 0; i < iMax-1; i++ )
934  {
935  angleSum += integral.Legendre96(this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
936  lim[i],lim[i+1]);
937  // result += integral.Legendre10(this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
938  // lim[i],lim[i+1]);
939  }
940  }
941  return angleSum;
942 }
943 
945 //
946 // for photon angle distribution tables
947 //
948 
950 {
951  G4double result = GetStackFactor(energy,fGamma,fVarAngle);
952  if(result < 0) result = 0.0;
953  return result;
954 }
955 
957 //
958 // The XTR angular distribution based on transparent regular radiator
959 
961 {
962  // G4cout<<"angle2 = "<<varAngle<<"; fGamma = "<<fGamma<<G4endl;
963 
964  G4double result;
965  G4double sum = 0., tmp1, tmp2, tmp=0., cof1, cof2, cofMin, cofPHC, energy1, energy2;
966  G4int k, kMax, kMin, i;
967 
968  cofPHC = twopi*hbarc;
969 
970  cof1 = (fPlateThick + fGasThick)*(1./fGamma/fGamma + varAngle);
972 
973  // G4cout<<"cof1 = "<<cof1<<"; cof2 = "<<cof2<<"; cofPHC = "<<cofPHC<<G4endl;
974 
975  cofMin = std::sqrt(cof1*cof2);
976  cofMin /= cofPHC;
977 
978  kMin = G4int(cofMin);
979  if (cofMin > kMin) kMin++;
980 
981  kMax = kMin + 9; // 9; // kMin + G4int(tmp);
982 
983  // G4cout<<"cofMin = "<<cofMin<<"; kMin = "<<kMin<<"; kMax = "<<kMax<<G4endl;
984 
985  for( k = kMin; k <= kMax; k++ )
986  {
987  tmp1 = cofPHC*k;
988  tmp2 = std::sqrt(tmp1*tmp1-cof1*cof2);
989  energy1 = (tmp1+tmp2)/cof1;
990  energy2 = (tmp1-tmp2)/cof1;
991 
992  for(i = 0; i < 2; i++)
993  {
994  if( i == 0 )
995  {
996  if (energy1 > fTheMaxEnergyTR || energy1 < fTheMinEnergyTR) continue;
997 
998  tmp1 = ( energy1*energy1*(1./fGamma/fGamma + varAngle) + fSigma1 )
999  * fPlateThick/(4*hbarc*energy1);
1000  tmp2 = std::sin(tmp1);
1001  tmp = energy1*tmp2*tmp2;
1002  tmp2 = fPlateThick/(4.*tmp1);
1003  tmp1 = hbarc*energy1/( energy1*energy1*(1./fGamma/fGamma + varAngle) + fSigma2 );
1004  tmp *= (tmp1-tmp2)*(tmp1-tmp2);
1005  tmp1 = cof1/(4.*hbarc) - cof2/(4.*hbarc*energy1*energy1);
1006  tmp2 = std::abs(tmp1);
1007 
1008  if(tmp2 > 0.) tmp /= tmp2;
1009  else continue;
1010  }
1011  else
1012  {
1013  if (energy2 > fTheMaxEnergyTR || energy2 < fTheMinEnergyTR) continue;
1014 
1015  tmp1 = ( energy2*energy2*(1./fGamma/fGamma + varAngle) + fSigma1 )
1016  * fPlateThick/(4.*hbarc*energy2);
1017  tmp2 = std::sin(tmp1);
1018  tmp = energy2*tmp2*tmp2;
1019  tmp2 = fPlateThick/(4.*tmp1);
1020  tmp1 = hbarc*energy2/( energy2*energy2*(1./fGamma/fGamma + varAngle) + fSigma2 );
1021  tmp *= (tmp1-tmp2)*(tmp1-tmp2);
1022  tmp1 = cof1/(4.*hbarc) - cof2/(4.*hbarc*energy2*energy2);
1023  tmp2 = std::abs(tmp1);
1024 
1025  if(tmp2 > 0.) tmp /= tmp2;
1026  else continue;
1027  }
1028  sum += tmp;
1029  }
1030  // G4cout<<"k = "<<k<<"; energy1 = "<<energy1/keV<<" keV; energy2 = "<<energy2/keV
1031  // <<" keV; tmp = "<<tmp<<"; sum = "<<sum<<G4endl;
1032  }
1033  result = 4.*pi*fPlateNumber*sum*varAngle;
1034  result /= hbarc*hbarc;
1035 
1036  // old code based on general numeric integration
1037  // fVarAngle = varAngle;
1038  // G4Integrator<G4VXTRenergyLoss,G4double(G4VXTRenergyLoss::*)(G4double)> integral;
1039  // result = integral.Legendre10(this,&G4VXTRenergyLoss::AngleSpectralXTRdEdx,
1040  // fMinEnergyTR,fMaxEnergyTR);
1041  return result;
1042 }
1043 
1044 
1048 //
1049 // Calculates formation zone for plates. Omega is energy !!!
1050 
1052  G4double gamma ,
1053  G4double varAngle )
1054 {
1055  G4double cof, lambda;
1056  lambda = 1.0/gamma/gamma + varAngle + fSigma1/omega/omega;
1057  cof = 2.0*hbarc/omega/lambda;
1058  return cof;
1059 }
1060 
1062 //
1063 // Calculates complex formation zone for plates. Omega is energy !!!
1064 
1066  G4double gamma ,
1067  G4double varAngle )
1068 {
1069  G4double cof, length,delta, real_v, image_v;
1070 
1071  length = 0.5*GetPlateFormationZone(omega,gamma,varAngle);
1072  delta = length*GetPlateLinearPhotoAbs(omega);
1073  cof = 1.0/(1.0 + delta*delta);
1074 
1075  real_v = length*cof;
1076  image_v = real_v*delta;
1077 
1078  G4complex zone(real_v,image_v);
1079  return zone;
1080 }
1081 
1083 //
1084 // Computes matrix of Sandia photo absorption cross section coefficients for
1085 // plate material
1086 
1088 {
1089  const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
1090  const G4Material* mat = (*theMaterialTable)[fMatIndex1];
1092 
1093  return;
1094 }
1095 
1096 
1097 
1099 //
1100 // Returns the value of linear photo absorption coefficient (in reciprocal
1101 // length) for plate for given energy of X-ray photon omega
1102 
1104 {
1105  // G4int i;
1106  G4double omega2, omega3, omega4;
1107 
1108  omega2 = omega*omega;
1109  omega3 = omega2*omega;
1110  omega4 = omega2*omega2;
1111 
1112  const G4double* SandiaCof = fPlatePhotoAbsCof->GetSandiaCofForMaterial(omega);
1113  G4double cross = SandiaCof[0]/omega + SandiaCof[1]/omega2 +
1114  SandiaCof[2]/omega3 + SandiaCof[3]/omega4;
1115  return cross;
1116 }
1117 
1118 
1120 //
1121 // Calculates formation zone for gas. Omega is energy !!!
1122 
1124  G4double gamma ,
1125  G4double varAngle )
1126 {
1127  G4double cof, lambda;
1128  lambda = 1.0/gamma/gamma + varAngle + fSigma2/omega/omega;
1129  cof = 2.0*hbarc/omega/lambda;
1130  return cof;
1131 }
1132 
1133 
1135 //
1136 // Calculates complex formation zone for gas gaps. Omega is energy !!!
1137 
1139  G4double gamma ,
1140  G4double varAngle )
1141 {
1142  G4double cof, length,delta, real_v, image_v;
1143 
1144  length = 0.5*GetGasFormationZone(omega,gamma,varAngle);
1145  delta = length*GetGasLinearPhotoAbs(omega);
1146  cof = 1.0/(1.0 + delta*delta);
1147 
1148  real_v = length*cof;
1149  image_v = real_v*delta;
1150 
1151  G4complex zone(real_v,image_v);
1152  return zone;
1153 }
1154 
1155 
1156 
1158 //
1159 // Computes matrix of Sandia photo absorption cross section coefficients for
1160 // gas material
1161 
1163 {
1164  const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
1165  const G4Material* mat = (*theMaterialTable)[fMatIndex2];
1167  return;
1168 }
1169 
1171 //
1172 // Returns the value of linear photo absorption coefficient (in reciprocal
1173 // length) for gas
1174 
1176 {
1177  G4double omega2, omega3, omega4;
1178 
1179  omega2 = omega*omega;
1180  omega3 = omega2*omega;
1181  omega4 = omega2*omega2;
1182 
1183  const G4double* SandiaCof = fGasPhotoAbsCof->GetSandiaCofForMaterial(omega);
1184  G4double cross = SandiaCof[0]/omega + SandiaCof[1]/omega2 +
1185  SandiaCof[2]/omega3 + SandiaCof[3]/omega4;
1186  return cross;
1187 
1188 }
1189 
1191 //
1192 // Calculates the product of linear cof by formation zone for plate.
1193 // Omega is energy !!!
1194 
1196  G4double gamma ,
1197  G4double varAngle )
1198 {
1199  return GetPlateFormationZone(omega,gamma,varAngle)
1200  * GetPlateLinearPhotoAbs(omega);
1201 }
1203 //
1204 // Calculates the product of linear cof by formation zone for plate.
1205 // G4cout and output in file in some energy range.
1206 
1208 {
1209  std::ofstream outPlate("plateZmu.dat", std::ios::out );
1210  outPlate.setf( std::ios::scientific, std::ios::floatfield );
1211 
1212  G4int i;
1213  G4double omega, varAngle, gamma;
1214  gamma = 10000.;
1215  varAngle = 1/gamma/gamma;
1216  if(verboseLevel > 0)
1217  G4cout<<"energy, keV"<<"\t"<<"Zmu for plate"<<G4endl;
1218  for(i=0;i<100;i++)
1219  {
1220  omega = (1.0 + i)*keV;
1221  if(verboseLevel > 1)
1222  G4cout<<omega/keV<<"\t"<<GetPlateZmuProduct(omega,gamma,varAngle)<<"\t";
1223  if(verboseLevel > 0)
1224  outPlate<<omega/keV<<"\t\t"<<GetPlateZmuProduct(omega,gamma,varAngle)<<G4endl;
1225  }
1226  return;
1227 }
1228 
1230 //
1231 // Calculates the product of linear cof by formation zone for gas.
1232 // Omega is energy !!!
1233 
1235  G4double gamma ,
1236  G4double varAngle )
1237 {
1238  return GetGasFormationZone(omega,gamma,varAngle)*GetGasLinearPhotoAbs(omega);
1239 }
1241 //
1242 // Calculates the product of linear cof byformation zone for gas.
1243 // G4cout and output in file in some energy range.
1244 
1246 {
1247  std::ofstream outGas("gasZmu.dat", std::ios::out );
1248  outGas.setf( std::ios::scientific, std::ios::floatfield );
1249  G4int i;
1250  G4double omega, varAngle, gamma;
1251  gamma = 10000.;
1252  varAngle = 1/gamma/gamma;
1253  if(verboseLevel > 0)
1254  G4cout<<"energy, keV"<<"\t"<<"Zmu for gas"<<G4endl;
1255  for(i=0;i<100;i++)
1256  {
1257  omega = (1.0 + i)*keV;
1258  if(verboseLevel > 1)
1259  G4cout<<omega/keV<<"\t"<<GetGasZmuProduct(omega,gamma,varAngle)<<"\t";
1260  if(verboseLevel > 0)
1261  outGas<<omega/keV<<"\t\t"<<GetGasZmuProduct(omega,gamma,varAngle)<<G4endl;
1262  }
1263  return;
1264 }
1265 
1267 //
1268 // Computes Compton cross section for plate material in 1/mm
1269 
1271 {
1272  G4int i, numberOfElements;
1273  G4double xSection = 0., nowZ, sumZ = 0.;
1274 
1275  const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
1276  numberOfElements = (*theMaterialTable)[fMatIndex1]->GetNumberOfElements();
1277 
1278  for( i = 0; i < numberOfElements; i++ )
1279  {
1280  nowZ = (*theMaterialTable)[fMatIndex1]->GetElement(i)->GetZ();
1281  sumZ += nowZ;
1282  xSection += GetComptonPerAtom(omega,nowZ); // *nowZ;
1283  }
1284  xSection /= sumZ;
1285  xSection *= (*theMaterialTable)[fMatIndex1]->GetElectronDensity();
1286  return xSection;
1287 }
1288 
1289 
1291 //
1292 // Computes Compton cross section for gas material in 1/mm
1293 
1295 {
1296  G4int i, numberOfElements;
1297  G4double xSection = 0., nowZ, sumZ = 0.;
1298 
1299  const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
1300  numberOfElements = (*theMaterialTable)[fMatIndex2]->GetNumberOfElements();
1301 
1302  for( i = 0; i < numberOfElements; i++ )
1303  {
1304  nowZ = (*theMaterialTable)[fMatIndex2]->GetElement(i)->GetZ();
1305  sumZ += nowZ;
1306  xSection += GetComptonPerAtom(omega,nowZ); // *nowZ;
1307  }
1308  xSection /= sumZ;
1309  xSection *= (*theMaterialTable)[fMatIndex2]->GetElectronDensity();
1310  return xSection;
1311 }
1312 
1314 //
1315 // Computes Compton cross section per atom with Z electrons for gamma with
1316 // the energy GammaEnergy
1317 
1319 {
1320  G4double CrossSection = 0.0;
1321  if ( Z < 0.9999 ) return CrossSection;
1322  if ( GammaEnergy < 0.1*keV ) return CrossSection;
1323  if ( GammaEnergy > (100.*GeV/Z) ) return CrossSection;
1324 
1325  static const G4double a = 20.0 , b = 230.0 , c = 440.0;
1326 
1327  static const G4double
1328  d1= 2.7965e-1*barn, d2=-1.8300e-1*barn, d3= 6.7527 *barn, d4=-1.9798e+1*barn,
1329  e1= 1.9756e-5*barn, e2=-1.0205e-2*barn, e3=-7.3913e-2*barn, e4= 2.7079e-2*barn,
1330  f1=-3.9178e-7*barn, f2= 6.8241e-5*barn, f3= 6.0480e-5*barn, f4= 3.0274e-4*barn;
1331 
1332  G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
1333  p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
1334 
1335  G4double T0 = 15.0*keV;
1336  if (Z < 1.5) T0 = 40.0*keV;
1337 
1338  G4double X = std::max(GammaEnergy, T0) / electron_mass_c2;
1339  CrossSection = p1Z*std::log(1.+2.*X)/X
1340  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
1341 
1342  // modification for low energy. (special case for Hydrogen)
1343 
1344  if (GammaEnergy < T0)
1345  {
1346  G4double dT0 = 1.*keV;
1347  X = (T0+dT0) / electron_mass_c2;
1348  G4double sigma = p1Z*std::log(1.+2.*X)/X
1349  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
1350  G4double c1 = -T0*(sigma-CrossSection)/(CrossSection*dT0);
1351  G4double c2 = 0.150;
1352  if (Z > 1.5) c2 = 0.375-0.0556*std::log(Z);
1353  G4double y = std::log(GammaEnergy/T0);
1354  CrossSection *= std::exp(-y*(c1+c2*y));
1355  }
1356  // G4cout << "e= " << GammaEnergy << " Z= " << Z << " cross= " << CrossSection << G4endl;
1357  return CrossSection;
1358 }
1359 
1360 
1362 //
1363 // This function returns the spectral and angle density of TR quanta
1364 // in X-ray energy region generated forward when a relativistic
1365 // charged particle crosses interface between two materials.
1366 // The high energy small theta approximation is applied.
1367 // (matter1 -> matter2, or 2->1)
1368 // varAngle =2* (1 - std::cos(theta)) or approximately = theta*theta
1369 //
1370 
1371 G4double
1373  G4double varAngle ) const
1374 {
1375  G4double formationLength1, formationLength2;
1376  formationLength1 = 1.0/
1377  (1.0/(gamma*gamma)
1378  + fSigma1/(energy*energy)
1379  + varAngle);
1380  formationLength2 = 1.0/
1381  (1.0/(gamma*gamma)
1382  + fSigma2/(energy*energy)
1383  + varAngle);
1384  return (varAngle/energy)*(formationLength1 - formationLength2)
1385  *(formationLength1 - formationLength2);
1386 
1387 }
1388 
1390  G4double varAngle )
1391 {
1392  // return stack factor corresponding to one interface
1393 
1394  return std::real( OneInterfaceXTRdEdx(energy,gamma,varAngle) );
1395 }
1396 
1398 //
1399 // For photon energy distribution tables. Integrate first over angle
1400 //
1401 
1403 {
1404  return OneBoundaryXTRNdensity(fEnergy,fGamma,varAngle)*
1405  GetStackFactor(fEnergy,fGamma,varAngle);
1406 }
1407 
1409 //
1410 // For second integration over energy
1411 
1413 {
1414  fEnergy = energy;
1417  0.0,0.2*fMaxThetaTR) +
1419  0.2*fMaxThetaTR,fMaxThetaTR);
1420 }
1421 
1423 //
1424 // for photon angle distribution tables
1425 //
1426 
1428 {
1429  return OneBoundaryXTRNdensity(energy,fGamma,fVarAngle)*
1431 }
1432 
1434 //
1435 //
1436 
1438 {
1439  fVarAngle = varAngle;
1443 }
1444 
1446 //
1447 // Check number of photons for a range of Lorentz factors from both energy
1448 // and angular tables
1449 
1451 {
1452  G4int iTkin;
1453  G4double gamma, numberE;
1454 
1455  std::ofstream outEn("numberE.dat", std::ios::out );
1456  outEn.setf( std::ios::scientific, std::ios::floatfield );
1457 
1458  std::ofstream outAng("numberAng.dat", std::ios::out );
1459  outAng.setf( std::ios::scientific, std::ios::floatfield );
1460 
1461  for(iTkin=0;iTkin<fTotBin;iTkin++) // Lorentz factor loop
1462  {
1463  gamma = 1.0 + (fProtonEnergyVector->
1464  GetLowEdgeEnergy(iTkin)/proton_mass_c2);
1465  numberE = (*(*fEnergyDistrTable)(iTkin))(0);
1466  // numberA = (*(*fAngleDistrTable)(iTkin))(0);
1467  if(verboseLevel > 1)
1468  G4cout<<gamma<<"\t\t"<<numberE<<"\t" // <<numberA
1469  <<G4endl;
1470  if(verboseLevel > 0)
1471  outEn<<gamma<<"\t\t"<<numberE<<G4endl;
1472  }
1473  return;
1474 }
1475 
1477 //
1478 // Returns randon energy of a X-ray TR photon for given scaled kinetic energy
1479 // of a charged particle
1480 
1482 {
1483  G4int iTransfer, iPlace;
1484  G4double transfer = 0.0, position, E1, E2, W1, W2, W;
1485 
1486  iPlace = iTkin - 1;
1487 
1488  // G4cout<<"iPlace = "<<iPlace<<endl;
1489 
1490  if(iTkin == fTotBin) // relativistic plato, try from left
1491  {
1492  position = (*(*fEnergyDistrTable)(iPlace))(0)*G4UniformRand();
1493 
1494  for(iTransfer=0;;iTransfer++)
1495  {
1496  if(position >= (*(*fEnergyDistrTable)(iPlace))(iTransfer)) break;
1497  }
1498  transfer = GetXTRenergy(iPlace,position,iTransfer);
1499  }
1500  else
1501  {
1502  E1 = fProtonEnergyVector->GetLowEdgeEnergy(iTkin - 1);
1504  W = 1.0/(E2 - E1);
1505  W1 = (E2 - scaledTkin)*W;
1506  W2 = (scaledTkin - E1)*W;
1507 
1508  position =( (*(*fEnergyDistrTable)(iPlace))(0)*W1 +
1509  (*(*fEnergyDistrTable)(iPlace+1))(0)*W2 )*G4UniformRand();
1510 
1511  // G4cout<<position<<"\t";
1512 
1513  for(iTransfer=0;;iTransfer++)
1514  {
1515  if( position >=
1516  ( (*(*fEnergyDistrTable)(iPlace))(iTransfer)*W1 +
1517  (*(*fEnergyDistrTable)(iPlace+1))(iTransfer)*W2) ) break;
1518  }
1519  transfer = GetXTRenergy(iPlace,position,iTransfer);
1520 
1521  }
1522  // G4cout<<"XTR transfer = "<<transfer/keV<<" keV"<<endl;
1523  if(transfer < 0.0 ) transfer = 0.0;
1524  return transfer;
1525 }
1526 
1528 //
1529 // Returns approximate position of X-ray photon energy during random sampling
1530 // over integral energy distribution
1531 
1533  G4double /* position */,
1534  G4int iTransfer )
1535 {
1536  G4double x1, x2, y1, y2, result;
1537 
1538  if(iTransfer == 0)
1539  {
1540  result = (*fEnergyDistrTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
1541  }
1542  else
1543  {
1544  y1 = (*(*fEnergyDistrTable)(iPlace))(iTransfer-1);
1545  y2 = (*(*fEnergyDistrTable)(iPlace))(iTransfer);
1546 
1547  x1 = (*fEnergyDistrTable)(iPlace)->GetLowEdgeEnergy(iTransfer-1);
1548  x2 = (*fEnergyDistrTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
1549 
1550  if ( x1 == x2 ) result = x2;
1551  else
1552  {
1553  if ( y1 == y2 ) result = x1 + (x2 - x1)*G4UniformRand();
1554  else
1555  {
1556  // result = x1 + (position - y1)*(x2 - x1)/(y2 - y1);
1557  result = x1 + (x2 - x1)*G4UniformRand();
1558  }
1559  }
1560  }
1561  return result;
1562 }
1563 
1565 //
1566 // Get XTR photon angle at given energy and Tkin
1567 
1569 {
1570  G4int iTR, iAngle;
1572 
1573  if (iTkin == fTotBin) iTkin--;
1574 
1576 
1577  for( iTR = 0; iTR < fBinTR; iTR++ )
1578  {
1579  if( energyXTR < fXTREnergyVector->GetLowEdgeEnergy(iTR) ) break;
1580  }
1581  if (iTR == fBinTR) iTR--;
1582 
1583  position = (*(*fAngleForEnergyTable)(iTR))(0)*G4UniformRand();
1584 
1585  for( iAngle = 0;; iAngle++)
1586  {
1587  if( position >= (*(*fAngleForEnergyTable)(iTR))(iAngle) ) break;
1588  }
1589  angle = GetAngleXTR(iTR,position,iAngle);
1590  return angle;
1591 }
1592 
1594 //
1595 // Returns approximate position of X-ray photon angle at given energy during random sampling
1596 // over integral energy distribution
1597 
1599  G4double position,
1600  G4int iTransfer )
1601 {
1602  G4double x1, x2, y1, y2, result;
1603 
1604  if( iTransfer == 0 )
1605  {
1606  result = (*fAngleForEnergyTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
1607  }
1608  else
1609  {
1610  y1 = (*(*fAngleForEnergyTable)(iPlace))(iTransfer-1);
1611  y2 = (*(*fAngleForEnergyTable)(iPlace))(iTransfer);
1612 
1613  x1 = (*fAngleForEnergyTable)(iPlace)->GetLowEdgeEnergy(iTransfer-1);
1614  x2 = (*fAngleForEnergyTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
1615 
1616  if ( x1 == x2 ) result = x2;
1617  else
1618  {
1619  if ( y1 == y2 ) result = x1 + (x2 - x1)*G4UniformRand();
1620  else
1621  {
1622  result = x1 + (position - y1)*(x2 - x1)/(y2 - y1);
1623  }
1624  }
1625  }
1626  return result;
1627 }
1628 
1629 
1630 //
1631 //
1633