Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4AdjointPhotoElectricModel.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//
26//
27/////////////////////////////////////////////////////////////////////////////////
28// Module: G4AdjointPhotoElectricModel
29// Author: L. Desorgher
30// Organisation: SpaceIT GmbH
31// Contract: ESA contract 21435/08/NL/AT
32// Customer: ESA/ESTEC
33/////////////////////////////////////////////////////////////////////////////////
34//
35// CHANGE HISTORY
36// --------------
37// ChangeHistory:
38// -1 September 2007 creation by L. Desorgher
39//
40// -January 2009. L. Desorgher
41// Put a higher limit on the CS to avoid a high rate of Inverse Photo e- effect at low energy. The very high adjoint CS of the reverse
42// photo electric reaction produce a high rate of reverse photo electric reaction in the inner side of a shielding for eaxmple, the correction of this occurrence
43// by weight correction in the StepDoIt method is not statistically sufficient at small energy. The problem is partially solved by setting an higher CS limit
44// and compensating it by an extra weight correction factor. However when coupling it with other reverse processes the reverse photo-electric is still
45// the source of very occasional high weight that decrease the efficiency of the computation. A way to solve this problemn is still needed but is difficult
46// to find as it happens in rarea case but does give a weighrt that is outside the noemal distribution. (Very Tricky!)
47//
48// -October 2009 Correction of Element sampling. L. Desorgher
49//
50//-------------------------------------------------------------
51// Documentation:
52// Model for the adjoint photo electric process
53//
54#ifndef G4AdjointPhotoElectricModel_h
55#define G4AdjointPhotoElectricModel_h 1
56
57
58#include "globals.hh"
59#include "G4VEmAdjointModel.hh"
62
63{
64public:
65
68
69
70
71 virtual void SampleSecondaries(const G4Track& aTrack,
72 G4bool IsScatProjToProjCase,
73 G4ParticleChange* fParticleChange);
75 G4double primEnergy,
76 G4bool IsScatProjToProjCase);
78 G4double primEnergy,
79 G4bool IsScatProjToProjCase);
80
81 G4double AdjointCrossSectionPerAtom(const G4Element* anElement,G4double electronEnergy);
82
83
84
85 inline void SetTheDirectPEEffectModel(G4PEEffectFluoModel* aModel){theDirectPEEffectModel = aModel;
86 DefineDirectEMModel(aModel);}
87
88 virtual void CorrectPostStepWeight(G4ParticleChange* fParticleChange,
89 G4double old_weight,
90 G4double adjointPrimKinEnergy,
91 G4double projectileKinEnergy,
92 G4bool IsScatProjToProjCase);
93
94
95private:
96 G4double xsec[40];
97 G4double totAdjointCS;
98 G4double totBiasedAdjointCS;
99 G4double factorCSBiasing;
100 G4double pre_step_AdjointCS;
101 G4double post_step_AdjointCS;
102
103
104 G4double shell_prob[40][40];
105
106
107 G4PEEffectFluoModel* theDirectPEEffectModel;
108 size_t index_element;
109 G4double current_eEnergy;
110
111
112private:
113 void DefineCurrentMaterialAndElectronEnergy(const G4MaterialCutsCouple* aCouple,
114 G4double eEnergy);
115
116};
117
118#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
virtual void SampleSecondaries(const G4Track &aTrack, G4bool IsScatProjToProjCase, G4ParticleChange *fParticleChange)
virtual void CorrectPostStepWeight(G4ParticleChange *fParticleChange, G4double old_weight, G4double adjointPrimKinEnergy, G4double projectileKinEnergy, G4bool IsScatProjToProjCase)
virtual G4double AdjointCrossSection(const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)
void SetTheDirectPEEffectModel(G4PEEffectFluoModel *aModel)
G4double AdjointCrossSectionPerAtom(const G4Element *anElement, G4double electronEnergy)
virtual G4double GetAdjointCrossSection(const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)
void DefineDirectEMModel(G4VEmModel *aModel)