Geant4-11
Public Member Functions | Private Member Functions | Private Attributes
G4EmExtraParameters Class Reference

#include <G4EmExtraParameters.hh>

Public Member Functions

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 AddPAIModel (const G4String &particle, const G4String &region, const G4String &type)
 
void AddPhysics (const G4String &region, const G4String &type)
 
void DefineRegParamForEM (G4VEmProcess *) const
 
void DefineRegParamForLoss (G4VEnergyLossProcess *) const
 
void FillStepFunction (const G4ParticleDefinition *, G4VEnergyLossProcess *) const
 
 G4EmExtraParameters ()
 
 G4EmExtraParameters (G4EmExtraParameters &)=delete
 
G4bool GetDirectionalSplitting ()
 
G4double GetDirectionalSplittingRadius ()
 
G4ThreeVector GetDirectionalSplittingTarget () const
 
G4double GetStepFunctionIonsP1 () const
 
G4double GetStepFunctionIonsP2 () const
 
G4double GetStepFunctionLightIonsP1 () const
 
G4double GetStepFunctionLightIonsP2 () const
 
G4double GetStepFunctionMuHadP1 () const
 
G4double GetStepFunctionMuHadP2 () const
 
G4double GetStepFunctionP1 () const
 
G4double GetStepFunctionP2 () const
 
void Initialise ()
 
G4EmExtraParametersoperator= (const G4EmExtraParameters &right)=delete
 
const std::vector< G4String > & ParticlesPAI () const
 
G4bool QuantumEntanglement ()
 
const std::vector< G4String > & RegionsPAI () const
 
const std::vector< G4String > & RegionsPhysics () const
 
void SetDirectionalSplitting (G4bool v)
 
void SetDirectionalSplittingRadius (G4double r)
 
void SetDirectionalSplittingTarget (const G4ThreeVector &v)
 
void SetProcessBiasingFactor (const G4String &procname, G4double val, G4bool wflag)
 
void SetQuantumEntanglement (G4bool v)
 
void SetStepFunction (G4double v1, G4double v2)
 
void SetStepFunctionIons (G4double v1, G4double v2)
 
void SetStepFunctionLightIons (G4double v1, G4double v2)
 
void SetStepFunctionMuHad (G4double v1, G4double v2)
 
void SetSubCutRegion (const G4String &region)
 
const std::vector< G4String > & TypesPAI () const
 
const std::vector< G4String > & TypesPhysics () const
 
 ~G4EmExtraParameters ()
 

Private Member Functions

G4String CheckRegion (const G4String &) const
 
void PrintWarning (G4ExceptionDescription &ed) const
 

Private Attributes

G4bool directionalSplitting
 
G4double directionalSplittingRadius
 
G4ThreeVector directionalSplittingTarget
 
G4double dRoverRange
 
G4double dRoverRangeIons
 
G4double dRoverRangeLIons
 
G4double dRoverRangeMuHad
 
G4double finalRange
 
G4double finalRangeIons
 
G4double finalRangeLIons
 
G4double finalRangeMuHad
 
std::vector< G4doublem_elimBiasedSec
 
std::vector< G4doublem_factBiasedSec
 
std::vector< G4doublem_factBiasedXS
 
std::vector< G4doublem_lengthForced
 
std::vector< G4Stringm_particlesPAI
 
std::vector< G4Stringm_procBiasedSec
 
std::vector< G4Stringm_procBiasedXS
 
std::vector< G4Stringm_procForced
 
std::vector< G4Stringm_regnamesBiasedSec
 
std::vector< G4Stringm_regnamesForced
 
std::vector< G4Stringm_regnamesPAI
 
std::vector< G4Stringm_regnamesPhys
 
std::vector< G4Stringm_regnamesSubCut
 
std::vector< G4Stringm_typesPAI
 
std::vector< G4Stringm_typesPhys
 
std::vector< G4boolm_weightBiasedXS
 
