Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmProcessOptions.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// $Id$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class file
31//
32//
33// File name: G4EmProcessOptions
34//
35// Author: Vladimir Ivanchenko
36//
37// Creation date: 27.02.2004
38//
39// Modifications:
40// 30-06-04 G4EmProcess is pure discrete (V.Ivanchenko)
41// 24-03-05 Add ApplyCuts and RandomStep (V.Ivanchenko)
42// 10-01-06 PreciseRange -> CSDARange (V.Ivantchenko)
43// 10-05-06 Add command MscStepLimit to G4LossTableManager (V.Ivantchenko)
44// 22-05-06 Add SetBremsstrahlungTh (V.Ivanchenko)
45// 12-02-07 Add SetSkin, SetLinearLossLimit (V.Ivanchenko)
46// 30-05-12 Add biasing for G4VEmProcess (D. Sawkey)
47//
48// Class Description:
49//
50// -------------------------------------------------------------------
51//
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54
55#include "G4EmProcessOptions.hh"
56#include "G4SystemOfUnits.hh"
57#include "G4LossTableManager.hh"
58#include "G4VEmProcess.hh"
61#include "G4Region.hh"
62#include "G4RegionStore.hh"
64
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66
68{
69 theManager = G4LossTableManager::Instance();
70}
71
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
73
75{
76}
77
78//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
79
81{
82 theManager->SetLossFluctuations(val);
83}
84
85//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
86
88{
89 theManager->SetSubCutoff(val, r);
90}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
93
95{
96 theManager->SetIntegral(val);
97}
98
99//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
100
102{
103 theManager->SetMinSubRange(val);
104}
105
106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107
109{
110 theManager->SetMinEnergy(val);
111}
112
113//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
114
116{
117 theManager->SetMaxEnergy(val);
118}
119
120//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
121
123{
124 theManager->SetMaxEnergyForCSDARange(val);
125}
126
127//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
128
130{
131 theManager->SetMaxEnergyForMuons(val);
132}
133
134//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
135
137{
138 theManager->SetDEDXBinning(val);
139}
140
141//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
142
144{
145 theManager->SetDEDXBinningForCSDARange(val);
146}
147
148//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
149
151{
152 theManager->SetLambdaBinning(val);
153}
154
155//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
156
158{
159 theManager->SetStepFunction(v1, v2);
160}
161
162//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
163
165{
166 theManager->SetRandomStep(val);
167}
168
169//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
170
172{
173 const std::vector<G4VEmProcess*>& w =
174 theManager->GetEmProcessVector();
175 std::vector<G4VEmProcess*>::const_iterator itp;
176 for(itp = w.begin(); itp != w.end(); itp++) {
177 G4VEmProcess* q = *itp;
178 if(q) { q->SetApplyCuts(val); }
179 }
180}
181
182//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
183
185{
186 theManager->SetBuildCSDARange(val);
187}
188
189//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
190
192{
193 G4bool all = false;
194 if("all" == name) { all = true; }
195 const std::vector<G4VEnergyLossProcess*>& v =
196 theManager->GetEnergyLossProcessVector();
197
198 if(all) {
199 theManager->SetVerbose(val);
200 return;
201 }
202
203 std::vector<G4VEnergyLossProcess*>::const_iterator itr;
204 for(itr = v.begin(); itr != v.end(); ++itr) {
205 G4VEnergyLossProcess* p = *itr;
206 if(p) {
207 if (p->GetProcessName() == name) { p->SetVerboseLevel(val); }
208 }
209 }
210 const std::vector<G4VEmProcess*>& w =
211 theManager->GetEmProcessVector();
212 std::vector<G4VEmProcess*>::const_iterator itp;
213 for(itp = w.begin(); itp != w.end(); itp++) {
214 G4VEmProcess* q = *itp;
215 if(q) {
216 if (q->GetProcessName() == name) { q->SetVerboseLevel(val); }
217 }
218 }
219 const std::vector<G4VMultipleScattering*>& u =
220 theManager->GetMultipleScatteringVector();
221 std::vector<G4VMultipleScattering*>::const_iterator itm;
222 for(itm = u.begin(); itm != u.end(); itm++) {
223 G4VMultipleScattering* msc = *itm;
224 if(s) {
225 if (msc->GetProcessName() == name) { msc->SetVerboseLevel(val); }
226 }
227 }
228}
229
230//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
231
233{
234 const std::vector<G4VEnergyLossProcess*>& v =
235 theManager->GetEnergyLossProcessVector();
236 std::vector<G4VEnergyLossProcess*>::const_iterator itr;
237 for(itr = v.begin(); itr != v.end(); itr++) {
238 G4VEnergyLossProcess* p = *itr;
239 if(p) { p->SetLambdaFactor(val); }
240 }
241}
242
243//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
244
245void
247 G4bool valDeexcitation,
248 G4bool valAuger,
249 G4bool valPIXE)
250{
251 G4VAtomDeexcitation* ad = theManager-> AtomDeexcitation();
252 if(ad) {
253 ad->SetDeexcitationActiveRegion(rname, valDeexcitation,
254 valAuger,valPIXE);
255 }
256}
257
258//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
259
261{
262 G4VAtomDeexcitation* ad = theManager-> AtomDeexcitation();
263 if(ad) { ad->SetFluo(val); }
264}
265
266//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
267
269{
270 G4VAtomDeexcitation* ad = theManager-> AtomDeexcitation();
271 if(ad) { ad->SetAuger(val); }
272}
273
274//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
275
277{
278 G4VAtomDeexcitation* ad = theManager-> AtomDeexcitation();
279 if(ad) { ad->SetPIXE(val); }
280}
281
282//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
283
285{
286 G4VAtomDeexcitation* ad = theManager-> AtomDeexcitation();
287 if(ad) { ad->SetPIXECrossSectionModel(mname); }
288}
289
290//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
291
292void
294{
295 G4VAtomDeexcitation* ad = theManager-> AtomDeexcitation();
296 if(ad) { ad->SetPIXEElectronCrossSectionModel(mname); }
297}
298
299//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
300
302{
303 const std::vector<G4VMultipleScattering*>& u =
304 theManager->GetMultipleScatteringVector();
305 std::vector<G4VMultipleScattering*>::const_iterator itm;
306 for(itm = u.begin(); itm != u.end(); itm++) {
307 if(*itm) (*itm)->SetStepLimitType(val);
308 }
309}
310
311//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
312
314{
315 const std::vector<G4VMultipleScattering*>& u =
316 theManager->GetMultipleScatteringVector();
317 std::vector<G4VMultipleScattering*>::const_iterator itm;
318 for(itm = u.begin(); itm != u.end(); itm++) {
319 if(*itm) { (*itm)->SetLateralDisplasmentFlag(val); }
320 }
321}
322
323//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
324
326{
327 if(val < 0.0) return;
328 const std::vector<G4VMultipleScattering*>& u =
329 theManager->GetMultipleScatteringVector();
330 std::vector<G4VMultipleScattering*>::const_iterator itm;
331 for(itm = u.begin(); itm != u.end(); itm++) {
332 if(*itm) { (*itm)->SetSkin(val); }
333 }
334}
335
336//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
337
339{
340 if(val < 0.0) return;
341 const std::vector<G4VMultipleScattering*>& u =
342 theManager->GetMultipleScatteringVector();
343 std::vector<G4VMultipleScattering*>::const_iterator itm;
344 for(itm = u.begin(); itm != u.end(); itm++) {
345 if(*itm) { (*itm)->SetRangeFactor(val); }
346 }
347}
348
349//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
350
352{
353 if(val < 0.0) { return; }
354 const std::vector<G4VMultipleScattering*>& u =
355 theManager->GetMultipleScatteringVector();
356 std::vector<G4VMultipleScattering*>::const_iterator itm;
357 for(itm = u.begin(); itm != u.end(); itm++) {
358 if(*itm) { (*itm)->SetGeomFactor(val); }
359 }
360}
361
362//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
363
365{
366 const std::vector<G4VMultipleScattering*>& u =
367 theManager->GetMultipleScatteringVector();
368 std::vector<G4VMultipleScattering*>::const_iterator itm;
369 for(itm = u.begin(); itm != u.end(); itm++) {
370 if(*itm) { (*itm)->SetPolarAngleLimit(val); }
371 }
372 const std::vector<G4VEmProcess*>& w =
373 theManager->GetEmProcessVector();
374 std::vector<G4VEmProcess*>::const_iterator itp;
375 for(itp = w.begin(); itp != w.end(); itp++) {
376 G4VEmProcess* q = *itp;
377 if(q) { q->SetPolarAngleLimit(val); }
378 }
379}
380
381//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
382
384{
385 theManager->SetFactorForAngleLimit(val);
386}
387
388//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
389
391{
392 theManager->SetLPMFlag(val);
393}
394
395//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
396
398{
399 theManager->SetSplineFlag(val);
400}
401
402//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
403
405{
406 theManager->SetLinearLossLimit(val);
407}
408
409//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
410
412{
413 theManager->SetBremsstrahlungTh(val);
414}
415
416//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
417
418void
420 G4bool flag)
421{
422 const std::vector<G4VEnergyLossProcess*>& v =
423 theManager->GetEnergyLossProcessVector();
424 std::vector<G4VEnergyLossProcess*>::const_iterator itr;
425 for(itr = v.begin(); itr != v.end(); ++itr) {
426 G4VEnergyLossProcess* p = *itr;
427 if(p) {
428 if (p->GetProcessName() == name) {
429 p->SetCrossSectionBiasingFactor(val, flag);
430 }
431 }
432 }
433 const std::vector<G4VEmProcess*>& w =
434 theManager->GetEmProcessVector();
435 std::vector<G4VEmProcess*>::const_iterator itp;
436 for(itp = w.begin(); itp != w.end(); itp++) {
437 G4VEmProcess* q = *itp;
438 if(q) {
439 if (q->GetProcessName() == name) {
440 q->SetCrossSectionBiasingFactor(val, flag);
441 }
442 }
443 }
444}
445
446//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
447
448void
450 G4double length,
451 const G4String& region,
452 G4bool flag)
453{
454 const std::vector<G4VEnergyLossProcess*>& v =
455 theManager->GetEnergyLossProcessVector();
456 std::vector<G4VEnergyLossProcess*>::const_iterator itr;
457 for(itr = v.begin(); itr != v.end(); ++itr) {
458 G4VEnergyLossProcess* p = *itr;
459 if(p) {
460 if (p->GetProcessName() == name) {
461 p->ActivateForcedInteraction(length,region,flag);
462 }
463 }
464 }
465 const std::vector<G4VEmProcess*>& w =
466 theManager->GetEmProcessVector();
467 std::vector<G4VEmProcess*>::const_iterator itp;
468 for(itp = w.begin(); itp != w.end(); itp++) {
469 G4VEmProcess* q = *itp;
470 if(q) {
471 if (q->GetProcessName() == name) {
472 q->ActivateForcedInteraction(length,region,flag);
473 }
474 }
475 }
476}
477
478//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
479
480void
482 const G4String& region,
483 G4double factor,
484 G4double energyLimit)
485{
486 if(0.0 > factor) { return; }
487 const std::vector<G4VEnergyLossProcess*>& v =
488 theManager->GetEnergyLossProcessVector();
489 std::vector<G4VEnergyLossProcess*>::const_iterator itr;
490 for(itr = v.begin(); itr != v.end(); ++itr) {
491 G4VEnergyLossProcess* p = *itr;
492 if(p) {
493 if (p->GetProcessName() == name) {
494 p->ActivateSecondaryBiasing(region, factor, energyLimit);
495 }
496 }
497 }
498}
499
500//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
501
502void
504 const G4String& region,
505 G4double factor,
506 G4double energyLimit)
507{
508 if(0.0 > factor) { return; }
509 const std::vector<G4VEmProcess*>& v =
510 theManager->GetEmProcessVector();
511 std::vector<G4VEmProcess*>::const_iterator itr;
512 for(itr = v.begin(); itr != v.end(); ++itr) {
513 G4VEmProcess* p = *itr;
514 if(p) {
515 if (p->GetProcessName() == name) {
516 p->ActivateSecondaryBiasing(region, factor, energyLimit);
517 }
518 }
519 }
520}
521
522//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
523
G4MscStepLimitType
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
void SetBuildCSDARange(G4bool val)
void ActivateForcedInteraction(const G4String &name, G4double length=0.0, const G4String &region="", G4bool flag=true)
void SetAuger(G4bool val)
void SetLinearLossLimit(G4double val)
void SetMaxEnergy(G4double val)
void SetMscRangeFactor(G4double val)
void SetDeexcitationActiveRegion(const G4String &rname="", G4bool valDeexcitation=true, G4bool valAuger=true, G4bool valPIXE=true)
void SetIntegral(G4bool val)
void SetDEDXBinning(G4int val)
void SetLambdaBinning(G4int val)
void SetPolarAngleLimit(G4double val)
void SetVerbose(G4int val, const G4String &name="all")
void SetPIXECrossSectionModel(const G4String &val)
void SetBremsstrahlungTh(G4double val)
void SetDEDXBinningForCSDARange(G4int val)
void SetMinSubRange(G4double val)
void SetSkin(G4double val)
void SetPIXEElectronCrossSectionModel(const G4String &val)
void SetLPMFlag(G4bool val)
void SetMscLateralDisplacement(G4bool val)
void ActivateSecondaryBiasing(const G4String &name, const G4String &region, G4double factor, G4double energyLimit)
void SetMaxEnergyForCSDARange(G4double val)
void SetMinEnergy(G4double val)
void SetSplineFlag(G4bool val)
void SetMscGeomFactor(G4double val)
void SetLossFluctuations(G4bool val)
void SetRandomStep(G4bool val)
void SetLambdaFactor(G4double val)
void SetPIXE(G4bool val)
void SetMscStepLimitation(G4MscStepLimitType val)
void SetMaxEnergyForMuons(G4double val)
void SetFactorForAngleLimit(G4double val)
void SetApplyCuts(G4bool val)
void SetSubCutoff(G4bool val, const G4Region *r=0)
void SetStepFunction(G4double v1, G4double v2)
void ActivateSecondaryBiasingForGamma(const G4String &name, const G4String &region, G4double factor, G4double energyLimit)
void SetFluo(G4bool val)
void SetProcessBiasingFactor(const G4String &name, G4double val, G4bool flag=true)
void SetLinearLossLimit(G4double val)
void SetFactorForAngleLimit(G4double val)
void SetSplineFlag(G4bool val)
void SetMaxEnergy(G4double val)
static G4LossTableManager * Instance()
const std::vector< G4VEmProcess * > & GetEmProcessVector()
void SetSubCutoff(G4bool val, const G4Region *r=0)
void SetIntegral(G4bool val)
void SetBremsstrahlungTh(G4double val)
void SetMaxEnergyForCSDARange(G4double val)
const std::vector< G4VMultipleScattering * > & GetMultipleScatteringVector()
void SetVerbose(G4int val)
void SetBuildCSDARange(G4bool val)
void SetMinEnergy(G4double val)
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector()
void SetRandomStep(G4bool val)
void SetDEDXBinning(G4int val)
void SetStepFunction(G4double v1, G4double v2)
void SetLossFluctuations(G4bool val)
void SetMaxEnergyForMuons(G4double val)
void SetLambdaBinning(G4int val)
void SetLPMFlag(G4bool val)
void SetMinSubRange(G4double val)
void SetDEDXBinningForCSDARange(G4int val)
void SetPIXECrossSectionModel(const G4String &)
void SetPIXEElectronCrossSectionModel(const G4String &)
void SetDeexcitationActiveRegion(const G4String &rname, G4bool valDeexcitation, G4bool valAuger, G4bool valPIXE)
void SetCrossSectionBiasingFactor(G4double f, G4bool flag=true)
void SetApplyCuts(G4bool val)
void ActivateForcedInteraction(G4double length=0.0, const G4String &r="", G4bool flag=true)
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
void SetPolarAngleLimit(G4double a)
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
void SetCrossSectionBiasingFactor(G4double f, G4bool flag=true)
void SetLambdaFactor(G4double val)
void ActivateForcedInteraction(G4double length=0.0, const G4String &region="", G4bool flag=true)
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:408
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379