Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4UIcmdWith3VectorAndUnit.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 //
27 // $Id: G4UIcmdWith3VectorAndUnit.cc 77308 2013-11-22 11:22:28Z gcosmo $
28 //
29 //
30 
32 #include "G4Tokenizer.hh"
33 #include "G4UnitsTable.hh"
34 #include "G4UIcommandStatus.hh"
35 #include <sstream>
36 
38 (const char * theCommandPath,G4UImessenger * theMessenger)
39 :G4UIcommand(theCommandPath,theMessenger)
40 {
41  G4UIparameter * dblParamX = new G4UIparameter('d');
42  SetParameter(dblParamX);
43  G4UIparameter * dblParamY = new G4UIparameter('d');
44  SetParameter(dblParamY);
45  G4UIparameter * dblParamZ = new G4UIparameter('d');
46  SetParameter(dblParamZ);
47  G4UIparameter * untParam = new G4UIparameter('s');
48  SetParameter(untParam);
49  untParam->SetParameterName("Unit");
50 }
51 
53 {
54  std::vector<G4String> token_vector;
55  G4Tokenizer tkn(parameterList);
56  G4String str;
57  while( (str = tkn()) != "" ) {
58  token_vector.push_back(str);
59  }
60 
61  // convert a value in default unit
62  G4String converted_parameter;
63  G4String default_unit = GetParameter(3)-> GetDefaultValue();
64  if (default_unit != "" && token_vector.size() >= 4) {
65  if(CategoryOf(token_vector[3])!=CategoryOf(default_unit))
66  { return fParameterOutOfCandidates+3; }
67  G4double value_given = ValueOf(token_vector[3]);
68  G4double value_default = ValueOf(default_unit);
69  G4double x = ConvertToDouble(token_vector[0]) * value_given / value_default;
70  G4double y = ConvertToDouble(token_vector[1]) * value_given / value_default;
71  G4double z = ConvertToDouble(token_vector[2]) * value_given / value_default;
72 
73  // reconstruct parameter list
74  converted_parameter += ConvertToString(x);
75  converted_parameter += " ";
76  converted_parameter += ConvertToString(y);
77  converted_parameter += " ";
78  converted_parameter += ConvertToString(z);
79  converted_parameter += " ";
80  converted_parameter += default_unit;
81  for ( size_t i=4 ; i< token_vector.size(); i++) {
82  converted_parameter += " ";
83  converted_parameter += token_vector[i];
84  }
85  } else {
86  converted_parameter = parameterList;
87  }
88 
89  return G4UIcommand::DoIt(converted_parameter);
90 }
91 
93 {
94  return ConvertToDimensioned3Vector(paramString);
95 }
96 
98 {
99  G4double vx;
100  G4double vy;
101  G4double vz;
102  char unts[30];
103  std::istringstream is(paramString);
104  is >> vx >> vy >> vz >> unts;
105  return G4ThreeVector(vx,vy,vz);
106 }
107 
109 {
110  G4double vx;
111  G4double vy;
112  G4double vz;
113  char unts[30];
114  std::istringstream is(paramString);
115  is >> vx >> vy >> vz >> unts;
116  G4String unt = unts;
117  return ValueOf(unt);
118 }
119 
121 {
122  G4UIparameter* unitParam = GetParameter(3);
123  G4String canList = unitParam->GetParameterCandidates();
124  G4Tokenizer candidateTokenizer(canList);
125  G4String aToken = candidateTokenizer();
126 
127  std::ostringstream os;
128  os << G4BestUnit(vec,CategoryOf(aToken));
129  G4String st = os.str();
130 
131  return st;
132 }
133 
135 {
136  G4UIparameter* unitParam = GetParameter(3);
137  G4String st;
138  if(unitParam->IsOmittable())
139  { st = ConvertToString(vec,unitParam->GetDefaultValue()); }
140  else
141  { st = ConvertToStringWithBestUnit(vec); }
142  return st;
143 }
144 
146 (const char * theNameX,const char * theNameY,const char * theNameZ,
147 G4bool omittable,G4bool currentAsDefault)
148 {
149  G4UIparameter * theParamX = GetParameter(0);
150  theParamX->SetParameterName(theNameX);
151  theParamX->SetOmittable(omittable);
152  theParamX->SetCurrentAsDefault(currentAsDefault);
153  G4UIparameter * theParamY = GetParameter(1);
154  theParamY->SetParameterName(theNameY);
155  theParamY->SetOmittable(omittable);
156  theParamY->SetCurrentAsDefault(currentAsDefault);
157  G4UIparameter * theParamZ = GetParameter(2);
158  theParamZ->SetParameterName(theNameZ);
159  theParamZ->SetOmittable(omittable);
160  theParamZ->SetCurrentAsDefault(currentAsDefault);
161 }
162 
164 {
165  G4UIparameter * theParamX = GetParameter(0);
166  theParamX->SetDefaultValue(vec.x());
167  G4UIparameter * theParamY = GetParameter(1);
168  theParamY->SetDefaultValue(vec.y());
169  G4UIparameter * theParamZ = GetParameter(2);
170  theParamZ->SetDefaultValue(vec.z());
171 }
172 
173 void G4UIcmdWith3VectorAndUnit::SetUnitCategory(const char * unitCategory)
174 {
175  SetUnitCandidates(UnitsList(unitCategory));
176 }
177 
178 void G4UIcmdWith3VectorAndUnit::SetUnitCandidates(const char * candidateList)
179 {
180  G4UIparameter * untParam = GetParameter(3);
181  G4String canList = candidateList;
182  untParam->SetParameterCandidates(canList);
183 }
184 
186 {
187  G4UIparameter * untParam = GetParameter(3);
188  untParam->SetOmittable(true);
189  untParam->SetDefaultValue(defUnit);
190  SetUnitCategory(CategoryOf(defUnit));
191 }
192 
G4String GetParameterCandidates() const
CLHEP::Hep3Vector G4ThreeVector
double x() const
G4double z
Definition: TRTMaterials.hh:39
void SetOmittable(G4bool om)
void SetParameterCandidates(const char *theString)
static G4double GetNewUnitValue(const char *paramString)
void SetDefaultUnit(const char *defUnit)
void SetUnitCategory(const char *unitCategory)
static G4ThreeVector GetNew3VectorValue(const char *paramString)
void SetDefaultValue(const char *theDefaultValue)
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:357
G4String ConvertToStringWithDefaultUnit(G4ThreeVector vec)
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
void SetParameterName(const char *theName)
int G4int
Definition: G4Types.hh:78
void SetParameterName(const char *theNameX, const char *theNameY, const char *theNameZ, G4bool omittable, G4bool currentAsDefault=false)
double z() const
void SetCurrentAsDefault(G4bool val)
G4bool IsOmittable() const
virtual G4int DoIt(G4String parameterList)
static G4String UnitsList(const char *unitCategory)
Definition: G4UIcommand.cc:306
bool G4bool
Definition: G4Types.hh:79
G4UIparameter * GetParameter(G4int i) const
Definition: G4UIcommand.hh:145
static G4double ConvertToDouble(const char *st)
Definition: G4UIcommand.cc:429
void SetDefaultValue(G4ThreeVector defVal)
G4UIcmdWith3VectorAndUnit(const char *theCommandPath, G4UImessenger *theMessenger)
G4String GetDefaultValue() const
static G4double ValueOf(const char *unitName)
Definition: G4UIcommand.cc:294
double y() const
double G4double
Definition: G4Types.hh:76
void SetUnitCandidates(const char *candidateList)
static G4String CategoryOf(const char *unitName)
Definition: G4UIcommand.cc:301
static G4ThreeVector GetNew3VectorRawValue(const char *paramString)
virtual G4int DoIt(G4String parameterList)
Definition: G4UIcommand.cc:108
G4String ConvertToStringWithBestUnit(G4ThreeVector vec)
static G4ThreeVector ConvertToDimensioned3Vector(const char *st)
Definition: G4UIcommand.cc:459