G4DNATransformElectronModel.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: G4DNATransformElectronModel.cc 64057 2012-10-30 15:04:49Z gcosmo $
00027 //
00028 #include "G4DNATransformElectronModel.hh"
00029 #include "G4SystemOfUnits.hh"
00030 #include "G4ParticleChangeForGamma.hh"
00031 #include "G4Electron.hh"
00032 #include "G4NistManager.hh"
00033 #include "G4DNAChemistryManager.hh"
00034 #include "G4DNAMolecularMaterial.hh"
00035 
00036 G4DNATransformElectronModel::G4DNATransformElectronModel(const G4ParticleDefinition*,
00037                                                          const G4String& nam):
00038     G4VEmModel(nam),fIsInitialised(false)
00039 {
00040     fVerboseLevel = 0 ;
00041     SetLowEnergyLimit(0.*eV);
00042     SetHighEnergyLimit(0.025*eV);
00043     fParticleChangeForGamma = 0;
00044   //  fNistWater  = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
00045     fpWaterDensity = 0;
00046     fpWaterDensity = 0;
00047     fEpsilon = 0.0001*eV;
00048 }
00049 
00050 //______________________________________________________________________
00051 G4DNATransformElectronModel::~G4DNATransformElectronModel()
00052 {}
00053 
00054 //______________________________________________________________________
00055 void G4DNATransformElectronModel::Initialise(const G4ParticleDefinition* particleDefinition,
00056                                              const G4DataVector&)
00057 {
00058 #ifdef G4VERBOSE
00059     if (fVerboseLevel)
00060         G4cout << "Calling G4DNATransformElectronModel::Initialise()" << G4endl;
00061 #endif
00062 
00063     if (particleDefinition != G4Electron::ElectronDefinition())
00064     {
00065         G4ExceptionDescription exceptionDescription ;
00066         exceptionDescription << "Attempting to calculate cross section for wrong particle";
00067         G4Exception("G4DNATransformElectronModel::CrossSectionPerVolume","G4DNATransformElectronModel001",
00068                     FatalErrorInArgument,exceptionDescription);
00069         return;
00070     }
00071 
00072     // Initialize water density pointer
00073     fpWaterDensity = G4DNAMolecularMaterial::Instance()->GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
00074 
00075     if(!fIsInitialised)
00076     {
00077         fIsInitialised = true;
00078         fParticleChangeForGamma = GetParticleChangeForGamma();
00079     }
00080 }
00081 
00082 //______________________________________________________________________
00083 G4double G4DNATransformElectronModel::CrossSectionPerVolume(const G4Material* material,
00084                                                             const G4ParticleDefinition*,
00085                                                             G4double ekin,
00086                                                             G4double,
00087                                                             G4double)
00088 {
00089 #if G4VERBOSE
00090     if (fVerboseLevel > 1)
00091         G4cout << "Calling CrossSectionPerVolume() of G4DNATransformElectronModel" << G4endl;
00092 #endif
00093 
00094     if(ekin - fEpsilon > HighEnergyLimit())
00095     {
00096         return 0.0;
00097     }
00098 
00099     G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
00100 
00101     if(waterDensity!= 0.0)
00102    //  if (material == nistwater || material->GetBaseMaterial() == nistwater)
00103     {
00104         if (ekin - fEpsilon <= HighEnergyLimit())
00105         {
00106             return DBL_MAX;
00107         }
00108     }
00109 
00110     return 0.0 ;
00111 }
00112 
00113 //______________________________________________________________________
00114 void G4DNATransformElectronModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
00115                                                     const G4MaterialCutsCouple*,
00116                                                     const G4DynamicParticle* particle,
00117                                                     G4double,
00118                                                     G4double)
00119 {
00120 #if G4VERBOSE
00121     if (fVerboseLevel)
00122         G4cout << "Calling SampleSecondaries() of G4DNATransformElectronModel" << G4endl;
00123 #endif
00124 
00125     G4double k = particle->GetKineticEnergy();
00126 
00127 //    if (k - fEpsilon <= HighEnergyLimit())
00128 //    {
00129         const G4Track * track = fParticleChangeForGamma->GetCurrentTrack();
00130         G4DNAChemistryManager::Instance()->CreateSolvatedElectron(track);
00131         fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
00132         fParticleChangeForGamma->ProposeLocalEnergyDeposit(k);
00133 //    }
00134 }

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