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 (  ) 

Definition at line 67 of file G4EmProcessOptions.cc.

References G4LossTableManager::Instance().

00068 {
00069   theManager = G4LossTableManager::Instance();
00070 }

G4EmProcessOptions::~G4EmProcessOptions (  ) 

Definition at line 74 of file G4EmProcessOptions.cc.

00075 {
00076 }


Member Function Documentation

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

Definition at line 449 of file G4EmProcessOptions.cc.

References G4LossTableManager::GetEmProcessVector(), and G4LossTableManager::GetEnergyLossProcessVector().

Referenced by G4EnergyLossMessenger::SetNewValue().

00453 {
00454   const std::vector<G4VEnergyLossProcess*>& v =
00455         theManager->GetEnergyLossProcessVector();
00456   std::vector<G4VEnergyLossProcess*>::const_iterator itr;
00457   for(itr = v.begin(); itr != v.end(); ++itr) {
00458     G4VEnergyLossProcess* p = *itr;
00459     if(p) {
00460       if (p->GetProcessName() == name) { 
00461         p->ActivateForcedInteraction(length,region,flag); 
00462       }
00463     }
00464   }
00465   const std::vector<G4VEmProcess*>& w =
00466         theManager->GetEmProcessVector();
00467   std::vector<G4VEmProcess*>::const_iterator itp;
00468   for(itp = w.begin(); itp != w.end(); itp++) {
00469     G4VEmProcess* q = *itp;
00470     if(q) {
00471       if (q->GetProcessName() == name) { 
00472         q->ActivateForcedInteraction(length,region,flag); 
00473       }
00474     }
00475   }
00476 }

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

Definition at line 481 of file G4EmProcessOptions.cc.

References G4LossTableManager::GetEnergyLossProcessVector().

Referenced by G4EnergyLossMessenger::SetNewValue().

00485 {
00486   if(0.0 > factor) { return; }
00487   const std::vector<G4VEnergyLossProcess*>& v =
00488         theManager->GetEnergyLossProcessVector();
00489   std::vector<G4VEnergyLossProcess*>::const_iterator itr;
00490   for(itr = v.begin(); itr != v.end(); ++itr) {
00491     G4VEnergyLossProcess* p = *itr;
00492     if(p) {
00493       if (p->GetProcessName() == name) { 
00494         p->ActivateSecondaryBiasing(region, factor, energyLimit); 
00495       }
00496     }
00497   }
00498 }

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

Definition at line 503 of file G4EmProcessOptions.cc.

References G4LossTableManager::GetEmProcessVector().

Referenced by G4EnergyLossMessenger::SetNewValue().

00507 {
00508   if(0.0 > factor) { return; }
00509   const std::vector<G4VEmProcess*>& v =
00510         theManager->GetEmProcessVector();
00511   std::vector<G4VEmProcess*>::const_iterator itr;
00512   for(itr = v.begin(); itr != v.end(); ++itr) {
00513     G4VEmProcess* p = *itr;
00514     if(p) {
00515       if (p->GetProcessName() == name) { 
00516         p->ActivateSecondaryBiasing(region, factor, energyLimit); 
00517       }
00518     }
00519   }
00520 }

void G4EmProcessOptions::SetApplyCuts ( G4bool  val  ) 

Definition at line 171 of file G4EmProcessOptions.cc.

References G4LossTableManager::GetEmProcessVector().

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

00172 {
00173   const std::vector<G4VEmProcess*>& w =
00174         theManager->GetEmProcessVector();
00175   std::vector<G4VEmProcess*>::const_iterator itp;
00176   for(itp = w.begin(); itp != w.end(); itp++) {
00177     G4VEmProcess* q = *itp;
00178     if(q) { q->SetApplyCuts(val); }
00179   }
00180 }

void G4EmProcessOptions::SetAuger ( G4bool  val  ) 

Definition at line 268 of file G4EmProcessOptions.cc.

References G4VAtomDeexcitation::SetAuger().

Referenced by G4EnergyLossMessenger::SetNewValue().

00269 {
00270   G4VAtomDeexcitation* ad = theManager-> AtomDeexcitation();
00271   if(ad) { ad->SetAuger(val); }
00272 }

