source: trunk/source/processes/electromagnetic/adjoint/include/G4AdjointAlongStepWeightCorrection.hh@ 1305

Last change on this file since 1305 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

File size: 5.1 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.hh,v 1.4 2009/11/20 10:31:20 ldesorgh Exp $
27// GEANT4 tag $Name: geant4-09-03 $
28//
29/////////////////////////////////////////////////////////////////////////////////
30// Class: G4AdjointAlongStepWeightCorrection
31// Author: L. Desorgher
32// Organisation: SpaceIT GmbH
33// Contract: ESA contract 21435/08/NL/AT
34// Customer: ESA/ESTEC
35/////////////////////////////////////////////////////////////////////////////////
36//
37// CHANGE HISTORY
38// --------------
39// ChangeHistory:
40// 10 May 2007 creation by L. Desorgher
41// October 2009 implementation of the mode where the total adjoint and forward cross sections are equivalent. L. Desorgher
42//
43//-------------------------------------------------------------
44// Documentation:
45// Continuous processes acting on adjoint particles to correct continuously their weight during the adjoint reverse tracking.
46// Thi process is needed whene the adjoint cross section are not scaled such that the total adjoint cross section match the total forward cross section.
47// By default the mode where the total adjoint cross section is equal to the total forward cross section is used an therefore this along step weight
48// correction factor is 1.
49// However in some cases (some energy ranges) the total forward cross section or the total adjoint cross section can be null, in this case the along step
50// weight correction is neede and is given by exp(-(Sigma_tot_adj-Sigma_tot_fwd).dx)
51//
52//
53//
54
55
56#ifndef G4AdjointAlongStepWeightCorrection_h
57#define G4AdjointAlongStepWeightCorrection_h 1
58
59#include "G4VContinuousProcess.hh"
60#include "globals.hh"
61#include "G4Material.hh"
62#include "G4MaterialCutsCouple.hh"
63#include "G4Track.hh"
64#include "G4ParticleChange.hh"
65
66class G4Step;
67class G4ParticleDefinition;
68
69
70
71class G4AdjointAlongStepWeightCorrection : public G4VContinuousProcess
72{
73public:
74
75 G4AdjointAlongStepWeightCorrection(const G4String& name = "ContinuousWeightCorrection",
76 G4ProcessType type = fElectromagnetic);
77
78 virtual ~G4AdjointAlongStepWeightCorrection();
79
80
81protected:
82 virtual G4double GetContinuousStepLimit(const G4Track& track,
83 G4double previousStepSize,
84 G4double currentMinimumStep,
85 G4double& currentSafety);
86
87
88 //------------------------------------------------------------------------
89 // Generic methods common to all processes
90 //------------------------------------------------------------------------
91public:
92
93 void PreparePhysicsTable(const G4ParticleDefinition&);
94
95 void BuildPhysicsTable(const G4ParticleDefinition&);
96
97 G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&);
98
99
100private:
101
102 void DefineMaterial(const G4MaterialCutsCouple* couple);
103
104
105
106 G4AdjointAlongStepWeightCorrection(G4AdjointAlongStepWeightCorrection &);
107 G4AdjointAlongStepWeightCorrection & operator=(const G4AdjointAlongStepWeightCorrection &right);
108
109
110
111protected:
112
113 G4ParticleChange* fParticleChange;
114
115
116private:
117
118 const G4Material* currentMaterial;
119 const G4MaterialCutsCouple* currentCouple;
120 size_t currentMaterialIndex;
121 G4double currentTcut;
122 G4double preStepKinEnergy;
123
124};
125
126inline void G4AdjointAlongStepWeightCorrection::DefineMaterial(
127 const G4MaterialCutsCouple* couple)
128{
129 if(couple != currentCouple) {
130 currentCouple = couple;
131 currentMaterial = couple->GetMaterial();
132 currentMaterialIndex = couple->GetIndex();
133
134 }
135}
136
137
138
139#endif
Note: See TracBrowser for help on using the repository browser.