Geant4 Cross Reference |
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 // G4UIcmdWith3VectorAndUnit 27 // 28 // Author: M.Asai, 1998 29 // -------------------------------------------------------------------- 30 31 #include "G4UIcmdWith3VectorAndUnit.hh" 32 33 #include "G4Tokenizer.hh" 34 #include "G4UIcommandStatus.hh" 35 #include "G4UnitsTable.hh" 36 37 #include <sstream> 38 39 // -------------------------------------------------------------------- 40 G4UIcmdWith3VectorAndUnit::G4UIcmdWith3VectorAndUnit(const char* theCommandPath, 41 G4UImessenger* theMessenger) 42 : G4UIcommand(theCommandPath, theMessenger) 43 { 44 auto* dblParamX = new G4UIparameter('d'); 45 SetParameter(dblParamX); 46 auto* dblParamY = new G4UIparameter('d'); 47 SetParameter(dblParamY); 48 auto* dblParamZ = new G4UIparameter('d'); 49 SetParameter(dblParamZ); 50 auto* untParam = new G4UIparameter('s'); 51 untParam->SetParameterName("Unit"); 52 SetParameter(untParam); 53 SetCommandType(With3VectorAndUnitCmd); 54 } 55 56 // -------------------------------------------------------------------- 57 G4int G4UIcmdWith3VectorAndUnit::DoIt(const G4String& parameterList) 58 { 59 std::vector<G4String> token_vector; 60 G4Tokenizer tkn(parameterList); 61 G4String str; 62 while (!(str = tkn()).empty()) { 63 token_vector.push_back(str); 64 } 65 66 // convert a value in default unit 67 G4String converted_parameter; 68 G4String default_unit = GetParameter(3)->GetDefaultValue(); 69 if (!default_unit.empty() && token_vector.size() >= 4) { 70 if (CategoryOf(token_vector[3]) != CategoryOf(default_unit)) { 71 return fParameterOutOfCandidates + 3; 72 } 73 G4double value_given = ValueOf(token_vector[3]); 74 G4double value_default = ValueOf(default_unit); 75 G4double x = ConvertToDouble(token_vector[0]) * value_given / value_default; 76 G4double y = ConvertToDouble(token_vector[1]) * value_given / value_default; 77 G4double z = ConvertToDouble(token_vector[2]) * value_given / value_default; 78 79 // reconstruct parameter list 80 converted_parameter += ConvertToString(x); 81 converted_parameter += " "; 82 converted_parameter += ConvertToString(y); 83 converted_parameter += " "; 84 converted_parameter += ConvertToString(z); 85 converted_parameter += " "; 86 converted_parameter += default_unit; 87 for (std::size_t i = 4; i < token_vector.size(); ++i) { 88 converted_parameter += " "; 89 converted_parameter += token_vector[i]; 90 } 91 } 92 else { 93 converted_parameter = parameterList; 94 } 95 96 return G4UIcommand::DoIt(converted_parameter); 97 } 98 99 // -------------------------------------------------------------------- 100 G4ThreeVector G4UIcmdWith3VectorAndUnit::GetNew3VectorValue(const char* paramString) 101 { 102 return ConvertToDimensioned3Vector(paramString); 103 } 104 105 // -------------------------------------------------------------------- 106 G4ThreeVector G4UIcmdWith3VectorAndUnit::GetNew3VectorRawValue(const char* paramString) 107 { 108 G4double vx; 109 G4double vy; 110 G4double vz; 111 char unts[30]; 112 std::istringstream is(paramString); 113 is >> vx >> vy >> vz >> unts; 114 return G4ThreeVector(vx, vy, vz); 115 } 116 117 // -------------------------------------------------------------------- 118 G4double G4UIcmdWith3VectorAndUnit::GetNewUnitValue(const char* paramString) 119 { 120 G4double vx; 121 G4double vy; 122 G4double vz; 123 char unts[30]; 124 std::istringstream is(paramString); 125 is >> vx >> vy >> vz >> unts; 126 G4String unt = unts; 127 return ValueOf(unt); 128 } 129 130 // -------------------------------------------------------------------- 131 G4String G4UIcmdWith3VectorAndUnit::ConvertToStringWithBestUnit(const G4ThreeVector& vec) 132 { 133 G4UIparameter* unitParam = GetParameter(3); 134 G4String canList = unitParam->GetParameterCandidates(); 135 G4Tokenizer candidateTokenizer(canList); 136 G4String aToken = candidateTokenizer(); 137 138 std::ostringstream os; 139 os << G4BestUnit(vec, CategoryOf(aToken)); 140 G4String st = os.str(); 141 142 return st; 143 } 144 145 // -------------------------------------------------------------------- 146 G4String G4UIcmdWith3VectorAndUnit::ConvertToStringWithDefaultUnit(const G4ThreeVector& vec) 147 { 148 G4UIparameter* unitParam = GetParameter(3); 149 G4String st; 150 if (unitParam->IsOmittable()) { 151 st = ConvertToString(vec, unitParam->GetDefaultValue()); 152 } 153 else { 154 st = ConvertToStringWithBestUnit(vec); 155 } 156 return st; 157 } 158 159 // -------------------------------------------------------------------- 160 void G4UIcmdWith3VectorAndUnit::SetParameterName(const char* theNameX, const char* theNameY, 161 const char* theNameZ, G4bool omittable, 162 G4bool currentAsDefault) 163 { 164 G4UIparameter* theParamX = GetParameter(0); 165 theParamX->SetParameterName(theNameX); 166 theParamX->SetOmittable(omittable); 167 theParamX->SetCurrentAsDefault(currentAsDefault); 168 G4UIparameter* theParamY = GetParameter(1); 169 theParamY->SetParameterName(theNameY); 170 theParamY->SetOmittable(omittable); 171 theParamY->SetCurrentAsDefault(currentAsDefault); 172 G4UIparameter* theParamZ = GetParameter(2); 173 theParamZ->SetParameterName(theNameZ); 174 theParamZ->SetOmittable(omittable); 175 theParamZ->SetCurrentAsDefault(currentAsDefault); 176 } 177 178 // -------------------------------------------------------------------- 179 void G4UIcmdWith3VectorAndUnit::SetDefaultValue(const G4ThreeVector& vec) 180 { 181 G4UIparameter* theParamX = GetParameter(0); 182 theParamX->SetDefaultValue(vec.x()); 183 G4UIparameter* theParamY = GetParameter(1); 184 theParamY->SetDefaultValue(vec.y()); 185 G4UIparameter* theParamZ = GetParameter(2); 186 theParamZ->SetDefaultValue(vec.z()); 187 } 188 189 // -------------------------------------------------------------------- 190 void G4UIcmdWith3VectorAndUnit::SetUnitCategory(const char* unitCategory) 191 { 192 SetUnitCandidates(UnitsList(unitCategory)); 193 } 194 195 // -------------------------------------------------------------------- 196 void G4UIcmdWith3VectorAndUnit::SetUnitCandidates(const char* candidateList) 197 { 198 G4UIparameter* untParam = GetParameter(3); 199 G4String canList = candidateList; 200 untParam->SetParameterCandidates(canList); 201 } 202 203 // -------------------------------------------------------------------- 204 void G4UIcmdWith3VectorAndUnit::SetDefaultUnit(const char* defUnit) 205 { 206 G4UIparameter* untParam = GetParameter(3); 207 untParam->SetOmittable(true); 208 untParam->SetDefaultValue(defUnit); 209 SetUnitCategory(CategoryOf(defUnit)); 210 } 211