Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmProcessOptions Class Reference

#include <G4EmProcessOptions.hh>

Public Member Functions

 G4EmProcessOptions ()
 
 ~G4EmProcessOptions ()
 
void SetLossFluctuations (G4bool val)
 
void SetSubCutoff (G4bool val, const G4Region *r=0)
 
void SetIntegral (G4bool val)
 
void SetMinSubRange (G4double val)
 
void SetMinEnergy (G4double val)
 
void SetMaxEnergy (G4double val)
 
void SetMaxEnergyForCSDARange (G4double val)
 
void SetMaxEnergyForMuons (G4double val)
 
void SetDEDXBinning (G4int val)
 
void SetDEDXBinningForCSDARange (G4int val)
 
void SetLambdaBinning (G4int val)
 
void SetStepFunction (G4double v1, G4double v2)
 
void SetRandomStep (G4bool val)
 
void SetApplyCuts (G4bool val)
 
void SetBuildCSDARange (G4bool val)
 
void SetVerbose (G4int val, const G4String &name="all")
 
void SetLambdaFactor (G4double val)
 
void SetLinearLossLimit (G4double val)
 
void SetDeexcitationActiveRegion (const G4String &rname="", G4bool valDeexcitation=true, G4bool valAuger=true, G4bool valPIXE=true)
 
void SetFluo (G4bool val)
 
void SetAuger (G4bool val)
 
void SetPIXE (G4bool val)
 
void SetPIXECrossSectionModel (const G4String &val)
 
void SetPIXEElectronCrossSectionModel (const G4String &val)
 
void SetMscStepLimitation (G4MscStepLimitType val)
 
void SetMscLateralDisplacement (G4bool val)
 
void SetSkin (G4double val)
 
void SetMscRangeFactor (G4double val)
 
void SetMscGeomFactor (G4double val)
 
void SetLPMFlag (G4bool val)
 
void SetSplineFlag (G4bool val)
 
void SetBremsstrahlungTh (G4double val)
 
void SetPolarAngleLimit (G4double val)
 
void SetFactorForAngleLimit (G4double val)
 
void SetProcessBiasingFactor (const G4String &name, G4double val, G4bool flag=true)
 
void ActivateForcedInteraction (const G4String &name, G4double length=0.0, const G4String &region="", G4bool flag=true)
 
void ActivateSecondaryBiasing (const G4String &name, const G4String &region, G4double factor, G4double energyLimit)
 
void ActivateSecondaryBiasingForGamma (const G4String &name, const G4String &region, G4double factor, G4double energyLimit)
 

Detailed Description

Definition at line 63 of file G4EmProcessOptions.hh.

Constructor & Destructor Documentation

◆ G4EmProcessOptions()

G4EmProcessOptions::G4EmProcessOptions ( )

Definition at line 67 of file G4EmProcessOptions.cc.

68{
69 theManager = G4LossTableManager::Instance();
70}
static G4LossTableManager * Instance()

◆ ~G4EmProcessOptions()

G4EmProcessOptions::~G4EmProcessOptions ( )

Definition at line 74 of file G4EmProcessOptions.cc.

75{
76}

Member Function Documentation

◆ ActivateForcedInteraction()

void G4EmProcessOptions::ActivateForcedInteraction ( const G4String name,
G4double  length = 0.0,
const G4String region = "",
G4bool  flag = true 
)

Definition at line 449 of file G4EmProcessOptions.cc.

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}
const std::vector< G4VEmProcess * > & GetEmProcessVector()
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector()
void ActivateForcedInteraction(G4double length=0.0, const G4String &r="", G4bool flag=true)
void ActivateForcedInteraction(G4double length=0.0, const G4String &region="", G4bool flag=true)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ ActivateSecondaryBiasing()

void G4EmProcessOptions::ActivateSecondaryBiasing ( const G4String name,
const G4String region,
G4double  factor,
G4double  energyLimit 
)

Definition at line 481 of file G4EmProcessOptions.cc.

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}
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ ActivateSecondaryBiasingForGamma()

void G4EmProcessOptions::ActivateSecondaryBiasingForGamma ( const G4String name,
const G4String region,
G4double  factor,
G4double  energyLimit 
)

Definition at line 503 of file G4EmProcessOptions.cc.

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}
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetApplyCuts()

void G4EmProcessOptions::SetApplyCuts ( G4bool  val)

Definition at line 171 of file G4EmProcessOptions.cc.

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}
void SetApplyCuts(G4bool val)

Referenced by G4EmStandardPhysics_option1::ConstructProcess(), G4EmStandardPhysics_option2::ConstructProcess(), and G4EnergyLossMessenger::SetNewValue().

◆ SetAuger()

void G4EmProcessOptions::SetAuger ( G4bool  val)

Definition at line 268 of file G4EmProcessOptions.cc.

269{
270 G4VAtomDeexcitation* ad = theManager-> AtomDeexcitation();
271 if(ad) { ad->SetAuger(val); }
272}

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetBremsstrahlungTh()

