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

#include <G4MatScanMessenger.hh>

+ Inheritance diagram for G4MatScanMessenger:

Public Member Functions

 G4MatScanMessenger (G4MaterialScanner *p1)
 
virtual ~G4MatScanMessenger ()
 
virtual G4String GetCurrentValue (G4UIcommand *command)
 
virtual void SetNewValue (G4UIcommand *command, G4String newValue)
 
- 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)
 
virtual void SetNewValue (G4UIcommand *command, G4String newValue)
 
G4bool operator== (const G4UImessenger &messenger) const
 
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 (const G4String &s)
 
G4long StoL (const G4String &s)
 
G4double StoD (const 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 = nullptr
 
G4String baseDirName = ""
 
G4bool commandsShouldBeInMaster = false
 

Detailed Description

Definition at line 48 of file G4MatScanMessenger.hh.

Constructor & Destructor Documentation

◆ G4MatScanMessenger()

G4MatScanMessenger::G4MatScanMessenger ( G4MaterialScanner p1)

Definition at line 47 of file G4MatScanMessenger.cc.

48{
49 theScanner = p1;
50 G4UIparameter* par;
51 msDirectory = new G4UIdirectory("/control/matScan/");
52 msDirectory->SetGuidance("Material scanner commands.");
53
54 scanCmd = new G4UIcmdWithoutParameter("/control/matScan/scan", this);
55 scanCmd->SetGuidance("Start material scanning.");
56 scanCmd->SetGuidance("Scanning range should be defined with");
57 scanCmd->SetGuidance(
58 "/control/matScan/theta and /control/matSca/phi commands.");
60
61 thetaCmd = new G4UIcommand("/control/matScan/theta", this);
62 thetaCmd->SetGuidance("Define theta range.");
63 thetaCmd->SetGuidance(
64 "Usage : /control/matScan/theta [nbin] [thetaMin] [thetaSpan] [unit]");
65 thetaCmd->SetGuidance("Notation of angles :");
66 thetaCmd->SetGuidance(
67 " theta --- +Z axis : +90 deg. / X-Y plane : 0 deg. / -Z axis : -90 deg.");
68 par = new G4UIparameter("nbin", 'i', false);
69 par->SetParameterRange("nbin>0");
70 thetaCmd->SetParameter(par);
71 par = new G4UIparameter("thetaMin", 'd', false);
72 thetaCmd->SetParameter(par);
73 par = new G4UIparameter("thetaSpan", 'd', true);
74 par->SetParameterRange("thetaSpan>=0.");
75 par->SetDefaultValue(0.);
76 thetaCmd->SetParameter(par);
77 par = new G4UIparameter("unit", 'c', true);
78 par->SetDefaultValue("deg");
79 par->SetParameterCandidates(thetaCmd->UnitsList(thetaCmd->CategoryOf("deg")));
80 thetaCmd->SetParameter(par);
81
82 phiCmd = new G4UIcommand("/control/matScan/phi", this);
83 phiCmd->SetGuidance("Define phi range.");
84 phiCmd->SetGuidance(
85 "Usage : /control/matScan/phi [nbin] [phiMin] [phiSpan] [unit]");
86 phiCmd->SetGuidance("Notation of angles :");
87 phiCmd->SetGuidance(
88 " phi --- +X axis : 0 deg. / +Y axis : 90 deg. / -X axis : 180 "
89 "deg. / -Y axis : 270 deg.");
90 par = new G4UIparameter("nbin", 'i', false);
91 par->SetParameterRange("nbin>0");
92 phiCmd->SetParameter(par);
93 par = new G4UIparameter("phiMin", 'd', false);
94 phiCmd->SetParameter(par);
95 par = new G4UIparameter("phiSpan", 'd', true);
96 par->SetParameterRange("phiSpan>=0.");
97 par->SetDefaultValue(0.);
98 phiCmd->SetParameter(par);
99 par = new G4UIparameter("unit", 'c', true);
100 par->SetDefaultValue("deg");
101 par->SetParameterCandidates(phiCmd->UnitsList(phiCmd->CategoryOf("deg")));
102 phiCmd->SetParameter(par);
103
104 singleCmd = new G4UIcommand("/control/matScan/singleMeasure", this);
105 singleCmd->SetGuidance("Measure thickness for one particular direction.");
106 singleCmd->SetGuidance("Notation of angles :");
107 singleCmd->SetGuidance(
108 " theta --- +Z axis : +90 deg. / X-Y plane : 0 deg. / -Z axis : -90 deg.");
109 singleCmd->SetGuidance(
110 " phi --- +X axis : 0 deg. / +Y axis : 90 deg. / -X axis : "
111 "180 deg. / -Y axis : 270 deg.");
113 par = new G4UIparameter("theta", 'd', false);
114 singleCmd->SetParameter(par);
115 par = new G4UIparameter("phi", 'd', false);
116 singleCmd->SetParameter(par);
117 par = new G4UIparameter("unit", 'c', true);
118 par->SetDefaultValue("deg");
120 singleCmd->UnitsList(singleCmd->CategoryOf("deg")));
121 singleCmd->SetParameter(par);
122
123 single2Cmd = new G4UIcmdWith3Vector("/control/matScan/singleTo", this);
124 single2Cmd->SetGuidance(
125 "Measure thickness for one direction defined by a unit vector.");
126 single2Cmd->SetParameterName("X", "Y", "Z", false);
127
128 eyePosCmd =
129 new G4UIcmdWith3VectorAndUnit("/control/matScan/eyePosition", this);
130 eyePosCmd->SetGuidance("Define the eye position.");
131 eyePosCmd->SetParameterName("X", "Y", "Z", true);
132 eyePosCmd->SetDefaultValue(G4ThreeVector(0., 0., 0.));
133 eyePosCmd->SetDefaultUnit("m");
134
135 regSenseCmd = new G4UIcmdWithABool("/control/matScan/regionSensitive", this);
136 regSenseCmd->SetGuidance("Set region sensitivity.");
137 regSenseCmd->SetGuidance("This command is automatically set to TRUE");
138 regSenseCmd->SetGuidance(" if /control/matScan/region command is issued.");
139 regSenseCmd->SetParameterName("senseFlag", true);
140 regSenseCmd->SetDefaultValue(false);
141
142 regionCmd = new G4UIcmdWithAString("/control/matScan/region", this);
143 regionCmd->SetGuidance("Define region name to be scanned.");
144 regionCmd->SetGuidance(
145 "/control/matScan/regionSensitive command is automatically");
146 regionCmd->SetGuidance("set to TRUE with this command.");
147 regionCmd->SetParameterName("region", true);
148 regionCmd->SetDefaultValue("DefaultRegionForTheWorld");
149}
@ G4State_Idle
CLHEP::Hep3Vector G4ThreeVector
void SetDefaultUnit(const char *defUnit)
void SetParameterName(const char *theNameX, const char *theNameY, const char *theNameZ, G4bool omittable, G4bool currentAsDefault=false)
void SetDefaultValue(const G4ThreeVector &defVal)
void SetParameterName(const char *theNameX, const char *theNameY, const char *theNameZ, G4bool omittable, G4bool currentAsDefault=false)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetDefaultValue(G4bool defVal)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetDefaultValue(const char *defVal)
static G4String CategoryOf(const char *unitName)
Definition: G4UIcommand.cc:370
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:147
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:157
static G4String UnitsList(const char *unitCategory)
Definition: G4UIcommand.cc:376
void AvailableForStates(G4ApplicationState s1)
Definition: G4UIcommand.cc:287
void SetDefaultValue(const char *theDefaultValue)
void SetParameterRange(const char *theRange)
void SetParameterCandidates(const char *theString)

