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: G4PolarizationHelper.cc,v 1.4 2007/11/01 17:30:25 schaelic Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-02 $ |
---|
28 | // |
---|
29 | // GEANT4 Class file |
---|
30 | // |
---|
31 | // |
---|
32 | // File name: G4PolarizationHelper |
---|
33 | // |
---|
34 | // Author: Andreas Schaelicke |
---|
35 | // |
---|
36 | // Creation date: 12.08.2006 |
---|
37 | // |
---|
38 | // Modifications: |
---|
39 | // |
---|
40 | // Class Description: |
---|
41 | // |
---|
42 | // Provides some basic polarization transformation routines. |
---|
43 | // |
---|
44 | #include "G4PolarizationHelper.hh" |
---|
45 | #include "G4StokesVector.hh" |
---|
46 | |
---|
47 | |
---|
48 | G4ThreeVector G4PolarizationHelper::GetFrame(const G4ThreeVector & mom1, const G4ThreeVector & mom2) |
---|
49 | { |
---|
50 | G4ThreeVector normal = mom1.cross(mom2); |
---|
51 | return 1./normal.mag()*normal; |
---|
52 | } |
---|
53 | |
---|
54 | G4ThreeVector G4PolarizationHelper::GetParticleFrameY(const G4ThreeVector &uZ) |
---|
55 | { |
---|
56 | // compare also G4ThreeVector::rotateUz() |
---|
57 | |
---|
58 | if (uZ.x()==0. && uZ.y()==0.) { |
---|
59 | return G4ThreeVector(0.,1.,0.); |
---|
60 | } |
---|
61 | |
---|
62 | G4double invPerp = 1./std::sqrt(sqr(uZ.x())+sqr(uZ.y())); |
---|
63 | return G4ThreeVector(-uZ.y()*invPerp,uZ.x()*invPerp,0); |
---|
64 | } |
---|
65 | |
---|
66 | G4ThreeVector G4PolarizationHelper::GetParticleFrameX(const G4ThreeVector &uZ) |
---|
67 | { |
---|
68 | // compare also G4ThreeVector::rotateUz() |
---|
69 | |
---|
70 | if (uZ.x()==0. && uZ.y()==0.) { |
---|
71 | if (uZ.z()>=0.) return G4ThreeVector(1.,0.,0.); |
---|
72 | return G4ThreeVector(-1.,0.,0.); |
---|
73 | } |
---|
74 | |
---|
75 | G4double perp = std::sqrt(sqr(uZ.x())+sqr(uZ.y())); |
---|
76 | G4double invPerp = uZ.z()/perp; |
---|
77 | return G4ThreeVector(uZ.x()*invPerp,uZ.y()*invPerp,-perp); |
---|
78 | } |
---|
79 | |
---|
80 | G4ThreeVector G4PolarizationHelper::GetRandomFrame(const G4ThreeVector & mom1) |
---|
81 | { |
---|
82 | G4double phi =2.*pi*G4UniformRand(); |
---|
83 | G4ThreeVector normal = std::cos(phi)*GetParticleFrameX(mom1) |
---|
84 | + std::sin(phi)*G4PolarizationHelper::GetParticleFrameY(mom1); |
---|
85 | return normal; |
---|
86 | } |
---|
87 | |
---|
88 | |
---|
89 | G4ThreeVector G4PolarizationHelper::GetSpinInPRF(const G4ThreeVector &uZ, const G4ThreeVector & spin) |
---|
90 | { |
---|
91 | // compare also G4ThreeVector::rotateUz() |
---|
92 | |
---|
93 | if (uZ.x()==0. && uZ.y()==0.) { |
---|
94 | if (uZ.z()>=0.) return spin; |
---|
95 | return G4ThreeVector(-spin.x(),spin.y(),-spin.z()); |
---|
96 | } |
---|
97 | |
---|
98 | G4double perp = std::sqrt(sqr(uZ.x())+sqr(uZ.y())); |
---|
99 | G4double invPerp = uZ.z()/perp; |
---|
100 | |
---|
101 | G4ThreeVector uX(uZ.x()*invPerp,uZ.y()*invPerp,-perp); |
---|
102 | G4ThreeVector uY(-uZ.y()*invPerp,uZ.x()*invPerp,0); |
---|
103 | |
---|
104 | return G4ThreeVector(spin*uX,spin*uY,spin*uZ); |
---|
105 | } |
---|
106 | |
---|
107 | void G4PolarizationHelper::TestPolarizationTransformations() |
---|
108 | { |
---|
109 | G4double theta=0.; |
---|
110 | G4cout<<"========================================\n\n"; |
---|
111 | for (G4int i=0; i<=10; ++i) { |
---|
112 | theta=pi*i/10.; |
---|
113 | G4ThreeVector zAxis = G4ThreeVector(std::sin(theta),0.,std::cos(theta)); |
---|
114 | if (i==5) zAxis = G4ThreeVector(1.,0.,0.); |
---|
115 | if (i==10) zAxis = G4ThreeVector(0.,0.,-1.); |
---|
116 | G4ThreeVector yAxis = GetParticleFrameY(zAxis); |
---|
117 | |
---|
118 | G4cout<<zAxis<<" "<<zAxis.mag()<<"\n"; |
---|
119 | G4cout<<yAxis<<" "<<yAxis.mag()<<"\n"; |
---|
120 | G4ThreeVector xAxis = yAxis.cross(zAxis); |
---|
121 | G4cout<<xAxis<<" "<<xAxis.mag()<<"\n\n"; |
---|
122 | } |
---|
123 | |
---|
124 | G4cout<<"========================================\n\n"; |
---|
125 | |
---|
126 | for (G4int i=0; i<=10; ++i) { |
---|
127 | theta=pi*i/10.; |
---|
128 | G4ThreeVector zAxis = G4ThreeVector(0.,std::sin(theta),std::cos(theta)); |
---|
129 | if (i==5) zAxis = G4ThreeVector(0.,1.,0.); |
---|
130 | if (i==10) zAxis = G4ThreeVector(0.,0.,-1.); |
---|
131 | G4ThreeVector yAxis = GetParticleFrameY(zAxis); |
---|
132 | |
---|
133 | G4cout<<zAxis<<" "<<zAxis.mag()<<"\n"; |
---|
134 | G4cout<<yAxis<<" "<<yAxis.mag()<<"\n"; |
---|
135 | G4ThreeVector xAxis = yAxis.cross(zAxis); |
---|
136 | G4cout<<xAxis<<" "<<xAxis.mag()<<"\n\n"; |
---|
137 | } |
---|
138 | G4cout<<"========================================\n\n"; |
---|
139 | } |
---|
140 | |
---|
141 | void G4PolarizationHelper::TestInteractionFrame() |
---|
142 | { |
---|
143 | // check transformation procedure for polarisation transfer |
---|
144 | // calculation in scattering processes |
---|
145 | // a) transfer target polarisation in beam particle reference frame (PRF) |
---|
146 | // b) calc correct asymmetry w.r.t. scattering plane |
---|
147 | // c) determine incomming polarisation in interaction frame (IF) |
---|
148 | // d) transfer outgoing polarisation from IF to PRF |
---|
149 | G4cout<<"========================================\n\n"; |
---|
150 | |
---|
151 | G4double theta=0.; |
---|
152 | |
---|
153 | G4ThreeVector dir0=G4ThreeVector(0.,0.,1.); |
---|
154 | G4ThreeVector dir2=G4ThreeVector(std::sin(theta),0.,std::cos(theta)); |
---|
155 | |
---|
156 | G4StokesVector pol0=G4ThreeVector(0.,0.,1.); |
---|
157 | G4StokesVector pol1=G4ThreeVector(0.,0.,1.); |
---|
158 | |
---|
159 | pol1.rotateUz(dir0); |
---|
160 | |
---|
161 | G4cout<<"========================================\n\n"; |
---|
162 | |
---|
163 | |
---|
164 | } |
---|