LDMX Software
simcore::GammaConversionToFCPs Class Reference

Discrete process for gamma -> fcp+ fcp- conversion in a nuclear field. More...

#include <GammaConversionToFCPs.h>

Public Member Functions

 GammaConversionToFCPs (const G4String &processName="GammaToFCPPair", G4ProcessType type=fElectromagnetic)
 Constructor.
 
 ~GammaConversionToFCPs ()=default
 Destructor Just the default one.
 
G4bool IsApplicable (const G4ParticleDefinition &particle) override
 Check if this process applies to the given particle.
 
void BuildPhysicsTable (const G4ParticleDefinition &p) override
 Build physics table (called during initialization).
 
void printInfoDefinition ()
 Print process information.
 
void setCrossSecFactor (G4double fac)
 Set a factor to artificially scale the cross section.
 
G4double getCrossSecFactor () const
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override
 Return the mean free path of the gamma in the current material.
 
G4double getCrossSectionPerAtom (const G4DynamicParticle *aDynamicGamma, const G4Element *anElement)
 Return the total cross section per atom.
 
G4VParticleChange * PostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 Generate the final state: gamma -> fcp+ fcp-.
 
G4double computeCrossSectionPerAtom (G4double GammaEnergy, G4int Z)
 Compute the cross section per atom for the given gamma energy and Z.
 
G4double computeMeanFreePath (G4double GammaEnergy, const G4Material *aMaterial)
 Compute the mean free path for the given energy and material.
 
 GammaConversionToFCPs (const GammaConversionToFCPs &)=delete
 
GammaConversionToFCPsoperator= (const GammaConversionToFCPs &)=delete
 

Static Public Attributes

static const std::string PROCESS_NAME = "GammaToFCPPair"
 The name of this process.
 

Private Member Functions

const G4Element * selectRandomAtom (const G4DynamicParticle *aDynamicGamma, const G4Material *aMaterial)
 Select a random atom from the material, weighted by cross section.
 
 enableLogging ("GammaConversionToFCPs")
 Enable logging for this class.
 

Private Attributes

G4double mfcp_
 fcp mass
 
G4double rc_
 Coupling-scaled reduced Compton wavelength: q^2 * elm_coupling / Mfcp.
 
G4double limit_energy_
 Energy limit for accurate cross section (5 * Mfcp)
 
G4double lowest_energy_limit_
 Absolute low energy limit of the model (2 * Mfcp)
 
G4double highest_energy_limit_
 High energy limit of the model.
 
G4double cross_sec_factor_ = 1.0
 Factor to artificially scale the cross section.
 
const G4ParticleDefinition * the_fcp_minus_
 Pointer to fcp- definition.
 
const G4ParticleDefinition * the_fcp_plus_
 Pointer to fcp+ definition.
 
std::vector< G4double > temp_
 Temporary vector for element selection.
 

Detailed Description

Discrete process for gamma -> fcp+ fcp- conversion in a nuclear field.

This process is the SM-photon analogue of G4GammaConversionToMuons (gamma -> mu+ mu-). An SM photon converts into a pair of massive fractionally charged particles (fcp+ fcp-) in the electromagnetic field of a nucleus.

Definition at line 31 of file GammaConversionToFCPs.h.

Constructor & Destructor Documentation

◆ GammaConversionToFCPs()

simcore::GammaConversionToFCPs::GammaConversionToFCPs ( const G4String & processName = "GammaToFCPPair",
G4ProcessType type = fElectromagnetic )
explicit

Constructor.

Parameters
[in]processNamename of this process (default: "GammaToFCPPair")
[in]typeGeant4 process type

Definition at line 22 of file GammaConversionToFCPs.cxx.

