G4VXTRenergyLoss.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 //
00027 // $Id$
00028 //
00029 // History:
00030 // 2001-2002 R&D by V.Grichine
00031 // 19.06.03 V. Grichine, modifications in BuildTable for the integration 
00032 //                       in respect of angle: range is increased, accuracy is
00033 //                       improved
00034 // 28.07.05, P.Gumplinger add G4ProcessType to constructor
00035 // 28.09.07, V.Ivanchenko general cleanup without change of algorithms
00036 //
00037 
00038 #include "G4VXTRenergyLoss.hh"
00039 
00040 #include "G4Timer.hh"
00041 #include "G4PhysicalConstants.hh"
00042 #include "G4SystemOfUnits.hh"
00043 
00044 #include "G4Poisson.hh"
00045 #include "G4MaterialTable.hh"
00046 #include "G4VDiscreteProcess.hh"
00047 #include "G4VParticleChange.hh"
00048 #include "G4VSolid.hh"
00049 
00050 #include "G4RotationMatrix.hh"
00051 #include "G4ThreeVector.hh"
00052 #include "G4AffineTransform.hh"
00053 #include "G4SandiaTable.hh"
00054 
00055 #include "G4PhysicsVector.hh"
00056 #include "G4PhysicsFreeVector.hh"
00057 #include "G4PhysicsLinearVector.hh"
00058 
00060 //
00061 // Constructor, destructor
00062 
00063 G4VXTRenergyLoss::G4VXTRenergyLoss(G4LogicalVolume *anEnvelope,
00064                                    G4Material* foilMat,G4Material* gasMat,
00065                                    G4double a, G4double b,
00066                                    G4int n,const G4String& processName,
00067                                    G4ProcessType type) :
00068   G4VDiscreteProcess(processName, type),
00069   fGammaCutInKineticEnergy(0),
00070   fGammaTkinCut(0),
00071   fAngleDistrTable(0),
00072   fEnergyDistrTable(0),
00073   fPlatePhotoAbsCof(0),
00074   fGasPhotoAbsCof(0),
00075   fAngleForEnergyTable(0)
00076 {
00077   verboseLevel = 1;
00078 
00079   fPtrGamma = 0;
00080   fMinEnergyTR = fMaxEnergyTR = fMaxThetaTR = fGamma = fEnergy = fVarAngle 
00081     = fLambda = fTotalDist = fPlateThick = fGasThick = fAlphaPlate = fAlphaGas = 0.0;
00082 
00083   // Initialization of local constants
00084   fTheMinEnergyTR = 1.0*keV;
00085   fTheMaxEnergyTR = 100.0*keV;
00086   fTheMaxAngle    = 1.0e-2;
00087   fTheMinAngle    = 5.0e-6;
00088   fBinTR          = 50;
00089 
00090   fMinProtonTkin  = 100.0*GeV;
00091   fMaxProtonTkin  = 100.0*TeV;
00092   fTotBin         = 50;
00093 
00094   // Proton energy vector initialization
00095 
00096   fProtonEnergyVector = new G4PhysicsLogVector(fMinProtonTkin,
00097                                                fMaxProtonTkin,
00098                                                fTotBin  );
00099 
00100   fXTREnergyVector = new G4PhysicsLogVector(fTheMinEnergyTR,
00101                                             fTheMaxEnergyTR,
00102                                             fBinTR  );
00103 
00104   fPlasmaCof = 4.0*pi*fine_structure_const*hbarc*hbarc*hbarc/electron_mass_c2;
00105 
00106   fCofTR     = fine_structure_const/pi;
00107 
00108   fEnvelope  = anEnvelope;
00109 
00110   fPlateNumber = n;
00111   if(verboseLevel > 0)
00112     G4cout<<"### G4VXTRenergyLoss: the number of TR radiator plates = "
00113           <<fPlateNumber<<G4endl;
00114   if(fPlateNumber == 0)
00115   {
00116     G4Exception("G4VXTRenergyLoss::G4VXTRenergyLoss()","VXTRELoss01",
00117     FatalException,"No plates in X-ray TR radiator");
00118   }
00119   // default is XTR dEdx, not flux after radiator
00120 
00121   fExitFlux      = false;
00122   fAngleRadDistr = false;
00123   fCompton       = false;
00124 
00125   fLambda = DBL_MAX;
00126   // Mean thicknesses of plates and gas gaps
00127 
00128   fPlateThick = a;
00129   fGasThick   = b;
00130   fTotalDist  = fPlateNumber*(fPlateThick+fGasThick);
00131   if(verboseLevel > 0)
00132     G4cout<<"total radiator thickness = "<<fTotalDist/cm<<" cm"<<G4endl;
00133 
00134   // index of plate material
00135   fMatIndex1 = foilMat->GetIndex();
00136   if(verboseLevel > 0)
00137     G4cout<<"plate material = "<<foilMat->GetName()<<G4endl;
00138 
00139   // index of gas material
00140   fMatIndex2 = gasMat->GetIndex();
00141   if(verboseLevel > 0)
00142     G4cout<<"gas material = "<<gasMat->GetName()<<G4endl;
00143 
00144   // plasma energy squared for plate material
00145 
00146   fSigma1 = fPlasmaCof*foilMat->GetElectronDensity();
00147   //  fSigma1 = (20.9*eV)*(20.9*eV);
00148   if(verboseLevel > 0)
00149     G4cout<<"plate plasma energy = "<<std::sqrt(fSigma1)/eV<<" eV"<<G4endl;
00150 
00151   // plasma energy squared for gas material
00152 
00153   fSigma2 = fPlasmaCof*gasMat->GetElectronDensity();
00154   if(verboseLevel > 0)
00155     G4cout<<"gas plasma energy = "<<std::sqrt(fSigma2)/eV<<" eV"<<G4endl;
00156 
00157   // Compute cofs for preparation of linear photo absorption
00158 
00159   ComputePlatePhotoAbsCof();
00160   ComputeGasPhotoAbsCof();
00161 
00162   pParticleChange = &fParticleChange;
00163 }
00164 
00166 
00167 G4VXTRenergyLoss::~G4VXTRenergyLoss()
00168 {
00169   if(fEnvelope) delete fEnvelope;
00170   delete fProtonEnergyVector;
00171   delete fXTREnergyVector;
00172   delete fEnergyDistrTable;
00173   if(fAngleRadDistr) delete fAngleDistrTable;
00174   delete fAngleForEnergyTable;
00175 }
00176 
00178 //
00179 // Returns condition for application of the model depending on particle type
00180 
00181 
00182 G4bool G4VXTRenergyLoss::IsApplicable(const G4ParticleDefinition& particle)
00183 {
00184   return  ( particle.GetPDGCharge() != 0.0 );
00185 }
00186 
00188 //
00189 // Calculate step size for XTR process inside raaditor
00190 
00191 G4double G4VXTRenergyLoss::GetMeanFreePath(const G4Track& aTrack,
00192                                            G4double, // previousStepSize,
00193                                            G4ForceCondition* condition)
00194 {
00195   G4int iTkin, iPlace;
00196   G4double lambda, sigma, kinEnergy, mass, gamma;
00197   G4double charge, chargeSq, massRatio, TkinScaled;
00198   G4double E1,E2,W,W1,W2;
00199 
00200   *condition = NotForced;
00201   
00202   if( aTrack.GetVolume()->GetLogicalVolume() != fEnvelope ) lambda = DBL_MAX;
00203   else
00204   {
00205     const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
00206     kinEnergy = aParticle->GetKineticEnergy();
00207     mass      = aParticle->GetDefinition()->GetPDGMass();
00208     gamma     = 1.0 + kinEnergy/mass;
00209     if(verboseLevel > 1)
00210     {
00211       G4cout<<" gamma = "<<gamma<<";   fGamma = "<<fGamma<<G4endl;
00212     }
00213 
00214     if ( std::fabs( gamma - fGamma ) < 0.05*gamma ) lambda = fLambda;
00215     else
00216     {
00217       charge = aParticle->GetDefinition()->GetPDGCharge();
00218       chargeSq  = charge*charge;
00219       massRatio = proton_mass_c2/mass;
00220       TkinScaled = kinEnergy*massRatio;
00221 
00222       for(iTkin = 0; iTkin < fTotBin; iTkin++)
00223       {
00224         if( TkinScaled < fProtonEnergyVector->GetLowEdgeEnergy(iTkin))  break;    
00225       }
00226       iPlace = iTkin - 1;
00227 
00228       if(iTkin == 0) lambda = DBL_MAX; // Tkin is too small, neglect of TR photon generation
00229       else          // general case: Tkin between two vectors of the material
00230       {
00231         if(iTkin == fTotBin) 
00232         {
00233           sigma = (*(*fEnergyDistrTable)(iPlace))(0)*chargeSq;
00234         }
00235         else
00236         {
00237           E1 = fProtonEnergyVector->GetLowEdgeEnergy(iTkin - 1); 
00238           E2 = fProtonEnergyVector->GetLowEdgeEnergy(iTkin);
00239           W = 1.0/(E2 - E1);
00240           W1 = (E2 - TkinScaled)*W;
00241           W2 = (TkinScaled - E1)*W;
00242           sigma = ( (*(*fEnergyDistrTable)(iPlace  ))(0)*W1 +
00243                 (*(*fEnergyDistrTable)(iPlace+1))(0)*W2   )*chargeSq;
00244       
00245         }
00246         if (sigma < DBL_MIN) lambda = DBL_MAX;
00247         else                 lambda = 1./sigma; 
00248         fLambda = lambda;
00249         fGamma  = gamma;   
00250         if(verboseLevel > 1)
00251         {
00252           G4cout<<" lambda = "<<lambda/mm<<" mm"<<G4endl;
00253         }
00254       }
00255     }
00256   }  
00257   return lambda;
00258 }
00259 
00261 //
00262 // Interface for build table from physics list
00263 
00264 void G4VXTRenergyLoss::BuildPhysicsTable(const G4ParticleDefinition& pd)
00265 {
00266   if(pd.GetPDGCharge()  == 0.) 
00267   {
00268     G4Exception("G4VXTRenergyLoss::BuildPhysicsTable", "Notification", JustWarning,
00269                  "XTR initialisation for neutral particle ?!" );   
00270   }
00271   BuildEnergyTable();
00272 
00273   if (fAngleRadDistr) 
00274   {
00275     if(verboseLevel > 0)
00276     {
00277       G4cout<<"Build angle for energy distribution according the current radiator"
00278             <<G4endl;
00279     }
00280     BuildAngleForEnergyBank();
00281   }
00282 }
00283 
00284 
00286 //
00287 // Build integral energy distribution of XTR photons
00288 
00289 void G4VXTRenergyLoss::BuildEnergyTable()
00290 {
00291   G4int iTkin, iTR, iPlace;
00292   G4double radiatorCof = 1.0;           // for tuning of XTR yield
00293   G4double energySum = 0.0;
00294 
00295   fEnergyDistrTable = new G4PhysicsTable(fTotBin);
00296   if(fAngleRadDistr) fAngleDistrTable = new G4PhysicsTable(fTotBin);
00297 
00298   fGammaTkinCut = 0.0;
00299   
00300   // setting of min/max TR energies 
00301   
00302   if(fGammaTkinCut > fTheMinEnergyTR)  fMinEnergyTR = fGammaTkinCut;
00303   else                                 fMinEnergyTR = fTheMinEnergyTR;
00304         
00305   if(fGammaTkinCut > fTheMaxEnergyTR) fMaxEnergyTR = 2.0*fGammaTkinCut;  
00306   else                                fMaxEnergyTR = fTheMaxEnergyTR;
00307     
00308   G4Integrator<G4VXTRenergyLoss,G4double(G4VXTRenergyLoss::*)(G4double)> integral;
00309 
00310   G4cout.precision(4);
00311   G4Timer timer;
00312   timer.Start();
00313 
00314   if(verboseLevel > 0) 
00315   {
00316     G4cout<<G4endl;
00317     G4cout<<"Lorentz Factor"<<"\t"<<"XTR photon number"<<G4endl;
00318     G4cout<<G4endl;
00319   }
00320   for( iTkin = 0; iTkin < fTotBin; iTkin++ )      // Lorentz factor loop
00321   {
00322     G4PhysicsLogVector* energyVector = new G4PhysicsLogVector( fMinEnergyTR,
00323                                                                fMaxEnergyTR,
00324                                                                fBinTR  );
00325 
00326     fGamma = 1.0 + (fProtonEnergyVector->
00327                     GetLowEdgeEnergy(iTkin)/proton_mass_c2);
00328 
00329     fMaxThetaTR = 2500.0/(fGamma*fGamma) ;  // theta^2
00330 
00331     fTheMinAngle = 1.0e-3; // was 5.e-6, e-6 !!!, e-5, e-4
00332  
00333     if(      fMaxThetaTR > fTheMaxAngle )  fMaxThetaTR = fTheMaxAngle; 
00334     else if( fMaxThetaTR < fTheMinAngle )  fMaxThetaTR = fTheMinAngle;
00335       
00336     energySum = 0.0;
00337 
00338     energyVector->PutValue(fBinTR-1,energySum);
00339 
00340     for( iTR = fBinTR - 2; iTR >= 0; iTR-- )
00341     {
00342         // Legendre96 or Legendre10
00343 
00344       energySum += radiatorCof*fCofTR*integral.Legendre10(
00345              this,&G4VXTRenergyLoss::SpectralXTRdEdx,
00346                    energyVector->GetLowEdgeEnergy(iTR),
00347                    energyVector->GetLowEdgeEnergy(iTR+1) ); 
00348       
00349       energyVector->PutValue(iTR,energySum/fTotalDist);
00350     }
00351     iPlace = iTkin;
00352     fEnergyDistrTable->insertAt(iPlace,energyVector);
00353 
00354     if(verboseLevel > 0)
00355     {
00356         G4cout
00357           // <<iTkin<<"\t"
00358           //   <<"fGamma = "
00359           <<fGamma<<"\t"  //  <<"  fMaxThetaTR = "<<fMaxThetaTR
00360           //  <<"sumN = "
00361           <<energySum      // <<"; sumA = "<<angleSum
00362           <<G4endl;
00363     }
00364   }     
00365   timer.Stop();
00366   G4cout.precision(6);
00367   if(verboseLevel > 0) 
00368   {
00369     G4cout<<G4endl;
00370     G4cout<<"total time for build X-ray TR energy loss tables = "
00371           <<timer.GetUserElapsed()<<" s"<<G4endl;
00372   }
00373   fGamma = 0.;
00374   return;
00375 }
00376 
00378 //
00379 // Bank of angle distributions for given energies (slow!)
00380 
00381 void G4VXTRenergyLoss::BuildAngleForEnergyBank()
00382 {
00383   if( this->GetProcessName() == "TranspRegXTRadiator" || 
00384       this->GetProcessName() == "TranspRegXTRmodel"   || 
00385       this->GetProcessName() == "RegularXTRadiator"   || 
00386       this->GetProcessName() == "RegularXTRmodel"       )
00387   {
00388     BuildAngleTable();
00389     return;
00390   }
00391   G4int i, iTkin, iTR;
00392   G4double angleSum  = 0.0;
00393 
00394 
00395   fGammaTkinCut = 0.0;
00396   
00397   // setting of min/max TR energies 
00398   
00399   if(fGammaTkinCut > fTheMinEnergyTR)  fMinEnergyTR = fGammaTkinCut;
00400   else                                 fMinEnergyTR = fTheMinEnergyTR;
00401         
00402   if(fGammaTkinCut > fTheMaxEnergyTR) fMaxEnergyTR = 2.0*fGammaTkinCut;  
00403   else                                fMaxEnergyTR = fTheMaxEnergyTR;
00404 
00405   G4PhysicsLogVector* energyVector = new G4PhysicsLogVector( fMinEnergyTR,
00406                                                                fMaxEnergyTR,
00407                                                                fBinTR  );
00408 
00409   G4Integrator<G4VXTRenergyLoss,G4double(G4VXTRenergyLoss::*)(G4double)> integral;
00410 
00411   G4cout.precision(4);
00412   G4Timer timer;
00413   timer.Start();
00414 
00415   for( iTkin = 0; iTkin < fTotBin; iTkin++ )      // Lorentz factor loop
00416   {
00417 
00418     fGamma = 1.0 + (fProtonEnergyVector->
00419                     GetLowEdgeEnergy(iTkin)/proton_mass_c2);
00420 
00421     fMaxThetaTR = 2500.0/(fGamma*fGamma) ;  // theta^2
00422 
00423     fTheMinAngle = 1.0e-3; // was 5.e-6, e-6 !!!, e-5, e-4
00424  
00425     if(      fMaxThetaTR > fTheMaxAngle )  fMaxThetaTR = fTheMaxAngle; 
00426     else if( fMaxThetaTR < fTheMinAngle )  fMaxThetaTR = fTheMinAngle;
00427       
00428     fAngleForEnergyTable = new G4PhysicsTable(fBinTR);
00429 
00430     for( iTR = 0; iTR < fBinTR; iTR++ )
00431     {
00432       angleSum     = 0.0;
00433       fEnergy      = energyVector->GetLowEdgeEnergy(iTR);    
00434       G4PhysicsLinearVector* angleVector = new G4PhysicsLinearVector(0.0,
00435                                                                    fMaxThetaTR,
00436                                                                    fBinTR  );
00437     
00438       angleVector ->PutValue(fBinTR - 1, angleSum);
00439 
00440       for( i = fBinTR - 2; i >= 0; i-- )
00441       {
00442           // Legendre96 or Legendre10
00443 
00444           angleSum  += integral.Legendre10(
00445                    this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
00446                    angleVector->GetLowEdgeEnergy(i),
00447                    angleVector->GetLowEdgeEnergy(i+1) );
00448 
00449           angleVector ->PutValue(i, angleSum);
00450       }
00451       fAngleForEnergyTable->insertAt(iTR, angleVector);
00452     }
00453     fAngleBank.push_back(fAngleForEnergyTable); 
00454   }    
00455   timer.Stop();
00456   G4cout.precision(6);
00457   if(verboseLevel > 0) 
00458   {
00459     G4cout<<G4endl;
00460     G4cout<<"total time for build X-ray TR angle for energy loss tables = "
00461           <<timer.GetUserElapsed()<<" s"<<G4endl;
00462   }
00463   fGamma = 0.;
00464   return;
00465 }
00466 
00468 //
00469 // Build XTR angular distribution at given energy based on the model 
00470 // of transparent regular radiator
00471 
00472 void G4VXTRenergyLoss::BuildAngleTable()
00473 {
00474   G4int iTkin, iTR;
00475   G4double  energy;
00476 
00477   fGammaTkinCut = 0.0;
00478                               
00479   // setting of min/max TR energies 
00480   
00481   if(fGammaTkinCut > fTheMinEnergyTR)  fMinEnergyTR = fGammaTkinCut;
00482   else                                 fMinEnergyTR = fTheMinEnergyTR;
00483         
00484   if(fGammaTkinCut > fTheMaxEnergyTR) fMaxEnergyTR = 2.0*fGammaTkinCut;  
00485   else                                fMaxEnergyTR = fTheMaxEnergyTR;
00486 
00487   G4cout.precision(4);
00488   G4Timer timer;
00489   timer.Start();
00490   if(verboseLevel > 0) 
00491   {
00492     G4cout<<G4endl;
00493     G4cout<<"Lorentz Factor"<<"\t"<<"XTR photon number"<<G4endl;
00494     G4cout<<G4endl;
00495   }
00496   for( iTkin = 0; iTkin < fTotBin; iTkin++ )      // Lorentz factor loop
00497   {
00498     
00499     fGamma = 1.0 + (fProtonEnergyVector->
00500                             GetLowEdgeEnergy(iTkin)/proton_mass_c2);
00501 
00502     fMaxThetaTR = 25.0/(fGamma*fGamma);  // theta^2
00503 
00504     fTheMinAngle = 1.0e-3; // was 5.e-6, e-6 !!!, e-5, e-4
00505  
00506     if( fMaxThetaTR > fTheMaxAngle )    fMaxThetaTR = fTheMaxAngle; 
00507     else
00508     {
00509        if( fMaxThetaTR < fTheMinAngle )  fMaxThetaTR = fTheMinAngle;
00510     }
00511 
00512     fAngleForEnergyTable = new G4PhysicsTable(fBinTR);
00513 
00514     for( iTR = 0; iTR < fBinTR; iTR++ )
00515     {
00516       // energy = fMinEnergyTR*(iTR+1);
00517 
00518       energy = fXTREnergyVector->GetLowEdgeEnergy(iTR);
00519 
00520       G4PhysicsFreeVector* angleVector = GetAngleVector(energy,fBinTR);
00521       // G4cout<<G4endl;
00522 
00523       fAngleForEnergyTable->insertAt(iTR,angleVector);
00524     }
00525     
00526     fAngleBank.push_back(fAngleForEnergyTable); 
00527   }     
00528   timer.Stop();
00529   G4cout.precision(6);
00530   if(verboseLevel > 0) 
00531   {
00532     G4cout<<G4endl;
00533     G4cout<<"total time for build XTR angle for given energy tables = "
00534           <<timer.GetUserElapsed()<<" s"<<G4endl;
00535   }
00536   fGamma = 0.;
00537   
00538   return;
00539 } 
00540 
00542 //
00543 // Vector of angles and angle integral distributions
00544 
00545 G4PhysicsFreeVector* G4VXTRenergyLoss::GetAngleVector(G4double energy, G4int n)
00546 {
00547   G4double theta=0., result, tmp=0., cof1, cof2, cofMin, cofPHC, angleSum  = 0.;
00548   G4int iTheta, k, /*kMax,*/ kMin;
00549 
00550   G4PhysicsFreeVector* angleVector = new G4PhysicsFreeVector(n);
00551   
00552   cofPHC  = 4*pi*hbarc;
00553   tmp     = (fSigma1 - fSigma2)/cofPHC/energy;
00554   cof1    = fPlateThick*tmp;
00555   cof2    = fGasThick*tmp;
00556 
00557   cofMin  =  energy*(fPlateThick + fGasThick)/fGamma/fGamma;
00558   cofMin += (fPlateThick*fSigma1 + fGasThick*fSigma2)/energy;
00559   cofMin /= cofPHC;
00560 
00561   kMin = G4int(cofMin);
00562   if (cofMin > kMin) kMin++;
00563 
00564   //kMax = kMin + fBinTR -1;
00565 
00566   if(verboseLevel > 2)
00567   {
00568     G4cout<<"n-1 = "<<n-1<<"; theta = "
00569           <<std::sqrt(fMaxThetaTR)*fGamma<<"; tmp = "
00570           <<0. 
00571           <<";    angleSum = "<<angleSum<<G4endl; 
00572   }
00573   //  angleVector->PutValue(n-1,fMaxThetaTR, angleSum);
00574 
00575   for( iTheta = n - 1; iTheta >= 1; iTheta-- )
00576   {
00577 
00578     k = iTheta- 1 + kMin;
00579 
00580     tmp    = pi*fPlateThick*(k + cof2)/(fPlateThick + fGasThick);
00581 
00582     result = (k - cof1)*(k - cof1)*(k + cof2)*(k + cof2);
00583 
00584     tmp = std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
00585 
00586     if( k == kMin && kMin == G4int(cofMin) )
00587     {
00588       angleSum   += 0.5*tmp; // 0.5*std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
00589     }
00590     else if(iTheta == n-1);
00591     else
00592     {
00593       angleSum   += tmp; // std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
00594     }
00595     theta = std::abs(k-cofMin)*cofPHC/energy/(fPlateThick + fGasThick);
00596 
00597     if(verboseLevel > 2)
00598     {
00599       G4cout<<"iTheta = "<<iTheta<<"; k = "<<k<<"; theta = "
00600             <<std::sqrt(theta)*fGamma<<"; tmp = "
00601             <<tmp // std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result
00602             <<";    angleSum = "<<angleSum<<G4endl;
00603     }
00604     angleVector->PutValue( iTheta, theta, angleSum );       
00605   }
00606   if (theta > 0.)
00607   {
00608     angleSum += 0.5*tmp;
00609     theta = 0.;
00610   }
00611   if(verboseLevel > 2)
00612   {
00613     G4cout<<"iTheta = "<<iTheta<<"; theta = "
00614           <<std::sqrt(theta)*fGamma<<"; tmp = "
00615           <<tmp 
00616           <<";    angleSum = "<<angleSum<<G4endl;
00617   }
00618   angleVector->PutValue( iTheta, theta, angleSum );
00619 
00620   return angleVector;
00621 }
00622 
00624 //
00625 // Build XTR angular distribution based on the model of transparent regular radiator
00626 
00627 void G4VXTRenergyLoss::BuildGlobalAngleTable()
00628 {
00629   G4int iTkin, iTR, iPlace;
00630   G4double radiatorCof = 1.0;           // for tuning of XTR yield
00631   G4double angleSum;
00632   fAngleDistrTable = new G4PhysicsTable(fTotBin);
00633 
00634   fGammaTkinCut = 0.0;
00635   
00636   // setting of min/max TR energies 
00637   
00638   if(fGammaTkinCut > fTheMinEnergyTR)  fMinEnergyTR = fGammaTkinCut;
00639   else                                 fMinEnergyTR = fTheMinEnergyTR;
00640         
00641   if(fGammaTkinCut > fTheMaxEnergyTR) fMaxEnergyTR = 2.0*fGammaTkinCut;  
00642   else                                fMaxEnergyTR = fTheMaxEnergyTR;
00643 
00644   G4cout.precision(4);
00645   G4Timer timer;
00646   timer.Start();
00647   if(verboseLevel > 0) {
00648     G4cout<<G4endl;
00649     G4cout<<"Lorentz Factor"<<"\t"<<"XTR photon number"<<G4endl;
00650     G4cout<<G4endl;
00651   }
00652   for( iTkin = 0; iTkin < fTotBin; iTkin++ )      // Lorentz factor loop
00653   {
00654     
00655     fGamma = 1.0 + (fProtonEnergyVector->
00656                             GetLowEdgeEnergy(iTkin)/proton_mass_c2);
00657 
00658     fMaxThetaTR = 25.0/(fGamma*fGamma);  // theta^2
00659 
00660     fTheMinAngle = 1.0e-3; // was 5.e-6, e-6 !!!, e-5, e-4
00661  
00662     if( fMaxThetaTR > fTheMaxAngle )    fMaxThetaTR = fTheMaxAngle; 
00663     else
00664     {
00665        if( fMaxThetaTR < fTheMinAngle )  fMaxThetaTR = fTheMinAngle;
00666     }
00667     G4PhysicsLinearVector* angleVector = new G4PhysicsLinearVector(0.0,
00668                                                                fMaxThetaTR,
00669                                                                fBinTR      );
00670 
00671     angleSum  = 0.0;
00672 
00673     G4Integrator<G4VXTRenergyLoss,G4double(G4VXTRenergyLoss::*)(G4double)> integral;
00674 
00675    
00676     angleVector->PutValue(fBinTR-1,angleSum);
00677 
00678     for( iTR = fBinTR - 2; iTR >= 0; iTR-- )
00679     {
00680 
00681       angleSum  += radiatorCof*fCofTR*integral.Legendre96(
00682                    this,&G4VXTRenergyLoss::AngleXTRdEdx,
00683                    angleVector->GetLowEdgeEnergy(iTR),
00684                    angleVector->GetLowEdgeEnergy(iTR+1) );
00685 
00686       angleVector ->PutValue(iTR,angleSum);
00687     }
00688     if(verboseLevel > 1) {
00689       G4cout
00690         // <<iTkin<<"\t"
00691         //   <<"fGamma = "
00692         <<fGamma<<"\t"  //  <<"  fMaxThetaTR = "<<fMaxThetaTR
00693         //  <<"sumN = "<<energySum      // <<"; sumA = "
00694         <<angleSum
00695         <<G4endl;
00696     }
00697     iPlace = iTkin;
00698     fAngleDistrTable->insertAt(iPlace,angleVector);
00699   }     
00700   timer.Stop();
00701   G4cout.precision(6);
00702   if(verboseLevel > 0) {
00703     G4cout<<G4endl;
00704     G4cout<<"total time for build X-ray TR angle tables = "
00705           <<timer.GetUserElapsed()<<" s"<<G4endl;
00706   }
00707   fGamma = 0.;
00708   
00709   return;
00710 } 
00711 
00712 
00714 //
00715 // The main function which is responsible for the treatment of a particle passage
00716 // trough G4Envelope with discrete generation of G4Gamma
00717 
00718 G4VParticleChange* G4VXTRenergyLoss::PostStepDoIt( const G4Track& aTrack, 
00719                                                   const G4Step&  aStep   )
00720 {
00721   G4int iTkin /*, iPlace*/;
00722   G4double energyTR, theta,theta2, phi, dirX, dirY, dirZ;
00723  
00724 
00725   fParticleChange.Initialize(aTrack);
00726 
00727   if(verboseLevel > 1)
00728   {
00729     G4cout<<"Start of G4VXTRenergyLoss::PostStepDoIt "<<G4endl;
00730     G4cout<<"name of current material =  "
00731           <<aTrack.GetVolume()->GetLogicalVolume()->GetMaterial()->GetName()<<G4endl;
00732   }
00733   if( aTrack.GetVolume()->GetLogicalVolume() != fEnvelope ) 
00734   {
00735     if(verboseLevel > 0)
00736     {
00737       G4cout<<"Go out from G4VXTRenergyLoss::PostStepDoIt: wrong volume "<<G4endl;
00738     }
00739     return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
00740   }
00741   else
00742   {
00743     G4StepPoint* pPostStepPoint        = aStep.GetPostStepPoint();
00744     const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
00745    
00746     // Now we are ready to Generate one TR photon
00747 
00748     G4double kinEnergy = aParticle->GetKineticEnergy();
00749     G4double mass      = aParticle->GetDefinition()->GetPDGMass();
00750     G4double gamma     = 1.0 + kinEnergy/mass;
00751 
00752     if(verboseLevel > 1 )
00753     {
00754       G4cout<<"gamma = "<<gamma<<G4endl;
00755     }
00756     G4double         massRatio   = proton_mass_c2/mass;
00757     G4double          TkinScaled = kinEnergy*massRatio;
00758     G4ThreeVector      position  = pPostStepPoint->GetPosition();
00759     G4ParticleMomentum direction = aParticle->GetMomentumDirection();
00760     G4double           startTime = pPostStepPoint->GetGlobalTime();
00761 
00762     for( iTkin = 0; iTkin < fTotBin; iTkin++ )
00763     {
00764       if(TkinScaled < fProtonEnergyVector->GetLowEdgeEnergy(iTkin))  break;    
00765     }
00766     //iPlace = iTkin - 1;
00767 
00768     if(iTkin == 0) // Tkin is too small, neglect of TR photon generation
00769     {
00770       if( verboseLevel > 0)
00771       {
00772         G4cout<<"Go out from G4VXTRenergyLoss::PostStepDoIt:iTkin = "<<iTkin<<G4endl;
00773       }
00774       return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
00775     } 
00776     else          // general case: Tkin between two vectors of the material
00777     {
00778       fParticleChange.SetNumberOfSecondaries(1);
00779 
00780       energyTR = GetXTRrandomEnergy(TkinScaled,iTkin);
00781 
00782       if( verboseLevel > 1)
00783       {
00784         G4cout<<"energyTR = "<<energyTR/keV<<" keV"<<G4endl;
00785       }
00786       if (fAngleRadDistr)
00787       {
00788         // theta = std::fabs(G4RandGauss::shoot(0.0,pi/gamma));
00789         theta2 = GetRandomAngle(energyTR,iTkin);
00790         if(theta2 > 0.) theta = std::sqrt(theta2);
00791         else            theta = 0.; // theta2;
00792       }
00793       else theta = std::fabs(G4RandGauss::shoot(0.0,pi/gamma));
00794 
00795       // theta = 0.;  // check no spread
00796 
00797       if( theta >= 0.1 ) theta = 0.1;
00798 
00799       // G4cout<<" : theta = "<<theta<<endl;
00800 
00801       phi = twopi*G4UniformRand();
00802 
00803       dirX = std::sin(theta)*std::cos(phi);
00804       dirY = std::sin(theta)*std::sin(phi);
00805       dirZ = std::cos(theta);
00806 
00807       G4ThreeVector directionTR(dirX,dirY,dirZ);
00808       directionTR.rotateUz(direction);
00809       directionTR.unit();
00810 
00811       G4DynamicParticle* aPhotonTR = new G4DynamicParticle(G4Gamma::Gamma(),
00812                                                            directionTR, energyTR);
00813 
00814       // A XTR photon is set on the particle track inside the radiator 
00815       // and is moved to the G4Envelope surface for standard X-ray TR models
00816       // only. The case of fExitFlux=true
00817 
00818       if( fExitFlux )
00819       {
00820         const G4RotationMatrix* rotM = pPostStepPoint->GetTouchable()->GetRotation();
00821         G4ThreeVector transl = pPostStepPoint->GetTouchable()->GetTranslation();
00822         G4AffineTransform transform = G4AffineTransform(rotM,transl);
00823         transform.Invert();
00824         G4ThreeVector localP = transform.TransformPoint(position);
00825         G4ThreeVector localV = transform.TransformAxis(directionTR);
00826 
00827         G4double distance = fEnvelope->GetSolid()->DistanceToOut(localP, localV);
00828         if(verboseLevel > 1)
00829         {
00830           G4cout<<"distance to exit = "<<distance/mm<<" mm"<<G4endl;
00831         }
00832         position         += distance*directionTR;
00833         startTime        += distance/c_light;
00834       }
00835       G4Track* aSecondaryTrack = new G4Track( aPhotonTR, 
00836                                                 startTime, position );
00837       aSecondaryTrack->SetTouchableHandle(
00838                          aStep.GetPostStepPoint()->GetTouchableHandle());
00839       aSecondaryTrack->SetParentID( aTrack.GetTrackID() );
00840 
00841       fParticleChange.AddSecondary(aSecondaryTrack);
00842       fParticleChange.ProposeEnergy(kinEnergy);     
00843     }
00844   }
00845   return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
00846 }
00847 
00849 //
00850 // This function returns the spectral and angle density of TR quanta
00851 // in X-ray energy region generated forward when a relativistic
00852 // charged particle crosses interface between two materials.
00853 // The high energy small theta approximation is applied.
00854 // (matter1 -> matter2, or 2->1)
00855 // varAngle =2* (1 - std::cos(theta)) or approximately = theta*theta
00856 //
00857 
00858 G4complex G4VXTRenergyLoss::OneInterfaceXTRdEdx( G4double energy,
00859                                            G4double gamma,
00860                                            G4double varAngle ) 
00861 {
00862   G4complex Z1    = GetPlateComplexFZ(energy,gamma,varAngle);
00863   G4complex Z2    = GetGasComplexFZ(energy,gamma,varAngle);
00864 
00865   G4complex zOut  = (Z1 - Z2)*(Z1 - Z2)
00866                     * (varAngle*energy/hbarc/hbarc);  
00867   return    zOut;
00868 
00869 }
00870 
00871 
00873 //
00874 // For photon energy distribution tables. Integrate first over angle
00875 //
00876 
00877 G4double G4VXTRenergyLoss::SpectralAngleXTRdEdx(G4double varAngle)
00878 {
00879   G4double result =  GetStackFactor(fEnergy,fGamma,varAngle);
00880   if(result < 0.0) result = 0.0;
00881   return result;
00882 }
00883 
00885 //
00886 // For second integration over energy
00887  
00888 G4double G4VXTRenergyLoss::SpectralXTRdEdx(G4double energy)
00889 {
00890   G4int i, iMax = 8;
00891   G4double angleSum = 0.0;
00892 
00893   G4double lim[8] = { 0.0, 0.01, 0.02, 0.05, 0.1, 0.2, 0.5, 1.0 };
00894 
00895   for( i = 0; i < iMax; i++ ) lim[i] *= fMaxThetaTR;
00896 
00897   G4Integrator<G4VXTRenergyLoss,G4double(G4VXTRenergyLoss::*)(G4double)> integral;
00898 
00899   fEnergy = energy;
00900   /*
00901   if( fAngleRadDistr && ( fEnergy == fEnergyForAngle ) )
00902   {
00903     fAngleVector ->PutValue(fBinTR - 1, angleSum);
00904 
00905     for( i = fBinTR - 2; i >= 0; i-- )
00906     {
00907 
00908       angleSum  += integral.Legendre10(
00909                    this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
00910                    fAngleVector->GetLowEdgeEnergy(i),
00911                    fAngleVector->GetLowEdgeEnergy(i+1) );
00912 
00913       fAngleVector ->PutValue(i, angleSum);
00914     }
00915   }
00916   else
00917   */
00918   {
00919     for( i = 0; i < iMax-1; i++ )
00920     {
00921       angleSum += integral.Legendre96(this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
00922                            lim[i],lim[i+1]);
00923       // result += integral.Legendre10(this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
00924       //                           lim[i],lim[i+1]);
00925     }
00926   }
00927   return angleSum;
00928 } 
00929  
00931 // 
00932 // for photon angle distribution tables
00933 //
00934 
00935 G4double G4VXTRenergyLoss::AngleSpectralXTRdEdx(G4double energy)
00936 {
00937   G4double result =  GetStackFactor(energy,fGamma,fVarAngle);
00938   if(result < 0) result = 0.0;
00939   return result;
00940 } 
00941 
00943 //
00944 // The XTR angular distribution based on transparent regular radiator
00945 
00946 G4double G4VXTRenergyLoss::AngleXTRdEdx(G4double varAngle) 
00947 {
00948   // G4cout<<"angle2 = "<<varAngle<<"; fGamma = "<<fGamma<<G4endl;
00949  
00950   G4double result;
00951   G4double sum = 0., tmp1, tmp2, tmp=0., cof1, cof2, cofMin, cofPHC, energy1, energy2;
00952   G4int k, kMax, kMin, i;
00953 
00954   cofPHC  = twopi*hbarc;
00955 
00956   cof1    = (fPlateThick + fGasThick)*(1./fGamma/fGamma + varAngle);
00957   cof2    = fPlateThick*fSigma1 + fGasThick*fSigma2;
00958 
00959   // G4cout<<"cof1 = "<<cof1<<"; cof2 = "<<cof2<<"; cofPHC = "<<cofPHC<<G4endl; 
00960 
00961   cofMin  =  std::sqrt(cof1*cof2); 
00962   cofMin /= cofPHC;
00963 
00964   kMin = G4int(cofMin);
00965   if (cofMin > kMin) kMin++;
00966 
00967   kMax = kMin + 9; //  9; // kMin + G4int(tmp);
00968 
00969   // G4cout<<"cofMin = "<<cofMin<<"; kMin = "<<kMin<<"; kMax = "<<kMax<<G4endl;
00970 
00971   for( k = kMin; k <= kMax; k++ )
00972   {
00973     tmp1 = cofPHC*k;
00974     tmp2 = std::sqrt(tmp1*tmp1-cof1*cof2);
00975     energy1 = (tmp1+tmp2)/cof1;
00976     energy2 = (tmp1-tmp2)/cof1;
00977 
00978     for(i = 0; i < 2; i++)
00979     {
00980       if( i == 0 )
00981       {
00982         if (energy1 > fTheMaxEnergyTR || energy1 < fTheMinEnergyTR) continue;
00983         tmp1 = ( energy1*energy1*(1./fGamma/fGamma + varAngle) + fSigma1 )
00984           * fPlateThick/(4*hbarc*energy1);
00985         tmp2 = std::sin(tmp1);
00986         tmp  = energy1*tmp2*tmp2;
00987         tmp2 = fPlateThick/(4*tmp1);
00988         tmp1 = hbarc*energy1/( energy1*energy1*(1./fGamma/fGamma + varAngle) + fSigma2 );
00989         tmp *= (tmp1-tmp2)*(tmp1-tmp2);
00990         tmp1 = cof1/(4*hbarc) - cof2/(4*hbarc*energy1*energy1);
00991         tmp2 = std::abs(tmp1);
00992         if(tmp2 > 0.) tmp /= tmp2;
00993         else continue;
00994       }
00995       else
00996       {
00997         if (energy2 > fTheMaxEnergyTR || energy2 < fTheMinEnergyTR) continue;
00998         tmp1 = ( energy2*energy2*(1./fGamma/fGamma + varAngle) + fSigma1 )
00999           * fPlateThick/(4*hbarc*energy2);
01000         tmp2 = std::sin(tmp1);
01001         tmp  = energy2*tmp2*tmp2;
01002         tmp2 = fPlateThick/(4*tmp1);
01003         tmp1 = hbarc*energy2/( energy2*energy2*(1./fGamma/fGamma + varAngle) + fSigma2 );
01004         tmp *= (tmp1-tmp2)*(tmp1-tmp2);
01005         tmp1 = cof1/(4*hbarc) - cof2/(4*hbarc*energy2*energy2);
01006         tmp2 = std::abs(tmp1);
01007         if(tmp2 > 0.) tmp /= tmp2;
01008         else continue;
01009       }
01010       sum += tmp;
01011     }
01012     // G4cout<<"k = "<<k<<"; energy1 = "<<energy1/keV<<" keV; energy2 = "<<energy2/keV
01013     //  <<" keV; tmp = "<<tmp<<"; sum = "<<sum<<G4endl;
01014   }
01015   result = 4.*pi*fPlateNumber*sum*varAngle;
01016   result /= hbarc*hbarc;
01017 
01018   // old code based on general numeric integration
01019   // fVarAngle = varAngle;
01020   // G4Integrator<G4VXTRenergyLoss,G4double(G4VXTRenergyLoss::*)(G4double)> integral;
01021   // result = integral.Legendre10(this,&G4VXTRenergyLoss::AngleSpectralXTRdEdx,
01022   //                         fMinEnergyTR,fMaxEnergyTR);
01023   return result;
01024 }
01025 
01026 
01030 //
01031 // Calculates formation zone for plates. Omega is energy !!!
01032 
01033 G4double G4VXTRenergyLoss::GetPlateFormationZone( G4double omega ,
01034                                                 G4double gamma ,
01035                                                 G4double varAngle    ) 
01036 {
01037   G4double cof, lambda;
01038   lambda = 1.0/gamma/gamma + varAngle + fSigma1/omega/omega;
01039   cof = 2.0*hbarc/omega/lambda;
01040   return cof;
01041 }
01042 
01044 //
01045 // Calculates complex formation zone for plates. Omega is energy !!!
01046 
01047 G4complex G4VXTRenergyLoss::GetPlateComplexFZ( G4double omega ,
01048                                              G4double gamma ,
01049                                              G4double varAngle    ) 
01050 {
01051   G4double cof, length,delta, real_v, image_v;
01052 
01053   length = 0.5*GetPlateFormationZone(omega,gamma,varAngle);
01054   delta  = length*GetPlateLinearPhotoAbs(omega);
01055   cof    = 1.0/(1.0 + delta*delta);
01056 
01057   real_v  = length*cof;
01058   image_v = real_v*delta;
01059 
01060   G4complex zone(real_v,image_v); 
01061   return zone;
01062 }
01063 
01065 //
01066 // Computes matrix of Sandia photo absorption cross section coefficients for
01067 // plate material
01068 
01069 void G4VXTRenergyLoss::ComputePlatePhotoAbsCof() 
01070 {
01071   const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
01072   const G4Material* mat = (*theMaterialTable)[fMatIndex1];
01073   fPlatePhotoAbsCof = mat->GetSandiaTable();
01074 
01075   return;
01076 }
01077 
01078 
01079 
01081 //
01082 // Returns the value of linear photo absorption coefficient (in reciprocal 
01083 // length) for plate for given energy of X-ray photon omega
01084 
01085 G4double G4VXTRenergyLoss::GetPlateLinearPhotoAbs(G4double omega) 
01086 {
01087   //  G4int i;
01088   G4double omega2, omega3, omega4; 
01089 
01090   omega2 = omega*omega;
01091   omega3 = omega2*omega;
01092   omega4 = omega2*omega2;
01093 
01094   G4double* SandiaCof = fPlatePhotoAbsCof->GetSandiaCofForMaterial(omega);
01095   G4double cross = SandiaCof[0]/omega  + SandiaCof[1]/omega2 +
01096                    SandiaCof[2]/omega3 + SandiaCof[3]/omega4;
01097   return cross;
01098 }
01099 
01100 
01102 //
01103 // Calculates formation zone for gas. Omega is energy !!!
01104 
01105 G4double G4VXTRenergyLoss::GetGasFormationZone( G4double omega ,
01106                                               G4double gamma ,
01107                                               G4double varAngle   ) 
01108 {
01109   G4double cof, lambda;
01110   lambda = 1.0/gamma/gamma + varAngle + fSigma2/omega/omega;
01111   cof = 2.0*hbarc/omega/lambda;
01112   return cof;
01113 }
01114 
01115 
01117 //
01118 // Calculates complex formation zone for gas gaps. Omega is energy !!!
01119 
01120 G4complex G4VXTRenergyLoss::GetGasComplexFZ( G4double omega ,
01121                                            G4double gamma ,
01122                                            G4double varAngle    ) 
01123 {
01124   G4double cof, length,delta, real_v, image_v;
01125 
01126   length = 0.5*GetGasFormationZone(omega,gamma,varAngle);
01127   delta  = length*GetGasLinearPhotoAbs(omega);
01128   cof    = 1.0/(1.0 + delta*delta);
01129 
01130   real_v   = length*cof;
01131   image_v  = real_v*delta;
01132 
01133   G4complex zone(real_v,image_v); 
01134   return zone;
01135 }
01136 
01137 
01138 
01140 //
01141 // Computes matrix of Sandia photo absorption cross section coefficients for
01142 // gas material
01143 
01144 void G4VXTRenergyLoss::ComputeGasPhotoAbsCof() 
01145 {
01146   const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
01147   const G4Material* mat = (*theMaterialTable)[fMatIndex2];
01148   fGasPhotoAbsCof = mat->GetSandiaTable();
01149   return;
01150 }
01151 
01153 //
01154 // Returns the value of linear photo absorption coefficient (in reciprocal 
01155 // length) for gas
01156 
01157 G4double G4VXTRenergyLoss::GetGasLinearPhotoAbs(G4double omega) 
01158 {
01159   G4double omega2, omega3, omega4; 
01160 
01161   omega2 = omega*omega;
01162   omega3 = omega2*omega;
01163   omega4 = omega2*omega2;
01164 
01165   G4double* SandiaCof = fGasPhotoAbsCof->GetSandiaCofForMaterial(omega);
01166   G4double cross = SandiaCof[0]/omega  + SandiaCof[1]/omega2 +
01167                    SandiaCof[2]/omega3 + SandiaCof[3]/omega4;
01168   return cross;
01169 
01170 }
01171 
01173 //
01174 // Calculates the product of linear cof by formation zone for plate. 
01175 // Omega is energy !!!
01176 
01177 G4double G4VXTRenergyLoss::GetPlateZmuProduct( G4double omega ,
01178                                              G4double gamma ,
01179                                              G4double varAngle   ) 
01180 {
01181   return GetPlateFormationZone(omega,gamma,varAngle)
01182     * GetPlateLinearPhotoAbs(omega);
01183 }
01185 //
01186 // Calculates the product of linear cof by formation zone for plate. 
01187 // G4cout and output in file in some energy range.
01188 
01189 void G4VXTRenergyLoss::GetPlateZmuProduct() 
01190 {
01191   std::ofstream outPlate("plateZmu.dat", std::ios::out );
01192   outPlate.setf( std::ios::scientific, std::ios::floatfield );
01193 
01194   G4int i;
01195   G4double omega, varAngle, gamma;
01196   gamma = 10000.;
01197   varAngle = 1/gamma/gamma;
01198   if(verboseLevel > 0)
01199     G4cout<<"energy, keV"<<"\t"<<"Zmu for plate"<<G4endl;
01200   for(i=0;i<100;i++)
01201   {
01202     omega = (1.0 + i)*keV;
01203     if(verboseLevel > 1)
01204       G4cout<<omega/keV<<"\t"<<GetPlateZmuProduct(omega,gamma,varAngle)<<"\t";
01205     if(verboseLevel > 0)
01206       outPlate<<omega/keV<<"\t\t"<<GetPlateZmuProduct(omega,gamma,varAngle)<<G4endl;
01207   }
01208   return;
01209 }
01210 
01212 //
01213 // Calculates the product of linear cof by formation zone for gas. 
01214 // Omega is energy !!!
01215 
01216 G4double G4VXTRenergyLoss::GetGasZmuProduct( G4double omega ,
01217                                              G4double gamma ,
01218                                              G4double varAngle   ) 
01219 {
01220   return GetGasFormationZone(omega,gamma,varAngle)*GetGasLinearPhotoAbs(omega);
01221 }
01223 //
01224 // Calculates the product of linear cof byformation zone for gas. 
01225 // G4cout and output in file in some energy range.
01226 
01227 void G4VXTRenergyLoss::GetGasZmuProduct() 
01228 {
01229   std::ofstream outGas("gasZmu.dat", std::ios::out );
01230   outGas.setf( std::ios::scientific, std::ios::floatfield );
01231   G4int i;
01232   G4double omega, varAngle, gamma;
01233   gamma = 10000.;
01234   varAngle = 1/gamma/gamma;
01235   if(verboseLevel > 0)
01236     G4cout<<"energy, keV"<<"\t"<<"Zmu for gas"<<G4endl;
01237   for(i=0;i<100;i++)
01238   {
01239     omega = (1.0 + i)*keV;
01240     if(verboseLevel > 1)
01241       G4cout<<omega/keV<<"\t"<<GetGasZmuProduct(omega,gamma,varAngle)<<"\t";
01242     if(verboseLevel > 0)
01243       outGas<<omega/keV<<"\t\t"<<GetGasZmuProduct(omega,gamma,varAngle)<<G4endl;
01244   }
01245   return;
01246 }
01247 
01249 //
01250 // Computes Compton cross section for plate material in 1/mm
01251 
01252 G4double G4VXTRenergyLoss::GetPlateCompton(G4double omega) 
01253 {
01254   G4int i, numberOfElements;
01255   G4double xSection = 0., nowZ, sumZ = 0.;
01256 
01257   const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
01258   numberOfElements = (*theMaterialTable)[fMatIndex1]->GetNumberOfElements();
01259 
01260   for( i = 0; i < numberOfElements; i++ )
01261   {
01262     nowZ      = (*theMaterialTable)[fMatIndex1]->GetElement(i)->GetZ();
01263     sumZ     += nowZ;
01264     xSection += GetComptonPerAtom(omega,nowZ); // *nowZ;
01265   }
01266   xSection /= sumZ;
01267   xSection *= (*theMaterialTable)[fMatIndex1]->GetElectronDensity();
01268   return xSection;
01269 }
01270 
01271 
01273 //
01274 // Computes Compton cross section for gas material in 1/mm
01275 
01276 G4double G4VXTRenergyLoss::GetGasCompton(G4double omega) 
01277 {
01278   G4int i, numberOfElements;
01279   G4double xSection = 0., nowZ, sumZ = 0.;
01280 
01281   const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
01282   numberOfElements = (*theMaterialTable)[fMatIndex2]->GetNumberOfElements();
01283 
01284   for( i = 0; i < numberOfElements; i++ )
01285   {
01286     nowZ      = (*theMaterialTable)[fMatIndex2]->GetElement(i)->GetZ();
01287     sumZ     += nowZ;
01288     xSection += GetComptonPerAtom(omega,nowZ); // *nowZ;
01289   }
01290   xSection /= sumZ;
01291   xSection *= (*theMaterialTable)[fMatIndex2]->GetElectronDensity();
01292   return xSection;
01293 }
01294 
01296 //
01297 // Computes Compton cross section per atom with Z electrons for gamma with
01298 // the energy GammaEnergy
01299 
01300 G4double G4VXTRenergyLoss::GetComptonPerAtom(G4double GammaEnergy, G4double Z) 
01301 {
01302   G4double CrossSection = 0.0;
01303   if ( Z < 0.9999 )                 return CrossSection;
01304   if ( GammaEnergy < 0.1*keV      ) return CrossSection;
01305   if ( GammaEnergy > (100.*GeV/Z) ) return CrossSection;
01306 
01307   static const G4double a = 20.0 , b = 230.0 , c = 440.0;
01308 
01309   static const G4double
01310   d1= 2.7965e-1*barn, d2=-1.8300e-1*barn, d3= 6.7527   *barn, d4=-1.9798e+1*barn,
01311   e1= 1.9756e-5*barn, e2=-1.0205e-2*barn, e3=-7.3913e-2*barn, e4= 2.7079e-2*barn,
01312   f1=-3.9178e-7*barn, f2= 6.8241e-5*barn, f3= 6.0480e-5*barn, f4= 3.0274e-4*barn;
01313 
01314   G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
01315            p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
01316 
01317   G4double T0  = 15.0*keV;
01318   if (Z < 1.5) T0 = 40.0*keV;
01319 
01320   G4double X   = std::max(GammaEnergy, T0) / electron_mass_c2;
01321   CrossSection = p1Z*std::log(1.+2.*X)/X
01322                + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
01323 
01324   //  modification for low energy. (special case for Hydrogen)
01325 
01326   if (GammaEnergy < T0) 
01327   {
01328     G4double dT0 = 1.*keV;
01329     X = (T0+dT0) / electron_mass_c2;
01330     G4double sigma = p1Z*std::log(1.+2*X)/X
01331                     + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
01332     G4double   c1 = -T0*(sigma-CrossSection)/(CrossSection*dT0);
01333     G4double   c2 = 0.150;
01334     if (Z > 1.5) c2 = 0.375-0.0556*std::log(Z);
01335     G4double    y = std::log(GammaEnergy/T0);
01336     CrossSection *= std::exp(-y*(c1+c2*y));
01337   }
01338   //  G4cout << "e= " << GammaEnergy << " Z= " << Z << " cross= " << CrossSection << G4endl;
01339   return CrossSection;  
01340 }
01341 
01342 
01344 //
01345 // This function returns the spectral and angle density of TR quanta
01346 // in X-ray energy region generated forward when a relativistic
01347 // charged particle crosses interface between two materials.
01348 // The high energy small theta approximation is applied.
01349 // (matter1 -> matter2, or 2->1)
01350 // varAngle =2* (1 - std::cos(theta)) or approximately = theta*theta
01351 //
01352 
01353 G4double
01354 G4VXTRenergyLoss::OneBoundaryXTRNdensity( G4double energy,G4double gamma,
01355                                          G4double varAngle ) const
01356 {
01357   G4double  formationLength1, formationLength2;
01358   formationLength1 = 1.0/
01359   (1.0/(gamma*gamma)
01360   + fSigma1/(energy*energy)
01361   + varAngle);
01362   formationLength2 = 1.0/
01363   (1.0/(gamma*gamma)
01364   + fSigma2/(energy*energy)
01365   + varAngle);
01366   return (varAngle/energy)*(formationLength1 - formationLength2)
01367               *(formationLength1 - formationLength2);
01368 
01369 }
01370 
01371 G4double G4VXTRenergyLoss::GetStackFactor( G4double energy, G4double gamma,
01372                                                      G4double varAngle )
01373 {
01374   // return stack factor corresponding to one interface
01375 
01376   return std::real( OneInterfaceXTRdEdx(energy,gamma,varAngle) );
01377 }
01378 
01380 //
01381 // For photon energy distribution tables. Integrate first over angle
01382 //
01383 
01384 G4double G4VXTRenergyLoss::XTRNSpectralAngleDensity(G4double varAngle)
01385 {
01386   return OneBoundaryXTRNdensity(fEnergy,fGamma,varAngle)*
01387          GetStackFactor(fEnergy,fGamma,varAngle);
01388 }
01389 
01391 //
01392 // For second integration over energy
01393  
01394 G4double G4VXTRenergyLoss::XTRNSpectralDensity(G4double energy)
01395 {
01396   fEnergy = energy;
01397   G4Integrator<G4VXTRenergyLoss,G4double(G4VXTRenergyLoss::*)(G4double)> integral;
01398   return integral.Legendre96(this,&G4VXTRenergyLoss::XTRNSpectralAngleDensity,
01399                              0.0,0.2*fMaxThetaTR) +
01400          integral.Legendre10(this,&G4VXTRenergyLoss::XTRNSpectralAngleDensity,
01401                              0.2*fMaxThetaTR,fMaxThetaTR);
01402 } 
01403  
01405 // 
01406 // for photon angle distribution tables
01407 //
01408 
01409 G4double G4VXTRenergyLoss::XTRNAngleSpectralDensity(G4double energy)
01410 {
01411   return OneBoundaryXTRNdensity(energy,fGamma,fVarAngle)*
01412          GetStackFactor(energy,fGamma,fVarAngle);
01413 } 
01414 
01416 //
01417 //
01418 
01419 G4double G4VXTRenergyLoss::XTRNAngleDensity(G4double varAngle) 
01420 {
01421   fVarAngle = varAngle;
01422   G4Integrator<G4VXTRenergyLoss,G4double(G4VXTRenergyLoss::*)(G4double)> integral;
01423   return integral.Legendre96(this,&G4VXTRenergyLoss::XTRNAngleSpectralDensity,
01424                              fMinEnergyTR,fMaxEnergyTR);
01425 }
01426 
01428 //
01429 // Check number of photons for a range of Lorentz factors from both energy 
01430 // and angular tables
01431 
01432 void G4VXTRenergyLoss::GetNumberOfPhotons()
01433 {
01434   G4int iTkin;
01435   G4double gamma, numberE;
01436 
01437   std::ofstream outEn("numberE.dat", std::ios::out );
01438   outEn.setf( std::ios::scientific, std::ios::floatfield );
01439 
01440   std::ofstream outAng("numberAng.dat", std::ios::out );
01441   outAng.setf( std::ios::scientific, std::ios::floatfield );
01442 
01443   for(iTkin=0;iTkin<fTotBin;iTkin++)      // Lorentz factor loop
01444   {
01445      gamma = 1.0 + (fProtonEnergyVector->
01446                             GetLowEdgeEnergy(iTkin)/proton_mass_c2);
01447      numberE = (*(*fEnergyDistrTable)(iTkin))(0);
01448      //  numberA = (*(*fAngleDistrTable)(iTkin))(0);
01449      if(verboseLevel > 1)
01450        G4cout<<gamma<<"\t\t"<<numberE<<"\t"    //  <<numberA
01451              <<G4endl; 
01452      if(verboseLevel > 0)
01453        outEn<<gamma<<"\t\t"<<numberE<<G4endl; 
01454   }
01455   return;
01456 }  
01457 
01459 //
01460 // Returns randon energy of a X-ray TR photon for given scaled kinetic energy
01461 // of a charged particle
01462 
01463 G4double G4VXTRenergyLoss::GetXTRrandomEnergy( G4double scaledTkin, G4int iTkin )
01464 {
01465   G4int iTransfer, iPlace;
01466   G4double transfer = 0.0, position, E1, E2, W1, W2, W;
01467 
01468   iPlace = iTkin - 1;
01469 
01470   //  G4cout<<"iPlace = "<<iPlace<<endl;
01471 
01472   if(iTkin == fTotBin) // relativistic plato, try from left
01473   {
01474       position = (*(*fEnergyDistrTable)(iPlace))(0)*G4UniformRand();
01475 
01476       for(iTransfer=0;;iTransfer++)
01477       {
01478         if(position >= (*(*fEnergyDistrTable)(iPlace))(iTransfer)) break;
01479       }
01480       transfer = GetXTRenergy(iPlace,position,iTransfer);
01481   }
01482   else
01483   {
01484     E1 = fProtonEnergyVector->GetLowEdgeEnergy(iTkin - 1); 
01485     E2 = fProtonEnergyVector->GetLowEdgeEnergy(iTkin);
01486     W  = 1.0/(E2 - E1);
01487     W1 = (E2 - scaledTkin)*W;
01488     W2 = (scaledTkin - E1)*W;
01489 
01490     position =( (*(*fEnergyDistrTable)(iPlace))(0)*W1 + 
01491                     (*(*fEnergyDistrTable)(iPlace+1))(0)*W2 )*G4UniformRand();
01492 
01493         // G4cout<<position<<"\t";
01494 
01495     for(iTransfer=0;;iTransfer++)
01496     {
01497           if( position >=
01498           ( (*(*fEnergyDistrTable)(iPlace))(iTransfer)*W1 + 
01499             (*(*fEnergyDistrTable)(iPlace+1))(iTransfer)*W2) ) break;
01500     }
01501     transfer = GetXTRenergy(iPlace,position,iTransfer);
01502     
01503   } 
01504   //  G4cout<<"XTR transfer = "<<transfer/keV<<" keV"<<endl; 
01505   if(transfer < 0.0 ) transfer = 0.0;
01506   return transfer;
01507 }
01508 
01510 //
01511 // Returns approximate position of X-ray photon energy during random sampling
01512 // over integral energy distribution
01513 
01514 G4double G4VXTRenergyLoss::GetXTRenergy( G4int    iPlace, 
01515                                        G4double  /* position */, 
01516                                        G4int    iTransfer )
01517 {
01518   G4double x1, x2, y1, y2, result;
01519 
01520   if(iTransfer == 0)
01521   {
01522     result = (*fEnergyDistrTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
01523   }  
01524   else
01525   {
01526     y1 = (*(*fEnergyDistrTable)(iPlace))(iTransfer-1);
01527     y2 = (*(*fEnergyDistrTable)(iPlace))(iTransfer);
01528 
01529     x1 = (*fEnergyDistrTable)(iPlace)->GetLowEdgeEnergy(iTransfer-1);
01530     x2 = (*fEnergyDistrTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
01531 
01532     if ( x1 == x2 )    result = x2;
01533     else
01534     {
01535       if ( y1 == y2  ) result = x1 + (x2 - x1)*G4UniformRand();
01536       else
01537       {
01538         // result = x1 + (position - y1)*(x2 - x1)/(y2 - y1);
01539         result = x1 + (x2 - x1)*G4UniformRand();
01540       }
01541     }
01542   }
01543   return result;
01544 }
01545 
01547 //
01548 //  Get XTR photon angle at given energy and Tkin
01549 
01550 G4double G4VXTRenergyLoss::GetRandomAngle( G4double energyXTR, G4int iTkin )
01551 {
01552   G4int iTR, iAngle;
01553   G4double position, angle;
01554 
01555   if (iTkin == fTotBin) iTkin--;
01556 
01557   fAngleForEnergyTable = fAngleBank[iTkin];
01558 
01559   for( iTR = 0; iTR < fBinTR; iTR++ )
01560   {
01561     if( energyXTR < fXTREnergyVector->GetLowEdgeEnergy(iTR) )  break;    
01562   }
01563   if (iTR == fBinTR) iTR--;
01564       
01565   position = (*(*fAngleForEnergyTable)(iTR))(0)*G4UniformRand();
01566 
01567   for(iAngle = 0;;iAngle++)
01568   {
01569     if(position >= (*(*fAngleForEnergyTable)(iTR))(iAngle)) break;
01570   }
01571   angle = GetAngleXTR(iTR,position,iAngle);
01572   return angle;
01573 }
01574 
01576 //
01577 // Returns approximate position of X-ray photon angle at given energy during random sampling
01578 // over integral energy distribution
01579 
01580 G4double G4VXTRenergyLoss::GetAngleXTR( G4int    iPlace, 
01581                                        G4double position, 
01582                                        G4int    iTransfer )
01583 {
01584   G4double x1, x2, y1, y2, result;
01585 
01586   if(iTransfer == 0)
01587   {
01588     result = (*fAngleForEnergyTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
01589   }  
01590   else
01591   {
01592     y1 = (*(*fAngleForEnergyTable)(iPlace))(iTransfer-1);
01593     y2 = (*(*fAngleForEnergyTable)(iPlace))(iTransfer);
01594 
01595     x1 = (*fAngleForEnergyTable)(iPlace)->GetLowEdgeEnergy(iTransfer-1);
01596     x2 = (*fAngleForEnergyTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
01597 
01598     if ( x1 == x2 )    result = x2;
01599     else
01600     {
01601       if ( y1 == y2  ) result = x1 + (x2 - x1)*G4UniformRand();
01602       else
01603       {
01604         result = x1 + (position - y1)*(x2 - x1)/(y2 - y1);
01605       }
01606     }
01607   }
01608   return result;
01609 }
01610 
01611 
01612 //
01613 //
01615 

Generated on Mon May 27 17:50:24 2013 for Geant4 by  doxygen 1.4.7