Geant4  10.02.p03
G4MatScanMessenger Class Reference

#include <G4MatScanMessenger.hh>

Inheritance diagram for G4MatScanMessenger:
Collaboration 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 ()
 
G4bool operator== (const G4UImessenger &messenger) const
 
G4bool CommandsShouldBeInMaster () const
 

Private Attributes

G4MaterialScannertheScanner
 
G4UIdirectorymsDirectory
 
G4UIcmdWithoutParameterscanCmd
 
G4UIcommandthetaCmd
 
G4UIcommandphiCmd
 
G4UIcommandsingleCmd
 
G4UIcmdWith3Vectorsingle2Cmd
 
G4UIcmdWithABoolregSenseCmd
 
G4UIcmdWithAStringregionCmd
 
G4UIcmdWith3VectorAndUniteyePosCmd
 

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)
 
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
 
G4String baseDirName
 
G4bool commandsShouldBeInMaster
 

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("/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");
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(" phi --- +X axis : 0 deg. / +Y axis : 90 deg. / -X axis : 180 deg. / -Y axis : 270 deg.");
84  par = new G4UIparameter("nbin",'i',false);
85  par->SetParameterRange("nbin>0");
86  phiCmd->SetParameter(par);
87  par = new G4UIparameter("phiMin",'d',false);
88  phiCmd->SetParameter(par);
89  par = new G4UIparameter("phiSpan",'d',true);
90  par->SetParameterRange("phiSpan>=0.");
91  par->SetDefaultValue(0.);
92  phiCmd->SetParameter(par);
93  par = new G4UIparameter("unit",'c',true);
94  par->SetDefaultValue("deg");
96  phiCmd->SetParameter(par);
97 
98  singleCmd = new G4UIcommand("/control/matScan/singleMeasure",this);
99  singleCmd->SetGuidance("Measure thickness for one particular direction.");
100  singleCmd->SetGuidance("Notation of angles :");
101  singleCmd->SetGuidance(" theta --- +Z axis : +90 deg. / X-Y plane : 0 deg. / -Z axis : -90 deg.");
102  singleCmd->SetGuidance(" phi --- +X axis : 0 deg. / +Y axis : 90 deg. / -X axis : 180 deg. / -Y axis : 270 deg.");
104  par = new G4UIparameter("theta",'d',false);
105  singleCmd->SetParameter(par);
106  par = new G4UIparameter("phi",'d',false);
107  singleCmd->SetParameter(par);
108  par = new G4UIparameter("unit",'c',true);
109  par->SetDefaultValue("deg");
111  singleCmd->SetParameter(par);
112 
113  single2Cmd = new G4UIcmdWith3Vector("/control/matScan/singleTo",this);
114  single2Cmd->SetGuidance("Measure thicknesss for one direction defined by a unit vector.");
115  single2Cmd->SetParameterName("X","Y","Z",false);
116 
117  eyePosCmd = new G4UIcmdWith3VectorAndUnit("/control/matScan/eyePosition",this);
118  eyePosCmd->SetGuidance("Define the eye position.");
119  eyePosCmd->SetParameterName("X","Y","Z",true);
122 
123  regSenseCmd = new G4UIcmdWithABool("/control/matScan/regionSensitive",this);
124  regSenseCmd->SetGuidance("Set region sensitivity.");
125  regSenseCmd->SetGuidance("This command is automatically set to TRUE");
126  regSenseCmd->SetGuidance(" if /control/matScan/region command is issued.");
127  regSenseCmd->SetParameterName("senseFlag",true);
129 
130  regionCmd = new G4UIcmdWithAString("/control/matScan/region",this);
131  regionCmd->SetGuidance("Define region name to be scanned.");
132  regionCmd->SetGuidance("/control/matScan/regionSensitive command is automatically");
133  regionCmd->SetGuidance("set to TRUE with this command.");
134  regionCmd->SetParameterName("region",true);
135  regionCmd->SetDefaultValue("DefaultRegionForTheWorld");
136 }
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
G4UIcmdWithoutParameter * scanCmd
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:152
G4UIcmdWith3VectorAndUnit * eyePosCmd
CLHEP::Hep3Vector G4ThreeVector
void SetParameterRange(const char *theRange)
void SetParameterCandidates(const char *theString)
void SetDefaultUnit(const char *defUnit)
void SetParameterName(const char *theNameX, const char *theNameY, const char *theNameZ, G4bool omittable, G4bool currentAsDefault=false)
G4UIdirectory * msDirectory
void SetDefaultValue(const char *theDefaultValue)
G4UIcmdWithAString * regionCmd
G4UIcmdWithABool * regSenseCmd
void SetParameterName(const char *theNameX, const char *theNameY, const char *theNameZ, G4bool omittable, G4bool currentAsDefault=false)
G4UIcmdWith3Vector * single2Cmd
void SetDefaultValue(G4bool defVal)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
static G4String UnitsList(const char *unitCategory)
Definition: G4UIcommand.cc:320
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:161
void AvailableForStates(G4ApplicationState s1)
Definition: G4UIcommand.cc:239
G4MaterialScanner * theScanner
void SetDefaultValue(G4ThreeVector defVal)
void SetDefaultValue(const char *defVal)
static G4String CategoryOf(const char *unitName)
Definition: G4UIcommand.cc:315
Here is the call graph for this function:

