Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes
G4MatScanMessenger Class Reference

#include <G4MatScanMessenger.hh>

Inheritance diagram for G4MatScanMessenger:
G4UImessenger

Public Member Functions

G4bool CommandsShouldBeInMaster () const
 
 G4MatScanMessenger (G4MaterialScanner *p1)
 
virtual G4String GetCurrentValue (G4UIcommand *command)
 
G4bool operator!= (const G4UImessenger &messenger) const
 
G4bool operator== (const G4UImessenger &messenger) const
 
virtual void SetNewValue (G4UIcommand *command, G4String newValue)
 
virtual ~G4MatScanMessenger ()
 

Protected Member Functions

void AddUIcommand (G4UIcommand *newCommand)
 
G4String BtoS (G4bool b)
 
template<typename T >
T * CreateCommand (const G4String &cname, const G4String &dsc)
 
void CreateDirectory (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
G4String DtoS (G4double a)
 
G4String ItoS (G4int i)
 
G4bool StoB (G4String s)
 
G4double StoD (G4String s)
 
G4int StoI (G4String s)
 
G4long StoL (G4String s)
 

Protected Attributes

G4UIdirectorybaseDir = nullptr
 
G4String baseDirName = ""
 
G4bool commandsShouldBeInMaster = false
 

Private Attributes

G4UIcmdWith3VectorAndUniteyePosCmd = nullptr
 
G4UIdirectorymsDirectory = nullptr
 
G4UIcommandphiCmd = nullptr
 
G4UIcmdWithAStringregionCmd = nullptr
 
G4UIcmdWithABoolregSenseCmd = nullptr
 
G4UIcmdWithoutParameterscanCmd = nullptr
 
G4UIcmdWith3Vectorsingle2Cmd = nullptr
 
G4UIcommandsingleCmd = nullptr
 
G4MaterialScannertheScanner = nullptr
 
G4UIcommandthetaCmd = nullptr
 

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");
58 "/control/matScan/theta and /control/matSca/phi commands.");
60
61 thetaCmd = new G4UIcommand("/control/matScan/theta", this);
62 thetaCmd->SetGuidance("Define theta range.");
64 "Usage : /control/matScan/theta [nbin] [thetaMin] [thetaSpan] [unit]");
65 thetaCmd->SetGuidance("Notation of angles :");
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");
71 par = new G4UIparameter("thetaMin", 'd', false);
73 par = new G4UIparameter("thetaSpan", 'd', true);
74 par->SetParameterRange("thetaSpan>=0.");
75 par->SetDefaultValue(0.);
77 par = new G4UIparameter("unit", 'c', true);
78 par->SetDefaultValue("deg");
81
82 phiCmd = new G4UIcommand("/control/matScan/phi", this);
83 phiCmd->SetGuidance("Define phi range.");
85 "Usage : /control/matScan/phi [nbin] [phiMin] [phiSpan] [unit]");
86 phiCmd->SetGuidance("Notation of angles :");
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");
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 :");
108 " theta --- +Z axis : +90 deg. / X-Y plane : 0 deg. / -Z axis : -90 deg.");
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);
115 par = new G4UIparameter("phi", 'd', false);
117 par = new G4UIparameter("unit", 'c', true);
118 par->SetDefaultValue("deg");
122
123 single2Cmd = new G4UIcmdWith3Vector("/control/matScan/singleTo", this);
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);
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);
141
142 regionCmd = new G4UIcmdWithAString("/control/matScan/region", this);
143 regionCmd->SetGuidance("Define region name to be scanned.");
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
G4UIcmdWithABool * regSenseCmd
G4UIcmdWithoutParameter * scanCmd
G4UIdirectory * msDirectory
G4UIcmdWith3Vector * single2Cmd
G4UIcmdWith3VectorAndUnit * eyePosCmd
G4UIcmdWithAString * regionCmd
G4MaterialScanner * theScanner
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:371
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:377
void AvailableForStates(G4ApplicationState s1)
Definition: G4UIcommand.cc:288
void SetDefaultValue(const char *theDefaultValue)
void SetParameterRange(const char *theRange)
void SetParameterCandidates(const char *theString)

