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

#include <G4QSSMessenger.hh>

Inheritance diagram for G4QSSMessenger:

Public Types

enum  StepperSelection { None = 1 , G4QSS2 = 2 , G4QSS3 = 3 , NumMethods }

Public Member Functions

 G4QSSMessenger ()
 ~G4QSSMessenger () override
void SetNewValue (G4UIcommand *command, G4String newValues) override
G4int GetQssOrder ()
G4double Get_dQRel ()
G4double Get_dQMin ()
G4int GetMaxSubsteps ()
void selectStepper (const std::string &)
StepperSelection selectedStepper ()
G4bool SetQssOrder (G4int order)
Public Member Functions inherited from G4UImessenger
 G4UImessenger ()=default
 G4UImessenger (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
virtual ~G4UImessenger ()
virtual G4String GetCurrentValue (G4UIcommand *command)
G4bool CommandsShouldBeInMaster () const

Static Public Member Functions

static G4QSSMessengerinstance ()

Additional Inherited Members

Protected Member Functions inherited from G4UImessenger
G4String ItoS (G4int i)
G4String LtoS (G4long l)
G4String DtoS (G4double a)
G4String BtoS (G4bool b)
G4int StoI (const G4String &s)
G4long StoL (const G4String &s)
G4double StoD (const G4String &s)
G4bool StoB (const 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 = nullptr
G4String baseDirName = ""
G4bool commandsShouldBeInMaster = false

Detailed Description

Definition at line 46 of file G4QSSMessenger.hh.

Member Enumeration Documentation

◆ StepperSelection

Enumerator
None 
G4QSS2 
G4QSS3 
NumMethods 

Definition at line 75 of file G4QSSMessenger.hh.

Constructor & Destructor Documentation

◆ G4QSSMessenger()

G4QSSMessenger::G4QSSMessenger ( )

Constructor and Destructor.

Definition at line 34 of file G4QSSMessenger.cc.

35{
37
38 qssCmdDir = new G4UIdirectory("/QSS/", false);
39 qssCmdDir->SetGuidance("G4QSStepper configuration.");
40
41 dQMinCmd = new G4UIcmdWithADoubleAndUnit("/QSS/dQMin",this);
42 dQMinCmd->SetDefaultUnit("mm");
43 dQMinCmd->SetParameterName("dQMinCmd",false);
44 dQMinCmd->SetUnitCategory("Length");
45
46 dQRelCmd = new G4UIcmdWithADouble("/QSS/dQRel",this);
47 dQRelCmd->SetGuidance("Default is 1e-5");
48 dQRelCmd->SetParameterName("dQRelCmd",false);
49
50 stepperSelectorCmd = new G4UIcmdWithAString("/QSS/selectStepper",this);
51 stepperSelectorCmd->SetGuidance("Select stepper.");
52 stepperSelectorCmd->SetParameterName("choice", false);
53 stepperSelectorCmd->SetCandidates("TemplatedDoPri OldRK45 G4QSS2");
54
55 maxSubstepsCmd = new G4UIcmdWithAnInteger("/QSS/maxSubsteps",this);
56 maxSubstepsCmd->SetGuidance("Default is 5000");
57 maxSubstepsCmd->SetDefaultValue(5000);
58 maxSubstepsCmd->SetParameterName("maxSubstepsCmd", false);
59
60}
G4bool commandsShouldBeInMaster

Referenced by instance().

◆ ~G4QSSMessenger()

G4QSSMessenger::~G4QSSMessenger ( )
override

Definition at line 62 of file G4QSSMessenger.cc.

63{
64 delete qssCmdDir;
65 delete dQMinCmd;
66 delete dQRelCmd;
67 delete stepperSelectorCmd;
68 delete maxSubstepsCmd;
69}

Member Function Documentation

◆ Get_dQMin()

G4double G4QSSMessenger::Get_dQMin ( )
inline

Definition at line 72 of file G4QSSMessenger.hh.

static G4QSSParameters * Instance()

◆ Get_dQRel()

G4double G4QSSMessenger::Get_dQRel ( )
inline

Definition at line 71 of file G4QSSMessenger.hh.

◆ GetMaxSubsteps()

G4int G4QSSMessenger::GetMaxSubsteps ( )
inline

Definition at line 73 of file G4QSSMessenger.hh.

Referenced by SetNewValue(), and G4QSStepper::Stepper().

◆ GetQssOrder()

G4int G4QSSMessenger::GetQssOrder ( )
inline

Accessors.

Definition at line 70 of file G4QSSMessenger.hh.

Referenced by G4ChordFinder::G4ChordFinder().

◆ instance()

G4QSSMessenger * G4QSSMessenger::instance ( )
static

Definition at line 71 of file G4QSSMessenger.cc.

72{
73 static G4QSSMessenger theSingleMessengerInstance;
74 return &theSingleMessengerInstance;
75}

Referenced by G4ChordFinder::G4ChordFinder(), G4QSStepper::G4QSStepper(), and G4QSStepper::Stepper().

◆ selectedStepper()

G4QSSMessenger::StepperSelection G4QSSMessenger::selectedStepper ( )

Definition at line 195 of file G4QSSMessenger.cc.

196{
197 return _selectedStepper;
198}

◆ selectStepper()

void G4QSSMessenger::selectStepper ( const std::string & newValue)

Stepper selection, G4QSS2 or G4QSS3.

Definition at line 136 of file G4QSSMessenger.cc.

137{
138 const std::map<std::string, StepperSelection> stepperMapping =
139 {
140 {"G4QSS2", G4QSS2}, {"G4QSS3", G4QSS3}
141 };
142 _selectedStepper = stepperMapping.at(newValue);
143
144 if( (StepperSelection::None <= _selectedStepper )
145 && (_selectedStepper < StepperSelection::NumMethods ) )
146 {
147 G4cout << "G4QSSMessenger: Selecting stepper " << newValue
148 << " ( which is method # " << _selectedStepper << " ) " << G4endl;
149
150 G4int new_order= 0;
151 if( _selectedStepper == G4QSS2 )
152 {
153 new_order= 2;
154 }
155 else
156 {
157 if( _selectedStepper == G4QSS3 )
158 {
159 new_order= 3;
160 }
161 }
162 SetQssOrder( new_order ); // Will check requested value
163 }
164 else
165 {
166 G4cerr << "G4QSSMessenger::selectStepper : Invalid type of stepper requested"
167 << " Cannot find the requested stepper type '" << newValue << G4endl;
168 }
169}
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4bool SetQssOrder(G4int order)

Referenced by SetNewValue().

◆ SetNewValue()

void G4QSSMessenger::SetNewValue ( G4UIcommand * command,
G4String newValues )
overridevirtual

Applies command to the associated object.

Reimplemented from G4UImessenger.

Definition at line 100 of file G4QSSMessenger.cc.

101{
102 G4bool good_value= true;
103
104 if ( command == dQMinCmd )
105 {
106 good_value= Set_dQMin( dQMinCmd->GetNewDoubleValue(newValue) );
107 }
108
109 if ( command == dQRelCmd )
110 {
111 good_value= Set_dQRel( dQRelCmd->GetNewDoubleValue(newValue) );
112 }
113
114 if (command == maxSubstepsCmd)
115 {
116 G4int oldMaxSubstp = GetMaxSubsteps();
117 good_value= SetMaxSubsteps ( maxSubstepsCmd->GetNewIntValue(newValue) );
118 if( good_value )
119 {
120 G4cout << "G4QSSMessenger: changed maxSubsteps = " << GetMaxSubsteps()
121 << " ( old value = " << oldMaxSubstp << " )." << G4endl;
122 }
123 }
124
125 if( !good_value )
126 {
127 G4cerr << " G4QSSMessenger: Command FAILED - parameter was not accepted." << G4endl;
128 }
129
130 if(command == stepperSelectorCmd)
131 {
132 this->selectStepper(newValue);
133 }
134}
bool G4bool
Definition G4Types.hh:86
void selectStepper(const std::string &)

◆ SetQssOrder()

G4bool G4QSSMessenger::SetQssOrder ( G4int order)

Sets QSS order. To be suppressed in favour of the method it calls in G4QSSParameters.

Definition at line 171 of file G4QSSMessenger.cc.

172{
173 G4bool good_value= false;
174
175 if( 2 <= value && value <= 3 )
176 {
177 if( value == 2)
178 {
179 _selectedStepper = G4QSS2;
180 }
181 else
182 {
183 _selectedStepper = G4QSS3;
184 }
185 good_value= G4QSSParameters::Instance()->SetQssOrder( value );
186 }
187 else
188 {
189 G4cerr << "G4QSSMessenger::SetQssOrder : requested invalid *order* of QSS-stepper"
190 << " Asked for " << value << " . Value must be either 2 or 3." << G4endl;
191 }
192 return good_value;
193}
G4bool SetQssOrder(G4int value, G4bool onlyWarn=false)

Referenced by G4ChordFinder::G4ChordFinder(), G4QSStepper::G4QSStepper(), and selectStepper().


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