◆ ~G4MatScanMessenger()

G4MatScanMessenger::~G4MatScanMessenger ( )
virtual

Definition at line 138 of file G4MatScanMessenger.cc.

139 {
140  delete scanCmd;
141  delete thetaCmd;
142  delete phiCmd;
143  delete singleCmd;
144  delete single2Cmd;
145  delete eyePosCmd;
146  delete regSenseCmd;
147  delete regionCmd;
148  delete msDirectory;
149 }
G4UIcmdWithoutParameter * scanCmd
G4UIcmdWith3VectorAndUnit * eyePosCmd
G4UIdirectory * msDirectory
G4UIcmdWithAString * regionCmd
G4UIcmdWithABool * regSenseCmd
G4UIcmdWith3Vector * single2Cmd

Member Function Documentation

◆ GetCurrentValue()

G4String G4MatScanMessenger::GetCurrentValue ( G4UIcommand command)
virtual

Reimplemented from G4UImessenger.

Definition at line 151 of file G4MatScanMessenger.cc.

152 {
153  G4String currentValue;
154  if(command==thetaCmd)
155  {
156  currentValue = thetaCmd->ConvertToString(theScanner->GetNTheta());
157  currentValue += " ";
158  currentValue += thetaCmd->ConvertToString((theScanner->GetThetaMin())/deg);
159  currentValue += " ";
160  currentValue += thetaCmd->ConvertToString((theScanner->GetThetaSpan())/deg);
161  }
162  else if(command==phiCmd)
163  {
164  currentValue = phiCmd->ConvertToString(theScanner->GetNPhi());
165  currentValue += " ";
166  currentValue += phiCmd->ConvertToString((theScanner->GetPhiMin())/deg);
167  currentValue += " ";
168  currentValue += phiCmd->ConvertToString((theScanner->GetPhiSpan())/deg);
169  }
170  else if(command==eyePosCmd)
171  { currentValue = eyePosCmd->ConvertToString(theScanner->GetEyePosition(),"m"); }
172  else if(command==regSenseCmd)
174  else if(command==regionCmd)
175  { currentValue = theScanner->GetRegionName(); }
176  return currentValue;
177 }
G4UIcmdWith3VectorAndUnit * eyePosCmd
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:371
G4bool GetRegionSensitive() const
G4UIcmdWithAString * regionCmd
G4UIcmdWithABool * regSenseCmd
G4double GetThetaSpan() const
static const double deg
Definition: G4SIunits.hh:151
G4int GetNTheta() const
G4MaterialScanner * theScanner
G4ThreeVector GetEyePosition() const
G4double GetPhiSpan() const
G4double GetPhiMin() const
G4double GetThetaMin() const
G4String GetRegionName() const
Here is the call graph for this function:

◆ SetNewValue()

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

Reimplemented from G4UImessenger.

Definition at line 179 of file G4MatScanMessenger.cc.