std::vector< G4boolm_weightForced
 
G4bool quantumEntanglement
 
G4EmExtraParametersMessengertheMessenger
 

Detailed Description

Definition at line 62 of file G4EmExtraParameters.hh.

Constructor & Destructor Documentation

◆ G4EmExtraParameters() [1/2]

G4EmExtraParameters::G4EmExtraParameters ( )
explicit

◆ ~G4EmExtraParameters()

G4EmExtraParameters::~G4EmExtraParameters ( )

Definition at line 63 of file G4EmExtraParameters.cc.

64{
65 delete theMessenger;
66}

References theMessenger.

◆ G4EmExtraParameters() [2/2]

G4EmExtraParameters::G4EmExtraParameters ( G4EmExtraParameters )
delete

Member Function Documentation

◆ ActivateForcedInteraction()

void G4EmExtraParameters::ActivateForcedInteraction ( const G4String procname,
const G4String region,
G4double  length,
G4bool  wflag 
)

Definition at line 321 of file G4EmExtraParameters.cc.

325{
326 const G4String& r = CheckRegion(region);
327 if(length >= 0.0) {
328 G4int n = m_procForced.size();
329 for(G4int 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}
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
int G4int
Definition: G4Types.hh:85
void PrintWarning(G4ExceptionDescription &ed) const
std::vector< G4double > m_lengthForced
std::vector< G4String > m_regnamesForced
std::vector< G4String > m_procForced
G4String CheckRegion(const G4String &) const
std::vector< G4bool > m_weightForced

References CheckRegion(), m_lengthForced, m_procForced, m_regnamesForced, m_weightForced, CLHEP::detail::n, and PrintWarning().

Referenced by G4EmParameters::ActivateForcedInteraction(), and G4EmExtraParametersMessenger::SetNewValue().

◆ ActivateSecondaryBiasing()

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

Definition at line 350 of file G4EmExtraParameters.cc.

354{
355 const G4String& r = CheckRegion(region);
356 if(factor >= 0.0 && energyLim >= 0.0) {
357 G4int n = m_procBiasedSec.size();
358 for(G4int 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}
std::vector< G4String > m_procBiasedSec
std::vector< G4double > m_factBiasedSec
std::vector< G4double > m_elimBiasedSec
std::vector< G4String > m_regnamesBiasedSec

References CheckRegion(), m_elimBiasedSec, m_factBiasedSec, m_procBiasedSec, m_regnamesBiasedSec, CLHEP::detail::n, and PrintWarning().

Referenced by G4EmParameters::ActivateSecondaryBiasing(), and G4EmExtraParametersMessenger::SetNewValue().

◆ AddPAIModel()

void G4EmExtraParameters::AddPAIModel ( const G4String particle,
const G4String region,
const G4String type 
)

Definition at line 222 of file G4EmExtraParameters.cc.

225{
226 G4String r = CheckRegion(region);
227 G4int nreg = m_regnamesPAI.size();
228 for(G4int 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}
std::vector< G4String > m_typesPAI
std::vector< G4String > m_particlesPAI
std::vector< G4String > m_regnamesPAI

References CheckRegion(), m_particlesPAI, m_regnamesPAI, and m_typesPAI.

Referenced by G4EmParameters::AddPAIModel(), and G4EmExtraParametersMessenger::SetNewValue().

◆ AddPhysics()

void G4EmExtraParameters::AddPhysics ( const G4String region,
const G4String type 
)

Definition at line 262 of file G4EmExtraParameters.cc.

264{
265 G4String r = CheckRegion(region);
266 G4int nreg = m_regnamesPhys.size();
267 for(G4int 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}
std::vector< G4String > m_typesPhys
std::vector< G4String > m_regnamesPhys

References CheckRegion(), m_regnamesPhys, and m_typesPhys.

Referenced by G4EmParameters::AddPhysics(), and G4EmExtraParametersMessenger::SetNewValue().

◆ CheckRegion()

G4String G4EmExtraParameters::CheckRegion ( const G4String reg) const
private

Definition at line 101 of file G4EmExtraParameters.cc.

102{
103 G4String r = reg;
104 if(r == "" || r == "world" || r == "World") {
105 r = "DefaultRegionForTheWorld";
106 }
107 return r;
108}
static const G4double reg

References reg.

Referenced by ActivateForcedInteraction(), ActivateSecondaryBiasing(), AddPAIModel(), AddPhysics(), and SetSubCutRegion().

◆ DefineRegParamForEM()

void G4EmExtraParameters::DefineRegParamForEM ( G4VEmProcess ptr) const

Definition at line 414 of file G4EmExtraParameters.cc.

415{
416 G4int n = m_procBiasedXS.size();
417 for(G4int i=0; i<n; ++i) {
418 if(ptr->GetProcessName() == m_procBiasedXS[i]) {
421 break;
422 }
423 }
424 n = m_procForced.size();
425 for(G4int i=0; i<n; ++i) {
426 if(ptr->GetProcessName() == m_procForced[i]) {
429 m_weightForced[i]);
430 break;
431 }
432 }
433 n = m_procBiasedSec.size();
434 for(G4int i=0; i<n; ++i) {
435 if(ptr->GetProcessName() == m_procBiasedSec[i]) {
437 m_factBiasedSec[i],
438 m_elimBiasedSec[i]);
439 break;
440 }
441 }
442}
std::vector< G4String > m_procBiasedXS
std::vector< G4double > m_factBiasedXS
std::vector< G4bool > m_weightBiasedXS
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)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382

References G4VEmProcess::ActivateForcedInteraction(), G4VEmProcess::ActivateSecondaryBiasing(), G4VProcess::GetProcessName(), m_elimBiasedSec, m_factBiasedSec, m_factBiasedXS, m_lengthForced, m_procBiasedSec, m_procBiasedXS, m_procForced, m_regnamesBiasedSec, m_regnamesForced, m_weightBiasedXS, m_weightForced, CLHEP::detail::n, and G4VEmProcess::SetCrossSectionBiasingFactor().

Referenced by G4EmParameters::DefineRegParamForEM().

◆ DefineRegParamForLoss()

void G4EmExtraParameters::DefineRegParamForLoss ( G4VEnergyLossProcess ptr) const

Definition at line 378 of file G4EmExtraParameters.cc.

379{
380 const G4RegionStore* regionStore = G4RegionStore::GetInstance();
381 G4int n = m_regnamesSubCut.size();
382 for(G4int 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(G4int i=0; i<n; ++i) {
388 if(ptr->GetProcessName() == m_procBiasedXS[i]) {
391 break;
392 }
393 }
394 n = m_procForced.size();
395 for(G4int i=0; i<n; ++i) {
396 if(ptr->GetProcessName() == m_procForced[i]) {
399 m_weightForced[i]);
400 break;
401 }
402 }
403 n = m_procBiasedSec.size();
404 for(G4int i=0; i<n; ++i) {
405 if(ptr->GetProcessName() == m_procBiasedSec[i]) {
407 m_factBiasedSec[i],
408 m_elimBiasedSec[i]);
409 break;
410 }
411 }
412}
std::vector< G4String > m_regnamesSubCut
static G4RegionStore * GetInstance()
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
void SetCrossSectionBiasingFactor(G4double f, G4bool flag=true)
void ActivateForcedInteraction(G4double length, const G4String &region, G4bool flag=true)
void ActivateSubCutoff(const G4Region *region)

References G4VEnergyLossProcess::ActivateForcedInteraction(), G4VEnergyLossProcess::ActivateSecondaryBiasing(), G4VEnergyLossProcess::ActivateSubCutoff(), G4RegionStore::GetInstance(), G4VProcess::GetProcessName(), G4RegionStore::GetRegion(), m_elimBiasedSec, m_factBiasedSec, m_factBiasedXS, m_lengthForced, m_procBiasedSec, m_procBiasedXS, m_procForced, m_regnamesBiasedSec, m_regnamesForced, m_regnamesSubCut, m_weightBiasedXS, m_weightForced, CLHEP::detail::n, reg, and G4VEnergyLossProcess::SetCrossSectionBiasingFactor().

Referenced by G4EmParameters::DefineRegParamForLoss().

◆ FillStepFunction()

void G4EmExtraParameters::FillStepFunction ( const G4ParticleDefinition part,
G4VEnergyLossProcess proc 
) const

Definition at line 202 of file G4EmExtraParameters.cc.

203{
204 // electron and positron
205 if (11 == std::abs(part->GetPDGEncoding())) {
207
208 // all heavy ions
209 } else if (part->IsGeneralIon()) {
211
212 // light nucleus and anti-nucleus
213 } else if (part->GetParticleType() == "nucleus" || part->GetParticleType() == "anti_nucleus") {
215
216 // other particles
217 } else {
219 }
220}
G4bool IsGeneralIon() const
const G4String & GetParticleType() const
void SetStepFunction(G4double v1, G4double v2)

References dRoverRange, dRoverRangeIons, dRoverRangeLIons, dRoverRangeMuHad, finalRange, finalRangeIons, finalRangeLIons, finalRangeMuHad, G4ParticleDefinition::GetParticleType(), G4ParticleDefinition::GetPDGEncoding(), G4ParticleDefinition::IsGeneralIon(), and G4VEnergyLossProcess::SetStepFunction().

Referenced by G4EmParameters::FillStepFunction().

◆ GetDirectionalSplitting()

G4bool G4EmExtraParameters::GetDirectionalSplitting ( )

Definition at line 454 of file G4EmExtraParameters.cc.

454 {
455 return directionalSplitting;
456}

References directionalSplitting.

Referenced by G4EmParameters::GetDirectionalSplitting().

◆ GetDirectionalSplittingRadius()

G4double G4EmExtraParameters::GetDirectionalSplittingRadius ( )

◆ GetDirectionalSplittingTarget()

G4ThreeVector G4EmExtraParameters::GetDirectionalSplittingTarget ( ) const

Definition at line 469 of file G4EmExtraParameters.cc.

470{
472}
G4ThreeVector directionalSplittingTarget

References directionalSplittingTarget.

Referenced by G4EmParameters::GetDirectionalSplittingTarget().

◆ GetStepFunctionIonsP1()

G4double G4EmExtraParameters::GetStepFunctionIonsP1 ( ) const

Definition at line 192 of file G4EmExtraParameters.cc.

193{
194 return dRoverRangeIons;
195}

References dRoverRangeIons.

Referenced by G4EmParameters::StreamInfo().

◆ GetStepFunctionIonsP2()

G4double G4EmExtraParameters::GetStepFunctionIonsP2 ( ) const

Definition at line 197 of file G4EmExtraParameters.cc.

198{
199 return finalRangeIons;
200}

References finalRangeIons.

Referenced by G4EmParameters::StreamInfo().

◆ GetStepFunctionLightIonsP1()

G4double G4EmExtraParameters::GetStepFunctionLightIonsP1 ( ) const

Definition at line 169 of file G4EmExtraParameters.cc.

170{
171 return dRoverRangeLIons;
172}

References dRoverRangeLIons.

Referenced by G4EmParameters::StreamInfo().

◆ GetStepFunctionLightIonsP2()

G4double G4EmExtraParameters::GetStepFunctionLightIonsP2 ( ) const

Definition at line 174 of file G4EmExtraParameters.cc.

175{
176 return finalRangeLIons;
177}

References finalRangeLIons.

Referenced by G4EmParameters::StreamInfo().

◆ GetStepFunctionMuHadP1()

G4double G4EmExtraParameters::GetStepFunctionMuHadP1 ( ) const

Definition at line 146 of file G4EmExtraParameters.cc.

147{
148 return dRoverRangeMuHad;
149}

References dRoverRangeMuHad.

Referenced by G4EmParameters::StreamInfo().

◆ GetStepFunctionMuHadP2()

G4double G4EmExtraParameters::GetStepFunctionMuHadP2 ( ) const

Definition at line 151 of file G4EmExtraParameters.cc.

152{
153 return finalRangeMuHad;
154}

References finalRangeMuHad.

Referenced by G4EmParameters::StreamInfo().

◆ GetStepFunctionP1()

G4double G4EmExtraParameters::GetStepFunctionP1 ( ) const

Definition at line 123 of file G4EmExtraParameters.cc.

124{
125 return dRoverRange;
126}

References dRoverRange.

Referenced by G4EmParameters::StreamInfo().

◆ GetStepFunctionP2()

G4double G4EmExtraParameters::GetStepFunctionP2 ( ) const

Definition at line 128 of file G4EmExtraParameters.cc.

129{
130 return finalRange;
131}

References finalRange.

Referenced by G4EmParameters::StreamInfo().

◆ Initialise()

void G4EmExtraParameters::Initialise ( )

◆ operator=()

G4EmExtraParameters & G4EmExtraParameters::operator= ( const G4EmExtraParameters right)
delete

◆ ParticlesPAI()

const std::vector< G4String > & G4EmExtraParameters::ParticlesPAI ( ) const

Definition at line 247 of file G4EmExtraParameters.cc.

248{
249 return m_particlesPAI;
250}

References m_particlesPAI.

Referenced by G4EmParameters::ParticlesPAI().

◆ PrintWarning()

void G4EmExtraParameters::PrintWarning ( G4ExceptionDescription ed) const
private

Definition at line 96 of file G4EmExtraParameters.cc.

97{
98 G4Exception("G4EmExtraParameters", "em0044", JustWarning, ed);
99}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35

References G4Exception(), and JustWarning.

Referenced by ActivateForcedInteraction(), ActivateSecondaryBiasing(), SetProcessBiasingFactor(), SetStepFunction(), SetStepFunctionIons(), SetStepFunctionLightIons(), and SetStepFunctionMuHad().

◆ QuantumEntanglement()

G4bool G4EmExtraParameters::QuantumEntanglement ( )

Definition at line 444 of file G4EmExtraParameters.cc.

445{
446 return quantumEntanglement;
447}

References quantumEntanglement.

Referenced by G4EmParameters::QuantumEntanglement(), and G4EmParameters::StreamInfo().

◆ RegionsPAI()

const std::vector< G4String > & G4EmExtraParameters::RegionsPAI ( ) const

Definition at line 252 of file G4EmExtraParameters.cc.

253{
254 return m_regnamesPAI;
255}

References m_regnamesPAI.

Referenced by G4EmParameters::RegionsPAI().

◆ RegionsPhysics()

const std::vector< G4String > & G4EmExtraParameters::RegionsPhysics ( ) const

Definition at line 274 of file G4EmExtraParameters.cc.

275{
276 return m_regnamesPhys;
277}

References m_regnamesPhys.

Referenced by G4EmParameters::RegionsPhysics().

◆ SetDirectionalSplitting()

void G4EmExtraParameters::SetDirectionalSplitting ( G4bool  v)

◆ SetDirectionalSplittingRadius()

void G4EmExtraParameters::SetDirectionalSplittingRadius ( G4double  r)

◆ SetDirectionalSplittingTarget()

void G4EmExtraParameters::SetDirectionalSplittingTarget ( const G4ThreeVector v)

◆ SetProcessBiasingFactor()

void G4EmExtraParameters::SetProcessBiasingFactor ( const G4String procname,
G4double  val,
G4bool  wflag 
)

Definition at line 297 of file G4EmExtraParameters.cc.

299{
300 if(val > 0.0) {
301 G4int n = m_procBiasedXS.size();
302 for(G4int 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}

References m_factBiasedXS, m_procBiasedXS, m_weightBiasedXS, CLHEP::detail::n, and PrintWarning().

Referenced by G4EmExtraParametersMessenger::SetNewValue(), and G4EmParameters::SetProcessBiasingFactor().

◆ SetQuantumEntanglement()

void G4EmExtraParameters::SetQuantumEntanglement ( G4bool  v)

◆ SetStepFunction()

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

Definition at line 110 of file G4EmExtraParameters.cc.

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}

References dRoverRange, finalRange, CLHEP::mm, and PrintWarning().

Referenced by G4EmExtraParametersMessenger::SetNewValue(), and G4EmParameters::SetStepFunction().

◆ SetStepFunctionIons()

void G4EmExtraParameters::SetStepFunctionIons ( G4double  v1,
G4double  v2 
)

Definition at line 179 of file G4EmExtraParameters.cc.

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}

References dRoverRangeIons, finalRangeIons, CLHEP::mm, and PrintWarning().

Referenced by G4EmExtraParametersMessenger::SetNewValue(), and G4EmParameters::SetStepFunctionIons().

◆ SetStepFunctionLightIons()

void G4EmExtraParameters::SetStepFunctionLightIons ( G4double  v1,
G4double  v2 
)

Definition at line 156 of file G4EmExtraParameters.cc.

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}

References dRoverRangeLIons, finalRangeLIons, CLHEP::mm, and PrintWarning().

Referenced by G4EmExtraParametersMessenger::SetNewValue(), and G4EmParameters::SetStepFunctionLightIons().

◆ SetStepFunctionMuHad()

void G4EmExtraParameters::SetStepFunctionMuHad ( G4double  v1,
G4double  v2 
)

Definition at line 133 of file G4EmExtraParameters.cc.

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}