◆ ~G4MatScanMessenger()

G4MatScanMessenger::~G4MatScanMessenger ( )
virtual

Definition at line 152 of file G4MatScanMessenger.cc.

153{
154 delete scanCmd;
155 delete thetaCmd;
156 delete phiCmd;
157 delete singleCmd;
158 delete single2Cmd;
159 delete eyePosCmd;
160 delete regSenseCmd;
161 delete regionCmd;
162 delete msDirectory;
163}

Member Function Documentation

◆ GetCurrentValue()

G4String G4MatScanMessenger::GetCurrentValue ( G4UIcommand command)
virtual

Reimplemented from G4UImessenger.

Definition at line 166 of file G4MatScanMessenger.cc.

167{
168 G4String currentValue;
169 if(command == thetaCmd)
170 {
171 currentValue = thetaCmd->ConvertToString(theScanner->GetNTheta());
172 currentValue += " ";
173 currentValue +=
174 thetaCmd->ConvertToString((theScanner->GetThetaMin()) / deg);
175 currentValue += " ";
176 currentValue +=
177 thetaCmd->ConvertToString((theScanner->GetThetaSpan()) / deg);
178 }
179 else if(command == phiCmd)
180 {
181 currentValue = phiCmd->ConvertToString(theScanner->GetNPhi());
182 currentValue += " ";
183 currentValue += phiCmd->ConvertToString((theScanner->GetPhiMin()) / deg);
184 currentValue += " ";
185 currentValue += phiCmd->ConvertToString((theScanner->GetPhiSpan()) / deg);
186 }
187 else if(command == eyePosCmd)
188 {
189 currentValue =
190 eyePosCmd->ConvertToString(theScanner->GetEyePosition(), "m");
191 }
192 else if(command == regSenseCmd)
193 {
194 currentValue =
195 regSenseCmd->ConvertToString(theScanner->GetRegionSensitive());
196 }
197 else if(command == regionCmd)
198 {
199 currentValue = theScanner->GetRegionName();
200 }
201 return currentValue;
202}
G4bool GetRegionSensitive() const
G4double GetThetaMin() const
G4double GetThetaSpan() const
G4int GetNPhi() const
G4ThreeVector GetEyePosition() const
G4double GetPhiMin() const
const G4String & GetRegionName() const
G4int GetNTheta() const
G4double GetPhiSpan() const
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:446