180 {
181  if(command==scanCmd)
182  { theScanner->Scan(); }
183  else if(command==thetaCmd)
184  {
185  G4Tokenizer next( newValue );
186  G4int nbin = StoI(next());
187  G4double thetaMin = StoD(next());
188  G4double thetaSpan = StoD(next());
189  G4String unit = next();
190  thetaMin *= thetaCmd->ValueOf(unit);
191  thetaSpan *= thetaCmd->ValueOf(unit);
192  theScanner->SetNTheta(nbin);
193  theScanner->SetThetaMin(thetaMin);
194  theScanner->SetThetaSpan(thetaSpan);
195  }
196  else if(command==phiCmd)
197  {
198  G4Tokenizer next( newValue );
199  G4int nbin = StoI(next());
200  G4double phiMin = StoD(next());
201  G4double phiSpan = StoD(next());
202  G4String unit = next();
203  phiMin *= phiCmd->ValueOf(unit);
204  phiSpan *= phiCmd->ValueOf(unit);
205  theScanner->SetNPhi(nbin);
206  theScanner->SetPhiMin(phiMin);
207  theScanner->SetPhiSpan(phiSpan);
208  }
209  else if(command==eyePosCmd)
211  else if(command==regSenseCmd)
213  else if(command==regionCmd)
214  { if(theScanner->SetRegionName(newValue)) theScanner->SetRegionSensitive(true); }
215  else if(command==singleCmd || command==single2Cmd)
216  {
217  G4int ntheta = theScanner->GetNTheta();
218  G4double thetaMin = theScanner->GetThetaMin();
219  G4double thetaSpan = theScanner->GetThetaSpan();
220  G4int nphi = theScanner->GetNPhi();
221  G4double phiMin = theScanner->GetPhiMin();
222  G4double phiSpan = theScanner->GetPhiSpan();
223 
224  G4double theta = 0.;
225  G4double phi = 0.;
226  if(command==singleCmd)
227  {
228  G4Tokenizer next( newValue );
229  theta = StoD(next());
230  phi = StoD(next());
231  G4String unit = next();
232  theta *= singleCmd->ValueOf(unit);
233  phi *= singleCmd->ValueOf(unit);
234  }
235  else if(command==single2Cmd)
236  {
238  theta = 90.*deg - v.theta();
239  phi = v.phi();
240  }
241  theScanner->SetNTheta(1);
242  theScanner->SetThetaMin(theta);
244  theScanner->SetNPhi(1);
245  theScanner->SetPhiMin(phi);
246  theScanner->SetPhiSpan(0.);
247  theScanner->Scan();
248 
249  theScanner->SetNTheta(ntheta);
250  theScanner->SetThetaMin(thetaMin);
251  theScanner->SetThetaSpan(thetaSpan);
252  theScanner->SetNPhi(nphi);
253  theScanner->SetPhiMin(phiMin);
254  theScanner->SetPhiSpan(phiSpan);
255  }
256 
257 }
void SetRegionSensitive(G4bool val=true)
G4UIcmdWithoutParameter * scanCmd
G4UIcmdWith3VectorAndUnit * eyePosCmd
double phi() const
void SetPhiMin(G4double val)
void SetPhiSpan(G4double val)
G4bool SetRegionName(const G4String &val)
static G4ThreeVector GetNew3VectorValue(const char *paramString)
G4UIcmdWithAString * regionCmd
G4UIcmdWithABool * regSenseCmd
void SetThetaMin(G4double val)
int G4int
Definition: G4Types.hh:78
static G4bool GetNewBoolValue(const char *paramString)
G4UIcmdWith3Vector * single2Cmd
void SetEyePosition(const G4ThreeVector &val)
G4double GetThetaSpan() const
void SetNPhi(G4int val)
static G4ThreeVector GetNew3VectorValue(const char *paramString)
static const double deg
Definition: G4SIunits.hh:151
double theta() const
G4int GetNTheta() const
void SetThetaSpan(G4double val)
G4MaterialScanner * theScanner
G4double GetPhiSpan() const
G4int StoI(G4String s)
G4double StoD(G4String s)
static G4double ValueOf(const char *unitName)
Definition: G4UIcommand.cc:308
G4double GetPhiMin() const
double G4double
Definition: G4Types.hh:76
G4double GetThetaMin() const
void SetNTheta(G4int val)
Here is the call graph for this function:

Member Data Documentation

◆ eyePosCmd

G4UIcmdWith3VectorAndUnit* G4MatScanMessenger::eyePosCmd
private

Definition at line 68 of file G4MatScanMessenger.hh.

◆ msDirectory

G4UIdirectory* G4MatScanMessenger::msDirectory
private

Definition at line 60 of file G4MatScanMessenger.hh.

◆ phiCmd

G4UIcommand* G4MatScanMessenger::phiCmd
private

Definition at line 63 of file G4MatScanMessenger.hh.

◆ regionCmd

G4UIcmdWithAString* G4MatScanMessenger::regionCmd
private

Definition at line 67 of file G4MatScanMessenger.hh.

◆ regSenseCmd

G4UIcmdWithABool* G4MatScanMessenger::regSenseCmd
private

Definition at line 66 of file G4MatScanMessenger.hh.

◆ scanCmd

G4UIcmdWithoutParameter* G4MatScanMessenger::scanCmd
private

Definition at line 61 of file G4MatScanMessenger.hh.

◆ single2Cmd

G4UIcmdWith3Vector* G4MatScanMessenger::single2Cmd
private

Definition at line 65 of file G4MatScanMessenger.hh.

◆ singleCmd

G4UIcommand* G4MatScanMessenger::singleCmd
private

Definition at line 64 of file G4MatScanMessenger.hh.

◆ theScanner

G4MaterialScanner* G4MatScanMessenger::theScanner
private

Definition at line 58 of file G4MatScanMessenger.hh.

◆ thetaCmd

G4UIcommand* G4MatScanMessenger::thetaCmd
private

Definition at line 62 of file G4MatScanMessenger.hh.


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