References dRoverRangeMuHad, finalRangeMuHad, CLHEP::mm, and PrintWarning().

Referenced by G4EmExtraParametersMessenger::SetNewValue(), and G4EmParameters::SetStepFunctionMuHad().

◆ SetSubCutRegion()

void G4EmExtraParameters::SetSubCutRegion ( const G4String region)

Definition at line 284 of file G4EmExtraParameters.cc.

285{
286 const G4String& r = CheckRegion(region);
287 G4int nreg = m_regnamesSubCut.size();
288 for(G4int i=0; i<nreg; ++i) {
289 if(r == m_regnamesSubCut[i]) {
290 return;
291 }
292 }
293 m_regnamesSubCut.push_back(r);
294}

References CheckRegion(), and m_regnamesSubCut.

Referenced by G4EmExtraParametersMessenger::SetNewValue(), and G4EmParameters::SetSubCutRegion().

◆ TypesPAI()

const std::vector< G4String > & G4EmExtraParameters::TypesPAI ( ) const

Definition at line 257 of file G4EmExtraParameters.cc.

258{
259 return m_typesPAI;
260}

References m_typesPAI.

Referenced by G4EmParameters::TypesPAI().

◆ TypesPhysics()

const std::vector< G4String > & G4EmExtraParameters::TypesPhysics ( ) const