24 : G4VDiscreteProcess(processName, type) {
25 // Get the fcp particle definitions
26 auto* fcp_minus = G4FractionallyCharged::FractionallyCharged();
27 mfcp_ = fcp_minus->GetPDGMass();
28
29 // Get fcp charge magnitude in units of e
30 G4double qfcp = std::abs(fcp_minus->GetPDGCharge() / CLHEP::eplus);
31
32 // Coupling-scaled reduced Compton wavelength
33 // rc_ = q^2 * elm_coupling / mfcp_
34 // This gives the correct q^4 scaling in sigfac = 4*alpha*Z^2*rc_^2
35 rc_ = qfcp * qfcp * CLHEP::elm_coupling / mfcp_;
36
37 limit_energy_ = 5. * mfcp_;
39 highest_energy_limit_ = 1e12 * CLHEP::GeV;
40 ldmx_log(debug) << "Energy limits: lowest="
41 << lowest_energy_limit_ / CLHEP::MeV
42 << " MeV, limit=" << limit_energy_ / CLHEP::MeV
43 << " MeV, highest=" << highest_energy_limit_ / CLHEP::GeV
44 << " GeV";
45
46 the_fcp_minus_ = fcp_minus;
47 // Look up the antiparticle (fcp+) from the particle table
48 the_fcp_plus_ = G4ParticleTable::GetParticleTable()->FindParticle("fcp+");
49
50 ldmx_log(debug) << "FCP+ particle: "
51 << (the_fcp_plus_ ? the_fcp_plus_->GetParticleName()
52 : "NULL");
53 ldmx_log(debug) << "FCP- particle: "
54 << (the_fcp_minus_ ? the_fcp_minus_->GetParticleName()
55 : "NULL");
56}
G4double limit_energy_
Energy limit for accurate cross section (5 * Mfcp)
const G4ParticleDefinition * the_fcp_plus_
Pointer to fcp+ definition.
G4double lowest_energy_limit_
Absolute low energy limit of the model (2 * Mfcp)
G4double rc_
Coupling-scaled reduced Compton wavelength: q^2 * elm_coupling / Mfcp.
G4double highest_energy_limit_
High energy limit of the model.
const G4ParticleDefinition * the_fcp_minus_
Pointer to fcp- definition.

References highest_energy_limit_, limit_energy_, lowest_energy_limit_, mfcp_, rc_, the_fcp_minus_, and the_fcp_plus_.

Member Function Documentation

◆ BuildPhysicsTable()

void simcore::GammaConversionToFCPs::BuildPhysicsTable ( const G4ParticleDefinition & p)
override

Build physics table (called during initialization).

Sets up the temp vector for element selection.

Definition at line 67 of file GammaConversionToFCPs.cxx.

67 {
68 ldmx_log(debug) << "=== BuildPhysicsTable for particle: "
69 << p.GetParticleName() << " ===";
70 auto table = G4Material::GetMaterialTable();
71 std::size_t nelm = 0;
72 for (auto const& mat : *table) {
73 std::size_t n = mat->GetNumberOfElements();
74 nelm = std::max(nelm, n);
75 ldmx_log(debug) << " Material: " << mat->GetName() << " with " << n
76 << " elements";
77 }
78 temp_.resize(nelm, 0);
79 ldmx_log(debug) << "Max elements in any material: " << nelm;
81}
std::vector< G4double > temp_
Temporary vector for element selection.
void printInfoDefinition()
Print process information.

References printInfoDefinition(), and temp_.

◆ computeCrossSectionPerAtom()

G4double simcore::GammaConversionToFCPs::computeCrossSectionPerAtom ( G4double GammaEnergy,
G4int Z )

Compute the cross section per atom for the given gamma energy and Z.

Definition at line 132 of file GammaConversionToFCPs.cxx.

