Geant4 11.2.2
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
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
229 // in previously defined region other particles may be already defined
230 // type should be overrided for the same region and particle
231 for(std::size_t i=0; i<nreg; ++i) {
232 if(m_regnamesPAI[i] == r) {
233 if (particle == "all") {
234 m_particlesPAI[i] = particle;
235 m_typesPAI[i] = type;
236 return;
237 } else if(m_particlesPAI[i] == particle || m_particlesPAI[i] == "all") {
238 m_typesPAI[i] = type;
239 return;
240 }
241 }
242 }
243 // new regions and/or particles
244 m_particlesPAI.push_back(particle);
245 m_regnamesPAI.push_back(r);
246 m_typesPAI.push_back(type);
247}
248
249const std::vector<G4String>& G4EmExtraParameters::ParticlesPAI() const
250{
251 return m_particlesPAI;
252}
253
254const std::vector<G4String>& G4EmExtraParameters::RegionsPAI() const
255{
256 return m_regnamesPAI;
257}
258
259const std::vector<G4String>& G4EmExtraParameters::TypesPAI() const
260{
261 return m_typesPAI;
262}
263
265 const G4String& type)
266{
267 G4String r = CheckRegion(region);
268 std::size_t nreg = m_regnamesPhys.size();
269 for(std::size_t i=0; i<nreg; ++i) {
270 if(r == m_regnamesPhys[i]) { return; }
271 }
272 m_regnamesPhys.push_back(r);
273 m_typesPhys.push_back(type);
274}
275
276const std::vector<G4String>& G4EmExtraParameters::RegionsPhysics() const
277{
278 return m_regnamesPhys;
279}
280
281const std::vector<G4String>& G4EmExtraParameters::TypesPhysics() const
282{
283 return m_typesPhys;
284}
285
287{
288 const G4String& r = CheckRegion(region);
289 std::size_t nreg = m_regnamesSubCut.size();
290 for(std::size_t i=0; i<nreg; ++i) {
291 if(r == m_regnamesSubCut[i]) {
292 return;
293 }
294 }
295 m_regnamesSubCut.push_back(r);
296}
297
298void
300 G4double val, G4bool wflag)
301{
302 if(val > 0.0) {
303 std::size_t n = m_procBiasedXS.size();
304 for(std::size_t i=0; i<n; ++i) {
305 if(procname == m_procBiasedXS[i]) {
306 m_factBiasedXS[i] = val;
307 m_weightBiasedXS[i]= wflag;
308 return;
309 }
310 }
311 m_procBiasedXS.push_back(procname);
312 m_factBiasedXS.push_back(val);
313 m_weightBiasedXS.push_back(wflag);
314 } else {
316 ed << "Process: " << procname << " XS biasing factor "
317 << val << " is negative - ignored";
318 PrintWarning(ed);
319 }
320}
321
322void
324 const G4String& region,
325 G4double length,
326 G4bool wflag)
327{
328 const G4String& r = CheckRegion(region);
329 if(length >= 0.0) {
330 std::size_t n = m_procForced.size();
331 for(std::size_t i=0; i<n; ++i) {
332 if(procname == m_procForced[i] && r == m_regnamesForced[i] ) {
333 m_lengthForced[i] = length;
334 m_weightForced[i] = wflag;
335 return;
336 }
337 }
338 m_regnamesForced.push_back(r);
339 m_procForced.push_back(procname);
340 m_lengthForced.push_back(length);
341 m_weightForced.push_back(wflag);
342 } else {
344 ed << "Process: " << procname << " in region " << r
345 << " : forced interacttion length= "
346 << length << " is negative - ignored";
347 PrintWarning(ed);
348 }
349}
350
351void
353 const G4String& region,
354 G4double factor,
355 G4double energyLim)
356{
357 const G4String& r = CheckRegion(region);
358 if(factor >= 0.0 && energyLim >= 0.0) {
359 std::size_t n = m_procBiasedSec.size();
360 for(std::size_t i=0; i<n; ++i) {
361 if(procname == m_procBiasedSec[i] && r == m_regnamesBiasedSec[i] ) {
362 m_factBiasedSec[i] = factor;
363 m_elimBiasedSec[i] = energyLim;
364 return;
365 }
366 }
367 m_regnamesBiasedSec.push_back(r);
368 m_procBiasedSec.push_back(procname);
369 m_factBiasedSec.push_back(factor);
370 m_elimBiasedSec.push_back(energyLim);
371 } else {
373 ed << "Process: " << procname << " in region " << r
374 << " : secondary bised factor= "
375 << factor << ", Elim= " << energyLim << " - ignored";
376 PrintWarning(ed);
377 }
378}
379
381{
382 const G4RegionStore* regionStore = G4RegionStore::GetInstance();
383 std::size_t n = m_regnamesSubCut.size();
384 for(std::size_t i=0; i<n; ++i) {
385 const G4Region* reg = regionStore->GetRegion(m_regnamesSubCut[i], false);
386 if(nullptr != reg) { ptr->ActivateSubCutoff(reg); }
387 }
388 n = m_procBiasedXS.size();
389 for(std::size_t i=0; i<n; ++i) {
390 if(ptr->GetProcessName() == m_procBiasedXS[i]) {
391 ptr->SetCrossSectionBiasingFactor(m_factBiasedXS[i],
392 m_weightBiasedXS[i]);
393 break;
394 }
395 }
396 n = m_procForced.size();
397 for(std::size_t i=0; i<n; ++i) {
398 if(ptr->GetProcessName() == m_procForced[i]) {
399 ptr->ActivateForcedInteraction(m_lengthForced[i],
400 m_regnamesForced[i],
401 m_weightForced[i]);
402 break;
403 }
404 }
405 n = m_procBiasedSec.size();
406 for(std::size_t i=0; i<n; ++i) {
407 if(ptr->GetProcessName() == m_procBiasedSec[i]) {
408 ptr->ActivateSecondaryBiasing(m_regnamesBiasedSec[i],
409 m_factBiasedSec[i],
410 m_elimBiasedSec[i]);
411 break;
412 }
413 }
414}
415
417{
418 std::size_t n = m_procBiasedXS.size();
419 for(std::size_t i=0; i<n; ++i) {
420 if(ptr->GetProcessName() == m_procBiasedXS[i]) {
421 ptr->SetCrossSectionBiasingFactor(m_factBiasedXS[i],
422 m_weightBiasedXS[i]);
423 break;
424 }
425 }
426 n = m_procForced.size();
427 for(std::size_t i=0; i<n; ++i) {
428 if(ptr->GetProcessName() == m_procForced[i]) {
429 ptr->ActivateForcedInteraction(m_lengthForced[i],
430 m_regnamesForced[i],
431 m_weightForced[i]);
432 break;
433 }
434 }
435 n = m_procBiasedSec.size();
436 for(std::size_t i=0; i<n; ++i) {
437 if(ptr->GetProcessName() == m_procBiasedSec[i]) {
438 ptr->ActivateSecondaryBiasing(m_regnamesBiasedSec[i],
439 m_factBiasedSec[i],
440 m_elimBiasedSec[i]);
441 break;
442 }
443 }
444}
445
447{
448 return quantumEntanglement;
449}
450
452{
453 quantumEntanglement = v;
454}
455
457 return directionalSplitting;
458}
459
461{
462 directionalSplitting = v;
463}
464
465void
467{
468 directionalSplittingTarget = v;
469}
470
472{
473 return directionalSplittingTarget;
474}
475
477{
478 directionalSplittingRadius = r;
479}
480
482{
483 return directionalSplittingRadius;
484}
485
486//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
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