void G4EmProcessOptions::SetBremsstrahlungTh ( G4double  val  ) 

Definition at line 411 of file G4EmProcessOptions.cc.

References G4LossTableManager::SetBremsstrahlungTh().

00412 {
00413   theManager->SetBremsstrahlungTh(val);
00414 }

void G4EmProcessOptions::SetBuildCSDARange ( G4bool  val  ) 

Definition at line 184 of file G4EmProcessOptions.cc.

References G4LossTableManager::SetBuildCSDARange().

Referenced by G4EnergyLossMessenger::SetNewValue().

00185 {
00186   theManager->SetBuildCSDARange(val);
00187 }

void G4EmProcessOptions::SetDEDXBinning ( G4int  val  ) 

Definition at line 136 of file G4EmProcessOptions.cc.

References G4LossTableManager::SetDEDXBinning().

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

00137 {
00138   theManager->SetDEDXBinning(val);
00139 }

void G4EmProcessOptions::SetDEDXBinningForCSDARange ( G4int  val  ) 

Definition at line 143 of file G4EmProcessOptions.cc.

References G4LossTableManager::SetDEDXBinningForCSDARange().

00144 {
00145   theManager->SetDEDXBinningForCSDARange(val);
00146 }

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

Definition at line 246 of file G4EmProcessOptions.cc.

References G4VAtomDeexcitation::SetDeexcitationActiveRegion().

Referenced by G4EnergyLossMessenger::SetNewValue().

00250 {
00251   G4VAtomDeexcitation* ad = theManager-> AtomDeexcitation();
00252   if(ad) { 
00253     ad->SetDeexcitationActiveRegion(rname, valDeexcitation,
00254                                     valAuger,valPIXE); 
00255   }
00256 }

void G4EmProcessOptions::SetFactorForAngleLimit ( G4double  val  ) 

Definition at line 383 of file G4EmProcessOptions.cc.

References G4LossTableManager::SetFactorForAngleLimit().

Referenced by G4EnergyLossMessenger::SetNewValue().

00384 {
00385   theManager->SetFactorForAngleLimit(val);
00386 }

void G4EmProcessOptions::SetFluo ( G4bool  val  ) 

Definition at line 260 of file G4EmProcessOptions.cc.

References G4VAtomDeexcitation::SetFluo().

Referenced by G4EnergyLossMessenger::SetNewValue().

00261 {
00262   G4VAtomDeexcitation* ad = theManager-> AtomDeexcitation();
00263   if(ad) { ad->SetFluo(val); }
00264 }

void G4EmProcessOptions::SetIntegral ( G4bool  val  ) 

Definition at line 94 of file G4EmProcessOptions.cc.

References G4LossTableManager::SetIntegral().

Referenced by G4EnergyLossMessenger::SetNewValue().

00095 {
00096   theManager->SetIntegral(val);
00097 }

void G4EmProcessOptions::SetLambdaBinning ( G4int  val  ) 

Definition at line 150 of file G4EmProcessOptions.cc.

References G4LossTableManager::SetLambdaBinning().

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

00151 {
00152   theManager->SetLambdaBinning(val);
00153 }

void G4EmProcessOptions::SetLambdaFactor ( G4double  val  ) 

Definition at line 232 of file G4EmProcessOptions.cc.

References G4LossTableManager::GetEnergyLossProcessVector().

Referenced by G4EnergyLossMessenger::SetNewValue().

00233 {
00234   const std::vector<G4VEnergyLossProcess*>& v =
00235         theManager->GetEnergyLossProcessVector();
00236   std::vector<G4VEnergyLossProcess*>::const_iterator itr;
00237   for(itr = v.begin(); itr != v.end(); itr++) {
00238     G4VEnergyLossProcess* p = *itr;
00239     if(p) { p->SetLambdaFactor(val); }
00240   }
00241 }

void G4EmProcessOptions::SetLinearLossLimit ( G4double  val  ) 

Definition at line 404 of file G4EmProcessOptions.cc.

References G4LossTableManager::SetLinearLossLimit().

Referenced by G4EnergyLossMessenger::SetNewValue().