Definition at line 279 of file G4EmExtraParameters.cc.

280{
281 return m_typesPhys;
282}

References m_typesPhys.

Referenced by G4EmParameters::TypesPhysics().

Field Documentation

◆ directionalSplitting

G4bool G4EmExtraParameters::directionalSplitting
private

◆ directionalSplittingRadius

G4double G4EmExtraParameters::directionalSplittingRadius
private

◆ directionalSplittingTarget

G4ThreeVector G4EmExtraParameters::directionalSplittingTarget
private

◆ dRoverRange

G4double G4EmExtraParameters::dRoverRange
private

◆ dRoverRangeIons

G4double G4EmExtraParameters::dRoverRangeIons
private

◆ dRoverRangeLIons

G4double G4EmExtraParameters::dRoverRangeLIons
private

◆ dRoverRangeMuHad

G4double G4EmExtraParameters::dRoverRangeMuHad
private

◆ finalRange

G4double G4EmExtraParameters::finalRange
private

◆ finalRangeIons

G4double G4EmExtraParameters::finalRangeIons
private

◆ finalRangeLIons

G4double G4EmExtraParameters::finalRangeLIons
private

◆ finalRangeMuHad

G4double G4EmExtraParameters::finalRangeMuHad
private

◆ m_elimBiasedSec

