Geant4-11
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//
27// Class: G4AdjointPhotoElectricModel
28// Author: L. Desorgher
29// Organisation: SpaceIT GmbH
30//
31// Model for the adjoint photo electric process.
32// Put a higher limit on the CS to avoid a high rate of Inverse Photo e-effect
33// at low energy. The very high adjoint CS of the reverse photo electric
34// reaction produce a high rate of reverse photo electric reaction in the inner
35// side of a shielding for eaxmple, the correction of this occurrence by weight
36// correction in the StepDoIt method is not statistically sufficient at small
37// energy. The problem is partially solved by setting a higher CS limit and
38// compensating it by an extra weight correction factor. However when coupling
39// it with other reverse processes the reverse photo-electric is still the
40// source of very occasional high weights that decrease the efficiency of the
41// computation. A way to solve this problemn is still needed but is difficult
42// to find as it happens in rare cases but does give a weight that is outside
43// the normal distribution. (Very Tricky!)
44//
46
47#ifndef G4AdjointPhotoElectricModel_h
48#define G4AdjointPhotoElectricModel_h 1
49
50#include "globals.hh"
51#include "G4VEmAdjointModel.hh"
52
54{
55 public:
58
59 void SampleSecondaries(const G4Track& aTrack, G4bool isScatProjToProj,
60 G4ParticleChange* fParticleChange) override;
61
63 G4double primEnergy,
64 G4bool isScatProjToProj) override;
65
67 G4double electronEnergy);
68
71 const G4AdjointPhotoElectricModel& right) = delete;
72
73 protected:
74 void CorrectPostStepWeight(G4ParticleChange* fParticleChange,
75 G4double old_weight, G4double adjointPrimKinEnergy,
76 G4double projectileKinEnergy,
77 G4bool isScatProjToProj) override;
78
79 private:
81 const G4MaterialCutsCouple* aCouple, G4double eEnergy);
82
90
91 size_t fIndexElement = 0;
92};
93
94#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
G4AdjointPhotoElectricModel(G4AdjointPhotoElectricModel &)=delete
void CorrectPostStepWeight(G4ParticleChange *fParticleChange, G4double old_weight, G4double adjointPrimKinEnergy, G4double projectileKinEnergy, G4bool isScatProjToProj) override
void DefineCurrentMaterialAndElectronEnergy(const G4MaterialCutsCouple *aCouple, G4double eEnergy)
G4double AdjointCrossSection(const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool isScatProjToProj) override
void SampleSecondaries(const G4Track &aTrack, G4bool isScatProjToProj, G4ParticleChange *fParticleChange) override
G4AdjointPhotoElectricModel & operator=(const G4AdjointPhotoElectricModel &right)=delete
G4double AdjointCrossSectionPerAtom(const G4Element *anElement, G4double electronEnergy)