Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/run/src/G4MatScanMessenger.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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"
 38 #include "G4UIcmdWith3VectorAndUnit.hh"
 39 #include "G4UIcmdWithABool.hh"
 40 #include "G4UIcmdWithAString.hh"
 41 #include "G4UIcmdWithoutParameter.hh"
 42 #include "G4UIcmdWithAnInteger.hh"
 43 #include "G4UIcommand.hh"
 44 #include "G4UIdirectory.hh"
 45 #include "G4UIparameter.hh"
 46 
 47 // --------------------------------------------------------------------
 48 G4MatScanMessenger::G4MatScanMessenger(G4MaterialScanner* p1)
 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 }
150 
151 // --------------------------------------------------------------------
152 G4MatScanMessenger::~G4MatScanMessenger()
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 }
164 
165 // --------------------------------------------------------------------
166 G4String G4MatScanMessenger::GetCurrentValue(G4UIcommand* command)
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 }
194 
195 // --------------------------------------------------------------------
196 void G4MatScanMessenger::SetNewValue(G4UIcommand* command, G4String newValue)
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 }
277