source: trunk/source/processes/electromagnetic/lowenergy/include/G4LivermorePolarizedGammaConversionModel.hh @ 1315

Last change on this file since 1315 was 1315, checked in by garnier, 14 years ago

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

File size: 5.0 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: G4LivermorePolarizedGammaConversionModel.hh,v 1.1 2009/10/30 14:52:05 flongo Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-cand-01 $
28//
29
30#ifndef G4LivermorePolarizedGammaConversionModel_h
31#define G4LivermorePolarizedGammaConversionModel_h 1
32
33#include "G4VEmModel.hh"
34#include "G4Electron.hh"
35#include "G4Positron.hh"
36#include "G4ParticleChangeForGamma.hh"
37#include "G4CrossSectionHandler.hh"
38#include "G4LogLogInterpolation.hh"
39#include "G4CompositeEMDataSet.hh"
40#include "G4ProductionCutsTable.hh"
41#include "G4ForceCondition.hh"
42
43
44class G4LivermorePolarizedGammaConversionModel : public G4VEmModel
45{
46
47public:
48
49  G4LivermorePolarizedGammaConversionModel(const G4ParticleDefinition* p = 0, 
50                                   const G4String& nam = "LivermorePolarizedGammaConversion");
51
52  virtual ~G4LivermorePolarizedGammaConversionModel();
53
54  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
55
56  virtual G4double ComputeCrossSectionPerAtom(
57                                const G4ParticleDefinition*,
58                                      G4double kinEnergy, 
59                                      G4double Z, 
60                                      G4double A=0, 
61                                      G4double cut=0,
62                                      G4double emax=DBL_MAX);
63
64  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
65                                 const G4MaterialCutsCouple*,
66                                 const G4DynamicParticle*,
67                                 G4double tmin,
68                                 G4double maxEnergy);
69
70protected:
71
72  G4ParticleChangeForGamma* fParticleChange;
73
74  G4double GetMeanFreePath(const G4Track& aTrack, 
75                           G4double previousStepSize, 
76                           G4ForceCondition* condition);
77private:
78
79  G4double lowEnergyLimit; 
80  G4double highEnergyLimit; 
81  G4bool isInitialised;
82  G4int verboseLevel;
83
84  G4VEMDataSet* meanFreePathTable;
85  G4VCrossSectionHandler* crossSectionHandler;
86
87
88  // specific methods for polarization
89 
90  G4ThreeVector GetRandomPolarization(G4ThreeVector& direction0); // Random Polarization
91  G4ThreeVector GetPerpendicularPolarization(const G4ThreeVector& direction0, const G4ThreeVector& polarization0) const;
92  G4ThreeVector SetPerpendicularVector(G4ThreeVector& a); // temporary
93  void SystemOfRefChange(G4ThreeVector& direction0, G4ThreeVector& direction1, 
94                         G4ThreeVector& polarization0);
95
96  // Gamma Conversion methods
97
98
99  G4double SetPhi(G4double);
100  G4double SetPsi(G4double, G4double); //
101
102  G4double ScreenFunction1(G4double screenVariable);
103  G4double ScreenFunction2(G4double screenVariable);
104  void SetTheta(G4double*, G4double*, G4double);
105
106  G4double Poli(G4double , G4double, G4double, G4double);
107  G4double Fln(G4double, G4double, G4double);
108
109  G4double Encu(G4double*, G4double*, G4double);
110
111  G4double Flor(G4double*, G4double);
112  G4double Glor(G4double*, G4double);
113
114  G4double Fdlor(G4double*, G4double);
115  G4double Fintlor(G4double*, G4double);
116  G4double Finvlor(G4double*, G4double,G4double);
117
118  G4double Ftan(G4double*, G4double);
119
120  //G4double Gtan(G4double*, G4double);
121
122  G4double Fdtan(G4double*, G4double);
123  G4double Finttan(G4double*, G4double);
124  G4double Finvtan(G4double*, G4double, G4double);
125
126  G4double smallEnergy;
127  G4double Psi, Phi;
128
129 
130  G4LivermorePolarizedGammaConversionModel & operator=(const  G4LivermorePolarizedGammaConversionModel &right);
131  G4LivermorePolarizedGammaConversionModel(const  G4LivermorePolarizedGammaConversionModel&);
132
133};
134
135//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
136
137#endif
Note: See TracBrowser for help on using the repository browser.