[814] | 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: G4ScoreQuantityMessenger.cc,v 1.8 2007/11/07 04:12:07 akimura Exp $ |
---|
[850] | 28 | // GEANT4 tag $Name: HEAD $ |
---|
[814] | 29 | // |
---|
| 30 | // --------------------------------------------------------------------- |
---|
| 31 | |
---|
| 32 | #include "G4ScoreQuantityMessenger.hh" |
---|
| 33 | #include "G4ScoringManager.hh" |
---|
| 34 | #include "G4VScoringMesh.hh" |
---|
| 35 | |
---|
| 36 | #include "G4PSCellCharge3D.hh" |
---|
| 37 | #include "G4PSCellFlux3D.hh" |
---|
| 38 | #include "G4PSPassageCellFlux3D.hh" |
---|
| 39 | #include "G4PSEnergyDeposit3D.hh" |
---|
| 40 | #include "G4PSDoseDeposit3D.hh" |
---|
| 41 | #include "G4PSNofStep3D.hh" |
---|
| 42 | #include "G4PSNofSecondary3D.hh" |
---|
| 43 | // |
---|
| 44 | #include "G4PSTrackLength3D.hh" |
---|
| 45 | #include "G4PSPassageCellCurrent3D.hh" |
---|
| 46 | #include "G4PSPassageTrackLength3D.hh" |
---|
| 47 | #include "G4PSFlatSurfaceCurrent3D.hh" |
---|
| 48 | #include "G4PSFlatSurfaceFlux3D.hh" |
---|
| 49 | #include "G4PSSphereSurfaceCurrent3D.hh" |
---|
| 50 | #include "G4PSSphereSurfaceFlux3D.hh" |
---|
| 51 | #include "G4PSCylinderSurfaceCurrent3D.hh" |
---|
| 52 | #include "G4PSCylinderSurfaceFlux3D.hh" |
---|
| 53 | #include "G4PSNofCollision3D.hh" |
---|
| 54 | #include "G4PSPopulation3D.hh" |
---|
| 55 | #include "G4PSTrackCounter3D.hh" |
---|
| 56 | #include "G4PSTermination3D.hh" |
---|
| 57 | |
---|
| 58 | #include "G4SDChargedFilter.hh" |
---|
| 59 | #include "G4SDNeutralFilter.hh" |
---|
| 60 | #include "G4SDKineticEnergyFilter.hh" |
---|
| 61 | #include "G4SDParticleFilter.hh" |
---|
| 62 | #include "G4SDParticleWithEnergyFilter.hh" |
---|
| 63 | |
---|
| 64 | #include "G4UIdirectory.hh" |
---|
| 65 | #include "G4UIcmdWithoutParameter.hh" |
---|
| 66 | #include "G4UIcmdWithAnInteger.hh" |
---|
| 67 | #include "G4UIcmdWithAString.hh" |
---|
| 68 | #include "G4UIcmdWithABool.hh" |
---|
| 69 | #include "G4UIcmdWithADoubleAndUnit.hh" |
---|
| 70 | #include "G4UIcmdWith3VectorAndUnit.hh" |
---|
| 71 | #include "G4UIcommand.hh" |
---|
| 72 | #include "G4Tokenizer.hh" |
---|
| 73 | #include "G4UnitsTable.hh" |
---|
| 74 | |
---|
| 75 | G4ScoreQuantityMessenger::G4ScoreQuantityMessenger(G4ScoringManager* SManager) |
---|
| 76 | :fSMan(SManager) |
---|
| 77 | { |
---|
| 78 | QuantityCommands(); |
---|
| 79 | FilterCommands(); |
---|
| 80 | } |
---|
| 81 | |
---|
| 82 | void G4ScoreQuantityMessenger::FilterCommands() |
---|
| 83 | { |
---|
| 84 | G4UIparameter* param; |
---|
| 85 | |
---|
| 86 | // |
---|
| 87 | // Filter commands |
---|
| 88 | filterDir = new G4UIdirectory("/score/filter/"); |
---|
| 89 | filterDir->SetGuidance(" Scoring filter commands."); |
---|
| 90 | // |
---|
| 91 | fchargedCmd = new G4UIcmdWithAString("/score/filter/charged",this); |
---|
| 92 | fchargedCmd->SetGuidance("Charged particle filter."); |
---|
| 93 | fchargedCmd->SetParameterName("fname",false); |
---|
| 94 | // |
---|
| 95 | fneutralCmd = new G4UIcmdWithAString("/score/filter/neutral",this); |
---|
| 96 | fneutralCmd->SetGuidance("Neutral particle filter."); |
---|
| 97 | fneutralCmd->SetParameterName("fname",false); |
---|
| 98 | // |
---|
| 99 | fkinECmd = new G4UIcommand("/score/filter/kineticEnergy",this); |
---|
| 100 | fkinECmd->SetGuidance("Kinetic energy filter."); |
---|
| 101 | fkinECmd->SetGuidance("[usage] /score/filter/kineticEnergy fname Elow Ehigh unit"); |
---|
| 102 | fkinECmd->SetGuidance(" fname :(String) Filter Name "); |
---|
| 103 | fkinECmd->SetGuidance(" Elow :(Double) Lower edge of kinetic energy"); |
---|
| 104 | fkinECmd->SetGuidance(" Ehigh :(Double) Higher edge of kinetic energy"); |
---|
| 105 | fkinECmd->SetGuidance(" unit :(String) unit of given kinetic energy"); |
---|
| 106 | param = new G4UIparameter("fname",'s',false); |
---|
| 107 | fkinECmd->SetParameter(param); |
---|
| 108 | param = new G4UIparameter("elow",'d',true); |
---|
| 109 | param->SetDefaultValue("0.0"); |
---|
| 110 | fkinECmd->SetParameter(param); |
---|
| 111 | param = new G4UIparameter("ehigh",'d',false); |
---|
| 112 | fkinECmd->SetParameter(param); |
---|
| 113 | G4String smax = DtoS(DBL_MAX); |
---|
| 114 | param->SetDefaultValue(smax); |
---|
| 115 | param = new G4UIparameter("unit",'s',false); |
---|
| 116 | param->SetDefaultValue("keV"); |
---|
| 117 | fkinECmd->SetParameter(param); |
---|
| 118 | // |
---|
| 119 | fparticleCmd = new G4UIcommand("/score/filter/particle",this); |
---|
| 120 | fparticleCmd->SetGuidance("Particle filter."); |
---|
| 121 | fparticleCmd->SetGuidance("[usage] /score/filter/particle fname p0 .. pn"); |
---|
| 122 | fparticleCmd->SetGuidance(" fname :(String) Filter Name "); |
---|
| 123 | fparticleCmd->SetGuidance(" p0 .. pn :(String) particle names"); |
---|
| 124 | param = new G4UIparameter("fname",'s',false); |
---|
| 125 | fparticleCmd->SetParameter(param); |
---|
| 126 | param = new G4UIparameter("particlelist",'s',false); |
---|
| 127 | param->SetDefaultValue(""); |
---|
| 128 | fparticleCmd->SetParameter(param); |
---|
| 129 | // |
---|
| 130 | // |
---|
| 131 | // |
---|
| 132 | fparticleKinECmd = new G4UIcommand("/score/filter/particleWithKineticEnergy",this); |
---|
| 133 | fparticleKinECmd->SetGuidance("Particle with kinetic energy filter."); |
---|
| 134 | fparticleKinECmd->SetGuidance("[usage] /score/filter/particleWithKineticEnergy fname Elow Ehigh unit p0 .. pn"); |
---|
| 135 | fparticleKinECmd->SetGuidance(" fname :(String) Filter Name "); |
---|
| 136 | fparticleKinECmd->SetGuidance(" Elow :(Double) Lower edge of kinetic energy"); |
---|
| 137 | fparticleKinECmd->SetGuidance(" Ehigh :(Double) Higher edge of kinetic energy"); |
---|
| 138 | fparticleKinECmd->SetGuidance(" unit :(String) unit of given kinetic energy"); |
---|
| 139 | fparticleKinECmd->SetGuidance(" p0 .. pn :(String) particle names"); |
---|
| 140 | param = new G4UIparameter("fname",'s',false); |
---|
| 141 | fparticleKinECmd->SetParameter(param); |
---|
| 142 | param = new G4UIparameter("elow",'d',false); |
---|
| 143 | param->SetDefaultValue("0.0"); |
---|
| 144 | fparticleKinECmd->SetParameter(param); |
---|
| 145 | param = new G4UIparameter("ehigh",'d',true); |
---|
| 146 | param->SetDefaultValue(smax); |
---|
| 147 | fparticleKinECmd->SetParameter(param); |
---|
| 148 | param = new G4UIparameter("unit",'s',true); |
---|
| 149 | param->SetDefaultValue("keV"); |
---|
| 150 | fparticleKinECmd->SetParameter(param); |
---|
| 151 | param = new G4UIparameter("particlelist",'s',false); |
---|
| 152 | param->SetDefaultValue(""); |
---|
| 153 | fparticleKinECmd->SetParameter(param); |
---|
| 154 | // |
---|
| 155 | // |
---|
| 156 | } |
---|
| 157 | |
---|
| 158 | G4ScoreQuantityMessenger::~G4ScoreQuantityMessenger() |
---|
| 159 | { |
---|
| 160 | delete quantityDir; |
---|
| 161 | delete qTouchCmd; |
---|
| 162 | // |
---|
| 163 | delete qCellChgCmd; |
---|
| 164 | delete qCellFluxCmd; |
---|
| 165 | delete qPassCellFluxCmd; |
---|
| 166 | delete qeDepCmd; |
---|
| 167 | delete qdoseDepCmd; |
---|
| 168 | delete qnOfStepCmd; |
---|
| 169 | delete qnOfSecondaryCmd; |
---|
| 170 | // |
---|
| 171 | delete qTrackLengthCmd; |
---|
| 172 | delete qPassCellCurrCmd; |
---|
| 173 | delete qPassTrackLengthCmd; |
---|
| 174 | delete qFlatSurfCurrCmd; |
---|
| 175 | delete qFlatSurfFluxCmd; |
---|
| 176 | // delete qSphereSurfCurrCmd; |
---|
| 177 | // delete qSphereSurfFluxCmd; |
---|
| 178 | // delete qCylSurfCurrCmd; |
---|
| 179 | // delete qCylSurfFluxCmd; |
---|
| 180 | delete qNofCollisionCmd; |
---|
| 181 | delete qPopulationCmd; |
---|
| 182 | delete qTrackCountCmd; |
---|
| 183 | delete qTerminationCmd; |
---|
| 184 | // |
---|
| 185 | delete filterDir; |
---|
| 186 | delete fchargedCmd; |
---|
| 187 | delete fneutralCmd; |
---|
| 188 | delete fkinECmd; |
---|
| 189 | delete fparticleCmd; |
---|
| 190 | delete fparticleKinECmd; |
---|
| 191 | } |
---|
| 192 | |
---|
| 193 | void G4ScoreQuantityMessenger::SetNewValue(G4UIcommand * command,G4String newVal) |
---|
| 194 | { |
---|
| 195 | // |
---|
| 196 | // Get Current Mesh |
---|
| 197 | // |
---|
| 198 | G4VScoringMesh* mesh = fSMan->GetCurrentMesh(); |
---|
| 199 | if(!mesh) |
---|
| 200 | { |
---|
| 201 | G4cerr << "ERROR : No mesh is currently open. Open/create a mesh first. Command ignored." << G4endl; |
---|
| 202 | return; |
---|
| 203 | } |
---|
| 204 | // Tokens |
---|
| 205 | G4TokenVec token; |
---|
| 206 | FillTokenVec(newVal,token); |
---|
| 207 | // |
---|
| 208 | // Commands for Current Mesh |
---|
| 209 | if(command==qTouchCmd) { |
---|
| 210 | mesh->SetCurrentPrimitiveScorer(newVal); |
---|
| 211 | } else if(command== qCellChgCmd) { |
---|
| 212 | if(!mesh->FindPrimitiveScorer(newVal)) { |
---|
| 213 | mesh->SetPrimitiveScorer(new G4PSCellCharge3D(newVal)); |
---|
| 214 | } else { |
---|
| 215 | G4cout << "WARNING[" << qTouchCmd->GetCommandPath() |
---|
| 216 | << "] : Quantity name, \"" << newVal << "\", is already existing." << G4endl; |
---|
| 217 | mesh->SetNullToCurrentPrimitiveScorer(); |
---|
| 218 | } |
---|
| 219 | } else if(command== qCellFluxCmd) { |
---|
| 220 | if(!mesh->FindPrimitiveScorer(newVal)) { |
---|
| 221 | mesh->SetPrimitiveScorer(new G4PSCellFlux3D(newVal)); |
---|
| 222 | } else { |
---|
| 223 | G4cout << "WARNING[" << qCellFluxCmd->GetCommandPath() |
---|
| 224 | << "] : Quantity name, \"" << newVal << "\", is already existing." << G4endl; |
---|
| 225 | mesh->SetNullToCurrentPrimitiveScorer(); |
---|
| 226 | } |
---|
| 227 | } else if(command== qPassCellFluxCmd) { |
---|
| 228 | if(!mesh->FindPrimitiveScorer(newVal)) { |
---|
| 229 | mesh->SetPrimitiveScorer(new G4PSPassageCellFlux3D(newVal)); |
---|
| 230 | } else { |
---|
| 231 | G4cout << "WARNING[" << qPassCellFluxCmd->GetCommandPath() |
---|
| 232 | << "] : Quantity name, \"" << newVal << "\", is already existing." << G4endl; |
---|
| 233 | mesh->SetNullToCurrentPrimitiveScorer(); |
---|
| 234 | } |
---|
| 235 | } else if(command==qeDepCmd) { |
---|
| 236 | if(!mesh->FindPrimitiveScorer(newVal)) { |
---|
| 237 | mesh->SetPrimitiveScorer(new G4PSEnergyDeposit3D(newVal)); |
---|
| 238 | } else { |
---|
| 239 | G4cout << "WARNING[" << qeDepCmd->GetCommandPath() |
---|
| 240 | << "] : Quantity name, \"" << newVal << "\", is already existing." << G4endl; |
---|
| 241 | mesh->SetNullToCurrentPrimitiveScorer(); |
---|
| 242 | } |
---|
| 243 | } else if(command== qdoseDepCmd) { |
---|
| 244 | if(!mesh->FindPrimitiveScorer(newVal)) { |
---|
| 245 | mesh->SetPrimitiveScorer(new G4PSDoseDeposit3D(newVal)); |
---|
| 246 | } else { |
---|
| 247 | G4cout << "WARNING[" << qdoseDepCmd->GetCommandPath() |
---|
| 248 | << "] : Quantity name, \"" << newVal << "\", is already existing." << G4endl; |
---|
| 249 | mesh->SetNullToCurrentPrimitiveScorer(); |
---|
| 250 | } |
---|
| 251 | } else if(command== qnOfStepCmd) { |
---|
| 252 | if(!mesh->FindPrimitiveScorer(newVal)) { |
---|
| 253 | mesh->SetPrimitiveScorer(new G4PSNofStep3D(newVal)); |
---|
| 254 | } else { |
---|
| 255 | G4cout << "WARNING[" << qnOfStepCmd->GetCommandPath() |
---|
| 256 | << "] : Quantity name, \"" << newVal << "\", is already existing." << G4endl; |
---|
| 257 | mesh->SetNullToCurrentPrimitiveScorer(); |
---|
| 258 | } |
---|
| 259 | } else if(command== qnOfSecondaryCmd) { |
---|
| 260 | if(!mesh->FindPrimitiveScorer(newVal)) { |
---|
| 261 | mesh->SetPrimitiveScorer(new G4PSNofSecondary3D(newVal)); |
---|
| 262 | } else { |
---|
| 263 | G4cout << "WARNING[" << qnOfSecondaryCmd->GetCommandPath() |
---|
| 264 | << "] : Quantity name, \"" << newVal << "\", is already existing." << G4endl; |
---|
| 265 | mesh->SetNullToCurrentPrimitiveScorer(); |
---|
| 266 | } |
---|
| 267 | } else if(command== qTrackLengthCmd) { |
---|
| 268 | if(!mesh->FindPrimitiveScorer(newVal)) { |
---|
| 269 | G4PSTrackLength3D* ps = new G4PSTrackLength3D(token[0]); |
---|
| 270 | ps->Weighted(StoB(token[1])); |
---|
| 271 | ps->MultiplyKineticEnergy(StoB(token[2])); |
---|
| 272 | ps->DivideByVelocity(StoB(token[3])); |
---|
| 273 | mesh->SetPrimitiveScorer(ps); |
---|
| 274 | } else { |
---|
| 275 | G4cout << "WARNING[" << qTrackLengthCmd->GetCommandPath() |
---|
| 276 | << "] : Quantity name, \"" << newVal << "\", is already existing." << G4endl; |
---|
| 277 | mesh->SetNullToCurrentPrimitiveScorer(); |
---|
| 278 | } |
---|
| 279 | } else if(command== qPassCellCurrCmd){ |
---|
| 280 | if(!mesh->FindPrimitiveScorer(newVal)) { |
---|
| 281 | G4PSPassageCellCurrent* ps = new G4PSPassageCellCurrent3D(token[0]); |
---|
| 282 | ps->Weighted(StoB(token[1])); |
---|
| 283 | mesh->SetPrimitiveScorer(ps); |
---|
| 284 | } else { |
---|
| 285 | G4cout << "WARNING[" << qPassCellCurrCmd->GetCommandPath() |
---|
| 286 | << "] : Quantity name, \"" << newVal << "\", is already existing." << G4endl; |
---|
| 287 | mesh->SetNullToCurrentPrimitiveScorer(); |
---|
| 288 | } |
---|
| 289 | } else if(command== qPassTrackLengthCmd){ |
---|
| 290 | if(!mesh->FindPrimitiveScorer(newVal)) { |
---|
| 291 | G4PSPassageTrackLength* ps = new G4PSPassageTrackLength3D(token[0]); |
---|
| 292 | ps->Weighted(StoB(token[1])); |
---|
| 293 | mesh->SetPrimitiveScorer(ps); |
---|
| 294 | } else { |
---|
| 295 | G4cout << "WARNING[" << qPassTrackLengthCmd->GetCommandPath() |
---|
| 296 | << "] : Quantity name, \"" << newVal << "\", is already existing." << G4endl; |
---|
| 297 | mesh->SetNullToCurrentPrimitiveScorer(); |
---|
| 298 | } |
---|
| 299 | } else if(command== qFlatSurfCurrCmd){ |
---|
| 300 | if(!mesh->FindPrimitiveScorer(newVal)) { |
---|
| 301 | G4PSFlatSurfaceCurrent3D* ps = |
---|
| 302 | new G4PSFlatSurfaceCurrent3D(token[0],StoI(token[1])); |
---|
| 303 | ps->Weighted(StoB(token[2])); |
---|
| 304 | ps->DivideByArea(StoB(token[3])); |
---|
| 305 | mesh->SetPrimitiveScorer(ps); |
---|
| 306 | } else { |
---|
| 307 | G4cout << "WARNING[" << qFlatSurfCurrCmd->GetCommandPath() |
---|
| 308 | << "] : Quantity name, \"" << newVal << "\", is already existing." << G4endl; |
---|
| 309 | mesh->SetNullToCurrentPrimitiveScorer(); |
---|
| 310 | } |
---|
| 311 | } else if(command== qFlatSurfFluxCmd){ |
---|
| 312 | if(!mesh->FindPrimitiveScorer(newVal)) { |
---|
| 313 | mesh->SetPrimitiveScorer( |
---|
| 314 | new G4PSFlatSurfaceFlux3D(token[0],StoI(token[1]))); |
---|
| 315 | } else { |
---|
| 316 | G4cout << "WARNING[" << qFlatSurfFluxCmd->GetCommandPath() |
---|
| 317 | << "] : Quantity name, \"" << newVal << "\", is already existing." << G4endl; |
---|
| 318 | mesh->SetNullToCurrentPrimitiveScorer(); |
---|
| 319 | } |
---|
| 320 | // } else if(command== qSphereSurfCurrCmd){ |
---|
| 321 | // if(!mesh->FindPrimitiveScorer(newVal)) { |
---|
| 322 | // G4PSSphereSurfaceCurrent3D* ps = |
---|
| 323 | // new G4PSSphereSurfaceCurrent3D(token[0],StoI(token[1])); |
---|
| 324 | // ps->Weighted(StoB(token[2])); |
---|
| 325 | // ps->DivideByArea(StoB(token[3])); |
---|
| 326 | // mesh->SetPrimitiveScorer(ps); |
---|
| 327 | // } else { |
---|
| 328 | // G4cout << "WARNING[" << qSphereSurfCurrCmd->GetCommandPath() |
---|
| 329 | // << "] : Quantity name, \"" << newVal << "\", is already existing." << G4endl; |
---|
| 330 | // mesh->SetNullToCurrentPrimitiveScorer(); |
---|
| 331 | // } |
---|
| 332 | // } else if(command== qSphereSurfFluxCmd){ |
---|
| 333 | // if(!mesh->FindPrimitiveScorer(newVal)) { |
---|
| 334 | // mesh->SetPrimitiveScorer( |
---|
| 335 | // new G4PSSphereSurfaceFlux3D(token[0], StoI(token[1]))); |
---|
| 336 | // } else { |
---|
| 337 | // G4cout << "WARNING[" << qSphereSurfFluxCmd->GetCommandPath() |
---|
| 338 | // << "] : Quantity name, \"" << newVal << "\", is already existing." << G4endl; |
---|
| 339 | // mesh->SetNullToCurrentPrimitiveScorer(); |
---|
| 340 | // } |
---|
| 341 | // } else if(command== qCylSurfCurrCmd){ |
---|
| 342 | // if(!mesh->FindPrimitiveScorer(newVal)) { |
---|
| 343 | // G4PSCylinderSurfaceCurrent3D* ps = |
---|
| 344 | // new G4PSCylinderSurfaceCurrent3D(token[0],StoI(token[1])); |
---|
| 345 | // ps->Weighted(StoB(token[2])); |
---|
| 346 | // ps->DivideByArea(StoB(token[3])); |
---|
| 347 | // mesh->SetPrimitiveScorer(ps); |
---|
| 348 | // } else { |
---|
| 349 | // G4cout << "WARNING[" << qCylSurfCurrCmd->GetCommandPath() |
---|
| 350 | // << "] : Quantity name, \"" << newVal << "\", is already existing." << G4endl; |
---|
| 351 | // mesh->SetNullToCurrentPrimitiveScorer(); |
---|
| 352 | // } |
---|
| 353 | // } else if(command== qCylSurfFluxCmd){ |
---|
| 354 | // if(!mesh->FindPrimitiveScorer(newVal)) { |
---|
| 355 | // mesh->SetPrimitiveScorer( |
---|
| 356 | // new G4PSCylinderSurfaceFlux3D(token[0], StoI(token[1]))); |
---|
| 357 | // } else { |
---|
| 358 | // G4cout << "WARNING[" << qCylSurfFluxCmd->GetCommandPath() |
---|
| 359 | // << "] : Quantity name, \"" << newVal << "\", is already existing." << G4endl; |
---|
| 360 | // mesh->SetNullToCurrentPrimitiveScorer(); |
---|
| 361 | // } |
---|
| 362 | } else if(command== qNofCollisionCmd){ |
---|
| 363 | if(!mesh->FindPrimitiveScorer(newVal)) { |
---|
| 364 | G4PSNofCollision3D* ps =new G4PSNofCollision3D(token[0]); |
---|
| 365 | ps->Weighted(StoB(token[1])); |
---|
| 366 | mesh->SetPrimitiveScorer(ps); |
---|
| 367 | } else { |
---|
| 368 | G4cout << "WARNING[" << qNofCollisionCmd->GetCommandPath() |
---|
| 369 | << "] : Quantity name, \"" << newVal << "\", is already existing." << G4endl; |
---|
| 370 | mesh->SetNullToCurrentPrimitiveScorer(); |
---|
| 371 | } |
---|
| 372 | } else if(command== qPopulationCmd){ |
---|
| 373 | if(!mesh->FindPrimitiveScorer(newVal)) { |
---|
| 374 | G4PSPopulation3D* ps =new G4PSPopulation3D(token[0]); |
---|
| 375 | ps->Weighted(StoB(token[1])); |
---|
| 376 | mesh->SetPrimitiveScorer(ps); |
---|
| 377 | } else { |
---|
| 378 | G4cout << "WARNING[" << qPopulationCmd->GetCommandPath() |
---|
| 379 | << "] : Quantity name, \"" << newVal << "\", is already existing." << G4endl; |
---|
| 380 | mesh->SetNullToCurrentPrimitiveScorer(); |
---|
| 381 | } |
---|
| 382 | } else if(command== qTrackCountCmd){ |
---|
| 383 | if(!mesh->FindPrimitiveScorer(newVal)) { |
---|
| 384 | G4PSTrackCounter3D* ps =new G4PSTrackCounter3D(token[0],StoI(token[1])); |
---|
| 385 | mesh->SetPrimitiveScorer(ps); |
---|
| 386 | } else { |
---|
| 387 | G4cout << "WARNING[" << qTrackCountCmd->GetCommandPath() |
---|
| 388 | << "] : Quantity name, \"" << newVal << "\", is already existing." << G4endl; |
---|
| 389 | mesh->SetNullToCurrentPrimitiveScorer(); |
---|
| 390 | } |
---|
| 391 | } else if(command== qTerminationCmd){ |
---|
| 392 | if(!mesh->FindPrimitiveScorer(newVal)) { |
---|
| 393 | G4PSTermination3D* ps =new G4PSTermination3D(token[0]); |
---|
| 394 | ps->Weighted(StoB(token[1])); |
---|
| 395 | mesh->SetPrimitiveScorer(ps); |
---|
| 396 | } else { |
---|
| 397 | G4cout << "WARNING[" << qTerminationCmd->GetCommandPath() |
---|
| 398 | << "] : Quantity name, \"" << newVal << "\", is already existing." << G4endl; |
---|
| 399 | mesh->SetNullToCurrentPrimitiveScorer(); |
---|
| 400 | } |
---|
| 401 | |
---|
| 402 | // |
---|
| 403 | // Filters |
---|
| 404 | // |
---|
| 405 | }else if(command== fchargedCmd){ |
---|
| 406 | if(!mesh->IsCurrentPrimitiveScorerNull()) { |
---|
| 407 | mesh->SetFilter(new G4SDChargedFilter(token[0])); |
---|
| 408 | } else { |
---|
| 409 | G4cout << "WARNING[" << fchargedCmd->GetCommandPath() |
---|
| 410 | << "] : Current quantity is not set. Set or touch a quantity first." << G4endl; |
---|
| 411 | } |
---|
| 412 | }else if(command== fneutralCmd){ |
---|
| 413 | if(!mesh->IsCurrentPrimitiveScorerNull()) { |
---|
| 414 | mesh->SetFilter(new G4SDNeutralFilter(token[0])); |
---|
| 415 | } else { |
---|
| 416 | G4cout << "WARNING[" << fneutralCmd->GetCommandPath() |
---|
| 417 | << "] : Current quantity is not set. Set or touch a quantity first." << G4endl; |
---|
| 418 | } |
---|
| 419 | }else if(command== fkinECmd){ |
---|
| 420 | if(!mesh->IsCurrentPrimitiveScorerNull()) { |
---|
| 421 | G4String& name = token[0]; |
---|
| 422 | G4double elow = StoD(token[1]); |
---|
| 423 | G4double ehigh = StoD(token[2]); |
---|
| 424 | mesh->SetFilter(new G4SDKineticEnergyFilter(name,elow,ehigh)); |
---|
| 425 | } else { |
---|
| 426 | G4cout << "WARNING[" << fkinECmd->GetCommandPath() |
---|
| 427 | << "] : Current quantity is not set. Set or touch a quantity first." << G4endl; |
---|
| 428 | } |
---|
| 429 | }else if(command== fparticleKinECmd){ |
---|
| 430 | if(!mesh->IsCurrentPrimitiveScorerNull()) { |
---|
| 431 | FParticleWithEnergyCommand(mesh,token); |
---|
| 432 | } else { |
---|
| 433 | G4cout << "WARNING[" << fparticleKinECmd->GetCommandPath() |
---|
| 434 | << "] : Current quantity is not set. Set or touch a quantity first." << G4endl; |
---|
| 435 | } |
---|
| 436 | } else if(command==fparticleCmd) { |
---|
| 437 | if(!mesh->IsCurrentPrimitiveScorerNull()) { |
---|
| 438 | FParticleCommand(mesh,token); |
---|
| 439 | } else { |
---|
| 440 | G4cout << "WARNING[" << fparticleCmd->GetCommandPath() |
---|
| 441 | << "] : Current quantity is not set. Set or touch a quantity first." << G4endl; |
---|
| 442 | } |
---|
| 443 | } |
---|
| 444 | } |
---|
| 445 | |
---|
| 446 | G4String G4ScoreQuantityMessenger::GetCurrentValue(G4UIcommand * /*command*/) |
---|
| 447 | { |
---|
| 448 | G4String val; |
---|
| 449 | |
---|
| 450 | return val; |
---|
| 451 | } |
---|
| 452 | |
---|
| 453 | void G4ScoreQuantityMessenger::FillTokenVec(G4String newValues, G4TokenVec& token){ |
---|
| 454 | |
---|
| 455 | G4Tokenizer next(newValues); |
---|
| 456 | G4String val; |
---|
| 457 | while ( !(val = next()).isNull() ) { |
---|
| 458 | token.push_back(val); |
---|
| 459 | } |
---|
| 460 | } |
---|
| 461 | |
---|
| 462 | |
---|
| 463 | void G4ScoreQuantityMessenger::FParticleCommand(G4VScoringMesh* mesh, G4TokenVec& token){ |
---|
| 464 | // |
---|
| 465 | // Filter name |
---|
| 466 | G4String name = token[0]; |
---|
| 467 | // |
---|
| 468 | // particle list |
---|
| 469 | std::vector<G4String> pnames; |
---|
| 470 | for ( G4int i = 1; i<(G4int)token.size(); i++){ |
---|
| 471 | pnames.push_back(token[i]); |
---|
| 472 | } |
---|
| 473 | // |
---|
| 474 | // Attach Filter |
---|
| 475 | mesh->SetFilter(new G4SDParticleFilter(name,pnames)); |
---|
| 476 | } |
---|
| 477 | |
---|
| 478 | void G4ScoreQuantityMessenger::FParticleWithEnergyCommand(G4VScoringMesh* mesh,G4TokenVec& token){ |
---|
| 479 | G4String& name = token[0]; |
---|
| 480 | G4double elow = StoD(token[1]); |
---|
| 481 | G4double ehigh= StoD(token[2]); |
---|
| 482 | G4double unitVal = G4UnitDefinition::GetValueOf(token[3]); |
---|
| 483 | G4SDParticleWithEnergyFilter* filter = |
---|
| 484 | new G4SDParticleWithEnergyFilter(name,elow*unitVal,ehigh*unitVal); |
---|
| 485 | for ( G4int i = 4; i < (G4int)token.size(); i++){ |
---|
| 486 | filter->add(token[i]); |
---|
| 487 | } |
---|
| 488 | mesh->SetFilter(filter); |
---|
| 489 | } |
---|
| 490 | |
---|