◆ SetNewValue()

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

Reimplemented from G4UImessenger.

Definition at line 205 of file G4MatScanMessenger.cc.

206{
207 if(command == scanCmd)
208 {
209 theScanner->Scan();
210 }
211 else if(command == thetaCmd)
212 {
213 G4Tokenizer next(newValue);
214 G4int nbin = StoI(next());
215 G4double thetaMin = StoD(next());
216 G4double thetaSpan = StoD(next());
217 G4String unit = next();
218 thetaMin *= thetaCmd->ValueOf(unit);
219 thetaSpan *= thetaCmd->ValueOf(unit);
220 theScanner->SetNTheta(nbin);
221 theScanner->SetThetaMin(thetaMin);
222 theScanner->SetThetaSpan(thetaSpan);
223 }
224 else if(command == phiCmd)
225 {
226 G4Tokenizer next(newValue);
227 G4int nbin = StoI(next());
228 G4double phiMin = StoD(next());
229 G4double phiSpan = StoD(next());
230 G4String unit = next();
231 phiMin *= phiCmd->ValueOf(unit);
232 phiSpan *= phiCmd->ValueOf(unit);
233 theScanner->SetNPhi(nbin);
234 theScanner->SetPhiMin(phiMin);
235 theScanner->SetPhiSpan(phiSpan);
236 }
237 else if(command == eyePosCmd)
238 {
239 theScanner->SetEyePosition(eyePosCmd->GetNew3VectorValue(newValue));
240 }
241 else if(command == regSenseCmd)
242 {
243 theScanner->SetRegionSensitive(regSenseCmd->GetNewBoolValue(newValue));
244 }
245 else if(command == regionCmd)
246 {
247 if(theScanner->SetRegionName(newValue))
248 theScanner->SetRegionSensitive(true);
249 }
250 else if(command == singleCmd || command == single2Cmd)
251 {
252 G4int ntheta = theScanner->GetNTheta();
253 G4double thetaMin = theScanner->GetThetaMin();
254 G4double thetaSpan = theScanner->GetThetaSpan();
255 G4int nphi = theScanner->GetNPhi();
256 G4double phiMin = theScanner->GetPhiMin();
257 G4double phiSpan = theScanner->GetPhiSpan();
258
259 G4double theta = 0.;
260 G4double phi = 0.;
261 if(command == singleCmd)
262 {
263 G4Tokenizer next(newValue);
264 theta = StoD(next());
265 phi = StoD(next());
266 G4String unit = next();
267 theta *= singleCmd->ValueOf(unit);
268 phi *= singleCmd->ValueOf(unit);
269 }
270 else if(command == single2Cmd)
271 {
272 G4ThreeVector v = single2Cmd->GetNew3VectorValue(newValue);
273 theta = 90. * deg - v.theta();
274 phi = v.phi();
275 }
276 theScanner->SetNTheta(1);
277 theScanner->SetThetaMin(theta);
278 theScanner->SetThetaSpan(0.);
279 theScanner->SetNPhi(1);
280 theScanner->SetPhiMin(phi);
281 theScanner->SetPhiSpan(0.);
282 theScanner->Scan();
283
284 theScanner->SetNTheta(ntheta);
285 theScanner->SetThetaMin(thetaMin);
286 theScanner->SetThetaSpan(thetaSpan);
287 theScanner->SetNPhi(nphi);
288 theScanner->SetPhiMin(phiMin);
289 theScanner->SetPhiSpan(phiSpan);
290 }
291}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
double phi() const
double theta() const
void SetNPhi(G4int val)
void SetThetaSpan(G4double val)
G4bool SetRegionName(const G4String &val)
void SetRegionSensitive(G4bool val=true)
void SetThetaMin(G4double val)
void SetPhiMin(G4double val)
void SetEyePosition(const G4ThreeVector &val)
void SetPhiSpan(G4double val)
void SetNTheta(G4int val)
static G4ThreeVector GetNew3VectorValue(const char *paramString)
static G4ThreeVector GetNew3VectorValue(const char *paramString)
static G4bool GetNewBoolValue(const char *paramString)
static G4double ValueOf(const char *unitName)
Definition: G4UIcommand.cc:362
G4double StoD(const G4String &s)
G4int StoI(const G4String &s)

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