source: trunk/source/processes/electromagnetic/adjoint/include/G4AdjointPhotoElectricModel.hh@ 1206

Last change on this file since 1206 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 4.8 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: G4AdjointPhotoElectricModel.hh,v 1.4 2009/11/20 10:31:20 ldesorgh Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
28//
29/////////////////////////////////////////////////////////////////////////////////
30// Module: G4AdjointPhotoElectricModel
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// -1 September 2007 creation by L. Desorgher
41//
42// -January 2009. L. Desorgher
43// 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
44// 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 occurence
45// 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
46// and compensating it by an extra weight correction factor. However when coupling it with other reverse processes the reverse photo-electric is still
47// 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
48// to find as it happens in rarea case but does give a weighrt that is outside the noemal distribution. (Very Tricky!)
49//
50// -October 2009 Correction of Element sampling. L. Desorgher
51//
52//-------------------------------------------------------------
53// Documentation:
54// Model for the adjoint photo electric process
55//
56#ifndef G4AdjointPhotoElectricModel_h
57#define G4AdjointPhotoElectricModel_h 1
58
59
60#include "globals.hh"
61#include "G4VEmAdjointModel.hh"
62#include "G4PEEffectModel.hh"
63class G4AdjointPhotoElectricModel: public G4VEmAdjointModel
64
65{
66public:
67
68 G4AdjointPhotoElectricModel();
69 ~G4AdjointPhotoElectricModel();
70
71
72
73 virtual void SampleSecondaries(const G4Track& aTrack,
74 G4bool IsScatProjToProjCase,
75 G4ParticleChange* fParticleChange);
76 virtual G4double AdjointCrossSection(const G4MaterialCutsCouple* aCouple,
77 G4double primEnergy,
78 G4bool IsScatProjToProjCase);
79
80 G4double AdjointCrossSectionPerAtom(const G4Element* anElement,G4double electronEnergy);
81
82
83
84 inline void SetTheDirectPEEffectModel(G4PEEffectModel* aModel){theDirectPEEffectModel = aModel;
85 DefineDirectEMModel(aModel);}
86
87 virtual void CorrectPostStepWeight(G4ParticleChange* fParticleChange,
88 G4double old_weight,
89 G4double adjointPrimKinEnergy,
90 G4double projectileKinEnergy,
91 G4bool IsScatProjToProjCase);
92
93
94private:
95 G4double xsec[40];
96 G4double totAdjointCS;
97 G4double totBiasedAdjointCS;
98 G4double factorCSBiasing;
99 G4double pre_step_AdjointCS;
100 G4double post_step_AdjointCS;
101
102
103 G4double shell_prob[40][40];
104
105
106 G4PEEffectModel* theDirectPEEffectModel;
107 size_t index_element;
108 G4double current_eEnergy;
109
110
111private:
112 void DefineCurrentMaterialAndElectronEnergy(const G4MaterialCutsCouple* aCouple,
113 G4double eEnergy);
114
115};
116
117#endif
Note: See TracBrowser for help on using the repository browser.