Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
G4AdjointPhysicsMessenger Class Reference

#include <Doxymodules_biasing.h>

Inheritance diagram for G4AdjointPhysicsMessenger:
G4UImessenger

Public Member Functions

 G4AdjointPhysicsMessenger (G4AdjointPhysicsList *)
 
virtual ~G4AdjointPhysicsMessenger ()
 
virtual void SetNewValue (G4UIcommand *, G4String)
 

Private Attributes

G4AdjointPhysicsListfPhysicsList
 
G4UIdirectoryfPhysicsDir
 
G4UIcmdWithABoolfUsepIonisationCmd
 
G4UIcmdWithABoolfUseBremCmd
 
G4UIcmdWithABoolfUseComptonCmd
 
G4UIcmdWithABoolfUseMSCmd
 
G4UIcmdWithABoolfUsePEEffectCmd
 
G4UIcmdWithABoolfUseGammaConversionCmd
 
G4UIcmdWithABoolfUseEgainFluctuationCmd
 
G4UIcmdWithADoubleAndUnitfSetEminAdjModelsCmd
 
G4UIcmdWithADoubleAndUnitfSetEmaxAdjModelsCmd
 

Detailed Description

Definition at line 168 of file Doxymodules_biasing.h.

Constructor & Destructor Documentation

◆ G4AdjointPhysicsMessenger()

G4AdjointPhysicsMessenger::G4AdjointPhysicsMessenger ( G4AdjointPhysicsList pPhysicsList)

Definition at line 55 of file G4AdjointPhysicsMessenger.cc.

58 fPhysicsList(pPhysicsList),
59 fPhysicsDir(0),
61 fUseBremCmd(0),
63 fUseMSCmd(0),
69{
70 fPhysicsDir = new G4UIdirectory("/adjoint_physics/");
71
72 fPhysicsDir->SetGuidance(
73 "Definition of adjoint and forward physics processes");
74 //-------
76 "/adjoint_physics/UseProtonIonisation",this);
77 fUsepIonisationCmd->SetGuidance(
78 "If true (false) the proton ionisation is (not) considered");
79 fUsepIonisationCmd->AvailableForStates(G4State_PreInit);
80
81 fUseBremCmd = new G4UIcmdWithABool("/adjoint_physics/UseBremsstrahlung",this);
82 fUseBremCmd->SetGuidance(
83 "If true (false) the bremsstrahlung process is (not) considered");
84 fUseBremCmd->AvailableForStates(G4State_PreInit);
85
86 fUseComptonCmd = new G4UIcmdWithABool("/adjoint_physics/UseCompton",this);
87 fUseComptonCmd->SetGuidance(
88 "If true (false) the Compton scattering is (not) considered");
89 fUseComptonCmd->AvailableForStates(G4State_PreInit);
90
91 fUseMSCmd = new G4UIcmdWithABool("/adjoint_physics/UseMS",this);
92 fUseMSCmd->SetGuidance(
93 "If true (false) the continuous multiple scattering is (not) considered");
94 fUseMSCmd->AvailableForStates(G4State_PreInit);
95
97 "/adjoint_physics/UseEgainElossFluctuation",this);
98 fUseEgainFluctuationCmd->SetGuidance(
99 "Switch on/off the fluctation for continuous energy gain/loss");
100 fUseEgainFluctuationCmd->AvailableForStates(G4State_PreInit);
101
102 fUsePEEffectCmd = new G4UIcmdWithABool("/adjoint_physics/UsePEEffect",this);
103 fUsePEEffectCmd->AvailableForStates(G4State_PreInit);
104 fUsePEEffectCmd->SetGuidance(
105 "If true (false) the photo electric effect is (not) considered");
106
108 "/adjoint_physics/UseGammaConversion",this);
109 fUseGammaConversionCmd->AvailableForStates(G4State_PreInit);
110 fUseGammaConversionCmd->SetGuidance(
111 "If true the fwd gamma pair conversion is considered");
112
114 "/adjoint_physics/SetEminForAdjointModels",this);
115 fSetEminAdjModelsCmd->SetGuidance(
116 "Set the minimum energy of the adjoint models");
117 fSetEminAdjModelsCmd->SetParameterName("Emin",false);
118 fSetEminAdjModelsCmd->SetUnitCategory("Energy");
119 fSetEminAdjModelsCmd->AvailableForStates(G4State_PreInit);
120
122 "/adjoint_physics/SetEmaxForAdjointModels",this);
123 fSetEmaxAdjModelsCmd->SetGuidance(
124 "Set the minimum energy of the adjoint models.");
125 fSetEmaxAdjModelsCmd->SetParameterName("Emax",false);
126 fSetEmaxAdjModelsCmd->SetUnitCategory("Energy");
127 fSetEmaxAdjModelsCmd->AvailableForStates(G4State_PreInit);
128}
G4UIcmdWithADoubleAndUnit * fSetEmaxAdjModelsCmd
G4UIcmdWithADoubleAndUnit * fSetEminAdjModelsCmd

