Geant4 11.4.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)
 ~G4MatScanMessenger () override
G4String GetCurrentValue (G4UIcommand *command) override
void SetNewValue (G4UIcommand *command, G4String newValue) override
Public Member Functions inherited from G4UImessenger
 G4UImessenger ()=default
 G4UImessenger (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
virtual ~G4UImessenger ()
G4bool CommandsShouldBeInMaster () const

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 49 of file G4MatScanMessenger.hh.

Constructor & Destructor Documentation

◆ G4MatScanMessenger()

G4MatScanMessenger::G4MatScanMessenger ( G4MaterialScanner * p1)

Definition at line 48 of file G4MatScanMessenger.cc.

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

◆ ~G4MatScanMessenger()

G4MatScanMessenger::~G4MatScanMessenger ( )
override

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)
overridevirtual

Reimplemented from G4UImessenger.

Definition at line 166 of file G4MatScanMessenger.cc.

167{
168 G4String currentValue;
169 if (command == thetaCmd) {
170 currentValue = thetaCmd->ConvertToString(theScanner->GetNTheta());
171 currentValue += " ";
172 currentValue += thetaCmd->ConvertToString((theScanner->GetThetaMin()) / deg);
173 currentValue += " ";
174 currentValue += thetaCmd->ConvertToString((theScanner->GetThetaSpan()) / deg);
175 }
176 else if (command == phiCmd) {
177 currentValue = phiCmd->ConvertToString(theScanner->GetNPhi());
178 currentValue += " ";
179 currentValue += phiCmd->ConvertToString((theScanner->GetPhiMin()) / deg);
180 currentValue += " ";
181 currentValue += phiCmd->ConvertToString((theScanner->GetPhiSpan()) / deg);
182 }
183 else if (command == eyePosCmd) {
184 currentValue = eyePosCmd->ConvertToString(theScanner->GetEyePosition(), "m");
185 }
186 else if (command == regSenseCmd) {
187 currentValue = regSenseCmd->ConvertToString(theScanner->GetRegionSensitive());
188 }
189 else if (command == regionCmd) {
190 currentValue = theScanner->GetRegionName();
191 }
192 return currentValue;
193}

◆ SetNewValue()

void G4MatScanMessenger::SetNewValue ( G4UIcommand * command,
G4String newValue )
overridevirtual

Reimplemented from G4UImessenger.

Definition at line 196 of file G4MatScanMessenger.cc.

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

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