source: trunk/source/processes/electromagnetic/adjoint/src/G4AdjointAlongStepWeightCorrection.cc@ 1211

Last change on this file since 1211 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 4.6 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// $Id: G4AdjointAlongStepWeightCorrection.cc,v 1.4 2009/11/20 10:31:20 ldesorgh Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
28//
29#include "G4AdjointAlongStepWeightCorrection.hh"
30#include "G4Step.hh"
31#include "G4ParticleDefinition.hh"
32#include "G4VParticleChange.hh"
33#include "G4AdjointCSManager.hh"
34
35#ifdef WIN32
36#include <G4float .h>
37#endif
38
39///////////////////////////////////////////////////////
40//
41
42G4AdjointAlongStepWeightCorrection::G4AdjointAlongStepWeightCorrection(const G4String& name,
43 G4ProcessType type): G4VContinuousProcess(name, type)
44{fParticleChange = new G4ParticleChange();
45}
46
47///////////////////////////////////////////////////////
48//
49
50G4AdjointAlongStepWeightCorrection::~G4AdjointAlongStepWeightCorrection()
51{;
52}
53///////////////////////////////////////////////////////
54//
55void G4AdjointAlongStepWeightCorrection::PreparePhysicsTable(
56 const G4ParticleDefinition& )
57{
58;
59}
60///////////////////////////////////////////////////////
61//
62
63void G4AdjointAlongStepWeightCorrection::BuildPhysicsTable(const G4ParticleDefinition& )
64{;
65}
66///////////////////////////////////////////////////////
67//
68G4VParticleChange* G4AdjointAlongStepWeightCorrection::AlongStepDoIt(const G4Track& track,
69 const G4Step& step)
70{
71
72 fParticleChange->Initialize(track);
73
74 // Get the actual (true) Step length
75 //----------------------------------
76 G4double length = step.GetStepLength();
77
78
79 G4double Tkin = step.GetPostStepPoint()->GetKineticEnergy();
80 G4ParticleDefinition* thePartDef= const_cast<G4ParticleDefinition*> (track.GetDynamicParticle()->GetDefinition());
81 G4double weight_correction=G4AdjointCSManager::GetAdjointCSManager()->GetContinuousWeightCorrection(thePartDef,
82 preStepKinEnergy,Tkin, currentCouple,length);
83
84
85
86 G4double new_weight=weight_correction*track.GetWeight();
87
88 //if (weight_correction >2.) new_weight=1.e-300;
89
90
91 //The following test check for zero weight.
92 //This happens after weight correction of gamma for photo electric effect.
93 //When the new weight is 0 it will be later on consider as nan by G4.
94 //Therefore we do put a lower limit of 1.e-300. for new_weight
95 //Correction by L.Desorgher on 15 July 2009
96#ifdef WIN32
97 if (!!_isnan(new_weight) || new_weight==0){
98#else
99 if (std::isnan(new_weight) || new_weight==0){
100#endif
101 //G4cout<<new_weight<<'\t'<<weight_correction<<'\t'<<track.GetWeight()<<G4endl;
102 new_weight=1.e-300;
103 }
104
105 //G4cout<<new_weight<<'\t'<<weight_correction<<'\t'<<track.GetWeight()<<G4endl;
106 fParticleChange->SetParentWeightByProcess(false);
107 fParticleChange->SetSecondaryWeightByProcess(false);
108 fParticleChange->ProposeParentWeight(new_weight);
109
110
111 return fParticleChange;
112
113}
114///////////////////////////////////////////////////////
115//
116G4double G4AdjointAlongStepWeightCorrection::GetContinuousStepLimit(const G4Track& track,
117 G4double , G4double , G4double& )
118{
119 G4double x = DBL_MAX;
120 DefineMaterial(track.GetMaterialCutsCouple());
121 preStepKinEnergy = track.GetKineticEnergy();
122 return x;
123}
Note: See TracBrowser for help on using the repository browser.