Geant4 10.7.0
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 ()
 
 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 (G4String s)
 
G4long StoL (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 = nullptr
 
G4String baseDirName = ""
 
G4bool commandsShouldBeInMaster = false
 

Detailed Description

Definition at line 46 of file G4MatScanMessenger.hh.

Constructor & Destructor Documentation

◆ G4MatScanMessenger()

G4MatScanMessenger::G4MatScanMessenger ( G4MaterialScanner p1)

Definition at line 46 of file G4MatScanMessenger.cc.

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

◆ ~G4MatScanMessenger()

G4MatScanMessenger::~G4MatScanMessenger ( )
virtual

Definition at line 150 of file G4MatScanMessenger.cc.

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

Member Function Documentation

◆ GetCurrentValue()

G4String G4MatScanMessenger::GetCurrentValue ( G4UIcommand command)
virtual

Reimplemented from G4UImessenger.

Definition at line 163 of file G4MatScanMessenger.cc.

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

◆ SetNewValue()

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

Reimplemented from G4UImessenger.

Definition at line 201 of file G4MatScanMessenger.cc.

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

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