133 {
134 // Cross section parametrization adapted from H. Burkhardt
135 // (G4GammaConversionToMuons) with Mmuon -> mfcp_ and
136 // coupling scaled by q_fcp^4 through rc_
137 if (Egam <= lowest_energy_limit_) {
138 return 0.0;
139 }
140
141 G4NistManager* nist = G4NistManager::Instance();
142
143 G4double pow_thres, ecor, b, dn, zthird, winfty, w_med_appr, wsatur, sigfac;
144
145 if (Z == 1) { // special case of Hydrogen
146 b = 202.4;
147 dn = 1.49;
148 } else {
149 b = 183.;
150 dn = 1.54 * nist->GetA27(Z);
151 }
152 zthird = 1. / nist->GetZ13(Z); // Z^(-1/3)
153 winfty = b * zthird * mfcp_ / (dn * electron_mass_c2);
154 w_med_appr = 1. / (4. * dn * SQRTE * mfcp_);
155 wsatur = winfty / w_med_appr;
156 sigfac = 4. * fine_structure_const * Z * Z * rc_ * rc_;
157 pow_thres = 1.479 + 0.00799 * dn;
158 ecor = -18. + 4347. / (b * zthird);
159
160 G4double cor_fuc = 1. + .04 * G4Log(1. + ecor / Egam);
161 G4double eg = G4Exp(G4Log(1. - 4. * mfcp_ / Egam) * pow_thres) *
162 G4Exp(G4Log(G4Exp(G4Log(wsatur) * POW_SAT) +
163 G4Exp(G4Log(Egam) * POW_SAT)) /
164 POW_SAT);
165
166 G4double cross_section =
167 7. / 9. * sigfac * G4Log(1. + w_med_appr * cor_fuc * eg);
168 cross_section *= cross_sec_factor_;
169 return cross_section;
170}
G4double cross_sec_factor_
Factor to artificially scale the cross section.

References cross_sec_factor_, lowest_energy_limit_, mfcp_, and rc_.

Referenced by computeMeanFreePath(), getCrossSectionPerAtom(), and selectRandomAtom().

◆ computeMeanFreePath()

G4double simcore::GammaConversionToFCPs::computeMeanFreePath ( G4double GammaEnergy,
const G4Material * aMaterial )

Compute the mean free path for the given energy and material.

Definition at line 96 of file GammaConversionToFCPs.cxx.

97 {
98 if (GammaEnergy <= lowest_energy_limit_) {
99 return DBL_MAX;
100 }
101
102 const G4ElementVector* the_element_vector = aMaterial->GetElementVector();
103 const G4double* nb_of_atoms_per_volume =
104 aMaterial->GetVecNbOfAtomsPerVolume();
105
106 G4double sigma = 0.0;
107 G4double fact = 1.0;
108 G4double e = GammaEnergy;
109
110 // Low energy approximation as in Bethe-Heitler model
111 if (e < limit_energy_) {
112 G4double y =
114 fact = y * y;
115 e = limit_energy_;
116 }
117
118 for (std::size_t i = 0; i < aMaterial->GetNumberOfElements(); ++i) {
119 sigma += nb_of_atoms_per_volume[i] * fact *
121 e, G4lrint((*the_element_vector)[i]->GetZ()));
122 }
123 return (sigma > 0.0) ? 1. / sigma : DBL_MAX;
124}
G4double computeCrossSectionPerAtom(G4double GammaEnergy, G4int Z)
Compute the cross section per atom for the given gamma energy and Z.

References computeCrossSectionPerAtom(), limit_energy_, and lowest_energy_limit_.

Referenced by GetMeanFreePath().

◆ getCrossSecFactor()

G4double simcore::GammaConversionToFCPs::getCrossSecFactor ( ) const
inline
Returns
the current cross section scaling factor

Definition at line 78 of file GammaConversionToFCPs.h.

78{ return cross_sec_factor_; }

References cross_sec_factor_.

◆ getCrossSectionPerAtom()

G4double simcore::GammaConversionToFCPs::getCrossSectionPerAtom ( const G4DynamicParticle * aDynamicGamma,
const G4Element * anElement )

Return the total cross section per atom.