std::vector<G4double> G4EmExtraParameters::m_elimBiasedSec
private

◆ m_factBiasedSec

std::vector<G4double> G4EmExtraParameters::m_factBiasedSec
private

◆ m_factBiasedXS

std::vector<G4double> G4EmExtraParameters::m_factBiasedXS
private

◆ m_lengthForced

std::vector<G4double> G4EmExtraParameters::m_lengthForced
private

◆ m_particlesPAI

std::vector<G4String> G4EmExtraParameters::m_particlesPAI
private

Definition at line 160 of file G4EmExtraParameters.hh.

Referenced by AddPAIModel(), and ParticlesPAI().

◆ m_procBiasedSec

std::vector<G4String> G4EmExtraParameters::m_procBiasedSec
private

◆ m_procBiasedXS

std::vector<G4String> G4EmExtraParameters::m_procBiasedXS
private

◆ m_procForced

std::vector<G4String> G4EmExtraParameters::m_procForced
private

◆ m_regnamesBiasedSec

std::vector<G4String> G4EmExtraParameters::m_regnamesBiasedSec
private

◆ m_regnamesForced

std::vector<G4String> G4EmExtraParameters::m_regnamesForced
private

◆ m_regnamesPAI

std::vector<G4String> G4EmExtraParameters::m_regnamesPAI
private

