Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PolarizedAnnihilationModel Class Reference

#include <G4PolarizedAnnihilationModel.hh>

+ Inheritance diagram for G4PolarizedAnnihilationModel:

Public Member Functions

 G4PolarizedAnnihilationModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="Polarized-Annihilation")
 
virtual ~G4PolarizedAnnihilationModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) final
 
virtual G4double ComputeCrossSectionPerElectron (G4double kinEnergy) final
 
void ComputeAsymmetriesPerElectron (G4double gammaEnergy, G4double &valueX, G4double &valueA, G4double &valueT)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) final
 
void SetTargetPolarization (const G4ThreeVector &pTarget)
 
void SetBeamPolarization (const G4ThreeVector &pBeam)
 
const G4ThreeVectorGetTargetPolarization () const
 
const G4ThreeVectorGetBeamPolarization () const
 
const G4ThreeVectorGetFinalGamma1Polarization () const
 
const G4ThreeVectorGetFinalGamma2Polarization () const
 
- Public Member Functions inherited from G4eeToTwoGammaModel
 G4eeToTwoGammaModel (const G4ParticleDefinition *p=0, const G4String &nam="eplus2gg")
 
virtual ~G4eeToTwoGammaModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual G4double ComputeCrossSectionPerElectron (G4double kinEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0., G4double maxEnergy=DBL_MAX) override
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4eeToTwoGammaModeloperator= (const G4eeToTwoGammaModel &right)=delete
 
 G4eeToTwoGammaModel (const G4eeToTwoGammaModel &)=delete
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

Detailed Description

Definition at line 62 of file G4PolarizedAnnihilationModel.hh.

Constructor & Destructor Documentation

◆ G4PolarizedAnnihilationModel()

G4PolarizedAnnihilationModel::G4PolarizedAnnihilationModel ( const G4ParticleDefinition p = nullptr,
const G4String nam = "Polarized-Annihilation" 
)
explicit

Definition at line 63 of file G4PolarizedAnnihilationModel.cc.

65 : G4eeToTwoGammaModel(p,nam),
66 crossSectionCalculator(nullptr),
67 verboseLevel(0),
68 gParticleChange(nullptr)
69{
70 crossSectionCalculator = new G4PolarizedAnnihilationCrossSection();
71}

◆ ~G4PolarizedAnnihilationModel()

G4PolarizedAnnihilationModel::~G4PolarizedAnnihilationModel ( )
virtual

Definition at line 73 of file G4PolarizedAnnihilationModel.cc.

74{
75 delete crossSectionCalculator;
76}

Member Function Documentation

◆ ComputeAsymmetriesPerElectron()

void G4PolarizedAnnihilationModel::ComputeAsymmetriesPerElectron ( G4double  gammaEnergy,
G4double valueX,
G4double valueA,
G4double valueT 
)

Definition at line 107 of file G4PolarizedAnnihilationModel.cc.

