Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4AdjointPhysicsMessenger Class Reference

#include <G4AdjointPhysicsMessenger.hh>

Inheritance diagram for G4AdjointPhysicsMessenger:
Collaboration diagram for G4AdjointPhysicsMessenger:

Public Member Functions

 G4AdjointPhysicsMessenger (G4AdjointPhysicsList *)
 
virtual ~G4AdjointPhysicsMessenger ()
 
virtual void SetNewValue (G4UIcommand *, G4String)
 
- Public Member Functions inherited from G4UImessenger
 G4UImessenger ()
 
 G4UImessenger (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
virtual ~G4UImessenger ()
 
virtual G4String GetCurrentValue (G4UIcommand *command)
 
G4bool operator== (const G4UImessenger &messenger) const
 
G4bool CommandsShouldBeInMaster () const
 

Additional Inherited Members

- Protected Member Functions inherited from G4UImessenger
G4String ItoS (G4int i)
 
G4String DtoS (G4double a)
 
G4String BtoS (G4bool b)
 
G4int StoI (G4String s)
 
G4double StoD (G4String s)
 
G4bool StoB (G4String s)
 
void AddUIcommand (G4UIcommand *newCommand)
 
void CreateDirectory (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
template<typename T >
T * CreateCommand (const G4String &cname, const G4String &dsc)
 
- Protected Attributes inherited from G4UImessenger
G4UIdirectorybaseDir
 
G4String baseDirName
 
G4bool commandsShouldBeInMaster
 

Detailed Description

Definition at line 64 of file G4AdjointPhysicsMessenger.hh.

Constructor & Destructor Documentation

G4AdjointPhysicsMessenger::G4AdjointPhysicsMessenger ( G4AdjointPhysicsList pPhysicsList)

Definition at line 56 of file G4AdjointPhysicsMessenger.cc.

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

Here is the call graph for this function:

G4AdjointPhysicsMessenger::~G4AdjointPhysicsMessenger ( )
virtual

Definition at line 133 of file G4AdjointPhysicsMessenger.cc.

134 {
135  delete fUsepIonisationCmd;
136  delete fUseBremCmd;
137  delete fUseComptonCmd;
138  delete fUseMSCmd;
139  delete fUsePEEffectCmd;
140  delete fUseGammaConversionCmd;
141  delete fUseEgainFluctuationCmd;
142  delete fSetEminAdjModelsCmd;
143  delete fSetEmaxAdjModelsCmd;
144 }

Member Function Documentation

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

Reimplemented from G4UImessenger.

Definition at line 148 of file G4AdjointPhysicsMessenger.cc.

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

Here is the call graph for this function:


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