G4PolarizedCompton.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: G4PolarizedCompton.cc 69847 2013-05-16 09:36:18Z gcosmo $
00028 // 
00029 //
00030 // File name:     G4PolarizedCompton
00031 //
00032 // Author:        Andreas Schaelicke
00033 //                based on code by Michel Maire / Vladimir IVANTCHENKO
00034 // Class description
00035 //
00036 // modified version respecting media and beam polarization
00037 //     using the stokes formalism
00038 //
00039 // Creation date: 01.05.2005
00040 //
00041 // Modifications:
00042 //
00043 // 01-01-05, include polarization description (A.Stahl)
00044 // 01-01-05, create asymmetry table and determine interactionlength (A.Stahl)
00045 // 01-05-05, update handling of media polarization (A.Schalicke)
00046 // 01-05-05, update polarized differential cross section (A.Schalicke)
00047 // 20-05-05, added polarization transfer (A.Schalicke)
00048 // 10-06-05, transformation between different reference frames (A.Schalicke)
00049 // 17-10-05, correct reference frame dependence in GetMeanFreePath (A.Schalicke)
00050 // 26-07-06, cross section recalculated (P.Starovoitov)
00051 // 09-08-06, make it work under current geant4 release (A.Schalicke)
00052 // 11-06-07, add PostStepGetPhysicalInteractionLength (A.Schalicke)
00053 // -----------------------------------------------------------------------------
00054 
00055 
00056 #include "G4PolarizedCompton.hh"
00057 #include "G4SystemOfUnits.hh"
00058 #include "G4Electron.hh"
00059 
00060 #include "G4StokesVector.hh"
00061 #include "G4PolarizationManager.hh"
00062 #include "G4PolarizedComptonModel.hh"
00063 #include "G4ProductionCutsTable.hh"
00064 #include "G4PhysicsTableHelper.hh"
00065 #include "G4KleinNishinaCompton.hh"
00066 #include "G4PolarizedComptonModel.hh"
00067 
00068 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00069 
00070 G4PolarizedCompton::G4PolarizedCompton(const G4String& processName,
00071   G4ProcessType type):
00072   G4VEmProcess (processName, type),
00073   buildAsymmetryTable(true),
00074   useAsymmetryTable(true),
00075   isInitialised(false),
00076   selectedModel(0),
00077   mType(10),
00078   theAsymmetryTable(NULL)
00079 {
00080   SetLambdaBinning(90);
00081   SetMinKinEnergy(0.1*keV);
00082   SetMaxKinEnergy(100.0*GeV);
00083   SetProcessSubType(fComptonScattering);
00084   emModel = 0;
00085 }
00086 
00087 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00088  
00089 G4PolarizedCompton::~G4PolarizedCompton()
00090 {
00091   if (theAsymmetryTable) {
00092     theAsymmetryTable->clearAndDestroy();
00093     delete theAsymmetryTable;
00094   }
00095 }
00096 
00097 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00098 
00099 void G4PolarizedCompton::InitialiseProcess(const G4ParticleDefinition*)
00100 {
00101   if(!isInitialised) {
00102     isInitialised = true;
00103     SetBuildTableFlag(true);
00104     SetSecondaryParticle(G4Electron::Electron());
00105     G4double emin = MinKinEnergy();
00106     G4double emax = MaxKinEnergy();
00107     emModel = new G4PolarizedComptonModel();
00108     if(0 == mType) selectedModel = new G4KleinNishinaCompton();
00109     else if(10 == mType) selectedModel = emModel;
00110     selectedModel->SetLowEnergyLimit(emin);
00111     selectedModel->SetHighEnergyLimit(emax);
00112     AddEmModel(1, selectedModel);
00113   } 
00114 }
00115 
00116 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00117 
00118 void G4PolarizedCompton::PrintInfo()
00119 {
00120   G4cout << " Total cross sections has a good parametrisation"
00121          << " from 10 KeV to (100/Z) GeV" 
00122          << "\n      Sampling according " << selectedModel->GetName() << " model" 
00123          << G4endl;
00124 }         
00125 
00126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00127 
00128 void G4PolarizedCompton::SetModel(const G4String& ss)
00129 {
00130   if(ss == "Klein-Nishina") mType = 0;
00131   if(ss == "Polarized-Compton") mType = 10;
00132 }
00133 
00134 
00135 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00136 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00137 
00138 
00139 
00140 G4double G4PolarizedCompton::GetMeanFreePath(
00141                                    const G4Track& aTrack,
00142                                    G4double   previousStepSize,
00143                                    G4ForceCondition* condition)
00144 {
00145   // *** get unploarised mean free path from lambda table ***
00146   G4double mfp = G4VEmProcess::GetMeanFreePath(aTrack, previousStepSize, condition);
00147 
00148 
00149    if (theAsymmetryTable && useAsymmetryTable) {
00150      // *** get asymmetry, if target is polarized ***
00151      const G4DynamicParticle* aDynamicGamma = aTrack.GetDynamicParticle();
00152      const G4double GammaEnergy = aDynamicGamma->GetKineticEnergy();
00153      const G4StokesVector GammaPolarization = aTrack.GetPolarization();
00154      const G4ParticleMomentum GammaDirection0 = aDynamicGamma->GetMomentumDirection();
00155 
00156      G4Material*         aMaterial = aTrack.GetMaterial();
00157      G4VPhysicalVolume*  aPVolume  = aTrack.GetVolume();
00158      G4LogicalVolume*    aLVolume  = aPVolume->GetLogicalVolume();
00159 
00160      //   G4Material* bMaterial = aLVolume->GetMaterial();
00161      G4PolarizationManager * polarizationManger = G4PolarizationManager::GetInstance();
00162 
00163      const G4bool VolumeIsPolarized = polarizationManger->IsPolarized(aLVolume);
00164      G4StokesVector ElectronPolarization = polarizationManger->GetVolumePolarization(aLVolume);
00165 
00166      if (!VolumeIsPolarized || mfp == DBL_MAX) return mfp;
00167      
00168      if (verboseLevel>=2) {
00169 
00170        G4cout << " Mom " << GammaDirection0  << G4endl;
00171        G4cout << " Polarization " << GammaPolarization  << G4endl;
00172        G4cout << " MaterialPol. " << ElectronPolarization  << G4endl;
00173        G4cout << " Phys. Volume " << aPVolume->GetName() << G4endl;
00174        G4cout << " Log. Volume  " << aLVolume->GetName() << G4endl;
00175        G4cout << " Material     " << aMaterial          << G4endl;
00176      }
00177 
00178      G4int midx= CurrentMaterialCutsCoupleIndex();
00179      G4PhysicsVector * aVector=(*theAsymmetryTable)(midx);
00180      
00181      G4double asymmetry=0;
00182      if (aVector) {
00183        G4bool isOutRange;
00184        asymmetry = aVector->GetValue(GammaEnergy, isOutRange);
00185      } else {
00186        G4cout << " MaterialIndex     " << midx << " is out of range \n";
00187        asymmetry=0;
00188      }
00189 
00190      //  we have to determine angle between particle motion 
00191      //  and target polarisation here  
00192      //      circ pol * Vec(ElectronPol)*Vec(PhotonMomentum)
00193      //  both vectors in global reference frame
00194      
00195      G4double pol=ElectronPolarization*GammaDirection0;
00196      
00197      G4double polProduct = GammaPolarization.p3() * pol;
00198      mfp *= 1. / ( 1. + polProduct * asymmetry );
00199 
00200      if (verboseLevel>=2) {
00201        G4cout << " MeanFreePath:  " << mfp / mm << " mm " << G4endl;
00202        G4cout << " Asymmetry:     " << asymmetry          << G4endl;
00203        G4cout << " PolProduct:    " << polProduct         << G4endl;
00204      }
00205    }
00206 
00207    return mfp;
00208 }
00209 
00210 G4double G4PolarizedCompton::PostStepGetPhysicalInteractionLength(
00211                                    const G4Track& aTrack,
00212                                    G4double   previousStepSize,
00213                                    G4ForceCondition* condition)
00214 {
00215   // *** get unploarised mean free path from lambda table ***
00216   G4double mfp = G4VEmProcess::PostStepGetPhysicalInteractionLength(aTrack, previousStepSize, condition);
00217 
00218 
00219    if (theAsymmetryTable && useAsymmetryTable) {
00220      // *** get asymmetry, if target is polarized ***
00221      const G4DynamicParticle* aDynamicGamma = aTrack.GetDynamicParticle();
00222      const G4double GammaEnergy = aDynamicGamma->GetKineticEnergy();
00223      const G4StokesVector GammaPolarization = aTrack.GetPolarization();
00224      const G4ParticleMomentum GammaDirection0 = aDynamicGamma->GetMomentumDirection();
00225 
00226      G4Material*         aMaterial = aTrack.GetMaterial();
00227      G4VPhysicalVolume*  aPVolume  = aTrack.GetVolume();
00228      G4LogicalVolume*    aLVolume  = aPVolume->GetLogicalVolume();
00229 
00230      //   G4Material* bMaterial = aLVolume->GetMaterial();
00231      G4PolarizationManager * polarizationManger = G4PolarizationManager::GetInstance();
00232 
00233      const G4bool VolumeIsPolarized = polarizationManger->IsPolarized(aLVolume);
00234      G4StokesVector ElectronPolarization = polarizationManger->GetVolumePolarization(aLVolume);
00235 
00236      if (!VolumeIsPolarized || mfp == DBL_MAX) return mfp;
00237      
00238      if (verboseLevel>=2) {
00239 
00240        G4cout << " Mom " << GammaDirection0  << G4endl;
00241        G4cout << " Polarization " << GammaPolarization  << G4endl;
00242        G4cout << " MaterialPol. " << ElectronPolarization  << G4endl;
00243        G4cout << " Phys. Volume " << aPVolume->GetName() << G4endl;
00244        G4cout << " Log. Volume  " << aLVolume->GetName() << G4endl;
00245        G4cout << " Material     " << aMaterial          << G4endl;
00246      }
00247 
00248      G4int midx= CurrentMaterialCutsCoupleIndex();
00249      G4PhysicsVector * aVector=(*theAsymmetryTable)(midx);
00250      
00251      G4double asymmetry=0;
00252      if (aVector) {
00253        G4bool isOutRange;
00254        asymmetry = aVector->GetValue(GammaEnergy, isOutRange);
00255      } else {
00256        G4cout << " MaterialIndex     " << midx << " is out of range \n";
00257        asymmetry=0;
00258      }
00259 
00260      //  we have to determine angle between particle motion 
00261      //  and target polarisation here  
00262      //      circ pol * Vec(ElectronPol)*Vec(PhotonMomentum)
00263      //  both vectors in global reference frame
00264      
00265      G4double pol=ElectronPolarization*GammaDirection0;
00266      
00267      G4double polProduct = GammaPolarization.p3() * pol;
00268      mfp *= 1. / ( 1. + polProduct * asymmetry );
00269 
00270      if (verboseLevel>=2) {
00271        G4cout << " MeanFreePath:  " << mfp / mm << " mm " << G4endl;
00272        G4cout << " Asymmetry:     " << asymmetry          << G4endl;
00273        G4cout << " PolProduct:    " << polProduct         << G4endl;
00274      }
00275    }
00276 
00277    return mfp;
00278 }
00279 
00280 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00281 
00282 void G4PolarizedCompton::PreparePhysicsTable(const G4ParticleDefinition& part)
00283 {
00284   G4VEmProcess::PreparePhysicsTable(part);
00285   if(buildAsymmetryTable)
00286     theAsymmetryTable = G4PhysicsTableHelper::PreparePhysicsTable(theAsymmetryTable);
00287 }
00288 
00289 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00290 
00291 
00292 void G4PolarizedCompton::BuildPhysicsTable(const G4ParticleDefinition& part)
00293 {
00294   // *** build (unpolarized) cross section tables (Lambda)
00295   G4VEmProcess::BuildPhysicsTable(part);
00296   if(buildAsymmetryTable)
00297     BuildAsymmetryTable(part);
00298 }
00299 
00300 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00301 
00302 
00303 void G4PolarizedCompton::BuildAsymmetryTable(const G4ParticleDefinition& part)
00304 {
00305   // Access to materials
00306   const G4ProductionCutsTable* theCoupleTable=
00307         G4ProductionCutsTable::GetProductionCutsTable();
00308   size_t numOfCouples = theCoupleTable->GetTableSize();
00309   for(size_t i=0; i<numOfCouples; ++i) {
00310     if (!theAsymmetryTable) break;
00311     if (theAsymmetryTable->GetFlag(i)) {
00312 
00313       // create physics vector and fill it
00314       const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(i);
00315       // use same parameters as for lambda
00316       G4PhysicsVector* aVector = LambdaPhysicsVector(couple);
00317       //      modelManager->FillLambdaVector(aVector, couple, startFromNull);
00318 
00319       for (G4int j = 0 ; j < LambdaBinning() ; ++j ) {
00320         G4double lowEdgeEnergy = aVector->GetLowEdgeEnergy(j);
00321         G4double tasm=0.;
00322         G4double asym = ComputeAsymmetry(lowEdgeEnergy, couple, part, 0., tasm);
00323         aVector->PutValue(j,asym);
00324       }
00325 
00326       G4PhysicsTableHelper::SetPhysicsVector(theAsymmetryTable, i, aVector);
00327     }
00328   }
00329 
00330 }
00331 
00332 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00333 
00334 
00335 G4double G4PolarizedCompton::ComputeAsymmetry(G4double energy,
00336                          const G4MaterialCutsCouple* couple,
00337                        const G4ParticleDefinition& aParticle,
00338                                G4double cut,
00339                                G4double & tAsymmetry)
00340 {
00341   G4double lAsymmetry = 0.0;
00342   tAsymmetry=0;
00343 
00344   //
00345   // calculate polarized cross section
00346   //
00347   G4ThreeVector thePolarization=G4ThreeVector(0.,0.,1.);
00348   emModel->SetTargetPolarization(thePolarization);
00349   emModel->SetBeamPolarization(thePolarization);
00350   G4double sigma2=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
00351 
00352   //
00353   // calculate unpolarized cross section
00354   //
00355   thePolarization=G4ThreeVector();
00356   emModel->SetTargetPolarization(thePolarization);
00357   emModel->SetBeamPolarization(thePolarization);
00358   G4double sigma0=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
00359 
00360   // determine assymmetries
00361   if (sigma0>0.) {
00362     lAsymmetry=sigma2/sigma0-1.;
00363   }
00364   return lAsymmetry;
00365 }
00366 
00367 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....

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