G4XTRRegularRadModel.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 
00028 #include <complex>
00029 
00030 #include "G4XTRRegularRadModel.hh"
00031 #include "G4PhysicalConstants.hh"
00032 #include "Randomize.hh"
00033 
00034 #include "G4Gamma.hh"
00035 using namespace std;
00036 
00038 //
00039 // Constructor, destructor
00040 
00041 G4XTRRegularRadModel::G4XTRRegularRadModel(G4LogicalVolume *anEnvelope,
00042                                          G4Material* foilMat,G4Material* gasMat,
00043                                          G4double a, G4double b, G4int n,
00044                                          const G4String& processName) :
00045   G4VXTRenergyLoss(anEnvelope,foilMat,gasMat,a,b,n,processName)  
00046 {
00047   G4cout<<" XTR Regular discrete radiator model is called"<<G4endl ;
00048 
00049   fExitFlux = true;
00050 
00051   // Build energy and angular integral spectra of X-ray TR photons from
00052   // a radiator
00053 
00054   // BuildTable() ;
00055 }
00056 
00058 
00059 G4XTRRegularRadModel::~G4XTRRegularRadModel()
00060 {
00061   ;
00062 }
00063 
00065 //
00066 //
00067 
00068 G4double G4XTRRegularRadModel::SpectralXTRdEdx(G4double energy)
00069 {
00070   G4double result, sum = 0., tmp, cof1, cof2, cofMin, cofPHC, theta2, theta2k; 
00071     G4double aMa, bMb ,sigma, dump;
00072   G4int k, kMax, kMin;
00073 
00074   aMa = fPlateThick*GetPlateLinearPhotoAbs(energy);
00075   bMb = fGasThick*GetGasLinearPhotoAbs(energy);
00076   sigma = 0.5*(aMa + bMb);
00077   dump = std::exp(-fPlateNumber*sigma);
00078   if(verboseLevel > 2)  G4cout<<" dump = "<<dump<<G4endl;  
00079   cofPHC  = 4*pi*hbarc;
00080   tmp     = (fSigma1 - fSigma2)/cofPHC/energy;  
00081   cof1    = fPlateThick*tmp;
00082   cof2    = fGasThick*tmp;
00083 
00084   cofMin  =  energy*(fPlateThick + fGasThick)/fGamma/fGamma;
00085   cofMin += (fPlateThick*fSigma1 + fGasThick*fSigma2)/energy;
00086   cofMin /= cofPHC;
00087 
00088   theta2 = cofPHC/(energy*(fPlateThick + fGasThick));
00089 
00090   //  if (fGamma < 1200) kMin = G4int(cofMin);  // 1200 ?
00091   // else               kMin = 1;
00092 
00093 
00094   kMin = G4int(cofMin);
00095   if (cofMin > kMin) kMin++;
00096 
00097   // tmp  = (fPlateThick + fGasThick)*energy*fMaxThetaTR;
00098   // tmp /= cofPHC;
00099   // kMax = G4int(tmp);
00100   // if(kMax < 0) kMax = 0;
00101   // kMax += kMin;
00102   
00103 
00104   kMax = kMin + 49; //  19; // kMin + G4int(tmp);
00105 
00106   // tmp /= fGamma;
00107   // if( G4int(tmp) < kMin ) kMin = G4int(tmp);
00108 
00109   if(verboseLevel > 2)
00110   {    
00111     G4cout<<cof1<<"     "<<cof2<<"        "<<cofMin<<G4endl;
00112     G4cout<<"kMin = "<<kMin<<";    kMax = "<<kMax<<G4endl;
00113   }
00114   for( k = kMin; k <= kMax; k++ )
00115   {
00116     tmp    = pi*fPlateThick*(k + cof2)/(fPlateThick + fGasThick);
00117     result = (k - cof1)*(k - cof1)*(k + cof2)*(k + cof2);
00118     // tmp = std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
00119     if( k == kMin && kMin == G4int(cofMin) )
00120     {
00121       sum   += 0.5*std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
00122     }
00123     else
00124     {
00125       sum   += std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
00126     }
00127     theta2k = std::sqrt(theta2*std::abs(k-cofMin));
00128 
00129     if(verboseLevel > 2)
00130     {    
00131       // G4cout<<"k = "<<k<<"; sqrt(theta2k) = "<<theta2k<<"; tmp = "<<std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result
00132       //     <<";    sum = "<<sum<<G4endl;  
00133       G4cout<<k<<"   "<<theta2k<<"     "<<std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result
00134               <<"      "<<sum<<G4endl;  
00135     }  
00136   }
00137   result = 2*( cof1 + cof2 )*( cof1 + cof2 )*sum/energy;
00138   // result *= ( 1 - std::exp(-0.5*fPlateNumber*sigma) )/( 1 - std::exp(-0.5*sigma) );  
00139   // fPlateNumber;
00140   result *= dump*( -1 + dump + 2*fPlateNumber ); 
00141   /*  
00142   fEnergy = energy;
00143   //  G4Integrator<G4VXTRenergyLoss,G4double(G4VXTRenergyLoss::*)(G4double)> integral;
00144   G4Integrator<G4TransparentRegXTRadiator,G4double(G4VXTRenergyLoss::*)(G4double)> integral;
00145  
00146   tmp = integral.Legendre96(this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
00147                              0.0,0.3*fMaxThetaTR) +
00148       integral.Legendre96(this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
00149                              0.3*fMaxThetaTR,0.6*fMaxThetaTR) +         
00150       integral.Legendre96(this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
00151                              0.6*fMaxThetaTR,fMaxThetaTR) ;
00152   result += tmp;
00153   */
00154   return result;
00155 }
00156 
00157 
00158 
00160 //
00161 // Approximation for radiator interference factor for the case of
00162 // fully Regular radiator. The plate and gas gap thicknesses are fixed .
00163 // The mean values of the plate and gas gap thicknesses 
00164 // are supposed to be about XTR formation zones but much less than 
00165 // mean absorption length of XTR photons in coresponding material.
00166 
00167 G4double 
00168 G4XTRRegularRadModel::GetStackFactor( G4double energy, 
00169                                          G4double gamma, G4double varAngle )
00170 {
00171   G4double result, Qa, Qb, Q, aZa, bZb, aMa, bMb, I2 ;
00172   
00173   aZa = fPlateThick/GetPlateFormationZone(energy,gamma,varAngle) ;
00174   bZb = fGasThick/GetGasFormationZone(energy,gamma,varAngle) ;
00175 
00176   aMa = fPlateThick*GetPlateLinearPhotoAbs(energy) ;
00177   bMb = fGasThick*GetGasLinearPhotoAbs(energy) ;
00178 
00179   Qa = std::exp(-aMa) ;
00180   Qb = std::exp(-bMb) ;
00181   Q  = Qa*Qb ;
00182 
00183   //  G4complex Ca(1.0+0.5*fPlateThick*Ma,fPlateThick/Za) ; 
00184   //  G4complex Cb(1.0+0.5*fGasThick*Mb,fGasThick/Zb) ; 
00185 
00186   G4complex Ha( std::exp(-0.5*aMa)*std::cos(aZa),
00187                -std::exp(-0.5*aMa)*std::sin(aZa)   ) ; 
00188  
00189   G4complex Hb( std::exp(-0.5*bMb)*std::cos(bZb),
00190                -std::exp(-0.5*bMb)*std::sin(bZb)    ) ;
00191 
00192   G4complex H  = Ha*Hb ;
00193 
00194   G4complex Hs = std::conj(H) ;
00195 
00196   //  G4complex F1 = ( 0.5*(1+Qa)*(1+H) - Ha - Qa*Hb )/(1-H) ;
00197 
00198   G4complex F2 = (1.0-Ha)*(Qa-Ha)*Hb*(1.0-Hs)*(Q-Hs) ;
00199 
00200   F2          *= std::pow(Q,G4double(fPlateNumber)) - std::pow(H,fPlateNumber) ;
00201 
00202   result       = ( 1 - std::pow(Q,G4double(fPlateNumber)) )/( 1 - Q ) ;
00203 
00204   result      *= (1 - Qa)*(1 + Qa - 2*std::sqrt(Qa)*std::cos(aZa)) ;
00205 
00206   result      /= (1 - std::sqrt(Q))*(1 - std::sqrt(Q)) + 
00207                   4*std::sqrt(Q)*std::sin(0.5*(aZa+bZb))*std::sin(0.5*(aZa+bZb)) ;
00208 
00209   I2           = 1.; // 2.0*std::real(F2) ;
00210 
00211   I2           /= (1 - std::sqrt(Q))*(1 - std::sqrt(Q)) + 
00212                   4*std::sqrt(Q)*std::sin(0.5*(aZa+bZb))*std::sin(0.5*(aZa+bZb)) ;
00213 
00214   I2           /= Q*( (std::sqrt(Q)-std::cos(aZa+bZb))*(std::sqrt(Q)-std::cos(aZa+bZb)) + 
00215                       std::sin(aZa+bZb)*std::sin(aZa+bZb)   ) ;
00216 
00217   G4complex stack  = 2.*I2*F2;
00218             stack += result;
00219             stack *= OneInterfaceXTRdEdx(energy,gamma,varAngle);
00220 
00221             // result       += I2 ;
00222   result = std::real(stack);
00223 
00224   return      result ;
00225 }
00226 
00227 
00228 //
00229 //
00231 
00232 
00233 
00234 
00235 
00236 
00237 
00238 

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