References G4UIcommand::AvailableForStates(), G4UIcommand::CategoryOf(), eyePosCmd, G4State_Idle, msDirectory, phiCmd, regionCmd, regSenseCmd, scanCmd, G4UIcmdWith3VectorAndUnit::SetDefaultUnit(), G4UIcmdWithAString::SetDefaultValue(), G4UIparameter::SetDefaultValue(), G4UIcmdWithABool::SetDefaultValue(), G4UIcmdWith3VectorAndUnit::SetDefaultValue(), G4UIcommand::SetGuidance(), G4UIcommand::SetParameter(), G4UIparameter::SetParameterCandidates(), G4UIcmdWithABool::SetParameterName(), G4UIcmdWithAString::SetParameterName(), G4UIcmdWith3Vector::SetParameterName(), G4UIcmdWith3VectorAndUnit::SetParameterName(), G4UIparameter::SetParameterRange(), single2Cmd, singleCmd, theScanner, thetaCmd, and G4UIcommand::UnitsList().

◆ ~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}

References eyePosCmd, msDirectory, phiCmd, regionCmd, regSenseCmd, scanCmd, single2Cmd, singleCmd, and thetaCmd.

Member Function Documentation

◆ AddUIcommand()

void G4UImessenger::AddUIcommand ( G4UIcommand newCommand)
protectedinherited

Definition at line 149 of file G4UImessenger.cc.

150{
151 G4cerr << "Warning : Old style definition of G4UIcommand <"
152 << newCommand->GetCommandPath() << ">." << G4endl;
153}
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
const G4String & GetCommandPath() const
Definition: G4UIcommand.hh:136

References G4cerr, G4endl, and G4UIcommand::GetCommandPath().

◆ BtoS()

G4String G4UImessenger::BtoS ( G4bool  b)
protectedinherited

Definition at line 98 of file G4UImessenger.cc.

99{
100 G4String vl = "0";
101 if(b)
102 vl = "true";
103 return vl;
104}

◆ CommandsShouldBeInMaster()

G4bool G4UImessenger::CommandsShouldBeInMaster ( ) const
inlineinherited

Definition at line 77 of file G4UImessenger.hh.

78 {
80 }
G4bool commandsShouldBeInMaster

References G4UImessenger::commandsShouldBeInMaster.

Referenced by G4UIcommand::G4UIcommandCommonConstructorCode().

◆ CreateCommand()

template<typename T >
T * G4UImessenger::CreateCommand ( const G4String cname,
const G4String dsc 
)
protectedinherited

Definition at line 110 of file G4UImessenger.hh.

111{
112 G4String path;
113 if(cname[0] != '/')
114 {
115 path = baseDirName + cname;
116 if(path[0] != '/')
117 path = "/" + path;
118 }
119
120 T* command = new T(path.c_str(), this);
121 command->SetGuidance(dsc.c_str());
122
123 return command;
124}
G4String baseDirName

References G4UImessenger::baseDirName.

◆ CreateDirectory()

void G4UImessenger::CreateDirectory ( const G4String path,
const G4String dsc,
G4bool  commandsToBeBroadcasted = true 
)
protectedinherited

Definition at line 156 of file G4UImessenger.cc.

158{
160
161 G4String fullpath = path;
162 if(fullpath.back() != '/')
163 fullpath.append("/");
164
165 G4UIcommandTree* tree = ui->GetTree()->FindCommandTree(fullpath.c_str());
166 if(tree != nullptr)
167 {
168 baseDirName = tree->GetPathName();
169 }
170 else
171 {
172 baseDir = new G4UIdirectory(fullpath.c_str(), commandsToBeBroadcasted);
173 baseDirName = fullpath;
174 baseDir->SetGuidance(dsc.c_str());
175 }
176}
const G4String & GetPathName() const
G4UIcommandTree * FindCommandTree(const char *commandPath)
G4UIcommandTree * GetTree() const
Definition: G4UImanager.hh:186
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77
G4UIdirectory * baseDir