void G4EmProcessOptions::SetBremsstrahlungTh ( G4double  val)

Definition at line 411 of file G4EmProcessOptions.cc.

412{
413 theManager->SetBremsstrahlungTh(val);
414}
void SetBremsstrahlungTh(G4double val)

◆ SetBuildCSDARange()

void G4EmProcessOptions::SetBuildCSDARange ( G4bool  val)

Definition at line 184 of file G4EmProcessOptions.cc.

185{
186 theManager->SetBuildCSDARange(val);
187}
void SetBuildCSDARange(G4bool val)

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetDEDXBinning()

◆ SetDEDXBinningForCSDARange()

void G4EmProcessOptions::SetDEDXBinningForCSDARange ( G4int  val)

Definition at line 143 of file G4EmProcessOptions.cc.

144{
145 theManager->SetDEDXBinningForCSDARange(val);
146}
void SetDEDXBinningForCSDARange(G4int val)

◆ SetDeexcitationActiveRegion()

void G4EmProcessOptions::SetDeexcitationActiveRegion ( const G4String rname = "",
G4bool  valDeexcitation = true,
G4bool  valAuger = true,
G4bool  valPIXE = true 
)

Definition at line 246 of file G4EmProcessOptions.cc.

250{
251 G4VAtomDeexcitation* ad = theManager-> AtomDeexcitation();
252 if(ad) {
253 ad->SetDeexcitationActiveRegion(rname, valDeexcitation,
254 valAuger,valPIXE);
255 }
256}
void SetDeexcitationActiveRegion(const G4String &rname, G4bool valDeexcitation, G4bool valAuger, G4bool valPIXE)

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetFactorForAngleLimit()

void G4EmProcessOptions::SetFactorForAngleLimit ( G4double  val)

Definition at line 383 of file G4EmProcessOptions.cc.

384{
385 theManager->SetFactorForAngleLimit(val);
386}
void SetFactorForAngleLimit(G4double val)

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetFluo()

void G4EmProcessOptions::SetFluo ( G4bool  val)

Definition at line 260 of file G4EmProcessOptions.cc.

261{
262 G4VAtomDeexcitation* ad = theManager-> AtomDeexcitation();
263 if(ad) { ad->SetFluo(val); }
264}

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetIntegral()

void G4EmProcessOptions::SetIntegral ( G4bool  val)

Definition at line 94 of file G4EmProcessOptions.cc.

95{
96 theManager->SetIntegral(val);
97}
void SetIntegral(G4bool val)

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetLambdaBinning()

◆ SetLambdaFactor()

void G4EmProcessOptions::SetLambdaFactor ( G4double  val)

Definition at line 232 of file G4EmProcessOptions.cc.

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}
void SetLambdaFactor(G4double val)

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetLinearLossLimit()

void G4EmProcessOptions::SetLinearLossLimit ( G4double  val)

Definition at line 404 of file G4EmProcessOptions.cc.

405{
406 theManager->SetLinearLossLimit(val);
407}
void SetLinearLossLimit(G4double val)

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetLossFluctuations()

void G4EmProcessOptions::SetLossFluctuations ( G4bool  val)

Definition at line 80 of file G4EmProcessOptions.cc.

81{
82 theManager->SetLossFluctuations(val);
83}
void SetLossFluctuations(G4bool val)

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetLPMFlag()

void G4EmProcessOptions::SetLPMFlag ( G4bool  val)

Definition at line 390 of file G4EmProcessOptions.cc.

391{
392 theManager->SetLPMFlag(val);
393}
void SetLPMFlag(G4bool val)

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetMaxEnergy()

◆ SetMaxEnergyForCSDARange()

void G4EmProcessOptions::SetMaxEnergyForCSDARange ( G4double  val)

Definition at line 122 of file G4EmProcessOptions.cc.

123{
124 theManager->SetMaxEnergyForCSDARange(val);
125}
void SetMaxEnergyForCSDARange(G4double val)

◆ SetMaxEnergyForMuons()

void G4EmProcessOptions::SetMaxEnergyForMuons ( G4double  val)

Definition at line 129 of file G4EmProcessOptions.cc.

130{
131 theManager->SetMaxEnergyForMuons(val);
132}
void SetMaxEnergyForMuons(G4double val)

◆ SetMinEnergy()

◆ SetMinSubRange()

void G4EmProcessOptions::SetMinSubRange ( G4double  val)

Definition at line 101 of file G4EmProcessOptions.cc.

102{
103 theManager->SetMinSubRange(val);
104}
void SetMinSubRange(G4double val)

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetMscGeomFactor()

void G4EmProcessOptions::SetMscGeomFactor ( G4double  val)

Definition at line 351 of file G4EmProcessOptions.cc.

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}
const std::vector< G4VMultipleScattering * > & GetMultipleScatteringVector()

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetMscLateralDisplacement()

void G4EmProcessOptions::SetMscLateralDisplacement ( G4bool  val)

Definition at line 313 of file G4EmProcessOptions.cc.

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}

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetMscRangeFactor()

void G4EmProcessOptions::SetMscRangeFactor ( G4double  val)