111{
112 // *** calculate asymmetries
113 G4double gam = 1. + ene/electron_mass_c2;
114 G4double xs0=crossSectionCalculator->TotalXSection(0.,1.,gam,
117 G4double xsA=crossSectionCalculator->TotalXSection(0.,1.,gam,
120 G4double xsT1=crossSectionCalculator->TotalXSection(0.,1.,gam,
123 G4double xsT2=crossSectionCalculator->TotalXSection(0.,1.,gam,
126 G4double xsT=0.5*(xsT1+xsT2);
127
128 valueX=xs0;
129 valueA=xsA/xs0-1.;
130 valueT=xsT/xs0-1.;
131 // G4cout<<valueX<<"\t"<<valueA<<"\t"<<valueT<<" energy = "<<gam<<G4endl;
132 if ( (valueA < -1) || (1 < valueA)) {
133 G4cout<< " ERROR PolarizedAnnihilationPS::ComputeAsymmetries \n";
134 G4cout<< " something wrong in total cross section calculation (valueA)\n";
135 G4cout<< " LONG: "<<valueX<<"\t"<<valueA<<"\t"<<valueT<<" energy = "<<gam<<G4endl;
136 }
137 if ( (valueT < -1) || (1 < valueT)) {
138 G4cout<< " ERROR PolarizedAnnihilationPS::ComputeAsymmetries \n";
139 G4cout<< " something wrong in total cross section calculation (valueT)\n";
140 G4cout<< " TRAN: "<<valueX<<"\t"<<valueA<<"\t"<<valueT<<" energy = "<<gam<<G4endl;
141 }
142}
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
virtual G4double TotalXSection(G4double xmin, G4double xmax, G4double y, const G4StokesVector &pol0, const G4StokesVector &pol1) override
static const G4StokesVector P3
static const G4StokesVector ZERO
static const G4StokesVector P2
static const G4StokesVector P1

Referenced by ComputeCrossSectionPerElectron().

◆ ComputeCrossSectionPerElectron()

G4double G4PolarizedAnnihilationModel::ComputeCrossSectionPerElectron ( G4double  kinEnergy)
finalvirtual

Reimplemented from G4eeToTwoGammaModel.

Definition at line 90 of file G4PolarizedAnnihilationModel.cc.

91{
92 // cross section from base model
94
95 G4double polzz = theBeamPolarization.z()*theTargetPolarization.z();
96 G4double poltt = theBeamPolarization.x()*theTargetPolarization.x()
97 + theBeamPolarization.y()*theTargetPolarization.y();
98 if (polzz!=0 || poltt!=0) {
99 G4double xval,lasym,tasym;
100 ComputeAsymmetriesPerElectron(kinEnergy,xval,lasym,tasym);
101 xs*=(1.+polzz*lasym+poltt*tasym);
102 }
103
104 return xs;
105}
double z() const
double x() const
double y() const
void ComputeAsymmetriesPerElectron(G4double gammaEnergy, G4double &valueX, G4double &valueA, G4double &valueT)
virtual G4double ComputeCrossSectionPerElectron(G4double kinEnergy)

◆ GetBeamPolarization()

const G4ThreeVector & G4PolarizedAnnihilationModel::GetBeamPolarization ( ) const
inline

Definition at line 131 of file G4PolarizedAnnihilationModel.hh.

132{
133 return theBeamPolarization;
134}

◆ GetFinalGamma1Polarization()

const G4ThreeVector & G4PolarizedAnnihilationModel::GetFinalGamma1Polarization ( ) const
inline

Definition at line 135 of file G4PolarizedAnnihilationModel.hh.

136{
137 return finalGamma1Polarization;
138}

◆ GetFinalGamma2Polarization()

const G4ThreeVector & G4PolarizedAnnihilationModel::GetFinalGamma2Polarization ( ) const
inline

Definition at line 139 of file G4PolarizedAnnihilationModel.hh.

140{
141 return finalGamma2Polarization;
142}

◆ GetTargetPolarization()

const G4ThreeVector & G4PolarizedAnnihilationModel::GetTargetPolarization ( ) const
inline

Definition at line 127 of file G4PolarizedAnnihilationModel.hh.

128{
129 return theTargetPolarization;
130}

◆ Initialise()

void G4PolarizedAnnihilationModel::Initialise ( const G4ParticleDefinition part,
const G4DataVector dv 
)
finalvirtual

Reimplemented from G4eeToTwoGammaModel.

Definition at line 81 of file G4PolarizedAnnihilationModel.cc.

83{
85 if(gParticleChange) { return; }
86 gParticleChange = GetParticleChangeForGamma();
87}
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override

◆ SampleSecondaries()

void G4PolarizedAnnihilationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple ,
const G4DynamicParticle dp,
G4double  tmin,
G4double  maxEnergy 
)
finalvirtual

Reimplemented from G4eeToTwoGammaModel.

Definition at line 144 of file G4PolarizedAnnihilationModel.cc.

148{
149 const G4Track * aTrack = gParticleChange->GetCurrentTrack();
150
151 // kill primary
152 gParticleChange->SetProposedKineticEnergy(0.);
153 gParticleChange->ProposeTrackStatus(fStopAndKill);
154
155 // V.Ivanchenko add protection against zero kin energy
156 G4double PositKinEnergy = dp->GetKineticEnergy();
157
158 if(PositKinEnergy == 0.0) {
159
160 G4double cosTeta = 2.*G4UniformRand()-1.;
161 G4double sinTeta = std::sqrt((1.0 - cosTeta)*(1.0 + cosTeta));
162 G4double phi = twopi * G4UniformRand();
163 G4ThreeVector dir(sinTeta*std::cos(phi), sinTeta*std::sin(phi), cosTeta);
164 fvect->push_back( new G4DynamicParticle(G4Gamma::Gamma(), dir, electron_mass_c2));
165 fvect->push_back( new G4DynamicParticle(G4Gamma::Gamma(),-dir, electron_mass_c2));
166 return;
167 }
168
169 // *** obtain and save target and beam polarization ***
171
172 // obtain polarization of the beam
173 theBeamPolarization = aTrack->GetPolarization();
174
175 // obtain polarization of the media
176 G4VPhysicalVolume* aPVolume = aTrack->GetVolume();
177 G4LogicalVolume* aLVolume = aPVolume->GetLogicalVolume();
178 const G4bool targetIsPolarized = polarizationManager->IsPolarized(aLVolume);
179 theTargetPolarization = polarizationManager->GetVolumePolarization(aLVolume);
180
181 if (verboseLevel >= 1) {
182 G4cout << "G4PolarizedComptonModel::SampleSecondaries in "
183 << aLVolume->GetName() << G4endl;
184 }
185
186 // transfer target electron polarization in frame of positron
187 if (targetIsPolarized)
188 theTargetPolarization.rotateUz(dp->GetMomentumDirection());
189
190 G4ParticleMomentum PositDirection = dp->GetMomentumDirection();
191
192 // polar asymmetry:
193 G4double polarization = theBeamPolarization.p3()*theTargetPolarization.p3();
194
195 G4double gamam1 = PositKinEnergy/electron_mass_c2;
196 G4double gama = gamam1+1. , gamap1 = gamam1+2.;
197 G4double sqgrate = std::sqrt(gamam1/gamap1)/2. , sqg2m1 = std::sqrt(gamam1*gamap1);
198
199 // limits of the energy sampling
200 G4double epsilmin = 0.5 - sqgrate , epsilmax = 0.5 + sqgrate;
201 G4double epsilqot = epsilmax/epsilmin;
202
203 //
204 // sample the energy rate of the created gammas
205 // note: for polarized partices, the actual dicing strategy
206 // will depend on the energy, and the degree of polarization !!
207 //
208 G4double epsil;
209 G4double gmax=1. + std::fabs(polarization); // crude estimate
210
211 //G4bool check_range=true;
212
213 crossSectionCalculator->Initialize(epsilmin, gama, 0., theBeamPolarization, theTargetPolarization);
214 if (crossSectionCalculator->DiceEpsilon()<0) {
215 G4cout<<"ERROR in PolarizedAnnihilationPS::PostStepDoIt\n"
216 <<"epsilmin DiceRoutine not appropriate ! "<<crossSectionCalculator->DiceEpsilon()<<G4endl;
217 //check_range=false;
218 }
219
220 crossSectionCalculator->Initialize(epsilmax, gama, 0., theBeamPolarization, theTargetPolarization);
221 if (crossSectionCalculator->DiceEpsilon()<0) {
222 G4cout<<"ERROR in PolarizedAnnihilationPS::PostStepDoIt\n"
223 <<"epsilmax DiceRoutine not appropriate ! "<<crossSectionCalculator->DiceEpsilon()<<G4endl;
224 //check_range=false;
225 }
226
227 G4int ncount=0;
228 G4double trejectmax=0.;
229 G4double treject;
230
231
232 do {
233 //
234 epsil = epsilmin*std::pow(epsilqot,G4UniformRand());
235
236 crossSectionCalculator->Initialize(epsil, gama, 0., theBeamPolarization, theTargetPolarization,1);
237
238 treject = crossSectionCalculator->DiceEpsilon();
239 treject*=epsil;
240
241 if (treject>gmax || treject<0.)
242 G4cout<<"ERROR in PolarizedAnnihilationPS::PostStepDoIt\n"
243 <<" eps ("<<epsil<<") rejection does not work properly: "<<treject<<G4endl;
244 ++ncount;
245 if (treject>trejectmax) trejectmax=treject;
246 if (ncount>1000) {
247 G4cout<<"WARNING in PolarizedAnnihilationPS::PostStepDoIt\n"
248 <<"eps dicing very inefficient ="<<trejectmax/gmax
249 <<", "<<treject/gmax<<". For secondary energy = "<<epsil<<" "<<ncount<<G4endl;
250 break;
251 }
252
253 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
254 } while( treject < gmax*G4UniformRand() );
255
256 //
257 // scattered Gamma angles. ( Z - axis along the parent positron)
258 //
259
260 G4double cost = (epsil*gamap1-1.)/(epsil*sqg2m1);
261 G4double sint = std::sqrt((1.+cost)*(1.-cost));
262 G4double phi = 0.;
263 G4double beamTrans = std::sqrt(sqr(theBeamPolarization.p1()) + sqr(theBeamPolarization.p2()));
264 G4double targetTrans = std::sqrt(sqr(theTargetPolarization.p1()) + sqr(theTargetPolarization.p2()));
265
266 // G4cout<<"phi dicing START"<<G4endl;
267 do{
268 phi = twopi * G4UniformRand();
269 crossSectionCalculator->Initialize(epsil, gama, 0., theBeamPolarization, theTargetPolarization,2);
270
271 G4double gdiced =crossSectionCalculator->getVar(0);
272 gdiced += crossSectionCalculator->getVar(3)*theBeamPolarization.p3()*theTargetPolarization.p3();
273 gdiced += 1.*(std::fabs(crossSectionCalculator->getVar(1))
274 + std::fabs(crossSectionCalculator->getVar(2)))*beamTrans*targetTrans;
275 gdiced += 1.*std::fabs(crossSectionCalculator->getVar(4))
276 *(std::fabs(theBeamPolarization.p3())*targetTrans + std::fabs(theTargetPolarization.p3())*beamTrans);
277
278 G4double gdist = crossSectionCalculator->getVar(0);
279 gdist += crossSectionCalculator->getVar(3)*theBeamPolarization.p3()*theTargetPolarization.p3();
280 gdist += crossSectionCalculator->getVar(1)*(std::cos(phi)*theBeamPolarization.p1()
281 + std::sin(phi)*theBeamPolarization.p2())
282 *(std::cos(phi)*theTargetPolarization.p1()
283 + std::sin(phi)*theTargetPolarization.p2());
284 gdist += crossSectionCalculator->getVar(2)*(std::cos(phi)*theBeamPolarization.p2()
285 - std::sin(phi)*theBeamPolarization.p1())
286 *(std::cos(phi)*theTargetPolarization.p2()
287 - std::sin(phi)*theTargetPolarization.p1());
288 gdist += crossSectionCalculator->getVar(4)
289 *(std::cos(phi)*theBeamPolarization.p3()*theTargetPolarization.p1()
290 + std::cos(phi)*theBeamPolarization.p1()*theTargetPolarization.p3()
291 + std::sin(phi)*theBeamPolarization.p3()*theTargetPolarization.p2()
292 + std::sin(phi)*theBeamPolarization.p2()*theTargetPolarization.p3());
293
294 treject = gdist/gdiced;
295 //G4cout<<" treject = "<<treject<<" at phi = "<<phi<<G4endl;
296 if (treject>1.+1.e-10 || treject<0){
297 G4cout<<"!!!ERROR in PolarizedAnnihilationPS::PostStepDoIt\n"
298 <<" phi rejection does not work properly: "<<treject<<G4endl;
299 G4cout<<" gdiced = "<<gdiced<<G4endl;
300 G4cout<<" gdist = "<<gdist<<G4endl;
301 G4cout<<" epsil = "<<epsil<<G4endl;
302 }
303
304 if (treject<1.e-3) {
305 G4cout<<"!!!ERROR in PolarizedAnnihilationPS::PostStepDoIt\n"
306 <<" phi rejection does not work properly: "<<treject<<"\n";
307 G4cout<<" gdiced="<<gdiced<<" gdist="<<gdist<<"\n";
308 G4cout<<" epsil = "<<epsil<<G4endl;
309 }
310
311 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
312 } while( treject < G4UniformRand() );
313 // G4cout<<"phi dicing END"<<G4endl;
314
315 G4double dirx = sint*std::cos(phi) , diry = sint*std::sin(phi) , dirz = cost;
316
317 //
318 // kinematic of the created pair
319 //
320 G4double TotalAvailableEnergy = PositKinEnergy + 2*electron_mass_c2;
321 G4double Phot1Energy = epsil*TotalAvailableEnergy;
322 G4double Phot2Energy =(1.-epsil)*TotalAvailableEnergy;
323
324 // *** prepare calculation of polarization transfer ***
325 G4ThreeVector Phot1Direction (dirx, diry, dirz);
326
327 // get interaction frame
328 G4ThreeVector nInteractionFrame =
329 G4PolarizationHelper::GetFrame(PositDirection,Phot1Direction);
330
331 // define proper in-plane and out-of-plane component of initial spins
332 theBeamPolarization.InvRotateAz(nInteractionFrame,PositDirection);
333 theTargetPolarization.InvRotateAz(nInteractionFrame,PositDirection);
334
335 // calculate spin transfere matrix
336
337 crossSectionCalculator->Initialize(epsil,gama,phi,theBeamPolarization,theTargetPolarization,2);
338
339 // **********************************************************************
340
341 Phot1Direction.rotateUz(PositDirection);
342 // create G4DynamicParticle object for the particle1
344 Phot1Direction, Phot1Energy);
345 finalGamma1Polarization=crossSectionCalculator->GetPol2();
346 G4double n1=finalGamma1Polarization.mag2();
347 if (n1>1) {
348 G4cout<<"ERROR: PolarizedAnnihilation Polarization Vector at epsil = "
349 <<epsil<<" is too large!!! \n"
350 <<"annihi pol1= "<<finalGamma1Polarization<<", ("<<n1<<")\n";
351 finalGamma1Polarization*=1./std::sqrt(n1);
352 }
353
354 // define polarization of first final state photon
355 finalGamma1Polarization.SetPhoton();
356 finalGamma1Polarization.RotateAz(nInteractionFrame,Phot1Direction);
357 aParticle1->SetPolarization(finalGamma1Polarization.p1(),
358 finalGamma1Polarization.p2(),
359 finalGamma1Polarization.p3());
360
361 fvect->push_back(aParticle1);
362
363
364 // **********************************************************************
365
366 G4double Eratio= Phot1Energy/Phot2Energy;
367 G4double PositP= std::sqrt(PositKinEnergy*(PositKinEnergy+2.*electron_mass_c2));
368 G4ThreeVector Phot2Direction (-dirx*Eratio, -diry*Eratio,
369 (PositP-dirz*Phot1Energy)/Phot2Energy);
370 Phot2Direction.rotateUz(PositDirection);
371 // create G4DynamicParticle object for the particle2
373 Phot2Direction, Phot2Energy);
374
375 // define polarization of second final state photon
376 finalGamma2Polarization=crossSectionCalculator->GetPol3();
377 G4double n2=finalGamma2Polarization.mag2();
378 if (n2>1) {
379 G4cout<<"ERROR: PolarizedAnnihilation Polarization Vector at epsil = "<<epsil<<" is too large!!! \n";
380 G4cout<<"annihi pol2= "<<finalGamma2Polarization<<", ("<<n2<<")\n";
381
382 finalGamma2Polarization*=1./std::sqrt(n2);
383 }
384 finalGamma2Polarization.SetPhoton();
385 finalGamma2Polarization.RotateAz(nInteractionFrame,Phot2Direction);
386 aParticle2->SetPolarization(finalGamma2Polarization.p1(),
387 finalGamma2Polarization.p2(),
388 finalGamma2Polarization.p3());
389
390 fvect->push_back(aParticle2);
391}
@ fStopAndKill
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
double mag2() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
void SetPolarization(const G4ThreeVector &)
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
const G4String & GetName() const
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
static G4ThreeVector GetFrame(const G4ThreeVector &, const G4ThreeVector &)
bool IsPolarized(G4LogicalVolume *lVol) const
static G4PolarizationManager * GetInstance()
const G4ThreeVector & GetVolumePolarization(G4LogicalVolume *lVol) const
virtual void Initialize(G4double eps, G4double gamma, G4double phi, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0) override
G4double p3() const
G4double p1() const
void InvRotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
G4double p2() const
void RotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPolarization() const
void ProposeTrackStatus(G4TrackStatus status)
G4LogicalVolume * GetLogicalVolume() const
T sqr(const T &x)
Definition: templates.hh:128

◆ SetBeamPolarization()

void G4PolarizedAnnihilationModel::SetBeamPolarization ( const G4ThreeVector pBeam)
inline

Definition at line 123 of file G4PolarizedAnnihilationModel.hh.

124{
125 theBeamPolarization = pBeam;
126}

◆ SetTargetPolarization()

void G4PolarizedAnnihilationModel::SetTargetPolarization ( const G4ThreeVector pTarget)
inline

Definition at line 119 of file G4PolarizedAnnihilationModel.hh.

120{
121 theTargetPolarization = pTarget;
122}

The documentation for this class was generated from the following files: