Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmExtraParameters.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//
28// GEANT4 Class file
29//
30//
31// File name: G4EmExtraParameters
32//
33// Author: Vladimir Ivanchenko
34//
35// Creation date: 07.05.2019
36//
37// -------------------------------------------------------------------
38//
39//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
40//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
41
45#include "G4UnitsTable.hh"
46#include "G4SystemOfUnits.hh"
47#include "G4VEmProcess.hh"
50#include "G4RegionStore.hh"
51#include "G4Region.hh"
52
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
54
56{
57 theMessenger = new G4EmExtraParametersMessenger(this);
58 Initialise();
59}
60
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
62
64{
65 delete theMessenger;
66}
67
68//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
69
71{
72 quantumEntanglement = false;
73 directionalSplitting = false;
74 directionalSplittingTarget.set(0.,0.,0.);
75 directionalSplittingRadius = 0.;
76
77 dRoverRange = 0.2;
78 finalRange = CLHEP::mm;
79 dRoverRangeMuHad = 0.2;
80 finalRangeMuHad = 0.1*CLHEP::mm;
81 dRoverRangeLIons = 0.2;
82 finalRangeLIons = 0.1*CLHEP::mm;
83 dRoverRangeIons = 0.2;
84 finalRangeIons = 0.1*CLHEP::mm;
85
86 m_regnamesForced.clear();
87 m_procForced.clear();
88 m_lengthForced.clear();
89 m_weightForced.clear();
90 m_regnamesSubCut.clear();
91 m_subCuts.clear();
92}
93
94//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
95
96
97void G4EmExtraParameters::PrintWarning(G4ExceptionDescription& ed) const
98{
99 G4Exception("G4EmExtraParameters", "em0044", JustWarning, ed);
100}
101
102G4String G4EmExtraParameters::CheckRegion(const G4String& reg) const
103{
104 G4String r = reg;
105 if(r == "" || r == "world" || r == "World") {
106 r = "DefaultRegionForTheWorld";
107 }
108 return r;
109}
110
112{
113 if(v1 > 0.0 && v1 <= 1.0 && v2 > 0.0) {
114 dRoverRange = v1;
115 finalRange = v2;
116 } else {
118 ed << "Values of step function are out of range: "
119 << v1 << ", " << v2/CLHEP::mm << " mm - are ignored";
120 PrintWarning(ed);
121 }
122}
123
125{
126 return dRoverRange;
127}
128
130{
131 return finalRange;
132}
133
135{
136 if(v1 > 0.0 && v1 <= 1.0 && v2 > 0.0) {
137 dRoverRangeMuHad = v1;
138 finalRangeMuHad = v2;
139 } else {
141 ed << "Values of step function are out of range: "
142 << v1 << ", " << v2/CLHEP::mm << " mm - are ignored";
143 PrintWarning(ed);
144 }
145}
146
148{
149 return dRoverRangeMuHad;
150}
151
153{
154 return finalRangeMuHad;
155}
156
158{
159 if(v1 > 0.0 && v1 <= 1.0 && v2 > 0.0) {
160 dRoverRangeLIons = v1;
161 finalRangeLIons = v2;
162 } else {
164 ed << "Values of step function are out of range: "
165 << v1 << ", " << v2/CLHEP::mm << " mm - are ignored";
166 PrintWarning(ed);
167 }
168}
169
171{
172 return dRoverRangeLIons;
173}
174
176{
177 return finalRangeLIons;
178}
179
181{
182 if(v1 > 0.0 && v1 <= 1.0 && v2 > 0.0) {
183 dRoverRangeIons = v1;
184 finalRangeIons = v2;
185 } else {
187 ed << "Values of step function are out of range: "
188 << v1 << ", " << v2/CLHEP::mm << " mm - are ignored";
189 PrintWarning(ed);
190 }
191}
192
194{
195 return dRoverRangeIons;
196}
197
199{
200 return finalRangeIons;
201}
202
204{
205 // electron and positron
206 if (11 == std::abs(part->GetPDGEncoding())) {
207 proc->SetStepFunction(dRoverRange, finalRange);
208
209 // all heavy ions
210 } else if (part->IsGeneralIon()) {
211 proc->SetStepFunction(dRoverRangeIons, finalRangeIons);
212
213 // light nucleus and anti-nucleus
214 } else if (part->GetParticleType() == "nucleus" || part->GetParticleType() == "anti_nucleus") {
215 proc->SetStepFunction(dRoverRangeLIons, finalRangeLIons);
216
217 // other particles
218 } else {
219 proc->SetStepFunction(dRoverRangeMuHad, finalRangeMuHad);
220 }
221}
222
224 const G4String& region,
225 const G4String& type)
226{
227 G4String r = CheckRegion(region);
228 G4int nreg = m_regnamesPAI.size();
229 for(G4int i=0; i<nreg; ++i) {
230 if((m_particlesPAI[i] == particle ||
231 m_particlesPAI[i] == "all" ||
232 particle == "all") &&
233 (m_regnamesPAI[i] == r ||
234 m_regnamesPAI[i] == "DefaultRegionForTheWorld" ||
235 r == "DefaultRegionForTheWorld") ) {
236
237 m_typesPAI[i] = type;
238 if(particle == "all") { m_particlesPAI[i] = particle; }
239 if(r == "DefaultRegionForTheWorld") { m_regnamesPAI[i] = r; }
240 return;
241 }
242 }
243 m_particlesPAI.push_back(particle);
244 m_regnamesPAI.push_back(r);
245 m_typesPAI.push_back(type);
246}
247
248const std::vector<G4String>& G4EmExtraParameters::ParticlesPAI() const
249{
250 return m_particlesPAI;
251}
252
253const std::vector<G4String>& G4EmExtraParameters::RegionsPAI() const
254{
255 return m_regnamesPAI;
256}
257
258const std::vector<G4String>& G4EmExtraParameters::TypesPAI() const
259{
260 return m_typesPAI;
261}
262
264 const G4String& type)
265{
266 G4String r = CheckRegion(region);
267 G4int nreg = m_regnamesPhys.size();
268 for(G4int i=0; i<nreg; ++i) {
269 if(r == m_regnamesPhys[i]) { return; }
270 }
271 m_regnamesPhys.push_back(r);
272 m_typesPhys.push_back(type);
273}
274
275const std::vector<G4String>& G4EmExtraParameters::RegionsPhysics() const
276{
277 return m_regnamesPhys;
278}
279
280const std::vector<G4String>& G4EmExtraParameters::TypesPhysics() const
281{
282 return m_typesPhys;
283}
284
286{
287 G4String r = CheckRegion(region);
288 G4int nreg = m_regnamesSubCut.size();
289 for(G4int i=0; i<nreg; ++i) {
290 if(r == m_regnamesSubCut[i]) {
291 m_subCuts[i] = val;
292 return;
293 }
294 }
295 m_regnamesSubCut.push_back(r);
296 m_subCuts.push_back(val);
297}
298
299void
301 G4double val, G4bool wflag)
302{
303 if(val > 0.0) {
304 G4int n = m_procBiasedXS.size();
305 for(G4int i=0; i<n; ++i) {
306 if(procname == m_procBiasedXS[i]) {
307 m_factBiasedXS[i] = val;
308 m_weightBiasedXS[i]= wflag;
309 return;
310 }
311 }
312 m_procBiasedXS.push_back(procname);
313 m_factBiasedXS.push_back(val);
314 m_weightBiasedXS.push_back(wflag);
315 } else {
317 ed << "Process: " << procname << " XS biasing factor "
318 << val << " is negative - ignored";
319 PrintWarning(ed);
320 }
321}
322
323void
325 const G4String& region,
326 G4double length,
327 G4bool wflag)
328{
329 G4String r = CheckRegion(region);
330 if(length >= 0.0) {
331 G4int n = m_procForced.size();
332 for(G4int i=0; i<n; ++i) {
333 if(procname == m_procForced[i] && r == m_regnamesForced[i] ) {
334 m_lengthForced[i] = length;
335 m_weightForced[i]= wflag;
336 return;
337 }
338 }
339 m_regnamesForced.push_back(r);
340 m_procForced.push_back(procname);
341 m_lengthForced.push_back(length);
342 m_weightForced.push_back(wflag);
343 } else {
345 ed << "Process: " << procname << " in region " << r
346 << " : forced interacttion length= "
347 << length << " is negative - ignored";
348 PrintWarning(ed);
349 }
350}
351
352void
354 const G4String& region,
355 G4double factor,
356 G4double energyLim)
357{
358 G4String r = CheckRegion(region);
359 if(factor >= 0.0 && energyLim >= 0.0) {
360 G4int n = m_procBiasedSec.size();
361 for(G4int i=0; i<n; ++i) {
362 if(procname == m_procBiasedSec[i] && r == m_regnamesBiasedSec[i] ) {
363 m_factBiasedSec[i] = factor;
364 m_elimBiasedSec[i] = energyLim;
365 return;
366 }
367 }
368 m_regnamesBiasedSec.push_back(r);
369 m_procBiasedSec.push_back(procname);
370 m_factBiasedSec.push_back(factor);
371 m_elimBiasedSec.push_back(energyLim);
372 } else {
374 ed << "Process: " << procname << " in region " << r
375 << " : secondary bised factor= "
376 << factor << ", Elim= " << energyLim << " - ignored";
377 PrintWarning(ed);
378 }
379}
380
382{
384 G4int n = m_regnamesSubCut.size();
385 for(G4int i=0; i<n; ++i) {
386 const G4Region* reg = regionStore->GetRegion(m_regnamesSubCut[i], false);
387 if(reg) { ptr->ActivateSubCutoff(m_subCuts[i], reg); }
388 }
389 n = m_procBiasedXS.size();
390 for(G4int i=0; i<n; ++i) {
391 if(ptr->GetProcessName() == m_procBiasedXS[i]) {
392 ptr->SetCrossSectionBiasingFactor(m_factBiasedXS[i],
393 m_weightBiasedXS[i]);
394 break;
395 }
396 }
397 n = m_procForced.size();
398 for(G4int i=0; i<n; ++i) {
399 if(ptr->GetProcessName() == m_procForced[i]) {
400 ptr->ActivateForcedInteraction(m_lengthForced[i],
401 m_regnamesForced[i],
402 m_weightForced[i]);
403 break;
404 }
405 }
406 n = m_procBiasedSec.size();
407 for(G4int i=0; i<n; ++i) {
408 if(ptr->GetProcessName() == m_procBiasedSec[i]) {
409 ptr->ActivateSecondaryBiasing(m_regnamesBiasedSec[i],
410 m_factBiasedSec[i],
411 m_elimBiasedSec[i]);
412 break;
413 }
414 }
415}
416
418{
419 G4int n = m_procBiasedXS.size();
420 for(G4int i=0; i<n; ++i) {
421 if(ptr->GetProcessName() == m_procBiasedXS[i]) {
422 ptr->SetCrossSectionBiasingFactor(m_factBiasedXS[i],
423 m_weightBiasedXS[i]);
424 break;
425 }
426 }
427 n = m_procForced.size();
428 for(G4int i=0; i<n; ++i) {
429 if(ptr->GetProcessName() == m_procForced[i]) {
430 ptr->ActivateForcedInteraction(m_lengthForced[i],
431 m_regnamesForced[i],
432 m_weightForced[i]);
433 break;
434 }
435 }
436 n = m_procBiasedSec.size();
437 for(G4int i=0; i<n; ++i) {
438 if(ptr->GetProcessName() == m_procBiasedSec[i]) {
439 ptr->ActivateSecondaryBiasing(m_regnamesBiasedSec[i],
440 m_factBiasedSec[i],
441 m_elimBiasedSec[i]);
442 break;
443 }
444 }
445}
446
448{
449 return quantumEntanglement;
450}
451
453{
454 quantumEntanglement = v;
455}
456
458 return directionalSplitting;
459}
460
462{
463 directionalSplitting = v;
464}
465
466void
468{
469 directionalSplittingTarget = v;
470}
471
473{
474 return directionalSplittingTarget;
475}
476
478{
479 directionalSplittingRadius = r;
480}
481
483{
484 return directionalSplittingRadius;
485}
486
487//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
void set(double x, double y, double z)
G4double GetStepFunctionP1() const
G4double GetStepFunctionP2() const
void SetQuantumEntanglement(G4bool v)
void SetStepFunctionMuHad(G4double v1, G4double v2)
void SetDirectionalSplittingRadius(G4double r)
G4double GetStepFunctionIonsP2() const
void SetSubCutoff(G4bool val, const G4String &region="")
void SetStepFunctionIons(G4double v1, G4double v2)
const std::vector< G4String > & ParticlesPAI() const
G4double GetStepFunctionMuHadP1() const
void ActivateForcedInteraction(const G4String &procname, const G4String &region, G4double length, G4bool wflag)
void ActivateSecondaryBiasing(const G4String &name, const G4String &region, G4double factor, G4double energyLimit)
void SetDirectionalSplitting(G4bool v)
const std::vector< G4String > & TypesPhysics() const
const std::vector< G4String > & TypesPAI() const
void DefineRegParamForEM(G4VEmProcess *) const
void SetDirectionalSplittingTarget(const G4ThreeVector &v)
const std::vector< G4String > & RegionsPAI() const
G4double GetStepFunctionLightIonsP1() const
void FillStepFunction(const G4ParticleDefinition *, G4VEnergyLossProcess *) const
void SetProcessBiasingFactor(const G4String &procname, G4double val, G4bool wflag)
G4ThreeVector GetDirectionalSplittingTarget() const
void DefineRegParamForLoss(G4VEnergyLossProcess *) const
const std::vector< G4String > & RegionsPhysics() const
G4double GetStepFunctionIonsP1() const
void SetStepFunction(G4double v1, G4double v2)
G4double GetStepFunctionLightIonsP2() const
G4double GetStepFunctionMuHadP2() const
G4double GetDirectionalSplittingRadius()
void SetStepFunctionLightIons(G4double v1, G4double v2)
void AddPAIModel(const G4String &particle, const G4String &region, const G4String &type)
void AddPhysics(const G4String &region, const G4String &type)
G4bool IsGeneralIon() const
const G4String & GetParticleType() const
static G4RegionStore * GetInstance()
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
void SetCrossSectionBiasingFactor(G4double f, G4bool flag=true)
void ActivateForcedInteraction(G4double length=0.0, const G4String &r="", G4bool flag=true)
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
void ActivateSubCutoff(G4bool val, const G4Region *region=nullptr)
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
void SetStepFunction(G4double v1, G4double v2)
void SetCrossSectionBiasingFactor(G4double f, G4bool flag=true)
void ActivateForcedInteraction(G4double length, const G4String &region, G4bool flag=true)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382