◆ ~G4AdjointPhysicsMessenger()

G4AdjointPhysicsMessenger::~G4AdjointPhysicsMessenger ( )
virtual

Definition at line 132 of file G4AdjointPhysicsMessenger.cc.

133{
134 delete fUsepIonisationCmd;
135 delete fUseBremCmd;
136 delete fUseComptonCmd;
137 delete fUseMSCmd;
138 delete fUsePEEffectCmd;
143}

Member Function Documentation

◆ SetNewValue()

void G4AdjointPhysicsMessenger::SetNewValue ( G4UIcommand command,
G4String  newValue 
)
virtual

Definition at line 147 of file G4AdjointPhysicsMessenger.cc.

149{
150 if ( command==fUsepIonisationCmd){
152 fUsepIonisationCmd->GetNewBoolValue(newValue));
153
154 }
155 else if ( command==fUseBremCmd){
156 fPhysicsList->SetUseBrem(fUseBremCmd->GetNewBoolValue(newValue));
157 }
158 else if ( command==fUseComptonCmd){
159 fPhysicsList->SetUseCompton(fUseComptonCmd->GetNewBoolValue(newValue));
160 }
161 else if ( command==fUseMSCmd){
162 fPhysicsList->SetUseMS(fUseMSCmd->GetNewBoolValue(newValue));
163 }
164 else if ( command==fUsePEEffectCmd){
165 fPhysicsList->SetUsePEEffect(fUsePEEffectCmd->GetNewBoolValue(newValue));
166 }
167 else if ( command==fUseGammaConversionCmd){
169 fUseGammaConversionCmd->GetNewBoolValue(newValue));
170 }
171 else if ( command==fUseEgainFluctuationCmd){
173 fUseEgainFluctuationCmd->GetNewBoolValue(newValue));
174 }
175
176 else if ( command== fSetEminAdjModelsCmd){
178 fSetEminAdjModelsCmd->GetNewDoubleValue(newValue));
179 }
180 else if ( command== fSetEmaxAdjModelsCmd){
182 fSetEmaxAdjModelsCmd->GetNewDoubleValue(newValue));
183 }
184}
void SetEminAdjModels(G4double aVal)
void SetEmaxAdjModels(G4double aVal)
void SetUsePEEffect(bool aBool)
void SetUseEgainFluctuation(bool aBool)
void SetUseCompton(bool aBool)
void SetUseProtonIonisation(bool aBool)
void SetUseGammaConversion(bool aBool)

Member Data Documentation

◆ fPhysicsList

G4AdjointPhysicsList* G4AdjointPhysicsMessenger::fPhysicsList
private

Definition at line 71 of file G4AdjointPhysicsMessenger.hh.

◆ fPhysicsDir

G4UIdirectory* G4AdjointPhysicsMessenger::fPhysicsDir
private

Definition at line 72 of file G4AdjointPhysicsMessenger.hh.

◆ fUsepIonisationCmd

G4UIcmdWithABool* G4AdjointPhysicsMessenger::fUsepIonisationCmd
private

Definition at line 75 of file G4AdjointPhysicsMessenger.hh.

◆ fUseBremCmd

G4UIcmdWithABool* G4AdjointPhysicsMessenger::fUseBremCmd
private

Definition at line 76 of file G4AdjointPhysicsMessenger.hh.

◆ fUseComptonCmd

G4UIcmdWithABool* G4AdjointPhysicsMessenger::fUseComptonCmd
private

Definition at line 77 of file G4AdjointPhysicsMessenger.hh.

◆ fUseMSCmd

G4UIcmdWithABool* G4AdjointPhysicsMessenger::fUseMSCmd
private

Definition at line 78 of file G4AdjointPhysicsMessenger.hh.

◆ fUsePEEffectCmd

G4UIcmdWithABool* G4AdjointPhysicsMessenger::fUsePEEffectCmd
private

Definition at line 79 of file G4AdjointPhysicsMessenger.hh.

◆ fUseGammaConversionCmd

G4UIcmdWithABool* G4AdjointPhysicsMessenger::fUseGammaConversionCmd
private

Definition at line 80 of file G4AdjointPhysicsMessenger.hh.

◆ fUseEgainFluctuationCmd

G4UIcmdWithABool* G4AdjointPhysicsMessenger::fUseEgainFluctuationCmd
private

Definition at line 81 of file G4AdjointPhysicsMessenger.hh.

◆ fSetEminAdjModelsCmd

G4UIcmdWithADoubleAndUnit* G4AdjointPhysicsMessenger::fSetEminAdjModelsCmd
private

Definition at line 82 of file G4AdjointPhysicsMessenger.hh.

◆ fSetEmaxAdjModelsCmd

G4UIcmdWithADoubleAndUnit* G4AdjointPhysicsMessenger::fSetEmaxAdjModelsCmd
private

Definition at line 83 of file G4AdjointPhysicsMessenger.hh.


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

Applications | User Support | Publications | Collaboration