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

Last change on this file since 966 was 966, checked in by garnier, 15 years ago

fichier ajoutes

File size: 4.7 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/////////////////////////////////////////////////////////////////////////////////
27//      Module:         G4AdjointAlongStepWeightCorrection.hh
28//      Author:         L. Desorgher
29//      Date:           10 May 2007
30//      Organisation:   SpaceIT GmbH
31//      Customer:       ESA/ESTEC
32/////////////////////////////////////////////////////////////////////////////////
33//
34// CHANGE HISTORY
35// --------------
36//      ChangeHistory:
37//              10 May 2007 creation by L. Desorgher           
38//
39//-------------------------------------------------------------
40//      Documentation:
41//              Continuous processes acting on adjoint particles to correct continuously their weight during the adjoint reverse tracking.
42//
43
44
45#ifndef G4AdjointAlongStepWeightCorrection_h
46#define G4AdjointAlongStepWeightCorrection_h 1
47
48#include "G4VContinuousProcess.hh"
49#include "globals.hh"
50#include "G4Material.hh"
51#include "G4MaterialCutsCouple.hh"
52#include "G4Track.hh"
53#include "G4ParticleChange.hh"
54
55class G4Step;
56class G4ParticleDefinition;
57
58
59
60class G4AdjointAlongStepWeightCorrection : public G4VContinuousProcess
61{
62public:
63
64  G4AdjointAlongStepWeightCorrection(const G4String& name = "ContinuousWeightCorrection",
65                         G4ProcessType type = fElectromagnetic);
66
67  virtual ~G4AdjointAlongStepWeightCorrection();
68
69
70protected:
71  virtual G4double GetContinuousStepLimit(const G4Track& track,
72                                                G4double previousStepSize,
73                                                G4double currentMinimumStep,
74                                                G4double& currentSafety);
75                                       
76
77  //------------------------------------------------------------------------
78  // Generic methods common to all processes
79  //------------------------------------------------------------------------
80public:
81
82  void PreparePhysicsTable(const G4ParticleDefinition&);
83
84  void BuildPhysicsTable(const G4ParticleDefinition&);
85
86  G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&);
87
88
89private:
90
91  void DefineMaterial(const G4MaterialCutsCouple* couple);
92 
93
94
95  G4AdjointAlongStepWeightCorrection(G4AdjointAlongStepWeightCorrection &);
96  G4AdjointAlongStepWeightCorrection & operator=(const G4AdjointAlongStepWeightCorrection &right);
97
98
99
100protected:
101
102  G4ParticleChange* fParticleChange;
103 
104 
105private:
106 
107  const G4Material*  currentMaterial;
108  const G4MaterialCutsCouple* currentCouple;
109  size_t   currentMaterialIndex; 
110  G4double currentTcut;
111  G4double preStepKinEnergy;
112 
113};
114
115inline void G4AdjointAlongStepWeightCorrection::DefineMaterial(
116            const G4MaterialCutsCouple* couple)
117{
118  if(couple != currentCouple) {
119    currentCouple   = couple;
120    currentMaterial = couple->GetMaterial();
121    currentMaterialIndex = couple->GetIndex();
122    //G4cout<<"Define Material"<<std::endl;
123    //if(!meanFreePath) ResetNumberOfInteractionLengthLeft();
124  }
125}
126
127
128///////////////////////////////////////////////////////
129//
130inline G4double G4AdjointAlongStepWeightCorrection::GetContinuousStepLimit(const G4Track& track,
131                G4double , G4double , G4double& )
132{ 
133  G4double x = DBL_MAX;
134  DefineMaterial(track.GetMaterialCutsCouple());
135  preStepKinEnergy = track.GetKineticEnergy();
136  return x;
137}
138#endif
Note: See TracBrowser for help on using the repository browser.