Definition at line 126 of file GammaConversionToFCPs.cxx.

127 {
128 return computeCrossSectionPerAtom(aDynamicGamma->GetKineticEnergy(),
129 G4lrint(anElement->GetZ()));
130}

References computeCrossSectionPerAtom().

◆ GetMeanFreePath()

G4double simcore::GammaConversionToFCPs::GetMeanFreePath ( const G4Track & aTrack,
G4double previousStepSize,
G4ForceCondition * condition )
override

Return the mean free path of the gamma in the current material.

Definition at line 83 of file GammaConversionToFCPs.cxx.

84 {
85 const G4DynamicParticle* a_dynamic_gamma = aTrack.GetDynamicParticle();
86 G4double gamma_energy = a_dynamic_gamma->GetKineticEnergy();
87 const G4Material* a_material = aTrack.GetMaterial();
88 G4double mfp = computeMeanFreePath(gamma_energy, a_material);
89 ldmx_log(trace) << "GetMeanFreePath: gamma KE=" << gamma_energy / CLHEP::MeV
90 << " MeV in " << a_material->GetName()
91 << " -> MFP=" << (mfp < DBL_MAX ? mfp / CLHEP::mm : -1)
92 << " mm";
93 return mfp;
94}
G4double computeMeanFreePath(G4double GammaEnergy, const G4Material *aMaterial)
Compute the mean free path for the given energy and material.

References computeMeanFreePath().

◆ IsApplicable()

G4bool simcore::GammaConversionToFCPs::IsApplicable ( const G4ParticleDefinition & particle)
override

Check if this process applies to the given particle.

Returns
true only for the SM photon (gamma)

Definition at line 58 of file GammaConversionToFCPs.cxx.

59 {
60 bool applicable = (&particle == G4Gamma::Gamma());
61 ldmx_log(debug) << "IsApplicable called for particle: "
62 << particle.GetParticleName() << " -> "
63 << (applicable ? "YES (this is gamma)" : "NO");
64 return applicable;
65}

◆ PostStepDoIt()

G4VParticleChange * simcore::GammaConversionToFCPs::PostStepDoIt ( const G4Track & aTrack,
const G4Step & aStep )
override

Generate the final state: gamma -> fcp+ fcp-.

Definition at line 179 of file GammaConversionToFCPs.cxx.