References G4UImessenger::baseDir, G4UImessenger::baseDirName, G4UIcommandTree::FindCommandTree(), G4UIcommandTree::GetPathName(), G4UImanager::GetTree(), G4UImanager::GetUIpointer(), and G4UIcommand::SetGuidance().

Referenced by G4MoleculeShootMessenger::G4MoleculeShootMessenger(), and G4UImessenger::G4UImessenger().

◆ DtoS()

G4String G4UImessenger::DtoS ( G4double  a)
protectedinherited

Definition at line 90 of file G4UImessenger.cc.

91{
92 std::ostringstream os;
93 os << a;
94 return G4String(os.str());
95}

Referenced by G4ScoreQuantityMessenger::FilterCommands(), and G4UIcontrolMessenger::SetNewValue().

◆ 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 +=
175 currentValue += " ";
176 currentValue +=
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 =
191 }
192 else if(command == regSenseCmd)
193 {
194 currentValue =
196 }
197 else if(command == regionCmd)
198 {
199 currentValue = theScanner->GetRegionName();
200 }
201 return currentValue;
202}
static constexpr double deg
Definition: G4SIunits.hh:132
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:445

References G4UIcommand::ConvertToString(), deg, eyePosCmd, G4MaterialScanner::GetEyePosition(), G4MaterialScanner::GetNPhi(), G4MaterialScanner::GetNTheta(), G4MaterialScanner::GetPhiMin(), G4MaterialScanner::GetPhiSpan(), G4MaterialScanner::GetRegionName(), G4MaterialScanner::GetRegionSensitive(), G4MaterialScanner::GetThetaMin(), G4MaterialScanner::GetThetaSpan(), phiCmd, regionCmd, regSenseCmd, theScanner, and thetaCmd.

◆ ItoS()

G4String G4UImessenger::ItoS ( G4int  i)
protectedinherited

Definition at line 82 of file G4UImessenger.cc.

83{
84 std::ostringstream os;
85 os << i;
86 return G4String(os.str());
87}

Referenced by G4GenericMessenger::DeclareMethod(), and G4ParticleGunMessenger::GetCurrentValue().

◆ operator!=()

G4bool G4UImessenger::operator!= ( const G4UImessenger messenger) const
inherited

Definition at line 76 of file G4UImessenger.cc.

77{
78 return this != &messenger;
79}

◆ operator==()

G4bool G4UImessenger::operator== ( const G4UImessenger messenger) const
inherited

Definition at line 70 of file G4UImessenger.cc.

71{
72 return this == &messenger;
73}

