Geant4-11
G4AdjointAlongStepWeightCorrection.hh
Go to the documentation of this file.
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//
27// Class: G4AdjointAlongStepWeightCorrection
28// Author: L. Desorgher
29// Organisation: SpaceIT GmbH
31//
32// Documentation:
33//
34// Continuous processes act on adjoint particles to continuously correct their
35// weight during the adjoint reverse tracking. This process is needed when
36// the adjoint cross sections are not scaled such that the total adjoint cross
37// section matches the total forward cross section. By default the mode where
38// the total adjoint cross section is equal to the total forward cross section
39// is used and therefore this along step weightcorrection factor is 1. However
40// in some cases (some energy ranges) the total forward cross section or the
41// total adjoint cross section can be zero. In this case the along step weight
42// correction is needed and is given by exp(-(Sigma_tot_adj-Sigma_tot_fwd).dx)
43//
44//-------------------------------------------------------------
45
46#ifndef G4AdjointAlongStepWeightCorrection_h
47#define G4AdjointAlongStepWeightCorrection_h 1
48
49#include "globals.hh"
51
56class G4Step;
57class G4Track;
58
60{
61 public:
63 const G4String& name = "ContinuousWeightCorrection",
65
67
68 G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&) override;
69
70 void ProcessDescription(std::ostream&) const override;
71 void DumpInfo() const override { ProcessDescription(G4cout); };
72
74 delete;
76 const G4AdjointAlongStepWeightCorrection& right) = delete;
77
78 protected:
80 G4double previousStepSize,
81 G4double currentMinimumStep,
82 G4double& currentSafety) override;
83
84 private:
85 void DefineMaterial(const G4MaterialCutsCouple* couple);
86
90
92};
93
95 const G4MaterialCutsCouple* couple)
96{
97 if(couple != fCurrentCouple)
98 {
99 fCurrentCouple = couple;
100 }
101}
102
103#endif
G4ProcessType
@ fElectromagnetic
double G4double
Definition: G4Types.hh:83
G4GLOB_DLL std::ostream G4cout
G4AdjointAlongStepWeightCorrection(G4AdjointAlongStepWeightCorrection &)=delete
G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety) override
void ProcessDescription(std::ostream &) const override
void DefineMaterial(const G4MaterialCutsCouple *couple)
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &) override
G4AdjointAlongStepWeightCorrection(const G4String &name="ContinuousWeightCorrection", G4ProcessType type=fElectromagnetic)
G4AdjointAlongStepWeightCorrection & operator=(const G4AdjointAlongStepWeightCorrection &right)=delete
Definition: G4Step.hh:62
const char * name(G4int ptype)