Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LivermorePolarizedComptonModel.cc
Go to the documentation of this file.
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// Authors: G.Depaola & F.Longo
28//
29// History:
30// -------
31//
32// 05 Apr 2021 J Allison added quantum entanglement of e+ annihilation.
33// If the photons have been "tagged" as "quantum-entangled", for example by
34// G4eplusAnnihilation for annihilation into 2 photons, they are "analysed"
35// here if - and only if - both photons suffer Compton scattering. Theoretical
36// predictions from Pryce and Ward, Nature No 4065 (1947) p.435, and Snyder et al,
37// Physical Review 73 (1948) p.440. Experimental validation in "Photon quantum
38// entanglement in the MeV regime and its application in PET imaging",
39// D. Watts, J. Allison et al., Nature Communications (2021)12:2646
40// https://doi.org/10.1038/s41467-021-22907-5.
41//
42// 02 May 2009 S Incerti as V. Ivanchenko proposed in G4LivermoreComptonModel.cc
43//
44// Cleanup initialisation and generation of secondaries:
45// - apply internal high-energy limit only in constructor
46// - do not apply low-energy limit (default is 0)
47// - remove GetMeanFreePath method and table
48// - added protection against numerical problem in energy sampling
49// - use G4ElementSelector
50
53#include "G4SystemOfUnits.hh"
54#include "G4AutoLock.hh"
55#include "G4Electron.hh"
57#include "G4LossTableManager.hh"
59#include "G4AtomicShell.hh"
60#include "G4Gamma.hh"
61#include "G4ShellData.hh"
62#include "G4DopplerProfile.hh"
63#include "G4Log.hh"
64#include "G4Exp.hh"
65#include "G4Pow.hh"
70
71//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
72
73using namespace std;
74namespace { G4Mutex LivermorePolarizedComptonModelMutex = G4MUTEX_INITIALIZER; }
75
76
77G4PhysicsFreeVector* G4LivermorePolarizedComptonModel::data[] = {nullptr};
78G4ShellData* G4LivermorePolarizedComptonModel::shellData = nullptr;
79G4DopplerProfile* G4LivermorePolarizedComptonModel::profileData = nullptr;
80G4CompositeEMDataSet* G4LivermorePolarizedComptonModel::scatterFunctionData = nullptr;
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
83
85 :G4VEmModel(nam),isInitialised(false)
86{
87 verboseLevel= 1;
88 // Verbosity scale:
89 // 0 = nothing
90 // 1 = warning for energy non-conservation
91 // 2 = details of energy budget
92 // 3 = calculation of cross sections, file openings, sampling of atoms
93 // 4 = entering in methods
94
95 if( verboseLevel>1 )
96 G4cout << "Livermore Polarized Compton is constructed " << G4endl;
97
98 //Mark this model as "applicable" for atomic deexcitation
100
101 fParticleChange = nullptr;
102 fAtomDeexcitation = nullptr;
103 fEntanglementModelID = G4PhysicsModelCatalog::GetModelID("model_GammaGammaEntanglement");
104}
105
106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107
109{
110 if(IsMaster()) {
111 delete shellData;
112 shellData = nullptr;
113 delete profileData;
114 profileData = nullptr;
115 delete scatterFunctionData;
116 scatterFunctionData = nullptr;
117 for(G4int i=0; i<maxZ; ++i) {
118 if(data[i]) {
119 delete data[i];
120 data[i] = nullptr;
121 }
122 }
123 }
124}
125
126//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
127
129 const G4DataVector& cuts)
130{
131 if (verboseLevel > 1)
132 G4cout << "Calling G4LivermorePolarizedComptonModel::Initialise()" << G4endl;
133
134 // Initialise element selector
135 if(IsMaster()) {
136 // Access to elements
137 const char* path = G4FindDataDir("G4LEDATA");
138
139 G4ProductionCutsTable* theCoupleTable =
141
142 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
143
144 for(G4int i=0; i<numOfCouples; ++i) {
145 const G4Material* material =
146 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
147 const G4ElementVector* theElementVector = material->GetElementVector();
148 std::size_t nelm = material->GetNumberOfElements();
149
150 for (std::size_t j=0; j<nelm; ++j) {
151 G4int Z = G4lrint((*theElementVector)[j]->GetZ());
152 if(Z < 1) { Z = 1; }
153 else if(Z > maxZ){ Z = maxZ; }
154
155 if( (!data[Z]) ) { ReadData(Z, path); }
156 }
157 }
158
159 // For Doppler broadening
160 if(!shellData) {
161 shellData = new G4ShellData();
162 shellData->SetOccupancyData();
163 G4String file = "/doppler/shell-doppler";
164 shellData->LoadData(file);
165 }
166 if(!profileData) { profileData = new G4DopplerProfile(); }
167
168 // Scattering Function
169 if(!scatterFunctionData)
170 {
171
172 G4VDataSetAlgorithm* scatterInterpolation = new G4LogLogInterpolation;
173 G4String scatterFile = "comp/ce-sf-";
174 scatterFunctionData = new G4CompositeEMDataSet(scatterInterpolation, 1., 1.);
175 scatterFunctionData->LoadData(scatterFile);
176 }
177
178 InitialiseElementSelectors(particle, cuts);
179 }
180
181 if (verboseLevel > 2) {
182 G4cout << "Loaded cross section files" << G4endl;
183 }
184
185 if( verboseLevel>1 ) {
186 G4cout << "G4LivermoreComptonModel is initialized " << G4endl
187 << "Energy range: "
188 << LowEnergyLimit() / eV << " eV - "
189 << HighEnergyLimit() / GeV << " GeV"
190 << G4endl;
191 }
192 //
193 if(isInitialised) { return; }
194
195 fParticleChange = GetParticleChangeForGamma();
196 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
197 isInitialised = true;
198}
199
200
206
207//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
208
209void G4LivermorePolarizedComptonModel::ReadData(std::size_t Z, const char* path)
210{
211 if (verboseLevel > 1)
212 {
213 G4cout << "G4LivermorePolarizedComptonModel::ReadData()"
214 << G4endl;
215 }
216 if(data[Z]) { return; }
217 const char* datadir = path;
218 if(!datadir)
219 {
220 datadir = G4FindDataDir("G4LEDATA");
221 if(!datadir)
222 {
223 G4Exception("G4LivermorePolarizedComptonModel::ReadData()",
224 "em0006",FatalException,
225 "Environment variable G4LEDATA not defined");
226 return;
227 }
228 }
229
230 data[Z] = new G4PhysicsFreeVector();
231
232 std::ostringstream ost;
233 ost << datadir << "/livermore/comp/ce-cs-" << Z <<".dat";
234 std::ifstream fin(ost.str().c_str());
235
236 if( !fin.is_open())
237 {
239 ed << "G4LivermorePolarizedComptonModel data file <" << ost.str().c_str()
240 << "> is not opened!" << G4endl;
241 G4Exception("G4LivermoreComptonModel::ReadData()",
242 "em0003",FatalException,
243 ed,"G4LEDATA version should be G4EMLOW8.0 or later");
244 return;
245 } else {
246 if(verboseLevel > 3) {
247 G4cout << "File " << ost.str()
248 << " is opened by G4LivermorePolarizedComptonModel" << G4endl;
249 }
250 data[Z]->Retrieve(fin, true);
251 data[Z]->ScaleVector(MeV, MeV*barn);
252 }
253 fin.close();
254}
255
256//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
257
260 G4double GammaEnergy,
263{
264 if (verboseLevel > 3)
265 G4cout << "Calling ComputeCrossSectionPerAtom() of G4LivermorePolarizedComptonModel" << G4endl;
266
267 G4double cs = 0.0;
268
269 if (GammaEnergy < LowEnergyLimit())
270 return 0.0;
271
272 G4int intZ = G4lrint(Z);
273 if(intZ < 1 || intZ > maxZ) { return cs; }
274
275 G4PhysicsFreeVector* pv = data[intZ];
276
277 // if element was not initialised
278 // do initialisation safely for MT mode
279 if(!pv)
280 {
281 InitialiseForElement(0, intZ);
282 pv = data[intZ];
283 if(!pv) { return cs; }
284 }
285
286 G4int n = G4int(pv->GetVectorLength() - 1);
287 G4double e1 = pv->Energy(0);
288 G4double e2 = pv->Energy(n);
289
290 if(GammaEnergy <= e1) { cs = GammaEnergy/(e1*e1)*pv->Value(e1); }
291 else if(GammaEnergy <= e2) { cs = pv->Value(GammaEnergy)/GammaEnergy; }
292 else if(GammaEnergy > e2) { cs = pv->Value(e2)/GammaEnergy; }
293
294 return cs;
295
296}
297
298//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
299
300void G4LivermorePolarizedComptonModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
301 const G4MaterialCutsCouple* couple,
302 const G4DynamicParticle* aDynamicGamma,
303 G4double,
304 G4double)
305{
306 // The scattered gamma energy is sampled according to Klein - Nishina formula.
307 // The random number techniques of Butcher & Messel are used (Nuc Phys 20(1960),15).
308 // GEANT4 internal units
309 //
310 // Note : Effects due to binding of atomic electrons are negliged.
311
312 if (verboseLevel > 3)
313 G4cout << "Calling SampleSecondaries() of G4LivermorePolarizedComptonModel" << G4endl;
314
315 G4double gammaEnergy0 = aDynamicGamma->GetKineticEnergy();
316
317 // do nothing below the threshold
318 // should never get here because the XS is zero below the limit
319 if (gammaEnergy0 < LowEnergyLimit())
320 return ;
321
322 G4ThreeVector gammaPolarization0 = aDynamicGamma->GetPolarization();
323
324 // Protection: a polarisation parallel to the
325 // direction causes problems;
326 // in that case find a random polarization
327 G4ThreeVector gammaDirection0 = aDynamicGamma->GetMomentumDirection();
328
329 // Make sure that the polarization vector is perpendicular to the
330 // gamma direction. If not
331 if(!(gammaPolarization0.isOrthogonal(gammaDirection0, 1e-6))||(gammaPolarization0.mag()==0))
332 { // only for testing now
333 gammaPolarization0 = GetRandomPolarization(gammaDirection0);
334 }
335 else
336 {
337 if ( gammaPolarization0.howOrthogonal(gammaDirection0) != 0)
338 {
339 gammaPolarization0 = GetPerpendicularPolarization(gammaDirection0, gammaPolarization0);
340 }
341 }
342 // End of Protection
343
344 G4double E0_m = gammaEnergy0 / electron_mass_c2 ;
345
346 // Select randomly one element in the current material
347 //G4int Z = crossSectionHandler->SelectRandomAtom(couple,gammaEnergy0);
348 const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
349 const G4Element* elm = SelectRandomAtom(couple,particle,gammaEnergy0);
350 G4int Z = (G4int)elm->GetZ();
351
352 // Sample the energy and the polarization of the scattered photon
353 G4double epsilon, epsilonSq, onecost, sinThetaSqr, greject ;
354
355 G4double epsilon0Local = 1./(1. + 2*E0_m);
356 G4double epsilon0Sq = epsilon0Local*epsilon0Local;
357 G4double alpha1 = - G4Log(epsilon0Local);
358 G4double alpha2 = 0.5*(1.- epsilon0Sq);
359
360 G4double wlGamma = h_Planck*c_light/gammaEnergy0;
361 G4double gammaEnergy1;
362 G4ThreeVector gammaDirection1;
363
364 do {
365 if ( alpha1/(alpha1+alpha2) > G4UniformRand() )
366 {
367 epsilon = G4Exp(-alpha1*G4UniformRand());
368 epsilonSq = epsilon*epsilon;
369 }
370 else
371 {
372 epsilonSq = epsilon0Sq + (1.- epsilon0Sq)*G4UniformRand();
373 epsilon = std::sqrt(epsilonSq);
374 }
375
376 onecost = (1.- epsilon)/(epsilon*E0_m);
377 sinThetaSqr = onecost*(2.-onecost);
378
379 // Protection
380 if (sinThetaSqr > 1.)
381 {
382 G4cout
383 << " -- Warning -- G4LivermorePolarizedComptonModel::SampleSecondaries "
384 << "sin(theta)**2 = "
385 << sinThetaSqr
386 << "; set to 1"
387 << G4endl;
388 sinThetaSqr = 1.;
389 }
390 if (sinThetaSqr < 0.)
391 {
392 G4cout
393 << " -- Warning -- G4LivermorePolarizedComptonModel::SampleSecondaries "
394 << "sin(theta)**2 = "
395 << sinThetaSqr
396 << "; set to 0"
397 << G4endl;
398 sinThetaSqr = 0.;
399 }
400 // End protection
401
402 G4double x = std::sqrt(onecost/2.) / (wlGamma/cm);;
403 G4double scatteringFunction = scatterFunctionData->FindValue(x,Z-1);
404 greject = (1. - epsilon*sinThetaSqr/(1.+ epsilonSq))*scatteringFunction;
405
406 } while(greject < G4UniformRand()*Z);
407
408
409 // ****************************************************
410 // Phi determination
411 // ****************************************************
412 G4double phi = SetPhi(epsilon,sinThetaSqr);
413
414 //
415 // scattered gamma angles. ( Z - axis along the parent gamma)
416 //
417 G4double cosTheta = 1. - onecost;
418
419 // Protection
420 if (cosTheta > 1.)
421 {
422 G4cout
423 << " -- Warning -- G4LivermorePolarizedComptonModel::SampleSecondaries "
424 << "cosTheta = "
425 << cosTheta
426 << "; set to 1"
427 << G4endl;
428 cosTheta = 1.;
429 }
430 if (cosTheta < -1.)
431 {
432 G4cout
433 << " -- Warning -- G4LivermorePolarizedComptonModel::SampleSecondaries "
434 << "cosTheta = "
435 << cosTheta
436 << "; set to -1"
437 << G4endl;
438 cosTheta = -1.;
439 }
440 // End protection
441
442 G4double sinTheta = std::sqrt (sinThetaSqr);
443
444 // Protection
445 if (sinTheta > 1.)
446 {
447 G4cout
448 << " -- Warning -- G4LivermorePolarizedComptonModel::SampleSecondaries "
449 << "sinTheta = "
450 << sinTheta
451 << "; set to 1"
452 << G4endl;
453 sinTheta = 1.;
454 }
455 if (sinTheta < -1.)
456 {
457 G4cout
458 << " -- Warning -- G4LivermorePolarizedComptonModel::SampleSecondaries "
459 << "sinTheta = "
460 << sinTheta
461 << "; set to -1"
462 << G4endl;
463 sinTheta = -1.;
464 }
465 // End protection
466
467 // Check for entanglement and re-sample phi if necessary
468
469 const auto* auxInfo
470 = fParticleChange->GetCurrentTrack()->GetAuxiliaryTrackInformation(fEntanglementModelID);
471 if (auxInfo) {
472 const auto* entanglementAuxInfo = dynamic_cast<const G4EntanglementAuxInfo*>(auxInfo);
473 if (entanglementAuxInfo) {
474 auto* clipBoard = dynamic_cast<G4eplusAnnihilationEntanglementClipBoard*>
475 (entanglementAuxInfo->GetEntanglementClipBoard());
476 if (clipBoard) {
477 // This is an entangled photon from eplus annihilation at rest.
478 // If this is the first scatter of the first photon, place theta and
479 // phi on the clipboard.
480 // If this is the first scatter of the second photon, use theta and
481 // phi of the first scatter of the first photon, together with the
482 // theta of the second photon, to sample phi.
483 if (clipBoard->IsTrack1Measurement()) {
484 // Check we have the relevant track. Not sure this is strictly
485 // necessary but I want to be sure tracks from, say, more than one
486 // entangled system are properly paired.
487 // Note: the tracking manager pops the tracks in the reverse order. We
488 // will rely on that. (If not, the logic here would have to be a bit
489 // more complicated to ensure we matched the right tracks.)
490 // So our track 1 is clipboard track B.
491 if (clipBoard->GetTrackB() == fParticleChange->GetCurrentTrack()) {
492 // This is the first scatter of the first photon. Reset flag.
493 // // Debug
494 // auto* track1 = fParticleChange->GetCurrentTrack();
495 // G4cout
496 // << "This is the first scatter of the first photon. Reset flag."
497 // << "\nTrack: " << track1->GetTrackID()
498 // << ", Parent: " << track1->GetParentID()
499 // << ", Name: " << clipBoard->GetParentParticleDefinition()->GetParticleName()
500 // << G4endl;
501 // // End debug
502 clipBoard->ResetTrack1Measurement();
503 // Store cos(theta),phi of first photon.
504 clipBoard->SetComptonCosTheta1(cosTheta);
505 clipBoard->SetComptonPhi1(phi);
506 }
507 } else if (clipBoard->IsTrack2Measurement()) {
508 // Check we have the relevant track.
509 // Remember our track 2 is clipboard track A.
510 if (clipBoard->GetTrackA() == fParticleChange->GetCurrentTrack()) {
511 // This is the first scatter of the second photon. Reset flag.
512 // // Debug
513 // auto* track2 = fParticleChange->GetCurrentTrack();
514 // G4cout
515 // << "This is the first scatter of the second photon. Reset flag."
516 // << "\nTrack: " << track2->GetTrackID()
517 // << ", Parent: " << track2->GetParentID()
518 // << ", Name: " << clipBoard->GetParentParticleDefinition()->GetParticleName()
519 // << G4endl;
520 // // End debug
521 clipBoard->ResetTrack2Measurement();
522
523 // Get cos(theta),phi of first photon.
524 const G4double& cosTheta1 = clipBoard->GetComptonCosTheta1();
525 const G4double& phi1 = clipBoard->GetComptonPhi1();
526 // For clarity make aliases for the current cos(theta),phi.
527 const G4double& cosTheta2 = cosTheta;
528 G4double& phi2 = phi;
529 // G4cout << "cosTheta1,phi1: " << cosTheta1 << ',' << phi1 << G4endl;
530 // G4cout << "cosTheta2,phi2: " << cosTheta2 << ',' << phi2 << G4endl;
531
532 // Re-sample phi
533 // Draw the difference of azimuthal angles, deltaPhi, from
534 // A + B * cos(2*deltaPhi), or rather C + D * cos(2*deltaPhi), where
535 // C = A / (A + |B|) and D = B / (A + |B|), so that maximum is 1.
536 const G4double sin2Theta1 = 1.-cosTheta1*cosTheta1;
537 const G4double sin2Theta2 = 1.-cosTheta2*cosTheta2;
538
539 // Pryce and Ward, Nature No 4065 (1947) p.435.
540 auto* g4Pow = G4Pow::GetInstance();
541 const G4double A =
542 ((g4Pow->powN(1.-cosTheta1,3))+2.)*(g4Pow->powN(1.-cosTheta2,3)+2.)/
543 ((g4Pow->powN(2.-cosTheta1,3)*g4Pow->powN(2.-cosTheta2,3)));
544 const G4double B = -(sin2Theta1*sin2Theta2)/
545 ((g4Pow->powN(2.-cosTheta1,2)*g4Pow->powN(2.-cosTheta2,2)));
546
547 // // Snyder et al, Physical Review 73 (1948) p.440.
548 // // (This is an alternative formulation but result is identical.)
549 // const G4double& k0 = gammaEnergy0;
550 // const G4double k1 = k0/(2.-cosTheta1);
551 // const G4double k2 = k0/(2.-cosTheta2);
552 // const G4double gamma1 = k1/k0+k0/k1;
553 // const G4double gamma2 = k2/k0+k0/k2;
554 // const G4double A1 = gamma1*gamma2-gamma1*sin2Theta2-gamma2*sin2Theta1;
555 // const G4double B1 = 2.*sin2Theta1*sin2Theta2;
556 // // That's A1 + B1*sin2(deltaPhi) = A1 + B1*(0.5*(1.-cos(2.*deltaPhi).
557 // const G4double A = A1 + 0.5*B1;
558 // const G4double B = -0.5*B1;
559
560 const G4double maxValue = A + std::abs(B);
561 const G4double C = A / maxValue;
562 const G4double D = B / maxValue;
563 // G4cout << "A,B,C,D: " << A << ',' << B << ',' << C << ',' << D << G4endl;
564
565 // Sample delta phi
566 G4double deltaPhi;
567 const G4int maxCount = 999999;
568 G4int iCount = 0;
569 for (; iCount < maxCount; ++iCount) {
570 deltaPhi = twopi * G4UniformRand();
571 if (G4UniformRand() < C + D * cos(2.*deltaPhi)) break;
572 }
573 if (iCount >= maxCount ) {
574 G4cout << "G4LivermorePolarizedComptonModel::SampleSecondaries: "
575 << "Re-sampled delta phi not found in " << maxCount
576 << " tries - carrying on anyway." << G4endl;
577 }
578
579 // Thus, the desired second photon azimuth
580 phi2 = deltaPhi - phi1 + halfpi;
581 // The minus sign is in above statement because, since the two
582 // annihilation photons are in opposite directions, their phi's
583 // are measured in the opposite direction.
584 // halfpi is added for the following reason:
585 // In this function phi is relative to the polarisation - see
586 // SystemOfRefChange below. We know from G4eplusAnnihilation that
587 // the polarisations of the two annihilation photons are perpendicular
588 // to each other, i.e., halfpi different.
589 // Furthermore, only sin(phi) and cos(phi) are used below so no
590 // need to place any range constraints.
591 // if (phi2 > pi) {
592 // phi2 -= twopi;
593 // }
594 // if (phi2 < -pi) {
595 // phi2 += twopi;
596 // }
597 }
598 }
599 }
600 }
601 }
602
603 // End of entanglement
604
605 G4double dirx = sinTheta*std::cos(phi);
606 G4double diry = sinTheta*std::sin(phi);
607 G4double dirz = cosTheta ;
608
609 // oneCosT , eom
610
611 // Doppler broadening - Method based on:
612 // Y. Namito, S. Ban and H. Hirayama,
613 // "Implementation of the Doppler Broadening of a Compton-Scattered Photon Into the EGS4 Code"
614 // NIM A 349, pp. 489-494, 1994
615
616 // Maximum number of sampling iterations
617 static G4int maxDopplerIterations = 1000;
618 G4double bindingE = 0.;
619 G4double photonEoriginal = epsilon * gammaEnergy0;
620 G4double photonE = -1.;
621 G4int iteration = 0;
622 G4double eMax = gammaEnergy0;
623
624 G4int shellIdx = 0;
625
626 if (verboseLevel > 3) {
627 G4cout << "Started loop to sample broading" << G4endl;
628 }
629
630 do
631 {
632 iteration++;
633 // Select shell based on shell occupancy
634 shellIdx = shellData->SelectRandomShell(Z);
635 bindingE = shellData->BindingEnergy(Z,shellIdx);
636
637 if (verboseLevel > 3) {
638 G4cout << "Shell ID= " << shellIdx
639 << " Ebind(keV)= " << bindingE/keV << G4endl;
640 }
641 eMax = gammaEnergy0 - bindingE;
642
643 // Randomly sample bound electron momentum (memento: the data set is in Atomic Units)
644 G4double pSample = profileData->RandomSelectMomentum(Z,shellIdx);
645
646 if (verboseLevel > 3) {
647 G4cout << "pSample= " << pSample << G4endl;
648 }
649 // Rescale from atomic units
650 G4double pDoppler = pSample * fine_structure_const;
651 G4double pDoppler2 = pDoppler * pDoppler;
652 G4double var2 = 1. + onecost * E0_m;
653 G4double var3 = var2*var2 - pDoppler2;
654 G4double var4 = var2 - pDoppler2 * cosTheta;
655 G4double var = var4*var4 - var3 + pDoppler2 * var3;
656 if (var > 0.)
657 {
658 G4double varSqrt = std::sqrt(var);
659 G4double scale = gammaEnergy0 / var3;
660 // Random select either root
661 if (G4UniformRand() < 0.5) photonE = (var4 - varSqrt) * scale;
662 else photonE = (var4 + varSqrt) * scale;
663 }
664 else
665 {
666 photonE = -1.;
667 }
668 } while ( iteration <= maxDopplerIterations &&
669 (photonE < 0. || photonE > eMax || photonE < eMax*G4UniformRand()) );
670
671 // End of recalculation of photon energy with Doppler broadening
672 // Revert to original if maximum number of iterations threshold has been reached
673 if (iteration >= maxDopplerIterations)
674 {
675 photonE = photonEoriginal;
676 bindingE = 0.;
677 }
678
679 gammaEnergy1 = photonE;
680
681 //
682 // update G4VParticleChange for the scattered photon
683 //
684 // New polarization
685 G4ThreeVector gammaPolarization1 = SetNewPolarization(epsilon,
686 sinThetaSqr,
687 phi,
688 cosTheta);
689
690 // Set new direction
691 G4ThreeVector tmpDirection1( dirx,diry,dirz );
692 gammaDirection1 = tmpDirection1;
693
694 // Change reference frame.
695 SystemOfRefChange(gammaDirection0,gammaDirection1,
696 gammaPolarization0,gammaPolarization1);
697
698 if (gammaEnergy1 > 0.)
699 {
700 fParticleChange->SetProposedKineticEnergy( gammaEnergy1 ) ;
701 fParticleChange->ProposeMomentumDirection( gammaDirection1 );
702 fParticleChange->ProposePolarization( gammaPolarization1 );
703 }
704 else
705 {
706 gammaEnergy1 = 0.;
707 fParticleChange->SetProposedKineticEnergy(0.) ;
708 fParticleChange->ProposeTrackStatus(fStopAndKill);
709 }
710
711 //
712 // kinematic of the scattered electron
713 //
714 G4double ElecKineEnergy = gammaEnergy0 - gammaEnergy1 -bindingE;
715
716 // SI -protection against negative final energy: no e- is created
717 // like in G4LivermoreComptonModel.cc
718 if(ElecKineEnergy < 0.0) {
719 fParticleChange->ProposeLocalEnergyDeposit(gammaEnergy0 - gammaEnergy1);
720 return;
721 }
722
723 G4double ElecMomentum = std::sqrt(ElecKineEnergy*(ElecKineEnergy+2.*electron_mass_c2));
724
725 G4ThreeVector ElecDirection((gammaEnergy0 * gammaDirection0 -
726 gammaEnergy1 * gammaDirection1) * (1./ElecMomentum));
727
728 G4DynamicParticle* dp =
729 new G4DynamicParticle (G4Electron::Electron(),ElecDirection.unit(),ElecKineEnergy) ;
730 fvect->push_back(dp);
731
732 // sample deexcitation
733 //
734 if (verboseLevel > 3) {
735 G4cout << "Started atomic de-excitation " << fAtomDeexcitation << G4endl;
736 }
737
738 if(fAtomDeexcitation && iteration < maxDopplerIterations) {
739 G4int index = couple->GetIndex();
740 if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
741 std::size_t nbefore = fvect->size();
743 const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
744 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
745 std::size_t nafter = fvect->size();
746 if(nafter > nbefore) {
747 for (std::size_t i=nbefore; i<nafter; ++i) {
748 //Check if there is enough residual energy
749 if (bindingE >= ((*fvect)[i])->GetKineticEnergy())
750 {
751 //Ok, this is a valid secondary: keep it
752 bindingE -= ((*fvect)[i])->GetKineticEnergy();
753 }
754 else
755 {
756 //Invalid secondary: not enough energy to create it!
757 //Keep its energy in the local deposit
758 delete (*fvect)[i];
759 (*fvect)[i]=0;
760 }
761 }
762 }
763 }
764 }
765 //This should never happen
766 if(bindingE < 0.0)
767 G4Exception("G4LivermoreComptonModel::SampleSecondaries()",
768 "em2050",FatalException,"Negative local energy deposit");
769
770 fParticleChange->ProposeLocalEnergyDeposit(bindingE);
771
772}
773
774//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
775
776G4double G4LivermorePolarizedComptonModel::SetPhi(G4double energyRate,
777 G4double sinSqrTh)
778{
779 G4double rand1;
780 G4double rand2;
781 G4double phiProbability;
782 G4double phi;
783 G4double a, b;
784
785 do
786 {
787 rand1 = G4UniformRand();
788 rand2 = G4UniformRand();
789 phiProbability=0.;
790 phi = twopi*rand1;
791
792 a = 2*sinSqrTh;
793 b = energyRate + 1/energyRate;
794
795 phiProbability = 1 - (a/b)*(std::cos(phi)*std::cos(phi));
796 }
797 while ( rand2 > phiProbability );
798 return phi;
799}
800
801
802//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
803
804G4ThreeVector G4LivermorePolarizedComptonModel::SetPerpendicularVector(G4ThreeVector& a)
805{
806 G4double dx = a.x();
807 G4double dy = a.y();
808 G4double dz = a.z();
809 G4double x = dx < 0.0 ? -dx : dx;
810 G4double y = dy < 0.0 ? -dy : dy;
811 G4double z = dz < 0.0 ? -dz : dz;
812 if (x < y) {
813 return x < z ? G4ThreeVector(-dy,dx,0) : G4ThreeVector(0,-dz,dy);
814 }else{
815 return y < z ? G4ThreeVector(dz,0,-dx) : G4ThreeVector(-dy,dx,0);
816 }
817}
818
819//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
820
821G4ThreeVector G4LivermorePolarizedComptonModel::GetRandomPolarization(G4ThreeVector& direction0)
822{
823 G4ThreeVector d0 = direction0.unit();
824 G4ThreeVector a1 = SetPerpendicularVector(d0); //different orthogonal
825 G4ThreeVector a0 = a1.unit(); // unit vector
826
827 G4double rand1 = G4UniformRand();
828
829 G4double angle = twopi*rand1; // random polar angle
830 G4ThreeVector b0 = d0.cross(a0); // cross product
831
833
834 c.setX(std::cos(angle)*(a0.x())+std::sin(angle)*b0.x());
835 c.setY(std::cos(angle)*(a0.y())+std::sin(angle)*b0.y());
836 c.setZ(std::cos(angle)*(a0.z())+std::sin(angle)*b0.z());
837
838 G4ThreeVector c0 = c.unit();
839
840 return c0;
841}
842
843//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
844
845G4ThreeVector G4LivermorePolarizedComptonModel::GetPerpendicularPolarization
846(const G4ThreeVector& gammaDirection, const G4ThreeVector& gammaPolarization) const
847{
848 //
849 // The polarization of a photon is always perpendicular to its momentum direction.
850 // Therefore this function removes those vector component of gammaPolarization, which
851 // points in direction of gammaDirection
852 //
853 // Mathematically we search the projection of the vector a on the plane E, where n is the
854 // plains normal vector.
855 // The basic equation can be found in each geometry book (e.g. Bronstein):
856 // p = a - (a o n)/(n o n)*n
857
858 return gammaPolarization - gammaPolarization.dot(gammaDirection)/gammaDirection.dot(gammaDirection) * gammaDirection;
859}
860
861//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
862
863G4ThreeVector G4LivermorePolarizedComptonModel::SetNewPolarization(G4double epsilon,
864 G4double sinSqrTh,
865 G4double phi,
866 G4double costheta)
867{
868 G4double rand1;
869 G4double rand2;
870 G4double cosPhi = std::cos(phi);
871 G4double sinPhi = std::sin(phi);
872 G4double sinTheta = std::sqrt(sinSqrTh);
873 G4double cosSqrPhi = cosPhi*cosPhi;
874 // G4double cossqrth = 1.-sinSqrTh;
875 // G4double sinsqrphi = sinPhi*sinPhi;
876 G4double normalisation = std::sqrt(1. - cosSqrPhi*sinSqrTh);
877
878 // Determination of Theta
879 G4double theta;
880
881 // Dan Xu method (IEEE TNS, 52, 1160 (2005))
882 rand1 = G4UniformRand();
883 rand2 = G4UniformRand();
884
885 if (rand1<(epsilon+1.0/epsilon-2)/(2.0*(epsilon+1.0/epsilon)-4.0*sinSqrTh*cosSqrPhi))
886 {
887 if (rand2<0.5)
888 theta = pi/2.0;
889 else
890 theta = 3.0*pi/2.0;
891 }
892 else
893 {
894 if (rand2<0.5)
895 theta = 0;
896 else
897 theta = pi;
898 }
899 G4double cosBeta = std::cos(theta);
900 G4double sinBeta = std::sqrt(1-cosBeta*cosBeta);
901
902 G4ThreeVector gammaPolarization1;
903
904 G4double xParallel = normalisation*cosBeta;
905 G4double yParallel = -(sinSqrTh*cosPhi*sinPhi)*cosBeta/normalisation;
906 G4double zParallel = -(costheta*sinTheta*cosPhi)*cosBeta/normalisation;
907 G4double xPerpendicular = 0.;
908 G4double yPerpendicular = (costheta)*sinBeta/normalisation;
909 G4double zPerpendicular = -(sinTheta*sinPhi)*sinBeta/normalisation;
910
911 G4double xTotal = (xParallel + xPerpendicular);
912 G4double yTotal = (yParallel + yPerpendicular);
913 G4double zTotal = (zParallel + zPerpendicular);
914
915 gammaPolarization1.setX(xTotal);
916 gammaPolarization1.setY(yTotal);
917 gammaPolarization1.setZ(zTotal);
918
919 return gammaPolarization1;
920}
921
922//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
923
924void G4LivermorePolarizedComptonModel::SystemOfRefChange(G4ThreeVector& direction0,
925 G4ThreeVector& direction1,
926 G4ThreeVector& polarization0,
927 G4ThreeVector& polarization1)
928{
929 // direction0 is the original photon direction ---> z
930 // polarization0 is the original photon polarization ---> x
931 // need to specify y axis in the real reference frame ---> y
932 G4ThreeVector Axis_Z0 = direction0.unit();
933 G4ThreeVector Axis_X0 = polarization0.unit();
934 G4ThreeVector Axis_Y0 = (Axis_Z0.cross(Axis_X0)).unit(); // to be confirmed;
935
936 G4double direction_x = direction1.getX();
937 G4double direction_y = direction1.getY();
938 G4double direction_z = direction1.getZ();
939
940 direction1 = (direction_x*Axis_X0 + direction_y*Axis_Y0 + direction_z*Axis_Z0).unit();
941 G4double polarization_x = polarization1.getX();
942 G4double polarization_y = polarization1.getY();
943 G4double polarization_z = polarization1.getZ();
944
945 polarization1 = (polarization_x*Axis_X0 + polarization_y*Axis_Y0 + polarization_z*Axis_Z0).unit();
946
947}
948
949
950//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
951//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
952
953void
955 G4int Z)
956{
957 G4AutoLock l(&LivermorePolarizedComptonModelMutex);
958 if(!data[Z]) { ReadData(Z); }
959 l.unlock();
960}
G4double epsilon(G4double density, G4double temperature)
G4double C(G4double temp)
G4double B(G4double temperature)
G4double D(G4double temp)
std::vector< const G4Element * > G4ElementVector
const char * G4FindDataDir(const char *)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:180
const G4double a0
G4double G4Log(G4double x)
Definition G4Log.hh:227
#define G4MUTEX_INITIALIZER
std::mutex G4Mutex
CLHEP::Hep3Vector G4ThreeVector
@ fStopAndKill
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
const G4double A[17]
const G4double alpha2
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
double z() const
Hep3Vector unit() const
double getZ() const
double x() const
void setY(double)
double y() const
Hep3Vector cross(const Hep3Vector &) const
double dot(const Hep3Vector &) const
void setZ(double)
bool isOrthogonal(const Hep3Vector &v, double epsilon=tolerance) const
double mag() const
double getX() const
double howOrthogonal(const Hep3Vector &v) const
void setX(double)
double getY() const
virtual G4double FindValue(G4double x, G4int componentId=0) const
virtual G4bool LoadData(const G4String &fileName)
G4double RandomSelectMomentum(G4int Z, G4int shellIndex) const
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4ThreeVector & GetPolarization() const
static G4Electron * Electron()
Definition G4Electron.cc:91
G4double GetZ() const
Definition G4Element.hh:119
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void InitialiseForElement(const G4ParticleDefinition *, G4int Z) override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
G4LivermorePolarizedComptonModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="LivermorePolarizedCompton")
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
std::size_t GetNumberOfElements() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposePolarization(const G4ThreeVector &dir)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
static G4int GetModelID(const G4int modelIndex)
void ScaleVector(const G4double factorE, const G4double factorV)
G4double Energy(const std::size_t index) const
G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
G4double Value(const G4double energy, std::size_t &lastidx) const
std::size_t GetVectorLength() const
static G4Pow * GetInstance()
Definition G4Pow.cc:41
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
void SetOccupancyData()
G4double BindingEnergy(G4int Z, G4int shellIndex) const
void LoadData(const G4String &fileName)
G4int SelectRandomShell(G4int Z) const
G4VAuxiliaryTrackInformation * GetAuxiliaryTrackInformation(G4int id) const
Definition G4Track.cc:229
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double LowEnergyLimit() const
std::vector< G4EmElementSelector * > * GetElementSelectors()
G4bool IsMaster() const
G4double HighEnergyLimit() const
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
void SetDeexcitationFlag(G4bool val)
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
void ProposeTrackStatus(G4TrackStatus status)
const G4Track * GetCurrentTrack() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
int G4lrint(double ad)
Definition templates.hh:134