G4ePolarizedIonisation.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: G4ePolarizedIonisation.cc 69847 2013-05-16 09:36:18Z gcosmo $
00028 // -------------------------------------------------------------------
00029 //
00030 // GEANT4 Class file
00031 //
00032 //
00033 // File name:     G4ePolarizedIonisation
00034 //
00035 // Author:        A.Schaelicke on base of Vladimir Ivanchenko code
00036 //
00037 // Creation date: 10.11.2005
00038 //
00039 // Modifications:
00040 //
00041 // 10-11-05, include polarization description (A.Schaelicke)
00042 // , create asymmetry table and determine interactionlength 
00043 // , update polarized differential cross section 
00044 //
00045 // 20-08-06, modified interface (A.Schaelicke)
00046 // 11-06-07, add PostStepGetPhysicalInteractionLength (A.Schalicke)
00047 //
00048 // Class Description:
00049 //
00050 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00051 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00052 
00053 #include "G4ePolarizedIonisation.hh"
00054 #include "G4Electron.hh"
00055 #include "G4UniversalFluctuation.hh"
00056 #include "G4BohrFluctuations.hh"
00057 #include "G4UnitsTable.hh"
00058 
00059 #include "G4PolarizedMollerBhabhaModel.hh"
00060 #include "G4ProductionCutsTable.hh"
00061 #include "G4PolarizationManager.hh"
00062 #include "G4PolarizationHelper.hh"
00063 #include "G4StokesVector.hh"
00064 
00065 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00066 
00067 G4ePolarizedIonisation::G4ePolarizedIonisation(const G4String& name)
00068   : G4VEnergyLossProcess(name),
00069     theElectron(G4Electron::Electron()),
00070     isElectron(true),
00071     isInitialised(false),
00072     theAsymmetryTable(NULL),
00073     theTransverseAsymmetryTable(NULL)
00074 {
00075   verboseLevel=0;
00076   //  SetDEDXBinning(120);
00077   //  SetLambdaBinning(120);
00078   //  numBinAsymmetryTable=78;
00079 
00080   //  SetMinKinEnergy(0.1*keV);
00081   //  SetMaxKinEnergy(100.0*TeV);
00082   //  PrintInfoDefinition();
00083   SetProcessSubType(fIonisation);
00084   flucModel = 0;
00085   emModel = 0; 
00086 }
00087 
00088 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00089 
00090 G4ePolarizedIonisation::~G4ePolarizedIonisation()
00091 {
00092   if (theAsymmetryTable) {
00093     theAsymmetryTable->clearAndDestroy();
00094     delete theAsymmetryTable;
00095   }
00096   if (theTransverseAsymmetryTable) {
00097     theTransverseAsymmetryTable->clearAndDestroy();
00098     delete theTransverseAsymmetryTable;
00099   }
00100 }
00101 
00102 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00103 
00104 void G4ePolarizedIonisation::InitialiseEnergyLossProcess(
00105                     const G4ParticleDefinition* part,
00106                     const G4ParticleDefinition* /*part2*/)
00107 {
00108   if(!isInitialised) {
00109 
00110     if(part == G4Positron::Positron()) isElectron = false;
00111     SetSecondaryParticle(theElectron);
00112 
00113 
00114 
00115     flucModel = new G4UniversalFluctuation();
00116     //flucModel = new G4BohrFluctuations();
00117 
00118     //    G4VEmModel* em = new G4MollerBhabhaModel();
00119     emModel = new G4PolarizedMollerBhabhaModel;
00120     emModel->SetLowEnergyLimit(MinKinEnergy());
00121     emModel->SetHighEnergyLimit(MaxKinEnergy());
00122     AddEmModel(1, emModel, flucModel);
00123 
00124     isInitialised = true;
00125   }
00126 }
00127 
00128 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00129 
00130 void G4ePolarizedIonisation::PrintInfo()
00131 {
00132   G4cout << "      Delta cross sections from Moller+Bhabha, "
00133          << "good description from 1 KeV to 100 GeV."
00134          << G4endl;
00135 }
00136 
00137 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00138 
00139 G4double G4ePolarizedIonisation::GetMeanFreePath(const G4Track& track,
00140                                                  G4double step,
00141                                                  G4ForceCondition* cond)
00142 {
00143   // *** get unploarised mean free path from lambda table ***
00144   G4double mfp = G4VEnergyLossProcess::GetMeanFreePath(track, step, cond);
00145 
00146 
00147   // *** get asymmetry, if target is polarized ***
00148   G4VPhysicalVolume*  aPVolume  = track.GetVolume();
00149   G4LogicalVolume*    aLVolume  = aPVolume->GetLogicalVolume();
00150 
00151   G4PolarizationManager * polarizationManger = G4PolarizationManager::GetInstance();
00152   G4bool volumeIsPolarized = polarizationManger->IsPolarized(aLVolume);
00153   const G4StokesVector ePolarization = track.GetPolarization();
00154 
00155   if (mfp != DBL_MAX &&  volumeIsPolarized && !ePolarization.IsZero()) {
00156     const G4DynamicParticle* aDynamicElectron = track.GetDynamicParticle();
00157     G4double eEnergy = aDynamicElectron->GetKineticEnergy();
00158     const G4ParticleMomentum eDirection0 = aDynamicElectron->GetMomentumDirection();
00159 
00160     G4StokesVector volumePolarization = polarizationManger->GetVolumePolarization(aLVolume);
00161 
00162     G4bool isOutRange;
00163     size_t idx = CurrentMaterialCutsCoupleIndex();
00164     G4double lAsymmetry = (*theAsymmetryTable)(idx)->
00165                                   GetValue(eEnergy, isOutRange);
00166     G4double tAsymmetry = (*theTransverseAsymmetryTable)(idx)->
00167                                   GetValue(eEnergy, isOutRange);
00168 
00169     // calculate longitudinal spin component
00170     G4double polZZ = ePolarization.z()*
00171                         volumePolarization*eDirection0;
00172     // calculate transvers spin components
00173     G4double polXX = ePolarization.x()*
00174                         volumePolarization*G4PolarizationHelper::GetParticleFrameX(eDirection0);
00175     G4double polYY = ePolarization.y()*
00176                         volumePolarization*G4PolarizationHelper::GetParticleFrameY(eDirection0);
00177 
00178 
00179     G4double impact = 1. + polZZ*lAsymmetry + (polXX + polYY)*tAsymmetry;
00180     // determine polarization dependent mean free path
00181     mfp /= impact;
00182     if (mfp <=0.) {
00183      G4cout <<"PV impact ( "<<polXX<<" , "<<polYY<<" , "<<polZZ<<" )"<<G4endl;
00184      G4cout << " impact on MFP is "<< impact <<G4endl;
00185      G4cout<<" lAsymmetry= "<<lAsymmetry<<" ("<<std::fabs(lAsymmetry)-1.<<")\n";
00186      G4cout<<" tAsymmetry= "<<tAsymmetry<<" ("<<std::fabs(tAsymmetry)-1.<<")\n";
00187     }
00188   }
00189 
00190   return mfp;
00191 }
00192 
00193 G4double G4ePolarizedIonisation::PostStepGetPhysicalInteractionLength(const G4Track& track,
00194                                               G4double step,
00195                                               G4ForceCondition* cond)
00196 {
00197   // *** get unploarised mean free path from lambda table ***
00198   G4double mfp = G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength(track, step, cond);
00199 
00200 
00201   // *** get asymmetry, if target is polarized ***
00202   G4VPhysicalVolume*  aPVolume  = track.GetVolume();
00203   G4LogicalVolume*    aLVolume  = aPVolume->GetLogicalVolume();
00204 
00205   G4PolarizationManager * polarizationManger = G4PolarizationManager::GetInstance();
00206   G4bool volumeIsPolarized = polarizationManger->IsPolarized(aLVolume);
00207   const G4StokesVector ePolarization = track.GetPolarization();
00208 
00209   if (mfp != DBL_MAX &&  volumeIsPolarized && !ePolarization.IsZero()) {
00210     const G4DynamicParticle* aDynamicElectron = track.GetDynamicParticle();
00211     G4double eEnergy = aDynamicElectron->GetKineticEnergy();
00212     const G4ParticleMomentum eDirection0 = aDynamicElectron->GetMomentumDirection();
00213 
00214     G4StokesVector volumePolarization = polarizationManger->GetVolumePolarization(aLVolume);
00215 
00216     size_t idx = CurrentMaterialCutsCoupleIndex();
00217     G4double lAsymmetry = (*theAsymmetryTable)(idx)->Value(eEnergy);
00218     G4double tAsymmetry = (*theTransverseAsymmetryTable)(idx)->Value(eEnergy);
00219 
00220     // calculate longitudinal spin component
00221     G4double polZZ = ePolarization.z()*
00222                         volumePolarization*eDirection0;
00223     // calculate transvers spin components
00224     G4double polXX = ePolarization.x()*
00225                         volumePolarization*G4PolarizationHelper::GetParticleFrameX(eDirection0);
00226     G4double polYY = ePolarization.y()*
00227                         volumePolarization*G4PolarizationHelper::GetParticleFrameY(eDirection0);
00228 
00229 
00230     G4double impact = 1. + polZZ*lAsymmetry + (polXX + polYY)*tAsymmetry;
00231     // determine polarization dependent mean free path
00232     mfp /= impact;
00233     if (mfp <=0.) {
00234      G4cout <<"PV impact ( "<<polXX<<" , "<<polYY<<" , "<<polZZ<<" )"<<G4endl;
00235      G4cout << " impact on MFP is "<< impact <<G4endl;
00236      G4cout<<" lAsymmetry= "<<lAsymmetry<<" ("<<std::fabs(lAsymmetry)-1.<<")\n";
00237      G4cout<<" tAsymmetry= "<<tAsymmetry<<" ("<<std::fabs(tAsymmetry)-1.<<")\n";
00238     }
00239   }
00240 
00241   return mfp;
00242 }
00243 
00244 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00245 void G4ePolarizedIonisation::BuildPhysicsTable(const G4ParticleDefinition& part)
00246 {
00247   // *** build DEDX and (unpolarized) cross section tables
00248   G4VEnergyLossProcess::BuildPhysicsTable(part);
00249   //  G4PhysicsTable* pt =
00250   //  BuildDEDXTable();
00251 
00252 
00253   // *** build asymmetry-table
00254   if (theAsymmetryTable) {
00255     theAsymmetryTable->clearAndDestroy(); delete theAsymmetryTable;}
00256   if (theTransverseAsymmetryTable) {
00257     theTransverseAsymmetryTable->clearAndDestroy(); delete theTransverseAsymmetryTable;}
00258 
00259   const G4ProductionCutsTable* theCoupleTable=
00260         G4ProductionCutsTable::GetProductionCutsTable();
00261   size_t numOfCouples = theCoupleTable->GetTableSize();
00262 
00263   theAsymmetryTable = new G4PhysicsTable(numOfCouples);
00264   theTransverseAsymmetryTable = new G4PhysicsTable(numOfCouples);
00265 
00266   for (size_t j=0 ; j < numOfCouples; j++ ) {
00267     // get cut value
00268     const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(j);
00269 
00270     G4double cut = (*theCoupleTable->GetEnergyCutsVector(1))[j];
00271 
00272     //create physics vectors then fill it (same parameters as lambda vector)
00273     G4PhysicsVector * ptrVectorA = LambdaPhysicsVector(couple,cut);
00274     G4PhysicsVector * ptrVectorB = LambdaPhysicsVector(couple,cut);
00275     size_t bins = ptrVectorA->GetVectorLength();
00276 
00277     for (size_t i = 0 ; i < bins ; i++ ) {
00278       G4double lowEdgeEnergy = ptrVectorA->Energy(i);
00279       G4double tasm=0.;
00280       G4double asym = ComputeAsymmetry(lowEdgeEnergy, couple, part, cut, tasm);
00281       ptrVectorA->PutValue(i,asym);
00282       ptrVectorB->PutValue(i,tasm);
00283     }
00284     theAsymmetryTable->insertAt( j , ptrVectorA ) ;
00285     theTransverseAsymmetryTable->insertAt( j , ptrVectorB ) ;
00286   }
00287 
00288 }
00289 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00290 
00291 G4double G4ePolarizedIonisation::ComputeAsymmetry(G4double energy,
00292                                          const G4MaterialCutsCouple* couple,
00293                                                const G4ParticleDefinition& aParticle,
00294                                                G4double cut,
00295                                                G4double & tAsymmetry)
00296 {
00297   G4double lAsymmetry = 0.0;
00298            tAsymmetry = 0.0;
00299   if (isElectron) {lAsymmetry = tAsymmetry = -1.0;}
00300 
00301   // calculate polarized cross section
00302   theTargetPolarization=G4ThreeVector(0.,0.,1.);
00303   emModel->SetTargetPolarization(theTargetPolarization);
00304   emModel->SetBeamPolarization(theTargetPolarization);
00305   G4double sigma2=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
00306 
00307   // calculate transversely polarized cross section
00308   theTargetPolarization=G4ThreeVector(1.,0.,0.);
00309   emModel->SetTargetPolarization(theTargetPolarization);
00310   emModel->SetBeamPolarization(theTargetPolarization);
00311   G4double sigma3=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
00312 
00313   // calculate unpolarized cross section
00314   theTargetPolarization=G4ThreeVector();
00315   emModel->SetTargetPolarization(theTargetPolarization);
00316   emModel->SetBeamPolarization(theTargetPolarization);
00317   G4double sigma0=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
00318   // determine assymmetries
00319   if (sigma0>0.) {
00320     lAsymmetry=sigma2/sigma0-1.;
00321     tAsymmetry=sigma3/sigma0-1.;
00322   }
00323   if (std::fabs(lAsymmetry)>1.) {
00324     G4cout<<" energy="<<energy<<"\n";
00325     G4cout<<"WARNING lAsymmetry= "<<lAsymmetry<<" ("<<std::fabs(lAsymmetry)-1.<<")\n";
00326   }
00327   if (std::fabs(tAsymmetry)>1.) {
00328     G4cout<<" energy="<<energy<<"\n";
00329     G4cout<<"WARNING tAsymmetry= "<<tAsymmetry<<" ("<<std::fabs(tAsymmetry)-1.<<")\n";
00330   }
00331 //   else {
00332 //     G4cout<<"        tAsymmetry= "<<tAsymmetry<<" ("<<std::fabs(tAsymmetry)-1.<<")\n";
00333 //   }
00334   return lAsymmetry;
00335 }
00336 
00337 

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