G4AdjointAlongStepWeightCorrection.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: G4AdjointAlongStepWeightCorrection.cc 69844 2013-05-16 09:19:33Z gcosmo $
00027 //
00028 #include "G4AdjointAlongStepWeightCorrection.hh"
00029 #include "G4Step.hh"
00030 #include "G4ParticleDefinition.hh"
00031 #include "G4VParticleChange.hh"
00032 #include "G4AdjointCSManager.hh"
00033 
00035 //
00036 
00037 G4AdjointAlongStepWeightCorrection::G4AdjointAlongStepWeightCorrection(const G4String& name, 
00038   G4ProcessType type): G4VContinuousProcess(name, type)
00039 {fParticleChange = new G4ParticleChange();
00040  currentMaterialIndex=0;
00041  preStepKinEnergy=1.;
00042  currentCouple=0;
00043 }
00044 
00046 //
00047 
00048 G4AdjointAlongStepWeightCorrection::~G4AdjointAlongStepWeightCorrection()
00049 {delete fParticleChange;
00050 }
00052 //
00053 void G4AdjointAlongStepWeightCorrection::PreparePhysicsTable(
00054      const G4ParticleDefinition& )
00055 {
00056 ; 
00057 }
00059 //
00060 
00061 void G4AdjointAlongStepWeightCorrection::BuildPhysicsTable(const G4ParticleDefinition& )
00062 {;
00063 }
00065 //
00066 G4VParticleChange* G4AdjointAlongStepWeightCorrection::AlongStepDoIt(const G4Track& track,
00067                                                        const G4Step& step)
00068 {
00069    
00070   fParticleChange->Initialize(track);
00071   
00072   // Get the actual (true) Step length
00073   //----------------------------------
00074   G4double length = step.GetStepLength();
00075  
00076 
00077   G4double Tkin = step.GetPostStepPoint()->GetKineticEnergy();
00078   G4ParticleDefinition* thePartDef= const_cast<G4ParticleDefinition*>  (track.GetDynamicParticle()->GetDefinition());
00079   G4double weight_correction=G4AdjointCSManager::GetAdjointCSManager()->GetContinuousWeightCorrection(thePartDef,
00080                                                                         preStepKinEnergy,Tkin, currentCouple,length);
00081         
00082   
00083   
00084 
00085   //Caution!!!
00086   // It is important  to select the weight of the post_step_point
00087   // as the current weight and not the weight of the track, as t
00088   // the  weight of the track is changed after having applied all
00089   // the along_step_do_it.
00090 
00091   // G4double new_weight=weight_correction*track.GetWeight(); //old
00092   G4double new_weight=weight_correction*step.GetPostStepPoint()->GetWeight();
00093 
00094 
00095 
00096   //if (weight_correction >2.) new_weight=1.e-300;
00097   
00098   
00099   //The following test check for zero weight.
00100   //This happens after weight correction of gamma for photo electric effect.
00101   //When the new weight is 0 it will be later on consider as nan by G4.
00102   //Therefore we do put a lower limit of 1.e-300. for new_weight 
00103   //Correction by L.Desorgher on 15 July 2009 
00104   if (new_weight==0 || (new_weight<=0 && new_weight>0)){
00105                 //G4cout<<new_weight<<'\t'<<weight_correction<<'\t'<<track.GetWeight()<<G4endl;
00106                 new_weight=1.e-300;
00107   }
00108   
00109   //G4cout<<new_weight<<'\t'<<weight_correction<<'\t'<<track.GetWeight()<<G4endl;
00110   fParticleChange->SetParentWeightByProcess(false);
00111   fParticleChange->SetSecondaryWeightByProcess(false);
00112   fParticleChange->ProposeParentWeight(new_weight);
00113   
00114 
00115   return fParticleChange;
00116 
00117 }
00119 //
00120 G4double G4AdjointAlongStepWeightCorrection::GetContinuousStepLimit(const G4Track& track,
00121                 G4double , G4double , G4double& )
00122 { 
00123   G4double x = DBL_MAX;
00124   DefineMaterial(track.GetMaterialCutsCouple());
00125   preStepKinEnergy = track.GetKineticEnergy();
00126   return x;
00127 }

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