Geant4 11.1.1
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}
92
93//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
94
95
96void G4EmExtraParameters::PrintWarning(G4ExceptionDescription& ed) const
97{
98 G4Exception("G4EmExtraParameters", "em0044", JustWarning, ed);
99}
100
101G4String G4EmExtraParameters::CheckRegion(const G4String& reg) const
102{
103 G4String r = reg;
104 if(r == "" || r == "world" || r == "World") {
105 r = "DefaultRegionForTheWorld";
106 }
107 return r;
108}
109
111{
112 if(v1 > 0.0 && v1 <= 1.0 && v2 > 0.0) {
113 dRoverRange = v1;
114 finalRange = v2;
115 } else {
117 ed << "Values of step function are out of range: "
118 << v1 << ", " << v2/CLHEP::mm << " mm - are ignored";
119 PrintWarning(ed);
120 }
121}
122
124{
125 return dRoverRange;
126}
127
129{
130 return finalRange;
131}
132
134{
135 if(v1 > 0.0 && v1 <= 1.0 && v2 > 0.0) {
136 dRoverRangeMuHad = v1;
137 finalRangeMuHad = v2;
138 } else {
140 ed << "Values of step function are out of range: "
141 << v1 << ", " << v2/CLHEP::mm << " mm - are ignored";
142 PrintWarning(ed);
143 }
144}
145
147{
148 return dRoverRangeMuHad;
149}
150
152{
153 return finalRangeMuHad;
154}
155
157{
158 if(v1 > 0.0 && v1 <= 1.0 && v2 > 0.0) {
159 dRoverRangeLIons = v1;
160 finalRangeLIons = v2;
161 } else {
163 ed << "Values of step function are out of range: "
164 << v1 << ", " << v2/CLHEP::mm << " mm - are ignored";
165 PrintWarning(ed);
166 }
167}
168
170{
171 return dRoverRangeLIons;
172}
173
175{
176 return finalRangeLIons;
177}
178
180{
181 if(v1 > 0.0 && v1 <= 1.0 && v2 > 0.0) {
182 dRoverRangeIons = v1;
183 finalRangeIons = v2;
184 } else {
186 ed << "Values of step function are out of range: "
187 << v1 << ", " << v2/CLHEP::mm << " mm - are ignored";
188 PrintWarning(ed);
189 }
190}
191
193{
194 return dRoverRangeIons;
195}
196
198{
199 return finalRangeIons;
200}
201
203{
204 // electron and positron
205 if (11 == std::abs(part->GetPDGEncoding())) {
206 proc->SetStepFunction(dRoverRange, finalRange);
207
208 // all heavy ions
209 } else if ("GenericIon" == part->GetParticleName()) {
210 proc->SetStepFunction(dRoverRangeIons, finalRangeIons);
211
212 // light nucleus and anti-nucleus
213 } else if (part->GetParticleType() == "nucleus" || part->GetParticleType() == "anti_nucleus") {
214 proc->SetStepFunction(dRoverRangeLIons, finalRangeLIons);
215
216 // other particles
217 } else {
218 proc->SetStepFunction(dRoverRangeMuHad, finalRangeMuHad);
219 }
220}
221
223 const G4String& region,
224 const G4String& type)
225{
226 G4String r = CheckRegion(region);
227 std::size_t nreg = m_regnamesPAI.size();
228 for(std::size_t i=0; i<nreg; ++i) {
229 if((m_particlesPAI[i] == particle ||
230 m_particlesPAI[i] == "all" ||
231 particle == "all") &&
232 (m_regnamesPAI[i] == r ||
233 m_regnamesPAI[i] == "DefaultRegionForTheWorld" ||
234 r == "DefaultRegionForTheWorld") ) {
235
236 m_typesPAI[i] = type;
237 if(particle == "all") { m_particlesPAI[i] = particle; }
238 if(r == "DefaultRegionForTheWorld") { m_regnamesPAI[i] = r; }
239 return;
240 }
241 }
242 m_particlesPAI.push_back(particle);
243 m_regnamesPAI.push_back(r);
244 m_typesPAI.push_back(type);
245}
246
247const std::vector<G4String>& G4EmExtraParameters::ParticlesPAI() const
248{
249 return m_particlesPAI;
250}
251
252const std::vector<G4String>& G4EmExtraParameters::RegionsPAI() const
253{
254 return m_regnamesPAI;
255}
256
257const std::vector<G4String>& G4EmExtraParameters::TypesPAI() const
258{
259 return m_typesPAI;
260}
261
263 const G4String& type)
264{
265 G4String r = CheckRegion(region);
266 std::size_t nreg = m_regnamesPhys.size();
267 for(std::size_t i=0; i<nreg; ++i) {
268 if(r == m_regnamesPhys[i]) { return; }
269 }
270 m_regnamesPhys.push_back(r);
271 m_typesPhys.push_back(type);
272}
273
274const std::vector<G4String>& G4EmExtraParameters::RegionsPhysics() const
275{
276 return m_regnamesPhys;
277}
278
279const std::vector<G4String>& G4EmExtraParameters::TypesPhysics() const
280{
281 return m_typesPhys;
282}
283
285{
286 const G4String& r = CheckRegion(region);
287 std::size_t nreg = m_regnamesSubCut.size();
288 for(std::size_t i=0; i<nreg; ++i) {
289 if(r == m_regnamesSubCut[i]) {
290 return;
291 }
292 }
293 m_regnamesSubCut.push_back(r);
294}
295
296void
298 G4double val, G4bool wflag)
299{
300 if(val > 0.0) {
301 std::size_t n = m_procBiasedXS.size();
302 for(std::size_t i=0; i<n; ++i) {
303 if(procname == m_procBiasedXS[i]) {
304 m_factBiasedXS[i] = val;
305 m_weightBiasedXS[i]= wflag;
306 return;
307 }
308 }
309 m_procBiasedXS.push_back(procname);
310 m_factBiasedXS.push_back(val);
311 m_weightBiasedXS.push_back(wflag);
312 } else {
314 ed << "Process: " << procname << " XS biasing factor "
315 << val << " is negative - ignored";
316 PrintWarning(ed);
317 }
318}
319
320void
322 const G4String& region,
323 G4double length,
324 G4bool wflag)
325{
326 const G4String& r = CheckRegion(region);
327 if(length >= 0.0) {
328 std::size_t n = m_procForced.size();
329 for(std::size_t i=0; i<n; ++i) {
330 if(procname == m_procForced[i] && r == m_regnamesForced[i] ) {
331 m_lengthForced[i] = length;
332 m_weightForced[i] = wflag;
333 return;
334 }
335 }
336 m_regnamesForced.push_back(r);
337 m_procForced.push_back(procname);
338 m_lengthForced.push_back(length);
339 m_weightForced.push_back(wflag);
340 } else {
342 ed << "Process: " << procname << " in region " << r
343 << " : forced interacttion length= "
344 << length << " is negative - ignored";
345 PrintWarning(ed);
346 }
347}
348
349void
351 const G4String& region,
352 G4double factor,
353 G4double energyLim)
354{
355 const G4String& r = CheckRegion(region);
356 if(factor >= 0.0 && energyLim >= 0.0) {
357 std::size_t n = m_procBiasedSec.size();
358 for(std::size_t i=0; i<n; ++i) {
359 if(procname == m_procBiasedSec[i] && r == m_regnamesBiasedSec[i] ) {
360 m_factBiasedSec[i] = factor;
361 m_elimBiasedSec[i] = energyLim;
362 return;
363 }
364 }
365 m_regnamesBiasedSec.push_back(r);
366 m_procBiasedSec.push_back(procname);
367 m_factBiasedSec.push_back(factor);
368 m_elimBiasedSec.push_back(energyLim);
369 } else {
371 ed << "Process: " << procname << " in region " << r
372 << " : secondary bised factor= "
373 << factor << ", Elim= " << energyLim << " - ignored";
374 PrintWarning(ed);
375 }
376}
377
379{
380 const G4RegionStore* regionStore = G4RegionStore::GetInstance();
381 std::size_t n = m_regnamesSubCut.size();
382 for(std::size_t i=0; i<n; ++i) {
383 const G4Region* reg = regionStore->GetRegion(m_regnamesSubCut[i], false);
384 if(nullptr != reg) { ptr->ActivateSubCutoff(reg); }
385 }
386 n = m_procBiasedXS.size();
387 for(std::size_t i=0; i<n; ++i) {
388 if(ptr->GetProcessName() == m_procBiasedXS[i]) {
389 ptr->SetCrossSectionBiasingFactor(m_factBiasedXS[i],
390 m_weightBiasedXS[i]);
391 break;
392 }
393 }
394 n = m_procForced.size();
395 for(std::size_t i=0; i<n; ++i) {
396 if(ptr->GetProcessName() == m_procForced[i]) {
397 ptr->ActivateForcedInteraction(m_lengthForced[i],
398 m_regnamesForced[i],
399 m_weightForced[i]);
400 break;
401 }
402 }
403 n = m_procBiasedSec.size();
404 for(std::size_t i=0; i<n; ++i) {
405 if(ptr->GetProcessName() == m_procBiasedSec[i]) {
406 ptr->ActivateSecondaryBiasing(m_regnamesBiasedSec[i],
407 m_factBiasedSec[i],
408 m_elimBiasedSec[i]);
409 break;
410 }
411 }
412}
413
415{
416 std::size_t n = m_procBiasedXS.size();
417 for(std::size_t i=0; i<n; ++i) {
418 if(ptr->GetProcessName() == m_procBiasedXS[i]) {
419 ptr->SetCrossSectionBiasingFactor(m_factBiasedXS[i],
420 m_weightBiasedXS[i]);
421 break;
422 }
423 }
424 n = m_procForced.size();
425 for(std::size_t i=0; i<n; ++i) {
426 if(ptr->GetProcessName() == m_procForced[i]) {
427 ptr->ActivateForcedInteraction(m_lengthForced[i],
428 m_regnamesForced[i],
429 m_weightForced[i]);
430 break;
431 }
432 }
433 n = m_procBiasedSec.size();
434 for(std::size_t i=0; i<n; ++i) {
435 if(ptr->GetProcessName() == m_procBiasedSec[i]) {
436 ptr->ActivateSecondaryBiasing(m_regnamesBiasedSec[i],
437 m_factBiasedSec[i],
438 m_elimBiasedSec[i]);
439 break;
440 }
441 }
442}
443
445{
446 return quantumEntanglement;
447}
448
450{
451 quantumEntanglement = v;
452}
453
455 return directionalSplitting;
456}
457
459{
460 directionalSplitting = v;
461}
462
463void
465{
466 directionalSplittingTarget = v;
467}
468
470{
471 return directionalSplittingTarget;
472}
473
475{
476 directionalSplittingRadius = r;
477}
478
480{
481 return directionalSplittingRadius;
482}
483
484//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
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 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
void SetSubCutRegion(const G4String &region)
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)
const G4String & GetParticleType() const
const G4String & GetParticleName() 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 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)
void ActivateSubCutoff(const G4Region *region)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386