180 {
181 ldmx_log(debug) << "gamma -> fcp+ fcp- conversion starting!";
182 aParticleChange.Initialize(aTrack);
183 const G4Material* a_material = aTrack.GetMaterial();
184
185 ldmx_log(debug) << "Conversion in material: " << a_material->GetName();
186 ldmx_log(debug) << "Gamma Track ID: " << aTrack.GetTrackID()
187 << ", Parent ID: " << aTrack.GetParentID();
188
189 // Current gamma energy and direction
190 const G4DynamicParticle* a_dynamic_gamma = aTrack.GetDynamicParticle();
191 G4double egam = a_dynamic_gamma->GetKineticEnergy();
192
193 G4ThreeVector pos = aTrack.GetPosition();
194 ldmx_log(debug) << "Gamma kinetic energy: " << egam / CLHEP::MeV << " MeV";
195 ldmx_log(debug) << "Gamma position: (" << pos.x() / CLHEP::mm << ", "
196 << pos.y() / CLHEP::mm << ", " << pos.z() / CLHEP::mm
197 << ") mm";
198 ldmx_log(debug) << "Gamma momentum direction: "
199 << a_dynamic_gamma->GetMomentumDirection();
200
201 if (egam <= lowest_energy_limit_) {
202 ldmx_log(warn) << "Gamma energy (" << egam / CLHEP::MeV
203 << " MeV) below threshold ("
204 << lowest_energy_limit_ / CLHEP::MeV
205 << " MeV), skipping conversion!";
206 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
207 }
208
209 // Kill the incident gamma
210 ldmx_log(debug) << "Killing incident gamma and creating fcp+/fcp- pair...";
211 aParticleChange.ProposeMomentumDirection(0., 0., 0.);
212 aParticleChange.ProposeEnergy(0.);
213 aParticleChange.ProposeTrackStatus(fStopAndKill);
214
215 G4ParticleMomentum gamma_direction = a_dynamic_gamma->GetMomentumDirection();
216
217 // Select randomly one element constituting the material
218 const G4Element* an_element = selectRandomAtom(a_dynamic_gamma, a_material);
219 G4int z = G4lrint(an_element->GetZ());
220
221 G4NistManager* nist = G4NistManager::Instance();
222
223 G4double b, dn;
224 G4double a027 = nist->GetA27(z);
225
226 if (z == 1) {
227 b = 202.4;
228 dn = 1.49;
229 } else {
230 b = 183.;
231 dn = 1.54 * a027;
232 }
233 G4double zthird = 1. / nist->GetZ13(z);
234 G4double winfty = b * zthird * mfcp_ / (dn * electron_mass_c2);
235
236 G4double c1_num = 0.138 * a027;
237 G4double c1_num2 = c1_num * c1_num;
238 G4double c2_term2 = electron_mass_c2 / (183. * zthird * mfcp_);
239
240 G4double gamma_fcp_inv = mfcp_ / egam;
241
242 // Generate xPlus according to the differential cross section by rejection
243 G4double xmin =
244 (egam <= limit_energy_) ? 0.5 : 0.5 - std::sqrt(0.25 - gamma_fcp_inv);
245 G4double xmax = 1. - xmin;
246
247 G4double ds2 = (dn * SQRTE - 2.);
248 G4double s_bz = SQRTE * b * zthird / electron_mass_c2;
249 G4double log_wmax_inv = 1. / G4Log(winfty * (1. + 2. * ds2 * gamma_fcp_inv) /
250 (1. + 2. * s_bz * mfcp_ * gamma_fcp_inv));
251 G4double x_plus = 0.5;
252 G4double x_minus = 0.5;
253 G4double x_pm = 0.25;
254
255 G4int nn = 0;
256 const G4int nmax = 1000;
257
258 // Sampling for Egam > limit_energy_
259 if (xmin < 0.5) {
260 G4double result, w;
261 do {
262 x_plus = xmin + G4UniformRand() * (xmax - xmin);
263 x_minus = 1. - x_plus;
264 x_pm = x_plus * x_minus;
265 G4double del = mfcp_ * mfcp_ / (2. * egam * x_pm);
266 w = winfty * (1. + ds2 * del / mfcp_) / (1. + s_bz * del);
267 G4double xxp = 1. - 4. / 3. * x_pm;
268 result = (xxp > 0.) ? xxp * G4Log(w) * log_wmax_inv : 0.0;
269 if (result > 1.0) {
270 ldmx_log(warn) << "GammaConversionToFCPs::PostStepDoIt WARNING:"
271 << " in dSigxPlusGen, result=" << result << " > 1";
272 }
273 ++nn;
274 if (nn >= nmax) {
275 break;
276 }
277 } while (G4UniformRand() > result);
278 }
279
280 // Generate angular variables via auxiliary variables t, psi, rho
281 G4double t;
282 G4double psi;
283 G4double rho;
284
285 G4double a3 = (gamma_fcp_inv / (2. * x_pm));
286 G4double a33 = a3 * a3;
287 G4double f1;
288 G4double b1 = 1. / (4. * c1_num2);
289 G4double b3 = b1 * b1 * b1;
290 G4double a21 = a33 + b1;
291
292 G4double f1_max =
293 -(1. - x_pm) * (2. * b1 + (a21 + a33) * G4Log(a33 / a21)) / (2 * b3);
294
295 G4double theta_plus, theta_minus, phi_half;
296 nn = 0;
297
298 // t, psi, rho generation (while angle < pi)
299 do {
300 // Generate t by rejection method
301 do {
302 ++nn;
303 t = G4UniformRand();
304 G4double a34 = a33 / (t * t);
305 G4double a22 = a34 + b1;
306 if (std::abs(b1) < 0.0001 * a34) {
307 f1 = (1. - 2. * x_pm + 4. * x_pm * t * (1. - t)) /
308 (12. * a34 * a34 * a34 * a34);
309 } else {
310 f1 = -(1. - 2. * x_pm + 4. * x_pm * t * (1. - t)) *
311 (2. * b1 + (a22 + a34) * G4Log(a34 / a22)) / (2 * b3);
312 }
313 if (f1 < 0.0 || f1 > f1_max) {
314 f1 = 0.0;
315 }
316 if (nn > nmax) {
317 break;
318 }
319 } while (G4UniformRand() * f1_max > f1);
320
321 // Generate psi by rejection method
322 G4double f2_max = 1. - 2. * x_pm * (1. - 4. * t * (1. - t));
323 G4double f2;
324 do {
325 ++nn;
326 psi = CLHEP::twopi * G4UniformRand();
327 f2 =
328 1. - 2. * x_pm + 4. * x_pm * t * (1. - t) * (1. + std::cos(2. * psi));
329 if (f2 < 0 || f2 > f2_max) {
330 f2 = 0.0;
331 }
332 if (nn >= nmax) {
333 break;
334 }
335 } while (G4UniformRand() * f2_max > f2);
336
337 // Generate rho by direct transformation
338 G4double c2_term1 = gamma_fcp_inv / (2. * x_pm * t);
339 G4double c22 = c2_term1 * c2_term1 + c2_term2 * c2_term2;
340 G4double c2 = 4. * c22 * c22 / std::sqrt(x_pm);
341 G4double rhomax = (1. / t - 1.) * 1.9 / a027;
342 G4double beta = G4Log((c2 + rhomax * rhomax * rhomax * rhomax) / c2);
343 rho = G4Exp(G4Log(c2 * (G4Exp(beta * G4UniformRand()) - 1.)) * 0.25);
344
345 // Get kinematical variables from t and psi
346 G4double u = std::sqrt(1. / t - 1.);
347 G4double xi_half = 0.5 * rho * std::cos(psi);
348 phi_half = 0.5 * rho / u * std::sin(psi);
349
350 theta_plus = gamma_fcp_inv * (u + xi_half) / x_plus;
351 theta_minus = gamma_fcp_inv * (u - xi_half) / x_minus;
352
353 // Protection against infinite loop
354 if (nn > nmax) {
355 if (std::abs(theta_plus) > CLHEP::pi) {
356 theta_plus = 0.0;
357 }
358 if (std::abs(theta_minus) > CLHEP::pi) {
359 theta_minus = 0.0;
360 }
361 }
362 } while (std::abs(theta_plus) > CLHEP::pi ||
363 std::abs(theta_minus) > CLHEP::pi);
364
365 // Construct the momentum vectors
366 // Azimuthal symmetry: take phi0 random between 0 and 2pi
367 G4double phi0 = CLHEP::twopi * G4UniformRand();
368 G4double e_plus = x_plus * egam;
369 G4double e_minus = x_minus * egam;
370
371 // fcp+ fcp- directions for gamma in z-direction
372 G4ThreeVector fcp_plus_direction(
373 std::sin(theta_plus) * std::cos(phi0 + phi_half),
374 std::sin(theta_plus) * std::sin(phi0 + phi_half), std::cos(theta_plus));
375 G4ThreeVector fcp_minus_direction(
376 -std::sin(theta_minus) * std::cos(phi0 - phi_half),
377 -std::sin(theta_minus) * std::sin(phi0 - phi_half),
378 std::cos(theta_minus));
379
380 // Rotate to actual gamma direction
381 fcp_plus_direction.rotateUz(gamma_direction);
382 fcp_minus_direction.rotateUz(gamma_direction);
383
384 // Create fcp+ secondary
385 auto a_particle1 =
386 new G4DynamicParticle(the_fcp_plus_, fcp_plus_direction, e_plus - mfcp_);
387 aParticleChange.AddSecondary(a_particle1);
388
389 ldmx_log(debug) << "Created FCP+ with:";
390 ldmx_log(debug) << " - Kinetic energy: " << (e_plus - mfcp_) / CLHEP::MeV
391 << " MeV";
392 ldmx_log(debug) << " - Total energy: " << e_plus / CLHEP::MeV << " MeV";
393 ldmx_log(debug) << " - Direction: " << fcp_plus_direction;
394
395 // Create fcp- secondary
396 auto a_particle2 = new G4DynamicParticle(the_fcp_minus_, fcp_minus_direction,
397 e_minus - mfcp_);
398 aParticleChange.AddSecondary(a_particle2);
399
400 ldmx_log(debug) << "Created FCP- with:";
401 ldmx_log(debug) << " - Kinetic energy: " << (e_minus - mfcp_) / CLHEP::MeV
402 << " MeV";
403 ldmx_log(debug) << " - Total energy: " << e_minus / CLHEP::MeV << " MeV";
404 ldmx_log(debug) << " - Direction: " << fcp_minus_direction;
405
406 ldmx_log(debug) << "=== gamma -> fcp+ fcp- CONVERSION COMPLETE ===";
407
408 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
409} // end of PostStepDoIt
const G4Element * selectRandomAtom(const G4DynamicParticle *aDynamicGamma, const G4Material *aMaterial)
Select a random atom from the material, weighted by cross section.