00405 {
00406   theManager->SetLinearLossLimit(val);
00407 }

void G4EmProcessOptions::SetLossFluctuations ( G4bool  val  ) 

Definition at line 80 of file G4EmProcessOptions.cc.

References G4LossTableManager::SetLossFluctuations().

Referenced by G4EnergyLossMessenger::SetNewValue().

00081 {
00082   theManager->SetLossFluctuations(val);
00083 }

void G4EmProcessOptions::SetLPMFlag ( G4bool  val  ) 

Definition at line 390 of file G4EmProcessOptions.cc.

References G4LossTableManager::SetLPMFlag().

Referenced by G4EnergyLossMessenger::SetNewValue().

00391 {
00392   theManager->SetLPMFlag(val);
00393 }

void G4EmProcessOptions::SetMaxEnergy ( G4double  val  ) 

Definition at line 115 of file G4EmProcessOptions.cc.

References G4LossTableManager::SetMaxEnergy().

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

00116 {
00117   theManager->SetMaxEnergy(val);
00118 }

void G4EmProcessOptions::SetMaxEnergyForCSDARange ( G4double  val  ) 

Definition at line 122 of file G4EmProcessOptions.cc.

References G4LossTableManager::SetMaxEnergyForCSDARange().

00123 {
00124   theManager->SetMaxEnergyForCSDARange(val);
00125 }

void G4EmProcessOptions::SetMaxEnergyForMuons ( G4double  val  ) 

Definition at line 129 of file G4EmProcessOptions.cc.

References G4LossTableManager::SetMaxEnergyForMuons().

00130 {
00131   theManager->SetMaxEnergyForMuons(val);
00132 }

void G4EmProcessOptions::SetMinEnergy ( G4double  val  ) 

Definition at line 108 of file G4EmProcessOptions.cc.

References G4LossTableManager::SetMinEnergy().

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

00109 {
00110   theManager->SetMinEnergy(val);
00111 }

void G4EmProcessOptions::SetMinSubRange ( G4double  val  ) 

Definition at line 101 of file G4EmProcessOptions.cc.

References G4LossTableManager::SetMinSubRange().

Referenced by G4EnergyLossMessenger::SetNewValue().

00102 {
00103   theManager->SetMinSubRange(val);
00104 }

void G4EmProcessOptions::SetMscGeomFactor ( G4double  val  ) 

Definition at line 351 of file G4EmProcessOptions.cc.

References G4LossTableManager::GetMultipleScatteringVector().

Referenced by G4EnergyLossMessenger::SetNewValue().

00352 {
00353   if(val < 0.0) { return; }
00354   const std::vector<G4VMultipleScattering*>& u =
00355         theManager->GetMultipleScatteringVector();
00356   std::vector<G4VMultipleScattering*>::const_iterator itm;
00357   for(itm = u.begin(); itm != u.end(); itm++) {
00358     if(*itm) { (*itm)->SetGeomFactor(val); }
00359   }
00360 }

void G4EmProcessOptions::SetMscLateralDisplacement ( G4bool  val  ) 

Definition at line 313 of file G4EmProcessOptions.cc.

References G4LossTableManager::GetMultipleScatteringVector().

Referenced by G4EnergyLossMessenger::SetNewValue().

00314 {
00315   const std::vector<G4VMultipleScattering*>& u =
00316         theManager->GetMultipleScatteringVector();
00317   std::vector<G4VMultipleScattering*>::const_iterator itm;
00318   for(itm = u.begin(); itm != u.end(); itm++) {
00319     if(*itm) { (*itm)->SetLateralDisplasmentFlag(val); }
00320   }
00321 }

void G4EmProcessOptions::SetMscRangeFactor ( G4double  val  ) 

Definition at line 338 of file G4EmProcessOptions.cc.

References G4LossTableManager::GetMultipleScatteringVector().

Referenced by G4EnergyLossMessenger::SetNewValue().

00339 {
00340   if(val < 0.0) return;
00341   const std::vector<G4VMultipleScattering*>& u =
00342         theManager->GetMultipleScatteringVector();
00343   std::vector<G4VMultipleScattering*>::const_iterator itm;
00344   for(itm = u.begin(); itm != u.end(); itm++) {
00345     if(*itm) { (*itm)->SetRangeFactor(val); }
00346   }
00347 }