◆ 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 {
240 }
241 else if(command == regSenseCmd)
242 {
244 }
245 else if(command == regionCmd)
246 {
247 if(theScanner->SetRegionName(newValue))
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 {
273 theta = 90. * deg - v.theta();
274 phi = v.phi();
275 }
277 theScanner->SetThetaMin(theta);
280 theScanner->SetPhiMin(phi);
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:363
G4int StoI(G4String s)
G4double StoD(G4String s)

References deg, eyePosCmd, G4UIcmdWith3Vector::GetNew3VectorValue(), G4UIcmdWith3VectorAndUnit::GetNew3VectorValue(), G4UIcmdWithABool::GetNewBoolValue(), G4MaterialScanner::GetNPhi(), G4MaterialScanner::GetNTheta(), G4MaterialScanner::GetPhiMin(), G4MaterialScanner::GetPhiSpan(), G4MaterialScanner::GetThetaMin(), G4MaterialScanner::GetThetaSpan(), CLHEP::Hep3Vector::phi(), phiCmd, regionCmd, regSenseCmd, G4MaterialScanner::Scan(), scanCmd, G4MaterialScanner::SetEyePosition(), G4MaterialScanner::SetNPhi(), G4MaterialScanner::SetNTheta(), G4MaterialScanner::SetPhiMin(), G4MaterialScanner::SetPhiSpan(), G4MaterialScanner::SetRegionName(), G4MaterialScanner::SetRegionSensitive(), G4MaterialScanner::SetThetaMin(), G4MaterialScanner::SetThetaSpan(), single2Cmd, singleCmd, G4UImessenger::StoD(), G4UImessenger::StoI(), theScanner, CLHEP::Hep3Vector::theta(), thetaCmd, and G4UIcommand::ValueOf().

◆ StoB()

G4bool G4UImessenger::StoB ( G4String  s)
protectedinherited

Definition at line 137 of file G4UImessenger.cc.

138{
140 G4bool vl = false;
141 if(v == "Y" || v == "YES" || v == "1" || v == "T" || v == "TRUE")
142 {
143 vl = true;
144 }
145 return vl;
146}
bool G4bool
Definition: G4Types.hh:86
G4String to_upper_copy(G4String str)
Return uppercase copy of string.

References G4StrUtil::to_upper_copy().

Referenced by G4LocalThreadCoutMessenger::SetNewValue(), G4CascadeParamMessenger::SetNewValue(), G4ScoreQuantityMessenger::SetNewValue(), and G4ScoringMessenger::SetNewValue().

◆ StoD()

G4double G4UImessenger::StoD ( G4String  s)
protectedinherited

◆ StoI()

G4int G4UImessenger::StoI ( G4String  s)
protectedinherited

◆ StoL()

G4long G4UImessenger::StoL ( G4String  s)
protectedinherited

Definition at line 117 of file G4UImessenger.cc.

118{
119 G4long vl;
120 const char* t = str;
121 std::istringstream is(t);
122 is >> vl;
123 return vl;
124}
long G4long
Definition: G4Types.hh:87

Referenced by G4RunMessenger::SetNewValue().

Field Documentation

◆ baseDir

G4UIdirectory* G4UImessenger::baseDir = nullptr
protectedinherited

◆ baseDirName

G4String G4UImessenger::baseDirName = ""
protectedinherited

◆ commandsShouldBeInMaster

G4bool G4UImessenger::commandsShouldBeInMaster = false
protectedinherited

◆ eyePosCmd

G4UIcmdWith3VectorAndUnit* G4MatScanMessenger::eyePosCmd = nullptr
private

◆ msDirectory

G4UIdirectory* G4MatScanMessenger::msDirectory = nullptr
private

Definition at line 62 of file G4MatScanMessenger.hh.

Referenced by G4MatScanMessenger(), and ~G4MatScanMessenger().

◆ phiCmd

G4UIcommand* G4MatScanMessenger::phiCmd = nullptr
private

◆ regionCmd

G4UIcmdWithAString* G4MatScanMessenger::regionCmd = nullptr
private

◆ regSenseCmd

G4UIcmdWithABool* G4MatScanMessenger::regSenseCmd = nullptr
private

◆ scanCmd

G4UIcmdWithoutParameter* G4MatScanMessenger::scanCmd = nullptr
private

Definition at line 63 of file G4MatScanMessenger.hh.

Referenced by G4MatScanMessenger(), SetNewValue(), and ~G4MatScanMessenger().

◆ single2Cmd

G4UIcmdWith3Vector* G4MatScanMessenger::single2Cmd = nullptr
private

Definition at line 67 of file G4MatScanMessenger.hh.

Referenced by G4MatScanMessenger(), SetNewValue(), and ~G4MatScanMessenger().

◆ singleCmd

G4UIcommand* G4MatScanMessenger::singleCmd = nullptr
private

Definition at line 66 of file G4MatScanMessenger.hh.

Referenced by G4MatScanMessenger(), SetNewValue(), and ~G4MatScanMessenger().

◆ theScanner

G4MaterialScanner* G4MatScanMessenger::theScanner = nullptr
private

Definition at line 60 of file G4MatScanMessenger.hh.

Referenced by G4MatScanMessenger(), GetCurrentValue(), and SetNewValue().

◆ thetaCmd

G4UIcommand* G4MatScanMessenger::thetaCmd = nullptr
private

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