References limit_energy_, lowest_energy_limit_, mfcp_, selectRandomAtom(), the_fcp_minus_, and the_fcp_plus_.

◆ printInfoDefinition()

void simcore::GammaConversionToFCPs::printInfoDefinition ( )

Print process information.

Definition at line 438 of file GammaConversionToFCPs.cxx.

438 {
439 G4cout << G4endl << GetProcessName() << ": "
440 << "gamma -> fcp+ fcp- Bethe-Heitler-like process" << G4endl;
441 G4cout << " good cross section parametrization from "
442 << G4BestUnit(lowest_energy_limit_, "Energy") << " to "
443 << highest_energy_limit_ / GeV << " GeV for all Z." << G4endl;
444 G4cout << " cross section factor: " << cross_sec_factor_ << G4endl;
445}

References cross_sec_factor_, highest_energy_limit_, and lowest_energy_limit_.

Referenced by BuildPhysicsTable().

◆ selectRandomAtom()

const G4Element * simcore::GammaConversionToFCPs::selectRandomAtom ( const G4DynamicParticle * aDynamicGamma,
const G4Material * aMaterial )
private

Select a random atom from the material, weighted by cross section.

Definition at line 411 of file GammaConversionToFCPs.cxx.

412 {
413 const std::size_t number_of_elements = aMaterial->GetNumberOfElements();
414 const G4ElementVector* the_element_vector = aMaterial->GetElementVector();
415 const G4Element* elm = (*the_element_vector)[0];
416
417 if (number_of_elements > 1) {
418 G4double e = std::max(aDynamicGamma->GetKineticEnergy(), limit_energy_);
419 const G4double* natom = aMaterial->GetVecNbOfAtomsPerVolume();
420
421 G4double sum = 0.;
422 for (std::size_t i = 0; i < number_of_elements; ++i) {
423 elm = (*the_element_vector)[i];
424 sum += natom[i] * computeCrossSectionPerAtom(e, G4lrint(elm->GetZ()));
425 temp_[i] = sum;
426 }
427 sum *= G4UniformRand();
428 for (std::size_t i = 0; i < number_of_elements; ++i) {
429 if (sum <= temp_[i]) {
430 elm = (*the_element_vector)[i];
431 break;
432 }
433 }
434 }
435 return elm;
436}

