Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QSSMessenger.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4QSSMessenger implementation
27//
28// Author: Leandro Gomez Vidal (Univ. Buenos Aires), October 2021
29// --------------------------------------------------------------------
30
31#include "G4QSSMessenger.hh"
32#include "G4QSSParameters.hh"
33
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}
61
63{
64 delete qssCmdDir;
65 delete dQMinCmd;
66 delete dQRelCmd;
67 delete stepperSelectorCmd;
68 delete maxSubstepsCmd;
69}
70
72{
73 static G4QSSMessenger theSingleMessengerInstance;
74 return &theSingleMessengerInstance;
75}
76
77G4bool G4QSSMessenger::Set_dQMin( G4double dvalue )
78{
79 G4bool good_value= true;
80 good_value= G4QSSParameters::Instance()->Set_dQMin( dvalue );
81 return good_value;
82}
83
84G4bool G4QSSMessenger::Set_dQRel( G4double dvalue )
85{
86 G4bool good_value= true;
87 good_value= G4QSSParameters::Instance()->Set_dQRel( dvalue );
88 return good_value;
89}
90
91G4bool G4QSSMessenger::SetMaxSubsteps( G4int value )
92{
93 G4bool good_value= true;
94
95 good_value = G4QSSParameters::Instance()->SetMaxSubsteps( value );
96
97 return good_value;
98}
99
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}
135
136void G4QSSMessenger::selectStepper(const std::string &newValue)
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}
170
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}
194
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
StepperSelection selectedStepper()
G4bool SetQssOrder(G4int order)
void SetNewValue(G4UIcommand *command, G4String newValues) override
void selectStepper(const std::string &)
static G4QSSMessenger * instance()
~G4QSSMessenger() override
static G4QSSParameters * Instance()
G4bool Set_dQMin(G4double dQMin)
G4bool SetQssOrder(G4int value, G4bool onlyWarn=false)
G4bool Set_dQRel(G4double dQRel)
G4bool SetMaxSubsteps(G4int maxSubsteps)
G4bool commandsShouldBeInMaster