source: trunk/source/processes/electromagnetic/lowenergy/include/G4LowEnergyPolarizedCompton.hh@ 899

Last change on this file since 899 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 4.9 KB
RevLine 
[819]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// $Id: G4LowEnergyPolarizedCompton.hh,v 1.9 2006/06/29 19:36:05 gunter Exp $
28// GEANT4 tag $Name: $
29//
30// ------------------------------------------------------------
31// GEANT 4 class header file
32// CERN Geneva Switzerland
33//
34
35// --------- G4LowEnergyPolarizedCompton class -----
36//
37// by G.Depaola & F.Longo (21 may 2001)
38// 24 May 2001 - MGP Modified to inherit from G4VDiscreteProcess
39// 25 May 2001 - MGP Added protections to avoid crashes
40//
41// 17 October 2001 - F.Longo Major revision according to design iteration
42//
43// 21 February 2002 - F.Longo Revisions with A.Zoglauer and G.Depaola
44// - better description of parallelism
45// - system of ref change method improved
46//
47//
48// Class description:
49// Low Energy electromagnetic process, Polarised Compton scattering
50// Further documentation available from http://www.ge.infn.it/geant4/lowE
51
52// ------------------------------------------------------------
53
54#ifndef G4LOWENERGYPOLARIZEDCOMPTON_H
55#define G4LOWENERGYPOLARIZEDCOMPTON_H 1
56
57#include "globals.hh"
58#include "G4VDiscreteProcess.hh"
59
60class G4Track;
61class G4Step;
62class G4ParticleDefinition;
63class G4VParticleChange;
64class G4VEMDataSet;
65class G4VCrossSectionHandler;
66class G4VRangeTest;
67
68class G4LowEnergyPolarizedCompton : public G4VDiscreteProcess
69{
70public:
71
72 G4LowEnergyPolarizedCompton(const G4String& processName = "polarLowEnCompt");
73
74 ~G4LowEnergyPolarizedCompton();
75
76 G4bool IsApplicable(const G4ParticleDefinition& definition);
77
78 void BuildPhysicsTable(const G4ParticleDefinition& photon);
79 G4VParticleChange* PostStepDoIt(const G4Track& aTrack, const G4Step& aStep);
80
81
82 // For testing purpose only
83 G4double DumpMeanFreePath(const G4Track& aTrack,
84 G4double previousStepSize,
85 G4ForceCondition* condition)
86 { return GetMeanFreePath(aTrack, previousStepSize, condition); }
87
88
89protected:
90
91 G4double GetMeanFreePath(const G4Track& aTrack,
92 G4double previousStepSize,
93 G4ForceCondition* condition);
94private:
95
96 // Hide copy constructor and assignment operator as private
97
98 G4LowEnergyPolarizedCompton& operator=(const G4LowEnergyPolarizedCompton&
99 right);
100 G4LowEnergyPolarizedCompton(const G4LowEnergyPolarizedCompton& );
101
102 G4double lowEnergyLimit; // low energy limit applied to the process
103 G4double highEnergyLimit; // high energy limit applied to the process
104
105 G4VEMDataSet* meanFreePathTable;
106 G4VEMDataSet* scatterFunctionData;
107
108 G4VCrossSectionHandler* crossSectionHandler;
109 G4VRangeTest* rangeTest;
110
111 const G4double intrinsicLowEnergyLimit; // intrinsic validity range
112 const G4double intrinsicHighEnergyLimit;
113
114 // specific methods for polarization
115
116 G4ThreeVector GetRandomPolarization(G4ThreeVector& direction0); // Random Polarization
117 G4ThreeVector GetPerpendicularPolarization(const G4ThreeVector& direction0, const G4ThreeVector& polarization0) const;
118
119 G4ThreeVector SetPerpendicularVector(G4ThreeVector& a); // temporary
120 G4ThreeVector SetNewPolarization(G4double epsilon, G4double sinSqrTheta,
121 G4double phi, G4double cosTheta);
122 G4double SetPhi(G4double, G4double);
123
124 void SystemOfRefChange(G4ThreeVector& direction0, G4ThreeVector& direction1,
125 G4ThreeVector& polarization0, G4ThreeVector& polarization1);
126
127};
128
129#endif
130
131
132
133
134
135
136
137
138
Note: See TracBrowser for help on using the repository browser.