Definition at line 161 of file G4EmExtraParameters.hh.

Referenced by AddPAIModel(), and RegionsPAI().

◆ m_regnamesPhys

std::vector<G4String> G4EmExtraParameters::m_regnamesPhys
private

Definition at line 164 of file G4EmExtraParameters.hh.

Referenced by AddPhysics(), and RegionsPhysics().

◆ m_regnamesSubCut

std::vector<G4String> G4EmExtraParameters::m_regnamesSubCut
private

Definition at line 167 of file G4EmExtraParameters.hh.

Referenced by DefineRegParamForLoss(), Initialise(), and SetSubCutRegion().

◆ m_typesPAI

std::vector<G4String> G4EmExtraParameters::m_typesPAI
private

Definition at line 162 of file G4EmExtraParameters.hh.

Referenced by AddPAIModel(), and TypesPAI().

◆ m_typesPhys

std::vector<G4String> G4EmExtraParameters::m_typesPhys
private

Definition at line 165 of file G4EmExtraParameters.hh.

Referenced by AddPhysics(), and TypesPhysics().

◆ m_weightBiasedXS

std::vector<G4bool> G4EmExtraParameters::m_weightBiasedXS
private

◆ m_weightForced

std::vector<G4bool> G4EmExtraParameters::m_weightForced
private

◆ quantumEntanglement

G4bool G4EmExtraParameters::quantumEntanglement
private

Definition at line 146 of file G4EmExtraParameters.hh.

Referenced by Initialise(), QuantumEntanglement(), and SetQuantumEntanglement().

◆ theMessenger

G4EmExtraParametersMessenger* G4EmExtraParameters::theMessenger
private

Definition at line 143 of file G4EmExtraParameters.hh.

Referenced by G4EmExtraParameters(), and ~G4EmExtraParameters().


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