Definition at line 338 of file G4EmProcessOptions.cc.

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}

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetMscStepLimitation()

void G4EmProcessOptions::SetMscStepLimitation ( G4MscStepLimitType  val)

Definition at line 301 of file G4EmProcessOptions.cc.

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}

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetPIXE()

void G4EmProcessOptions::SetPIXE ( G4bool  val)

Definition at line 276 of file G4EmProcessOptions.cc.

277{
278 G4VAtomDeexcitation* ad = theManager-> AtomDeexcitation();
279 if(ad) { ad->SetPIXE(val); }
280}

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetPIXECrossSectionModel()

void G4EmProcessOptions::SetPIXECrossSectionModel ( const G4String val)

Definition at line 284 of file G4EmProcessOptions.cc.

285{
286 G4VAtomDeexcitation* ad = theManager-> AtomDeexcitation();
287 if(ad) { ad->SetPIXECrossSectionModel(mname); }
288}
void SetPIXECrossSectionModel(const G4String &)

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetPIXEElectronCrossSectionModel()

void G4EmProcessOptions::SetPIXEElectronCrossSectionModel ( const G4String val)

Definition at line 293 of file G4EmProcessOptions.cc.

294{
295 G4VAtomDeexcitation* ad = theManager-> AtomDeexcitation();
296 if(ad) { ad->SetPIXEElectronCrossSectionModel(mname); }
297}
void SetPIXEElectronCrossSectionModel(const G4String &)

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetPolarAngleLimit()

void G4EmProcessOptions::SetPolarAngleLimit ( G4double  val)

Definition at line 364 of file G4EmProcessOptions.cc.

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}
void SetPolarAngleLimit(G4double a)

Referenced by G4EmLivermorePhysics::ConstructProcess(), G4EmLivermorePolarizedPhysics::ConstructProcess(), G4EmLowEPPhysics::ConstructProcess(), G4EmPenelopePhysics::ConstructProcess(), G4EmStandardPhysics::ConstructProcess(), G4EmStandardPhysics_option1::ConstructProcess(), G4EmStandardPhysics_option2::ConstructProcess(), G4EmStandardPhysics_option3::ConstructProcess(), G4EmStandardPhysics_option4::ConstructProcess(), and G4EnergyLossMessenger::SetNewValue().

◆ SetProcessBiasingFactor()

void G4EmProcessOptions::SetProcessBiasingFactor ( const G4String name,
G4double  val,
G4bool  flag = true 
)

Definition at line 419 of file G4EmProcessOptions.cc.

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}
void SetCrossSectionBiasingFactor(G4double f, G4bool flag=true)
void SetCrossSectionBiasingFactor(G4double f, G4bool flag=true)

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetRandomStep()

void G4EmProcessOptions::SetRandomStep ( G4bool  val)

Definition at line 164 of file G4EmProcessOptions.cc.

165{
166 theManager->SetRandomStep(val);
167}
void SetRandomStep(G4bool val)

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetSkin()

void G4EmProcessOptions::SetSkin ( G4double  val)

Definition at line 325 of file G4EmProcessOptions.cc.

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}

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetSplineFlag()

void G4EmProcessOptions::SetSplineFlag ( G4bool  val)

Definition at line 397 of file G4EmProcessOptions.cc.

398{
399 theManager->SetSplineFlag(val);
400}
void SetSplineFlag(G4bool val)

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetStepFunction()

void G4EmProcessOptions::SetStepFunction ( G4double  v1,
G4double  v2 
)

Definition at line 157 of file G4EmProcessOptions.cc.

158{
159 theManager->SetStepFunction(v1, v2);
160}
void SetStepFunction(G4double v1, G4double v2)

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetSubCutoff()

void G4EmProcessOptions::SetSubCutoff ( G4bool  val,
const G4Region r = 0 
)

Definition at line 87 of file G4EmProcessOptions.cc.

88{
89 theManager->SetSubCutoff(val, r);
90}
void SetSubCutoff(G4bool val, const G4Region *r=0)

Referenced by G4EnergyLossMessenger::SetNewValue().

◆ SetVerbose()

void G4EmProcessOptions::SetVerbose ( G4int  val,
const G4String name = "all" 
)

Definition at line 191 of file G4EmProcessOptions.cc.

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}
bool G4bool
Definition: G4Types.hh:67
void SetVerbose(G4int val)
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:408

Referenced by G4EmLivermorePhysics::ConstructProcess(), G4EmLivermorePolarizedPhysics::ConstructProcess(), G4EmLowEPPhysics::ConstructProcess(), G4EmPenelopePhysics::ConstructProcess(), G4EmStandardPhysics::ConstructProcess(), G4EmStandardPhysics_option1::ConstructProcess(), G4EmStandardPhysics_option2::ConstructProcess(), G4EmStandardPhysics_option3::ConstructProcess(), G4EmStandardPhysics_option4::ConstructProcess(), G4QAtomicPhysics::ConstructProcess(), and G4EnergyLossMessenger::SetNewValue().


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