void G4EmProcessOptions::SetMscStepLimitation ( G4MscStepLimitType  val  ) 

Definition at line 301 of file G4EmProcessOptions.cc.

References G4LossTableManager::GetMultipleScatteringVector().

Referenced by G4EnergyLossMessenger::SetNewValue().

00302 {
00303   const std::vector<G4VMultipleScattering*>& u =
00304         theManager->GetMultipleScatteringVector();
00305   std::vector<G4VMultipleScattering*>::const_iterator itm;
00306   for(itm = u.begin(); itm != u.end(); itm++) {
00307     if(*itm) (*itm)->SetStepLimitType(val);
00308   }
00309 }

void G4EmProcessOptions::SetPIXE ( G4bool  val  ) 

Definition at line 276 of file G4EmProcessOptions.cc.

References G4VAtomDeexcitation::SetPIXE().

Referenced by G4EnergyLossMessenger::SetNewValue().

00277 {
00278   G4VAtomDeexcitation* ad = theManager-> AtomDeexcitation();
00279   if(ad) { ad->SetPIXE(val); }
00280 }

void G4EmProcessOptions::SetPIXECrossSectionModel ( const G4String val  ) 

Definition at line 284 of file G4EmProcessOptions.cc.

References G4VAtomDeexcitation::SetPIXECrossSectionModel().

Referenced by G4EnergyLossMessenger::SetNewValue().

00285 {
00286   G4VAtomDeexcitation* ad = theManager-> AtomDeexcitation();
00287   if(ad) { ad->SetPIXECrossSectionModel(mname); }
00288 }

void G4EmProcessOptions::SetPIXEElectronCrossSectionModel ( const G4String val  ) 

Definition at line 293 of file G4EmProcessOptions.cc.

References G4VAtomDeexcitation::SetPIXEElectronCrossSectionModel().

Referenced by G4EnergyLossMessenger::SetNewValue().

00294 {
00295   G4VAtomDeexcitation* ad = theManager-> AtomDeexcitation();
00296   if(ad) { ad->SetPIXEElectronCrossSectionModel(mname); }
00297 }

void G4EmProcessOptions::SetPolarAngleLimit ( G4double  val  ) 

Definition at line 364 of file G4EmProcessOptions.cc.

References G4LossTableManager::GetEmProcessVector(), and G4LossTableManager::GetMultipleScatteringVector().

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

00365 {
00366   const std::vector<G4VMultipleScattering*>& u =
00367         theManager->GetMultipleScatteringVector();
00368   std::vector<G4VMultipleScattering*>::const_iterator itm;
00369   for(itm = u.begin(); itm != u.end(); itm++) {
00370     if(*itm) { (*itm)->SetPolarAngleLimit(val); }
00371   }
00372   const std::vector<G4VEmProcess*>& w =
00373         theManager->GetEmProcessVector();
00374   std::vector<G4VEmProcess*>::const_iterator itp;
00375   for(itp = w.begin(); itp != w.end(); itp++) {
00376     G4VEmProcess* q = *itp;
00377     if(q) { q->SetPolarAngleLimit(val); }
00378   }
00379 }

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

Definition at line 419 of file G4EmProcessOptions.cc.

References G4LossTableManager::GetEmProcessVector(), and G4LossTableManager::GetEnergyLossProcessVector().

Referenced by G4EnergyLossMessenger::SetNewValue().

00421 {
00422   const std::vector<G4VEnergyLossProcess*>& v =
00423         theManager->GetEnergyLossProcessVector();
00424   std::vector<G4VEnergyLossProcess*>::const_iterator itr;
00425   for(itr = v.begin(); itr != v.end(); ++itr) {
00426     G4VEnergyLossProcess* p = *itr;
00427     if(p) {
00428       if (p->GetProcessName() == name) { 
00429         p->SetCrossSectionBiasingFactor(val, flag); 
00430       }
00431     }
00432   }
00433   const std::vector<G4VEmProcess*>& w =
00434         theManager->GetEmProcessVector();
00435   std::vector<G4VEmProcess*>::const_iterator itp;
00436   for(itp = w.begin(); itp != w.end(); itp++) {
00437     G4VEmProcess* q = *itp;
00438     if(q) {
00439       if (q->GetProcessName() == name) { 
00440         q->SetCrossSectionBiasingFactor(val, flag); 
00441       }
00442     }
00443   }
00444 }

