Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MatScanMessenger.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// G4MatScanMessenger implementation
27//
28// Author: M.Asai, May 2006
29// --------------------------------------------------------------------
30
31#include "G4MatScanMessenger.hh"
32
33#include "G4MaterialScanner.hh"
34#include "G4SystemOfUnits.hh"
35#include "G4ThreeVector.hh"
36#include "G4Tokenizer.hh"
37#include "G4UIcmdWith3Vector.hh"
39#include "G4UIcmdWithABool.hh"
40#include "G4UIcmdWithAString.hh"
42#include "G4UIcommand.hh"
43#include "G4UIdirectory.hh"
44#include "G4UIparameter.hh"
45
46// --------------------------------------------------------------------
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("/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("Usage : /control/matScan/theta [nbin] [thetaMin] [thetaSpan] [unit]");
63 thetaCmd->SetGuidance("Notation of angles :");
64 thetaCmd->SetGuidance(" theta --- +Z axis : +90 deg. / X-Y plane : 0 deg. / -Z axis : -90 deg.");
65 par = new G4UIparameter("nbin", 'i', false);
66 par->SetParameterRange("nbin>0");
67 thetaCmd->SetParameter(par);
68 par = new G4UIparameter("thetaMin", 'd', false);
69 thetaCmd->SetParameter(par);
70 par = new G4UIparameter("thetaSpan", 'd', true);
71 par->SetParameterRange("thetaSpan>=0.");
72 par->SetDefaultValue(0.);
73 thetaCmd->SetParameter(par);
74 par = new G4UIparameter("unit", 'c', true);
75 par->SetDefaultValue("deg");
76 par->SetParameterCandidates(thetaCmd->UnitsList(thetaCmd->CategoryOf("deg")));
77 thetaCmd->SetParameter(par);
78
79 phiCmd = new G4UIcommand("/control/matScan/phi", this);
80 phiCmd->SetGuidance("Define phi range.");
81 phiCmd->SetGuidance("Usage : /control/matScan/phi [nbin] [phiMin] [phiSpan] [unit]");
82 phiCmd->SetGuidance("Notation of angles :");
83 phiCmd->SetGuidance(
84 " phi --- +X axis : 0 deg. / +Y axis : 90 deg. / -X axis : 180 "
85 "deg. / -Y axis : 270 deg.");
86 par = new G4UIparameter("nbin", 'i', false);
87 par->SetParameterRange("nbin>0");
88 phiCmd->SetParameter(par);
89 par = new G4UIparameter("phiMin", 'd', false);
90 phiCmd->SetParameter(par);
91 par = new G4UIparameter("phiSpan", 'd', true);
92 par->SetParameterRange("phiSpan>=0.");
93 par->SetDefaultValue(0.);
94 phiCmd->SetParameter(par);
95 par = new G4UIparameter("unit", 'c', true);
96 par->SetDefaultValue("deg");
97 par->SetParameterCandidates(phiCmd->UnitsList(phiCmd->CategoryOf("deg")));
98 phiCmd->SetParameter(par);
99
100 singleCmd = new G4UIcommand("/control/matScan/singleMeasure", this);
101 singleCmd->SetGuidance("Measure thickness for one particular direction.");
102 singleCmd->SetGuidance("Notation of angles :");
103 singleCmd->SetGuidance(" theta --- +Z axis : +90 deg. / X-Y plane : 0 deg. / -Z axis : -90 deg.");
104 singleCmd->SetGuidance(
105 " phi --- +X axis : 0 deg. / +Y axis : 90 deg. / -X axis : "
106 "180 deg. / -Y axis : 270 deg.");
108 par = new G4UIparameter("theta", 'd', false);
109 singleCmd->SetParameter(par);
110 par = new G4UIparameter("phi", 'd', false);
111 singleCmd->SetParameter(par);
112 par = new G4UIparameter("unit", 'c', true);
113 par->SetDefaultValue("deg");
114 par->SetParameterCandidates(singleCmd->UnitsList(singleCmd->CategoryOf("deg")));
115 singleCmd->SetParameter(par);
116
117 single2Cmd = new G4UIcmdWith3Vector("/control/matScan/singleTo", this);
118 single2Cmd->SetGuidance("Measure thickness for one direction defined by a unit vector.");
119 single2Cmd->SetParameterName("X", "Y", "Z", false);
120
121 eyePosCmd = new G4UIcmdWith3VectorAndUnit("/control/matScan/eyePosition", this);
122 eyePosCmd->SetGuidance("Define the eye position.");
123 eyePosCmd->SetParameterName("X", "Y", "Z", true);
124 eyePosCmd->SetDefaultValue(G4ThreeVector(0., 0., 0.));
125 eyePosCmd->SetDefaultUnit("m");
126
127 regSenseCmd = new G4UIcmdWithABool("/control/matScan/regionSensitive", this);
128 regSenseCmd->SetGuidance("Set region sensitivity.");
129 regSenseCmd->SetGuidance("This command is automatically set to TRUE");
130 regSenseCmd->SetGuidance(" if /control/matScan/region command is issued.");
131 regSenseCmd->SetParameterName("senseFlag", true);
132 regSenseCmd->SetDefaultValue(false);
133
134 regionCmd = new G4UIcmdWithAString("/control/matScan/region", this);
135 regionCmd->SetGuidance("Define region name to be scanned.");
136 regionCmd->SetGuidance("/control/matScan/regionSensitive command is automatically");
137 regionCmd->SetGuidance("set to TRUE with this command.");
138 regionCmd->SetParameterName("region", true);
139 regionCmd->SetDefaultValue("DefaultRegionForTheWorld");
140}
141
142// --------------------------------------------------------------------
144{
145 delete scanCmd;
146 delete thetaCmd;
147 delete phiCmd;
148 delete singleCmd;
149 delete single2Cmd;
150 delete eyePosCmd;
151 delete regSenseCmd;
152 delete regionCmd;
153 delete msDirectory;
154}
155
156// --------------------------------------------------------------------
158{
159 G4String currentValue;
160 if (command == thetaCmd) {
161 currentValue = thetaCmd->ConvertToString(theScanner->GetNTheta());
162 currentValue += " ";
163 currentValue += thetaCmd->ConvertToString((theScanner->GetThetaMin()) / deg);
164 currentValue += " ";
165 currentValue += thetaCmd->ConvertToString((theScanner->GetThetaSpan()) / deg);
166 }
167 else if (command == phiCmd) {
168 currentValue = phiCmd->ConvertToString(theScanner->GetNPhi());
169 currentValue += " ";
170 currentValue += phiCmd->ConvertToString((theScanner->GetPhiMin()) / deg);
171 currentValue += " ";
172 currentValue += phiCmd->ConvertToString((theScanner->GetPhiSpan()) / deg);
173 }
174 else if (command == eyePosCmd) {
175 currentValue = eyePosCmd->ConvertToString(theScanner->GetEyePosition(), "m");
176 }
177 else if (command == regSenseCmd) {
178 currentValue = regSenseCmd->ConvertToString(theScanner->GetRegionSensitive());
179 }
180 else if (command == regionCmd) {
181 currentValue = theScanner->GetRegionName();
182 }
183 return currentValue;
184}
185
186// --------------------------------------------------------------------
188{
189 if (command == scanCmd) {
190 theScanner->Scan();
191 }
192 else if (command == thetaCmd) {
193 G4Tokenizer next(newValue);
194 G4int nbin = StoI(next());
195 G4double thetaMin = StoD(next());
196 G4double thetaSpan = StoD(next());
197 G4String unit = next();
198 thetaMin *= thetaCmd->ValueOf(unit);
199 thetaSpan *= thetaCmd->ValueOf(unit);
200 theScanner->SetNTheta(nbin);
201 theScanner->SetThetaMin(thetaMin);
202 theScanner->SetThetaSpan(thetaSpan);
203 }
204 else if (command == phiCmd) {
205 G4Tokenizer next(newValue);
206 G4int nbin = StoI(next());
207 G4double phiMin = StoD(next());
208 G4double phiSpan = StoD(next());
209 G4String unit = next();
210 phiMin *= phiCmd->ValueOf(unit);
211 phiSpan *= phiCmd->ValueOf(unit);
212 theScanner->SetNPhi(nbin);
213 theScanner->SetPhiMin(phiMin);
214 theScanner->SetPhiSpan(phiSpan);
215 }
216 else if (command == eyePosCmd) {
217 theScanner->SetEyePosition(eyePosCmd->GetNew3VectorValue(newValue));
218 }
219 else if (command == regSenseCmd) {
220 theScanner->SetRegionSensitive(regSenseCmd->GetNewBoolValue(newValue));
221 }
222 else if (command == regionCmd) {
223 if (theScanner->SetRegionName(newValue)) theScanner->SetRegionSensitive(true);
224 }
225 else if (command == singleCmd || command == single2Cmd) {
226 G4int ntheta = theScanner->GetNTheta();
227 G4double thetaMin = theScanner->GetThetaMin();
228 G4double thetaSpan = theScanner->GetThetaSpan();
229 G4int nphi = theScanner->GetNPhi();
230 G4double phiMin = theScanner->GetPhiMin();
231 G4double phiSpan = theScanner->GetPhiSpan();
232
233 G4double theta = 0.;
234 G4double phi = 0.;
235 if (command == singleCmd) {
236 G4Tokenizer next(newValue);
237 theta = StoD(next());
238 phi = StoD(next());
239 G4String unit = next();
240 theta *= singleCmd->ValueOf(unit);
241 phi *= singleCmd->ValueOf(unit);
242 }
243 else if (command == single2Cmd) {
244 G4ThreeVector v = single2Cmd->GetNew3VectorValue(newValue);
245 theta = 90. * deg - v.theta();
246 phi = v.phi();
247 }
248 theScanner->SetNTheta(1);
249 theScanner->SetThetaMin(theta);
250 theScanner->SetThetaSpan(0.);
251 theScanner->SetNPhi(1);
252 theScanner->SetPhiMin(phi);
253 theScanner->SetPhiSpan(0.);
254 theScanner->Scan();
255
256 theScanner->SetNTheta(ntheta);
257 theScanner->SetThetaMin(thetaMin);
258 theScanner->SetThetaSpan(thetaSpan);
259 theScanner->SetNPhi(nphi);
260 theScanner->SetPhiMin(phiMin);
261 theScanner->SetPhiSpan(phiSpan);
262 }
263}
@ G4State_Idle
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
double phi() const
double theta() const
void SetNewValue(G4UIcommand *command, G4String newValue) override
G4MatScanMessenger(G4MaterialScanner *p1)
G4String GetCurrentValue(G4UIcommand *command) override
void SetNPhi(G4int val)
G4bool GetRegionSensitive() const
void SetThetaSpan(G4double val)
G4bool SetRegionName(const G4String &val)
G4double GetThetaMin() const
void SetRegionSensitive(G4bool val=true)
void SetThetaMin(G4double val)
G4double GetThetaSpan() const
void SetPhiMin(G4double val)
void SetEyePosition(const G4ThreeVector &val)
void SetPhiSpan(G4double val)
G4ThreeVector GetEyePosition() const
G4double GetPhiMin() const
const G4String & GetRegionName() const
G4int GetNTheta() const
G4double GetPhiSpan() const
void SetNTheta(G4int val)
void SetDefaultUnit(const char *defUnit)
static G4ThreeVector GetNew3VectorValue(const char *paramString)
void SetParameterName(const char *theNameX, const char *theNameY, const char *theNameZ, G4bool omittable, G4bool currentAsDefault=false)
void SetDefaultValue(const G4ThreeVector &defVal)
static G4ThreeVector GetNew3VectorValue(const char *paramString)
void SetParameterName(const char *theNameX, const char *theNameY, const char *theNameZ, G4bool omittable, G4bool currentAsDefault=false)
static G4bool GetNewBoolValue(const char *paramString)
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)
static G4double ValueOf(const char *unitName)
static G4String ConvertToString(G4bool boolVal)
void SetParameter(G4UIparameter *const newParameter)
void SetGuidance(const char *aGuidance)
static G4String UnitsList(const char *unitCategory)
void AvailableForStates(G4ApplicationState s1)
G4double StoD(const G4String &s)
G4int StoI(const G4String &s)
void SetDefaultValue(const char *theDefaultValue)
void SetParameterRange(const char *theRange)
void SetParameterCandidates(const char *theString)