G4AdjointPhotoElectricModel.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 // $Id: G4AdjointPhotoElectricModel.cc 69844 2013-05-16 09:19:33Z gcosmo $
00027 //
00028 #include "G4AdjointPhotoElectricModel.hh"
00029 #include "G4AdjointCSManager.hh"
00030 
00031 #include "G4PhysicalConstants.hh"
00032 #include "G4Integrator.hh"
00033 #include "G4TrackStatus.hh"
00034 #include "G4ParticleChange.hh"
00035 #include "G4AdjointElectron.hh"
00036 #include  "G4Gamma.hh"
00037 #include "G4AdjointGamma.hh"
00038 
00039 
00041 //
00042 G4AdjointPhotoElectricModel::G4AdjointPhotoElectricModel():
00043  G4VEmAdjointModel("AdjointPEEffect")
00044 
00045 { SetUseMatrix(false);
00046   SetApplyCutInRange(false);
00047 
00048   //Initialization
00049   current_eEnergy =0.;
00050   totAdjointCS=0.;
00051   factorCSBiasing =1.;
00052   post_step_AdjointCS =0.;
00053   pre_step_AdjointCS =0.;
00054   totBiasedAdjointCS =0.;
00055 
00056   index_element=0;
00057 
00058   theAdjEquivOfDirectPrimPartDef =G4AdjointGamma::AdjointGamma();
00059   theAdjEquivOfDirectSecondPartDef=G4AdjointElectron::AdjointElectron();
00060   theDirectPrimaryPartDef=G4Gamma::Gamma();
00061   second_part_of_same_type=false;
00062   theDirectPEEffectModel = new G4PEEffectModel();
00063 }
00065 //
00066 G4AdjointPhotoElectricModel::~G4AdjointPhotoElectricModel()
00067 {;}
00068 
00070 //
00071 void G4AdjointPhotoElectricModel::SampleSecondaries(const G4Track& aTrack,
00072                                 G4bool IsScatProjToProjCase,
00073                                 G4ParticleChange* fParticleChange)
00074 { if (IsScatProjToProjCase) return ;
00075 
00076   //Compute the totAdjointCS vectors if not already done for the current couple and electron energy
00077   //-----------------------------------------------------------------------------------------------
00078   const G4MaterialCutsCouple* aCouple = aTrack.GetMaterialCutsCouple();
00079   const G4DynamicParticle* aDynPart =  aTrack.GetDynamicParticle() ;
00080   G4double electronEnergy = aDynPart->GetKineticEnergy();
00081   G4ThreeVector electronDirection= aDynPart->GetMomentumDirection() ;
00082   pre_step_AdjointCS = totAdjointCS; //The last computed CS was  at pre step point
00083   post_step_AdjointCS = AdjointCrossSection(aCouple, electronEnergy,IsScatProjToProjCase);
00084   post_step_AdjointCS = totAdjointCS; 
00085                                 
00086 
00087 
00088 
00089   //Sample element
00090   //-------------
00091    const G4ElementVector* theElementVector = currentMaterial->GetElementVector();
00092    size_t nelm =  currentMaterial->GetNumberOfElements();
00093    G4double rand_CS= G4UniformRand()*xsec[nelm-1];
00094    for (index_element=0; index_element<nelm-1; index_element++){
00095         if (rand_CS<xsec[index_element]) break;
00096    }
00097         
00098    //Sample shell and binding energy
00099    //-------------
00100    G4int nShells = (*theElementVector)[index_element]->GetNbOfAtomicShells();
00101    rand_CS= shell_prob[index_element][nShells-1]*G4UniformRand();
00102    G4int i  = 0;  
00103    for (i=0; i<nShells-1; i++){
00104         if (rand_CS<shell_prob[index_element][i]) break;
00105    }
00106    G4double gammaEnergy= electronEnergy+(*theElementVector)[index_element]->GetAtomicShell(i);
00107         
00108   //Sample cos theta
00109   //Copy of the G4PEEffectModel cos theta sampling method ElecCosThetaDistribution.   
00110   //This method cannot be used directly from G4PEEffectModel because it is a friend method. I should ask Vladimir to change that  
00111   //------------------------------------------------------------------------------------------------    
00112   //G4double cos_theta = theDirectPEEffectModel->ElecCosThetaDistribution(electronEnergy);
00113         
00114    G4double  cos_theta = 1.;
00115    G4double gamma   = 1. + electronEnergy/electron_mass_c2;
00116    if (gamma <= 5.) {
00117         G4double beta  = std::sqrt(gamma*gamma-1.)/gamma;
00118         G4double b     = 0.5*gamma*(gamma-1.)*(gamma-2);
00119 
00120         G4double rndm,term,greject,grejsup;
00121         if (gamma < 2.) grejsup = gamma*gamma*(1.+b-beta*b);
00122         else            grejsup = gamma*gamma*(1.+b+beta*b);
00123 
00124         do {    rndm = 1.-2*G4UniformRand();
00125                 cos_theta = (rndm+beta)/(rndm*beta+1.);
00126                 term = 1.-beta*cos_theta;
00127                 greject = (1.-cos_theta*cos_theta)*(1.+b*term)/(term*term);
00128         } while(greject < G4UniformRand()*grejsup);
00129   }
00130         
00131   // direction of the adjoint gamma electron
00132   //---------------------------------------
00133   
00134      
00135   G4double sin_theta = std::sqrt(1.-cos_theta*cos_theta);
00136   G4double Phi     = twopi * G4UniformRand();
00137   G4double dirx = sin_theta*std::cos(Phi),diry = sin_theta*std::sin(Phi),dirz = cos_theta;
00138   G4ThreeVector adjoint_gammaDirection(dirx,diry,dirz);
00139   adjoint_gammaDirection.rotateUz(electronDirection);
00140   
00141   
00142   
00143   //Weight correction
00144  //-----------------------                                         
00145   CorrectPostStepWeight(fParticleChange, aTrack.GetWeight(), electronEnergy,gammaEnergy,IsScatProjToProjCase);  
00146  
00147   
00148   
00149   //Create secondary and modify fParticleChange 
00150   //--------------------------------------------
00151   G4DynamicParticle* anAdjointGamma = new G4DynamicParticle (
00152                        G4AdjointGamma::AdjointGamma(),adjoint_gammaDirection, gammaEnergy);
00153   
00154   
00155   
00156  
00157   
00158   fParticleChange->ProposeTrackStatus(fStopAndKill);
00159   fParticleChange->AddSecondary(anAdjointGamma);    
00160         
00161 
00162   
00163 
00164 }
00165 
00167 //
00168 void G4AdjointPhotoElectricModel::CorrectPostStepWeight(G4ParticleChange* fParticleChange, 
00169                                                             G4double old_weight,  
00170                                                             G4double adjointPrimKinEnergy, 
00171                                                             G4double projectileKinEnergy ,
00172                                                             G4bool  ) 
00173 {
00174  G4double new_weight=old_weight;
00175 
00176  G4double w_corr =G4AdjointCSManager::GetAdjointCSManager()->GetPostStepWeightCorrection()/factorCSBiasing;
00177  w_corr*=post_step_AdjointCS/pre_step_AdjointCS; 
00178  new_weight*=w_corr; 
00179  new_weight*=projectileKinEnergy/adjointPrimKinEnergy;
00180  fParticleChange->SetParentWeightByProcess(false);
00181  fParticleChange->SetSecondaryWeightByProcess(false);
00182  fParticleChange->ProposeParentWeight(new_weight);      
00183 }       
00184 
00186 //                      
00187 
00188 G4double G4AdjointPhotoElectricModel::AdjointCrossSection(const G4MaterialCutsCouple* aCouple,
00189                                 G4double electronEnergy,
00190                                 G4bool IsScatProjToProjCase)
00191 { 
00192   
00193 
00194   if (IsScatProjToProjCase)  return 0.;
00195 
00196         
00197   if (aCouple !=currentCouple || current_eEnergy !=electronEnergy) {
00198         totAdjointCS = 0.;
00199         DefineCurrentMaterialAndElectronEnergy(aCouple, electronEnergy);
00200         const G4ElementVector* theElementVector = currentMaterial->GetElementVector();
00201         const double* theAtomNumDensityVector = currentMaterial->GetVecNbOfAtomsPerVolume();
00202         size_t nelm =  currentMaterial->GetNumberOfElements();
00203         for (index_element=0;index_element<nelm;index_element++){
00204                 
00205                 totAdjointCS +=AdjointCrossSectionPerAtom((*theElementVector)[index_element],electronEnergy)*theAtomNumDensityVector[index_element];
00206                 xsec[index_element] = totAdjointCS;
00207         } 
00208 
00209         totBiasedAdjointCS=std::min(totAdjointCS,0.01);
00210 //      totBiasedAdjointCS=totAdjointCS;
00211         factorCSBiasing = totBiasedAdjointCS/totAdjointCS;
00212         lastCS=totBiasedAdjointCS;
00213         
00214         
00215   }
00216   return totBiasedAdjointCS;
00217 
00218   
00219 }
00221 //                      
00222 
00223 G4double G4AdjointPhotoElectricModel::GetAdjointCrossSection(const G4MaterialCutsCouple* aCouple,
00224                                 G4double electronEnergy,
00225                                 G4bool IsScatProjToProjCase)
00226 {       return AdjointCrossSection(aCouple,electronEnergy,IsScatProjToProjCase);
00227 }
00229 //                      
00230 
00231 G4double G4AdjointPhotoElectricModel::AdjointCrossSectionPerAtom(const G4Element* anElement,G4double electronEnergy)
00232 { 
00233   G4int nShells = anElement->GetNbOfAtomicShells();
00234   G4double Z= anElement->GetZ();
00235   G4int i  = 0;  
00236   G4double B0=anElement->GetAtomicShell(0);
00237   G4double gammaEnergy = electronEnergy+B0;
00238   G4double CS= theDirectPEEffectModel->ComputeCrossSectionPerAtom(G4Gamma::Gamma(),gammaEnergy,Z,0.,0.,0.);
00239   G4double adjointCS =0.;
00240   if (CS >0) adjointCS += CS/gammaEnergy; 
00241   shell_prob[index_element][0] = adjointCS;                                          
00242   for (i=1;i<nShells;i++){
00243         //G4cout<<i<<G4endl;
00244         G4double Bi_= anElement->GetAtomicShell(i-1);
00245         G4double Bi = anElement->GetAtomicShell(i);
00246         //G4cout<<Bi_<<'\t'<<Bi<<G4endl;
00247         if (electronEnergy <Bi_-Bi) {
00248                 gammaEnergy = electronEnergy+Bi;
00249                 
00250                 CS=theDirectPEEffectModel->ComputeCrossSectionPerAtom(G4Gamma::Gamma(),gammaEnergy,Z,0.,0.,0.);
00251                 if (CS>0) adjointCS +=CS/gammaEnergy;
00252         }
00253         shell_prob[index_element][i] = adjointCS;       
00254   
00255   }
00256   adjointCS*=electronEnergy;
00257   return adjointCS;
00258   
00259 }                               
00261 //                      
00262 
00263 void G4AdjointPhotoElectricModel::DefineCurrentMaterialAndElectronEnergy(const G4MaterialCutsCouple* couple, G4double anEnergy)
00264 { currentCouple   = const_cast<G4MaterialCutsCouple*> (couple);
00265   currentMaterial = const_cast<G4Material*> (couple->GetMaterial());
00266   currentCoupleIndex = couple->GetIndex();
00267   currentMaterialIndex = currentMaterial->GetIndex();
00268   current_eEnergy = anEnergy;   
00269 }

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