References computeCrossSectionPerAtom(), limit_energy_, and temp_.

Referenced by PostStepDoIt().

◆ setCrossSecFactor()

void simcore::GammaConversionToFCPs::setCrossSecFactor ( G4double fac)

Set a factor to artificially scale the cross section.

Parameters
[in]facmultiplicative factor (default 1.0)

Definition at line 172 of file GammaConversionToFCPs.cxx.

172 {
173 if (fac < 0.0) return;
174 cross_sec_factor_ = fac;
175 ldmx_log(info) << "GammaConversionToFCPs cross section factor set to "
177}

References cross_sec_factor_.

Member Data Documentation

◆ cross_sec_factor_

G4double simcore::GammaConversionToFCPs::cross_sec_factor_ = 1.0
private

Factor to artificially scale the cross section.

Definition at line 131 of file GammaConversionToFCPs.h.

Referenced by computeCrossSectionPerAtom(), getCrossSecFactor(), printInfoDefinition(), and setCrossSecFactor().

◆ highest_energy_limit_

G4double simcore::GammaConversionToFCPs::highest_energy_limit_
private

High energy limit of the model.

Definition at line 129 of file GammaConversionToFCPs.h.

Referenced by GammaConversionToFCPs(), and printInfoDefinition().

◆ limit_energy_