void G4EmProcessOptions::SetRandomStep ( G4bool  val  ) 

Definition at line 164 of file G4EmProcessOptions.cc.

References G4LossTableManager::SetRandomStep().

Referenced by G4EnergyLossMessenger::SetNewValue().

00165 {
00166   theManager->SetRandomStep(val);
00167 }

void G4EmProcessOptions::SetSkin ( G4double  val  ) 

Definition at line 325 of file G4EmProcessOptions.cc.

References G4LossTableManager::GetMultipleScatteringVector().

Referenced by G4EnergyLossMessenger::SetNewValue().

00326 {
00327   if(val < 0.0) return;
00328   const std::vector<G4VMultipleScattering*>& u =
00329         theManager->GetMultipleScatteringVector();
00330   std::vector<G4VMultipleScattering*>::const_iterator itm;
00331   for(itm = u.begin(); itm != u.end(); itm++) {
00332     if(*itm) { (*itm)->SetSkin(val); }
00333   }
00334 }

void G4EmProcessOptions::SetSplineFlag ( G4bool  val  ) 

Definition at line 397 of file G4EmProcessOptions.cc.

References G4LossTableManager::SetSplineFlag().

Referenced by G4EnergyLossMessenger::SetNewValue().

00398 {
00399   theManager->SetSplineFlag(val);
00400 }

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

Definition at line 157 of file G4EmProcessOptions.cc.

References G4LossTableManager::SetStepFunction().

Referenced by G4EnergyLossMessenger::SetNewValue().

00158 {
00159   theManager->SetStepFunction(v1, v2);
00160 }

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

Definition at line 87 of file G4EmProcessOptions.cc.

References G4LossTableManager::SetSubCutoff().

Referenced by G4EnergyLossMessenger::SetNewValue().

00088 {
00089   theManager->SetSubCutoff(val, r);
00090 }

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

Definition at line 191 of file G4EmProcessOptions.cc.

References G4LossTableManager::GetEmProcessVector(), G4LossTableManager::GetEnergyLossProcessVector(), G4LossTableManager::GetMultipleScatteringVector(), and G4LossTableManager::SetVerbose().

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

00192 {
00193   G4bool all = false;
00194   if("all" == name) { all = true; }
00195   const std::vector<G4VEnergyLossProcess*>& v =
00196         theManager->GetEnergyLossProcessVector();
00197 
00198   if(all) { 
00199     theManager->SetVerbose(val);
00200     return;
00201   }
00202 
00203   std::vector<G4VEnergyLossProcess*>::const_iterator itr;
00204   for(itr = v.begin(); itr != v.end(); ++itr) {
00205     G4VEnergyLossProcess* p = *itr;
00206     if(p) {
00207       if (p->GetProcessName() == name) { p->SetVerboseLevel(val); }
00208     }
00209   }
00210   const std::vector<G4VEmProcess*>& w =
00211         theManager->GetEmProcessVector();
00212   std::vector<G4VEmProcess*>::const_iterator itp;
00213   for(itp = w.begin(); itp != w.end(); itp++) {
00214     G4VEmProcess* q = *itp;
00215     if(q) {
00216       if (q->GetProcessName() == name) { q->SetVerboseLevel(val); }
00217     }
00218   }
00219   const std::vector<G4VMultipleScattering*>& u =
00220         theManager->GetMultipleScatteringVector();
00221   std::vector<G4VMultipleScattering*>::const_iterator itm;
00222   for(itm = u.begin(); itm != u.end(); itm++) {
00223     G4VMultipleScattering* msc = *itm;
00224     if(s) {
00225       if (msc->GetProcessName() == name) { msc->SetVerboseLevel(val); }
00226     }
00227   }
00228 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:51:53 2013 for Geant4 by  doxygen 1.4.7