G4double simcore::GammaConversionToFCPs::limit_energy_
private

Energy limit for accurate cross section (5 * Mfcp)

Definition at line 125 of file GammaConversionToFCPs.h.

Referenced by computeMeanFreePath(), GammaConversionToFCPs(), PostStepDoIt(), and selectRandomAtom().

◆ lowest_energy_limit_

G4double simcore::GammaConversionToFCPs::lowest_energy_limit_
private

Absolute low energy limit of the model (2 * Mfcp)

Definition at line 127 of file GammaConversionToFCPs.h.

Referenced by computeCrossSectionPerAtom(), computeMeanFreePath(), GammaConversionToFCPs(), PostStepDoIt(), and printInfoDefinition().

◆ mfcp_

G4double simcore::GammaConversionToFCPs::mfcp_
private

fcp mass

Definition at line 121 of file GammaConversionToFCPs.h.

Referenced by computeCrossSectionPerAtom(), GammaConversionToFCPs(), and PostStepDoIt().

◆ PROCESS_NAME

const std::string simcore::GammaConversionToFCPs::PROCESS_NAME = "GammaToFCPPair"
static

The name of this process.

Used to attach biasing operators to this process.

Definition at line 37 of file GammaConversionToFCPs.h.

◆ rc_

G4double simcore::GammaConversionToFCPs::rc_
private

Coupling-scaled reduced Compton wavelength: q^2 * elm_coupling / Mfcp.

Definition at line 123 of file GammaConversionToFCPs.h.

Referenced by computeCrossSectionPerAtom(), and GammaConversionToFCPs().

◆ temp_

std::vector<G4double> simcore::GammaConversionToFCPs::temp_
private

Temporary vector for element selection.

Definition at line 139 of file GammaConversionToFCPs.h.

Referenced by BuildPhysicsTable(), and selectRandomAtom().

◆ the_fcp_minus_

const G4ParticleDefinition* simcore::GammaConversionToFCPs::the_fcp_minus_
private

Pointer to fcp- definition.

Definition at line 134 of file GammaConversionToFCPs.h.

Referenced by GammaConversionToFCPs(), and PostStepDoIt().

◆ the_fcp_plus_

const G4ParticleDefinition* simcore::GammaConversionToFCPs::the_fcp_plus_
private

Pointer to fcp+ definition.

Definition at line 136 of file GammaConversionToFCPs.h.

Referenced by GammaConversionToFCPs(), and PostStepDoIt().


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