Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ScoreQuantityMessenger.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//
28// ---------------------------------------------------------------------
29// Modifications
30// 08-Oct-2010 T.Aso remove unit of G4PSPassageCellCurrent.
31// 24-Mar-2011 T.Aso Add StepChecker for debugging.
32// 24-Mar-2011 T.Aso Size and segmentation for replicated cylinder.
33// 01-Jun-2012 T.Aso Support weighted/dividedByArea options
34// in flatCurrent and flatFulx commands.
35// 27-Mar-2013 T.Aso Unit option in the kineticEnergy filter was
36// supported.
37//
38// ---------------------------------------------------------------------
39
41#include "G4ScoringManager.hh"
42#include "G4VScoringMesh.hh"
43#include "G4VPrimitiveScorer.hh"
44
45#include "G4PSCellCharge3D.hh"
46#include "G4PSCellFlux3D.hh"
51#include "G4PSDoseDeposit3D.hh"
53#include "G4PSNofStep3D.hh"
54#include "G4PSNofSecondary3D.hh"
55//
56#include "G4PSTrackLength3D.hh"
65#include "G4PSVolumeFlux3D.hh"
66#include "G4PSNofCollision3D.hh"
67#include "G4PSPopulation3D.hh"
68#include "G4PSTrackCounter3D.hh"
69#include "G4PSTermination3D.hh"
71
72#include "G4PSCellCharge.hh"
73#include "G4PSCellFlux.hh"
75#include "G4PSEnergyDeposit.hh"
76#include "G4PSDoseDeposit.hh"
77#include "G4PSNofStep.hh"
78#include "G4PSNofSecondary.hh"
79//
80#include "G4PSTrackLength.hh"
89#include "G4PSNofCollision.hh"
90#include "G4PSPopulation.hh"
91#include "G4PSTrackCounter.hh"
92#include "G4PSTermination.hh"
94
95//
96// For debug purpose
97#include "G4PSStepChecker3D.hh"
98
99#include "G4SDChargedFilter.hh"
100#include "G4SDNeutralFilter.hh"
102#include "G4SDParticleFilter.hh"
104
105#include "G4UIdirectory.hh"
108#include "G4UIcmdWithAString.hh"
109#include "G4UIcmdWithABool.hh"
112#include "G4UIcommand.hh"
113#include "G4UIparameter.hh"
114#include "G4Tokenizer.hh"
115#include "G4UnitsTable.hh"
116
118 : fSMan(SManager)
119{
120 QuantityCommands();
121 FilterCommands();
122}
123
124void G4ScoreQuantityMessenger::QuantityCommands()
125{
126 G4UIparameter* param;
127
128 //
129 // Quantity commands
130 quantityDir = new G4UIdirectory("/score/quantity/");
131 quantityDir->SetGuidance("Scoring quantity of the mesh.");
132 //
133 qTouchCmd = new G4UIcmdWithAString("/score/quantity/touch", this);
134 qTouchCmd->SetGuidance(
135 "Assign previously defined quantity to the current quantity.");
136 qTouchCmd->SetParameterName("qname", false);
137 //
138 qGetUnitCmd = new G4UIcmdWithoutParameter("/score/quantity/get/unit", this);
139 qGetUnitCmd->SetGuidance("Print output unit of the current quantity.");
140 //
141 qSetUnitCmd = new G4UIcmdWithAString("/score/quantity/set/unit", this);
142 qSetUnitCmd->SetGuidance("Set output unit of the current quantity.");
143 qSetUnitCmd->SetParameterName("unit", false);
144
145 // Primitive Scorers
146 qeDepCmd = new G4UIcommand("/score/quantity/energyDeposit", this);
147 qeDepCmd->SetGuidance("Energy deposit scorer.");
148 qeDepCmd->SetGuidance("[usage] /score/quantity/energyDeposit qname unit");
149 qeDepCmd->SetGuidance(" qname :(String) scorer name");
150 qeDepCmd->SetGuidance(" unit :(String) unit");
151 param = new G4UIparameter("qname", 's', false);
152 qeDepCmd->SetParameter(param);
153 param = new G4UIparameter("unit", 's', true);
154 param->SetDefaultUnit("MeV");
155 qeDepCmd->SetParameter(param);
156 //
157 qCellChgCmd = new G4UIcommand("/score/quantity/cellCharge", this);
158 qCellChgCmd->SetGuidance("Cell charge scorer.");
159 qCellChgCmd->SetGuidance("[usage] /score/quantity/cellCharge qname unit");
160 qCellChgCmd->SetGuidance(" qname :(String) scorer name");
161 qCellChgCmd->SetGuidance(" unit :(String) unit");
162 param = new G4UIparameter("qname", 's', false);
163 qCellChgCmd->SetParameter(param);
164 param = new G4UIparameter("unit", 's', true);
165 param->SetDefaultUnit("e+");
166 qCellChgCmd->SetParameter(param);
167 //
168 qCellFluxCmd = new G4UIcommand("/score/quantity/cellFlux", this);
169 qCellFluxCmd->SetGuidance("Cell flux scorer.");
170 qCellFluxCmd->SetGuidance("[usage] /score/quantity/cellFlux qname unit");
171 qCellFluxCmd->SetGuidance(" qname :(String) scorer name");
172 qCellFluxCmd->SetGuidance(" unit :(String) unit");
173 param = new G4UIparameter("qname", 's', false);
174 qCellFluxCmd->SetParameter(param);
175 param = new G4UIparameter("unit", 's', true);
176 param->SetDefaultValue("percm2");
177 qCellFluxCmd->SetParameter(param);
178 //
179 qPassCellFluxCmd = new G4UIcommand("/score/quantity/passageCellFlux", this);
180 qPassCellFluxCmd->SetGuidance("Passage cell flux scorer");
181 qPassCellFluxCmd->SetGuidance(
182 "[usage] /score/quantity/passageCellFlux qname unit");
183 qPassCellFluxCmd->SetGuidance(" qname :(String) scorer name");
184 qPassCellFluxCmd->SetGuidance(" unit :(String) unit");
185 param = new G4UIparameter("qname", 's', false);
186 qPassCellFluxCmd->SetParameter(param);
187 param = new G4UIparameter("unit", 's', true);
188 param->SetDefaultValue("percm2");
189 qPassCellFluxCmd->SetParameter(param);
190 //
191 qdoseDepCmd = new G4UIcommand("/score/quantity/doseDeposit", this);
192 qdoseDepCmd->SetGuidance("Dose deposit scorer.");
193 qdoseDepCmd->SetGuidance("[usage] /score/quantity/doseDeposit qname unit");
194 qdoseDepCmd->SetGuidance(" qname :(String) scorer name");
195 qdoseDepCmd->SetGuidance(" unit :(String) unit");
196 param = new G4UIparameter("qname", 's', false);
197 qdoseDepCmd->SetParameter(param);
198 param = new G4UIparameter("unit", 's', true);
199 param->SetDefaultUnit("Gy");
200 qdoseDepCmd->SetParameter(param);
201 //
202 qnOfStepCmd = new G4UIcommand("/score/quantity/nOfStep", this);
203 qnOfStepCmd->SetGuidance("Number of step scorer.");
204 qnOfStepCmd->SetGuidance("[usage] /score/quantity/nOfStep qname");
205 qnOfStepCmd->SetGuidance("[usage] /score/quantity/nOfStep qname bflag");
206 qnOfStepCmd->SetGuidance(" qname :(String) scorer name");
207 qnOfStepCmd->SetGuidance(" bflag :(Bool) Skip zero step ");
208 qnOfStepCmd->SetGuidance(" at geometry boundary if true");
209 param = new G4UIparameter("qname", 's', false);
210 qnOfStepCmd->SetParameter(param);
211 param = new G4UIparameter("bflag", 'b', true);
212 param->SetDefaultValue("false");
213 qnOfStepCmd->SetParameter(param);
214 //
215 qnOfSecondaryCmd = new G4UIcommand("/score/quantity/nOfSecondary", this);
216 qnOfSecondaryCmd->SetGuidance("Number of secondary scorer.");
217 qnOfSecondaryCmd->SetGuidance("[usage] /score/quantity/nOfSecondary qname");
218 qnOfSecondaryCmd->SetGuidance(" qname :(String) scorer name");
219 param = new G4UIparameter("qname", 's', false);
220 qnOfSecondaryCmd->SetParameter(param);
221 //
222 qTrackLengthCmd = new G4UIcommand("/score/quantity/trackLength", this);
223 qTrackLengthCmd->SetGuidance("Track length scorer.");
224 qTrackLengthCmd->SetGuidance(
225 "[usage] /score/quantity/trackLength qname wflag kflag vflag unit");
226 qTrackLengthCmd->SetGuidance(" qname :(String) scorer name");
227 qTrackLengthCmd->SetGuidance(" wflag :(Bool) weighted");
228 qTrackLengthCmd->SetGuidance(" kflag :(Bool) multiply kinetic energy");
229 qTrackLengthCmd->SetGuidance(" vflag :(Bool) divide by velocity");
230 qTrackLengthCmd->SetGuidance(" unit :(String) unit");
231 param = new G4UIparameter("qname", 's', false);
232 qTrackLengthCmd->SetParameter(param);
233 param = new G4UIparameter("wflag", 'b', true);
234 param->SetDefaultValue("false");
235 qTrackLengthCmd->SetParameter(param);
236 param = new G4UIparameter("kflag", 'b', true);
237 param->SetDefaultValue("false");
238 qTrackLengthCmd->SetParameter(param);
239 param = new G4UIparameter("vflag", 'b', true);
240 param->SetDefaultValue("false");
241 qTrackLengthCmd->SetParameter(param);
242 param = new G4UIparameter("unit", 's', true);
243 param->SetDefaultValue("mm");
244 qTrackLengthCmd->SetParameter(param);
245 //
246 qPassCellCurrCmd =
247 new G4UIcommand("/score/quantity/passageCellCurrent", this);
248 qPassCellCurrCmd->SetGuidance("Passage cell current scorer.");
249 qPassCellCurrCmd->SetGuidance(
250 "[usage] /score/quantity/passageCellCurrent qname wflag");
251 qPassCellCurrCmd->SetGuidance(" qname :(String) scorer name");
252 qPassCellCurrCmd->SetGuidance(" wflag :(Bool) weighted");
253 param = new G4UIparameter("qname", 's', false);
254 qPassCellCurrCmd->SetParameter(param);
255 param = new G4UIparameter("wflag", 'b', true);
256 param->SetDefaultValue("true");
257 qPassCellCurrCmd->SetParameter(param);
258 //
259 qPassTrackLengthCmd =
260 new G4UIcommand("/score/quantity/passageTrackLength", this);
261 qPassTrackLengthCmd->SetGuidance("Passage track length scorer.");
262 qPassTrackLengthCmd->SetGuidance(
263 "[usage] /score/quantity/passageTrackLength qname wflag unit");
264 qPassTrackLengthCmd->SetGuidance(" qname :(String) scorer name");
265 qPassTrackLengthCmd->SetGuidance(" wflag :(Bool) weighted");
266 qPassTrackLengthCmd->SetGuidance(" unit :(Bool) unit");
267 param = new G4UIparameter("qname", 's', false);
268 qPassTrackLengthCmd->SetParameter(param);
269 param = new G4UIparameter("wflag", 'b', true);
270 param->SetDefaultValue("true");
271 qPassTrackLengthCmd->SetParameter(param);
272 param = new G4UIparameter("unit", 's', true);
273 param->SetDefaultUnit("mm");
274 qPassTrackLengthCmd->SetParameter(param);
275 //
276 qFlatSurfCurrCmd =
277 new G4UIcommand("/score/quantity/flatSurfaceCurrent", this);
278 qFlatSurfCurrCmd->SetGuidance("Flat surface current Scorer.");
279 qFlatSurfCurrCmd->SetGuidance(
280 "[usage] /score/quantity/flatSurfaceCurrent qname dflag wflag aflag unit");
281 qFlatSurfCurrCmd->SetGuidance(" qname :(String) scorer name");
282 qFlatSurfCurrCmd->SetGuidance(" dflag :(Int) direction flag");
283 qFlatSurfCurrCmd->SetGuidance(" : 0 = Both In and Out");
284 qFlatSurfCurrCmd->SetGuidance(" : 1 = In only");
285 qFlatSurfCurrCmd->SetGuidance(" : 2 = Out only");
286 qFlatSurfCurrCmd->SetGuidance(" wflag :(Bool) weighted");
287 qFlatSurfCurrCmd->SetGuidance(" aflag :(Bool) divide by area");
288 qFlatSurfCurrCmd->SetGuidance(" unit :(String) unit");
289 param = new G4UIparameter("qname", 's', false);
290 qFlatSurfCurrCmd->SetParameter(param);
291 param = new G4UIparameter("dflag", 'i', true);
292 param->SetDefaultValue("0");
293 qFlatSurfCurrCmd->SetParameter(param);
294 param = new G4UIparameter("wflag", 'b', true);
295 param->SetDefaultValue("true");
296 qFlatSurfCurrCmd->SetParameter(param);
297 param = new G4UIparameter("aflag", 'b', true);
298 param->SetDefaultValue("true");
299 qFlatSurfCurrCmd->SetParameter(param);
300 param = new G4UIparameter("unit", 's', true);
301 param->SetDefaultValue("percm2");
302 qFlatSurfCurrCmd->SetParameter(param);
303 //
304 qFlatSurfFluxCmd = new G4UIcommand("/score/quantity/flatSurfaceFlux", this);
305 qFlatSurfFluxCmd->SetGuidance("Flat surface flux scorer.");
306 qFlatSurfFluxCmd->SetGuidance(
307 "[usage] /score/quantity/flatSurfaceFlux qname dflag unit");
308 qFlatSurfFluxCmd->SetGuidance(" qname :(String) scorer name");
309 qFlatSurfFluxCmd->SetGuidance(" dflag :(Int) direction flag");
310 qFlatSurfFluxCmd->SetGuidance(" : 0 = Both In and Out");
311 qFlatSurfFluxCmd->SetGuidance(" : 1 = In only");
312 qFlatSurfFluxCmd->SetGuidance(" : 2 = Out only");
313 qFlatSurfFluxCmd->SetGuidance(" wflag :(Bool) weighted");
314 qFlatSurfFluxCmd->SetGuidance(" aflag :(Bool) divide by area");
315 qFlatSurfFluxCmd->SetGuidance(" unit :(String) unit");
316 param = new G4UIparameter("qname", 's', false);
317 qFlatSurfFluxCmd->SetParameter(param);
318 param = new G4UIparameter("dflag", 'i', true);
319 param->SetDefaultValue("0");
320 qFlatSurfFluxCmd->SetParameter(param);
321 param = new G4UIparameter("wflag", 'b', true);
322 param->SetDefaultValue("true");
323 qFlatSurfFluxCmd->SetParameter(param);
324 param = new G4UIparameter("aflag", 'b', true);
325 param->SetDefaultValue("true");
326 qFlatSurfFluxCmd->SetParameter(param);
327 param = new G4UIparameter("unit", 's', true);
328 param->SetDefaultValue("percm2");
329 qFlatSurfFluxCmd->SetParameter(param);
330 //
331
332 qVolFluxCmd = new G4UIcommand("/score/quantity/volumeFlux", this);
333 qVolFluxCmd->SetGuidance("Volume flux scorer.");
334 qVolFluxCmd->SetGuidance(
335 "This scorer scores the number of particles getting into the volume "
336 "without normalized by the surface area.");
337 qVolFluxCmd->SetGuidance(
338 "[usage] /score/quantity/volumeFlux qname divcos dflag");
339 qVolFluxCmd->SetGuidance(" qname :(String) scorer name");
340 qVolFluxCmd->SetGuidance(" divcos :(Bool) divide by cos(theta), where theta "
341 "is the incident angle (default : false)");
342 qVolFluxCmd->SetGuidance(
343 " dflag :(Int) direction, 1 : inward (default), 2 : outward");
344 param = new G4UIparameter("qname", 's', false);
345 qVolFluxCmd->SetParameter(param);
346 param = new G4UIparameter("divcos", 'b', true);
347 param->SetDefaultValue(0);
348 qVolFluxCmd->SetParameter(param);
349 param = new G4UIparameter("dflag", 'i', true);
350 param->SetParameterRange("dflag>=1 && dflag<=2");
351 param->SetDefaultValue(1);
352 qVolFluxCmd->SetParameter(param);
353
354 qNofCollisionCmd = new G4UIcommand("/score/quantity/nOfCollision", this);
355 qNofCollisionCmd->SetGuidance("Number of collision scorer.");
356 qNofCollisionCmd->SetGuidance(
357 "[usage] /score/quantity/nOfCollision qname wflag");
358 qNofCollisionCmd->SetGuidance(" qname :(String) scorer name");
359 param = new G4UIparameter("qname", 's', false);
360 qNofCollisionCmd->SetParameter(param);
361 param = new G4UIparameter("wflag", 'b', true);
362 param->SetDefaultValue("false");
363 qNofCollisionCmd->SetParameter(param);
364 //
365 qPopulationCmd = new G4UIcommand("/score/quantity/population", this);
366 qPopulationCmd->SetGuidance("Population scorer.");
367 qPopulationCmd->SetGuidance("[usage] /score/quantity/population qname wflag");
368 qPopulationCmd->SetGuidance(" qname :(String) scorer name");
369 qPopulationCmd->SetGuidance(" wflag :(Bool) weighted");
370 param = new G4UIparameter("qname", 's', false);
371 qPopulationCmd->SetParameter(param);
372 param = new G4UIparameter("wflag", 'b', true);
373 param->SetDefaultValue("false");
374 qPopulationCmd->SetParameter(param);
375
376 //
377 qTrackCountCmd = new G4UIcommand("/score/quantity/nOfTrack", this);
378 qTrackCountCmd->SetGuidance("Number of track scorer.");
379 qTrackCountCmd->SetGuidance(
380 "[usage] /score/quantity/nOfTrack qname dflag wflag");
381 qTrackCountCmd->SetGuidance(" qname :(String) scorer name");
382 qTrackCountCmd->SetGuidance(" dflag :(Int) direction");
383 qTrackCountCmd->SetGuidance(" : 0 = Both In and Out");
384 qTrackCountCmd->SetGuidance(" : 1 = In only");
385 qTrackCountCmd->SetGuidance(" : 2 = Out only");
386 qTrackCountCmd->SetGuidance(" wflag :(Bool) weighted");
387 param = new G4UIparameter("qname", 's', false);
388 qTrackCountCmd->SetParameter(param);
389 param = new G4UIparameter("dflag", 'i', true);
390 param->SetDefaultValue("0");
391 qTrackCountCmd->SetParameter(param);
392 param = new G4UIparameter("wflag", 'b', true);
393 param->SetDefaultValue("false");
394 qTrackCountCmd->SetParameter(param);
395
396 //
397 qTerminationCmd = new G4UIcommand("/score/quantity/nOfTerminatedTrack", this);
398 qTerminationCmd->SetGuidance("Number of terminated tracks scorer.");
399 qTerminationCmd->SetGuidance(
400 "[usage] /score/quantity/nOfTerminatedTrack qname wflag");
401 qTerminationCmd->SetGuidance(" qname :(String) scorer name");
402 qTerminationCmd->SetGuidance(" wflag :(Bool) weighted");
403 param = new G4UIparameter("qname", 's', false);
404 qTerminationCmd->SetParameter(param);
405 param = new G4UIparameter("wflag", 'b', true);
406 param->SetDefaultValue("false");
407 qTerminationCmd->SetParameter(param);
408
409 //
410 qMinKinEAtGeneCmd =
411 new G4UIcommand("/score/quantity/minKinEAtGeneration", this);
412 qMinKinEAtGeneCmd->SetGuidance("Min Kinetic Energy at Generation");
413 qMinKinEAtGeneCmd->SetGuidance(
414 "[usage] /score/quantity/minKinEAtGeneration qname unit");
415 qMinKinEAtGeneCmd->SetGuidance(" qname :(String) scorer name");
416 qMinKinEAtGeneCmd->SetGuidance(" unit :(String) unit name");
417 param = new G4UIparameter("qname", 's', false);
418 qMinKinEAtGeneCmd->SetParameter(param);
419 param = new G4UIparameter("unit", 's', true);
420 param->SetDefaultUnit("MeV");
421 qMinKinEAtGeneCmd->SetParameter(param);
422 //
423 qStepCheckerCmd = new G4UIcommand("/score/quantity/stepChecker", this);
424 qStepCheckerCmd->SetGuidance("Display a comment when this PS is invoked");
425 qStepCheckerCmd->SetGuidance("[usage] /score/quantity/stepChecker qname");
426 qStepCheckerCmd->SetGuidance(" qname :(String) scorer name");
427 param = new G4UIparameter("qname", 's', false);
428 qStepCheckerCmd->SetParameter(param);
429}
430
431void G4ScoreQuantityMessenger::FilterCommands()
432{
433 G4UIparameter* param;
434
435 //
436 // Filter commands
437 filterDir = new G4UIdirectory("/score/filter/");
438 filterDir->SetGuidance(" Scoring filter commands.");
439 //
440 fchargedCmd = new G4UIcmdWithAString("/score/filter/charged", this);
441 fchargedCmd->SetGuidance("Charged particle filter.");
442 fchargedCmd->SetParameterName("fname", false);
443 //
444 fneutralCmd = new G4UIcmdWithAString("/score/filter/neutral", this);
445 fneutralCmd->SetGuidance("Neutral particle filter.");
446 fneutralCmd->SetParameterName("fname", false);
447 //
448 fkinECmd = new G4UIcommand("/score/filter/kineticEnergy", this);
449 fkinECmd->SetGuidance("Kinetic energy filter.");
450 fkinECmd->SetGuidance(
451 "[usage] /score/filter/kineticEnergy fname Elow Ehigh unit");
452 fkinECmd->SetGuidance(" fname :(String) Filter Name ");
453 fkinECmd->SetGuidance(" Elow :(Double) Lower edge of kinetic energy");
454 fkinECmd->SetGuidance(" Ehigh :(Double) Higher edge of kinetic energy");
455 fkinECmd->SetGuidance(" unit :(String) unit of given kinetic energy");
456 param = new G4UIparameter("fname", 's', false);
457 fkinECmd->SetParameter(param);
458 param = new G4UIparameter("elow", 'd', true);
459 param->SetDefaultValue("0.0");
460 fkinECmd->SetParameter(param);
461 param = new G4UIparameter("ehigh", 'd', true);
462 fkinECmd->SetParameter(param);
463 G4String smax = DtoS(DBL_MAX);
464 param->SetDefaultValue(smax);
465 param = new G4UIparameter("unit", 's', true);
466 param->SetDefaultUnit("keV");
467 fkinECmd->SetParameter(param);
468 //
469 fparticleCmd = new G4UIcommand("/score/filter/particle", this);
470 fparticleCmd->SetGuidance("Particle filter.");
471 fparticleCmd->SetGuidance("[usage] /score/filter/particle fname p0 .. pn");
472 fparticleCmd->SetGuidance(" fname :(String) Filter Name ");
473 fparticleCmd->SetGuidance(" p0 .. pn :(String) particle names");
474 param = new G4UIparameter("fname", 's', false);
475 fparticleCmd->SetParameter(param);
476 param = new G4UIparameter("particlelist", 's', false);
477 param->SetDefaultValue("");
478 fparticleCmd->SetParameter(param);
479 //
480 //
481 //
482 fparticleKinECmd =
483 new G4UIcommand("/score/filter/particleWithKineticEnergy", this);
484 fparticleKinECmd->SetGuidance("Particle with kinetic energy filter.");
485 fparticleKinECmd->SetGuidance(
486 "[usage] /score/filter/particleWithKineticEnergy fname Elow Ehigh unit p0 "
487 ".. pn");
488 fparticleKinECmd->SetGuidance(" fname :(String) Filter Name ");
489 fparticleKinECmd->SetGuidance(
490 " Elow :(Double) Lower edge of kinetic energy");
491 fparticleKinECmd->SetGuidance(
492 " Ehigh :(Double) Higher edge of kinetic energy");
493 fparticleKinECmd->SetGuidance(
494 " unit :(String) unit of given kinetic energy");
495 fparticleKinECmd->SetGuidance(" p0 .. pn :(String) particle names");
496 param = new G4UIparameter("fname", 's', false);
497 fparticleKinECmd->SetParameter(param);
498 param = new G4UIparameter("elow", 'd', false);
499 param->SetDefaultValue("0.0");
500 fparticleKinECmd->SetParameter(param);
501 param = new G4UIparameter("ehigh", 'd', true);
502 param->SetDefaultValue(smax);
503 fparticleKinECmd->SetParameter(param);
504 param = new G4UIparameter("unit", 's', true);
505 param->SetDefaultUnit("keV");
506 fparticleKinECmd->SetParameter(param);
507 param = new G4UIparameter("particlelist", 's', false);
508 param->SetDefaultValue("");
509 fparticleKinECmd->SetParameter(param);
510 //
511 //
512}
513
515{
516 delete quantityDir;
517 delete qTouchCmd;
518 delete qGetUnitCmd;
519 delete qSetUnitCmd;
520
521 //
522 delete qCellChgCmd;
523 delete qCellFluxCmd;
524 delete qPassCellFluxCmd;
525 delete qeDepCmd;
526 delete qdoseDepCmd;
527 delete qnOfStepCmd;
528 delete qnOfSecondaryCmd;
529 //
530 delete qTrackLengthCmd;
531 delete qPassCellCurrCmd;
532 delete qPassTrackLengthCmd;
533 delete qFlatSurfCurrCmd;
534 delete qFlatSurfFluxCmd;
535 delete qVolFluxCmd;
536
537 delete qNofCollisionCmd;
538 delete qPopulationCmd;
539 delete qTrackCountCmd;
540 delete qTerminationCmd;
541 delete qMinKinEAtGeneCmd;
542 //
543 delete qStepCheckerCmd;
544 //
545 delete filterDir;
546 delete fchargedCmd;
547 delete fneutralCmd;
548 delete fkinECmd;
549 delete fparticleCmd;
550 delete fparticleKinECmd;
551}
552
554 G4String newVal)
555{
556 using MeshShape = G4VScoringMesh::MeshShape;
557
559
560 //
561 // Get Current Mesh
562 //
563 G4VScoringMesh* mesh = fSMan->GetCurrentMesh();
564 if(mesh == nullptr)
565 {
566 ed << "ERROR : No mesh is currently open. Open/create a mesh first. "
567 "Command ignored.";
568 command->CommandFailed(ed);
569 return;
570 }
571 // Mesh type
572 auto shape = mesh->GetShape();
573 // Tokens
574 G4TokenVec token;
575 FillTokenVec(newVal, token);
576 //
577 // Commands for Current Mesh
578 if(command == qTouchCmd)
579 {
580 mesh->SetCurrentPrimitiveScorer(newVal);
581 }
582 else if(command == qGetUnitCmd)
583 {
584 G4cout << "Unit: " << mesh->GetCurrentPSUnit() << G4endl;
585 }
586 else if(command == qSetUnitCmd)
587 {
588 mesh->SetCurrentPSUnit(newVal);
589 }
590 else if(command == qCellChgCmd)
591 {
592 if(CheckMeshPS(mesh, token[0], command))
593 {
594 G4PSCellCharge* ps = nullptr;
595 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
596 {
597 ps = new G4PSCellCharge(token[0], mesh->GetCopyNumberLevel());
598 }
599 else
600 {
601 ps = new G4PSCellCharge3D(token[0]);
602 }
603 ps->SetUnit(token[1]);
604 mesh->SetPrimitiveScorer(ps);
605 }
606 }
607 else if(command == qCellFluxCmd)
608 {
609 if(CheckMeshPS(mesh, token[0], command))
610 {
611 G4PSCellFlux* ps = nullptr;
612 if(shape == MeshShape::box)
613 {
614 ps = new G4PSCellFlux3D(token[0]);
615 }
616 else if(shape == MeshShape::cylinder)
617 {
618 auto pps =
619 new G4PSCellFluxForCylinder3D(token[0]);
620 pps->SetCylinderSize(mesh->GetSize(),mesh->GetStartAngle(),mesh->GetAngleSpan());
621 G4int nSeg[3];
622 mesh->GetNumberOfSegments(nSeg);
623 pps->SetNumberOfSegments(nSeg);
624 ps = pps;
625 }
626 else if(shape == MeshShape::realWorldLogVol)
627 {
628 ed << "Cell flux for real world volume is not yet supported. Command "
629 "ignored.";
630 command->CommandFailed(ed);
631 return;
632 }
633 else if(shape == MeshShape::probe)
634 {
635 ps = new G4PSCellFlux(token[0]);
636 }
637 ps->SetUnit(token[1]);
638 mesh->SetPrimitiveScorer(ps);
639 }
640 }
641 else if(command == qPassCellFluxCmd)
642 {
643 if(CheckMeshPS(mesh, token[0], command))
644 {
645 G4PSPassageCellFlux* ps = nullptr;
646 if(shape == MeshShape::box)
647 {
648 ps = new G4PSPassageCellFlux3D(token[0]);
649 }
650 else if(shape == MeshShape::cylinder)
651 {
652 auto pps =
654 pps->SetCylinderSize(mesh->GetSize(),mesh->GetStartAngle(),mesh->GetAngleSpan());
655 G4int nSeg[3];
656 mesh->GetNumberOfSegments(nSeg);
657 pps->SetNumberOfSegments(nSeg);
658 ps = pps;
659 }
660 else if(shape == MeshShape::realWorldLogVol)
661 {
662 ed << "Passing cell flux for real world volume is not yet supported. "
663 "Command ignored.";
664 command->CommandFailed(ed);
665 return;
666 }
667 else if(shape == MeshShape::probe)
668 {
669 ps = new G4PSPassageCellFlux(token[0]);
670 }
671 ps->SetUnit(token[1]);
672 mesh->SetPrimitiveScorer(ps);
673 }
674 }
675 else if(command == qeDepCmd)
676 {
677 if(CheckMeshPS(mesh, token[0], command))
678 {
679 G4PSEnergyDeposit* ps = nullptr;
680 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
681 {
682 ps = new G4PSEnergyDeposit(token[0], mesh->GetCopyNumberLevel());
683 }
684 else
685 {
686 ps = new G4PSEnergyDeposit3D(token[0]);
687 }
688 ps->SetUnit(token[1]);
689 mesh->SetPrimitiveScorer(ps);
690 }
691 }
692 else if(command == qdoseDepCmd)
693 {
694 if(CheckMeshPS(mesh, token[0], command))
695 {
696 G4PSDoseDeposit* ps = nullptr;
697 if(shape == MeshShape::box)
698 {
699 ps = new G4PSDoseDeposit3D(token[0]);
700 }
701 else if(shape == MeshShape::cylinder)
702 {
703 auto pps =
704 new G4PSDoseDepositForCylinder3D(token[0]);
705 pps->SetUnit(token[1]);
706 pps->SetCylinderSize(mesh->GetSize(),mesh->GetStartAngle(),mesh->GetAngleSpan());
707 G4int nSeg[3];
708 mesh->GetNumberOfSegments(nSeg);
709 pps->SetNumberOfSegments(nSeg);
710 ps = pps;
711 }
712 else if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
713 {
714 ps = new G4PSDoseDeposit(token[0], mesh->GetCopyNumberLevel());
715 }
716 ps->SetUnit(token[1]);
717 mesh->SetPrimitiveScorer(ps);
718 }
719 }
720 else if(command == qnOfStepCmd)
721 {
722 if(CheckMeshPS(mesh, token[0], command))
723 {
724 G4PSNofStep* ps = nullptr;
725 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
726 {
727 ps = new G4PSNofStep(token[0], mesh->GetCopyNumberLevel());
728 }
729 else
730 {
731 ps = new G4PSNofStep3D(token[0]);
732 }
733 ps->SetBoundaryFlag(StoB(token[1]));
734 mesh->SetPrimitiveScorer(ps);
735 }
736 }
737 else if(command == qnOfSecondaryCmd)
738 {
739 if(CheckMeshPS(mesh, token[0], command))
740 {
741 G4PSNofSecondary* ps = nullptr;
742 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
743 {
744 ps = new G4PSNofSecondary(token[0], mesh->GetCopyNumberLevel());
745 }
746 else
747 {
748 ps = new G4PSNofSecondary3D(token[0]);
749 }
750 mesh->SetPrimitiveScorer(ps);
751 }
752 }
753 else if(command == qTrackLengthCmd)
754 {
755 if(CheckMeshPS(mesh, token[0], command))
756 {
757 G4PSTrackLength* ps = nullptr;
758 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
759 {
760 ps = new G4PSTrackLength(token[0], mesh->GetCopyNumberLevel());
761 }
762 else
763 {
764 ps = new G4PSTrackLength3D(token[0]);
765 }
766 ps->Weighted(StoB(token[1]));
767 ps->MultiplyKineticEnergy(StoB(token[2]));
768 ps->DivideByVelocity(StoB(token[3]));
769 ps->SetUnit(token[4]);
770 mesh->SetPrimitiveScorer(ps);
771 }
772 }
773 else if(command == qPassCellCurrCmd)
774 {
775 if(CheckMeshPS(mesh, token[0], command))
776 {
777 G4PSPassageCellCurrent* ps = nullptr;
778 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
779 {
780 ps = new G4PSPassageCellCurrent(token[0], mesh->GetCopyNumberLevel());
781 }
782 else
783 {
784 ps = new G4PSPassageCellCurrent3D(token[0]);
785 }
786 ps->Weighted(StoB(token[1]));
787 mesh->SetPrimitiveScorer(ps);
788 }
789 }
790 else if(command == qPassTrackLengthCmd)
791 {
792 if(CheckMeshPS(mesh, token[0], command))
793 {
794 G4PSPassageTrackLength* ps = nullptr;
795 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
796 {
797 ps = new G4PSPassageTrackLength(token[0], mesh->GetCopyNumberLevel());
798 }
799 else
800 {
801 ps = new G4PSPassageTrackLength3D(token[0]);
802 }
803 ps->Weighted(StoB(token[1]));
804 ps->SetUnit(token[2]);
805 mesh->SetPrimitiveScorer(ps);
806 }
807 }
808 else if(command == qFlatSurfCurrCmd)
809 {
810 if(CheckMeshPS(mesh, token[0], command))
811 {
812 G4PSFlatSurfaceCurrent* ps = nullptr;
813 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
814 {
815 ps = new G4PSFlatSurfaceCurrent(token[0], StoI(token[1]),
816 mesh->GetCopyNumberLevel());
817 }
818 else
819 {
820 ps = new G4PSFlatSurfaceCurrent3D(token[0], StoI(token[1]));
821 }
822 ps->Weighted(StoB(token[2]));
823 ps->DivideByArea(StoB(token[3]));
824 if(StoB(token[3]))
825 {
826 ps->SetUnit(token[4]);
827 }
828 else
829 {
830 ps->SetUnit("");
831 }
832 mesh->SetPrimitiveScorer(ps);
833 }
834 }
835 else if(command == qFlatSurfFluxCmd)
836 {
837 if(CheckMeshPS(mesh, token[0], command))
838 {
839 G4PSFlatSurfaceFlux* ps = nullptr;
840 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
841 {
842 ps = new G4PSFlatSurfaceFlux(token[0], StoI(token[1]),
843 mesh->GetCopyNumberLevel());
844 }
845 else
846 {
847 ps = new G4PSFlatSurfaceFlux3D(token[0], StoI(token[1]));
848 }
849 ps->Weighted(StoB(token[2]));
850 ps->DivideByArea(StoB(token[3]));
851 if(StoB(token[3]))
852 {
853 ps->SetUnit(token[4]);
854 }
855 else
856 {
857 ps->SetUnit("");
858 }
859 mesh->SetPrimitiveScorer(ps);
860 }
861 }
862 else if(command == qVolFluxCmd)
863 {
864 if(CheckMeshPS(mesh, token[0], command))
865 {
866 G4PSVolumeFlux* ps = nullptr;
867 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
868 {
869 ps = new G4PSVolumeFlux(token[0], StoI(token[2]),
870 mesh->GetCopyNumberLevel());
871 }
872 else
873 {
874 ps = new G4PSVolumeFlux3D(token[0], StoI(token[2]));
875 }
876 ps->SetDivCos(StoI(token[1]) != 0);
877 mesh->SetPrimitiveScorer(ps);
878 }
879 }
880 else if(command == qNofCollisionCmd)
881 {
882 if(CheckMeshPS(mesh, token[0], command))
883 {
884 G4PSNofCollision* ps = nullptr;
885 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
886 {
887 ps = new G4PSNofCollision3D(token[0], mesh->GetCopyNumberLevel());
888 }
889 else
890 {
891 ps = new G4PSNofCollision3D(token[0]);
892 }
893 ps->Weighted(StoB(token[1]));
894 mesh->SetPrimitiveScorer(ps);
895 }
896 }
897 else if(command == qPopulationCmd)
898 {
899 if(CheckMeshPS(mesh, token[0], command))
900 {
901 G4PSPopulation* ps = nullptr;
902 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
903 {
904 ps = new G4PSPopulation(token[0], mesh->GetCopyNumberLevel());
905 }
906 else
907 {
908 ps = new G4PSPopulation3D(token[0]);
909 }
910 ps->Weighted(StoB(token[1]));
911 mesh->SetPrimitiveScorer(ps);
912 }
913 }
914 else if(command == qTrackCountCmd)
915 {
916 if(CheckMeshPS(mesh, token[0], command))
917 {
918 G4PSTrackCounter* ps = nullptr;
919 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
920 {
921 ps = new G4PSTrackCounter(token[0], StoI(token[1]),
922 mesh->GetCopyNumberLevel());
923 }
924 else
925 {
926 ps = new G4PSTrackCounter3D(token[0], StoI(token[1]));
927 }
928 ps->Weighted(StoB(token[2]));
929 mesh->SetPrimitiveScorer(ps);
930 }
931 }
932 else if(command == qTerminationCmd)
933 {
934 if(CheckMeshPS(mesh, token[0], command))
935 {
936 G4PSTermination* ps = nullptr;
937 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
938 {
939 ps = new G4PSTermination(token[0], mesh->GetCopyNumberLevel());
940 }
941 else
942 {
943 ps = new G4PSTermination3D(token[0]);
944 }
945 ps->Weighted(StoB(token[1]));
946 mesh->SetPrimitiveScorer(ps);
947 }
948 }
949 else if(command == qMinKinEAtGeneCmd)
950 {
951 if(CheckMeshPS(mesh, token[0], command))
952 {
953 G4PSMinKinEAtGeneration* ps = nullptr;
954 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
955 {
956 ps = new G4PSMinKinEAtGeneration(token[0], mesh->GetCopyNumberLevel());
957 }
958 else
959 {
960 ps = new G4PSMinKinEAtGeneration3D(token[0]);
961 }
962 ps->SetUnit(token[1]);
963 mesh->SetPrimitiveScorer(ps);
964 }
965 }
966 else if(command == qStepCheckerCmd)
967 {
968 if(CheckMeshPS(mesh, token[0], command))
969 {
970 G4PSStepChecker* ps = nullptr;
971 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
972 {
973 ps = new G4PSStepChecker(token[0], mesh->GetCopyNumberLevel());
974 }
975 else
976 {
977 ps = new G4PSStepChecker3D(token[0]);
978 }
979 mesh->SetPrimitiveScorer(ps);
980 }
981
982 //
983 // Filters
984 //
985 }
986 else if(command == fchargedCmd)
987 {
989 {
990 mesh->SetFilter(new G4SDChargedFilter(token[0]));
991 }
992 else
993 {
994 ed << "WARNING[" << fchargedCmd->GetCommandPath()
995 << "] : Current quantity is not set. Set or touch a quantity first.";
996 command->CommandFailed(ed);
997 }
998 }
999 else if(command == fneutralCmd)
1000 {
1001 if(!mesh->IsCurrentPrimitiveScorerNull())
1002 {
1003 mesh->SetFilter(new G4SDNeutralFilter(token[0]));
1004 }
1005 else
1006 {
1007 ed << "WARNING[" << fneutralCmd->GetCommandPath()
1008 << "] : Current quantity is not set. Set or touch a quantity first.";
1009 command->CommandFailed(ed);
1010 }
1011 }
1012 else if(command == fkinECmd)
1013 {
1014 if(!mesh->IsCurrentPrimitiveScorerNull())
1015 {
1016 G4String& name = token[0];
1017 G4double elow = StoD(token[1]);
1018 G4double ehigh = StoD(token[2]);
1019 G4double unitVal = G4UnitDefinition::GetValueOf(token[3]);
1020 mesh->SetFilter(
1021 new G4SDKineticEnergyFilter(name, elow * unitVal, ehigh * unitVal));
1022 }
1023 else
1024 {
1025 ed << "WARNING[" << fkinECmd->GetCommandPath()
1026 << "] : Current quantity is not set. Set or touch a quantity first.";
1027 command->CommandFailed(ed);
1028 }
1029 }
1030 else if(command == fparticleKinECmd)
1031 {
1032 if(!mesh->IsCurrentPrimitiveScorerNull())
1033 {
1034 FParticleWithEnergyCommand(mesh, token);
1035 }
1036 else
1037 {
1038 ed << "WARNING[" << fparticleKinECmd->GetCommandPath()
1039 << "] : Current quantity is not set. Set or touch a quantity first.";
1040 command->CommandFailed(ed);
1041 }
1042 }
1043 else if(command == fparticleCmd)
1044 {
1045 if(!mesh->IsCurrentPrimitiveScorerNull())
1046 {
1047 FParticleCommand(mesh, token);
1048 }
1049 else
1050 {
1051 ed << "WARNING[" << fparticleCmd->GetCommandPath()
1052 << "] : Current quantity is not set. Set or touch a quantity first.";
1053 command->CommandFailed(ed);
1054 }
1055 }
1056}
1057
1059{
1060 G4String val;
1061
1062 return val;
1063}
1064
1066 G4TokenVec& token)
1067{
1068 G4Tokenizer next(newValues);
1069 G4String val;
1070 while(!(val = next()).empty())
1071 { // Loop checking 12.18.2015 M.Asai
1072 token.push_back(val);
1073 }
1074}
1075
1077 G4TokenVec& token)
1078{
1079 //
1080 // Filter name
1081 G4String name = token[0];
1082 //
1083 // particle list
1084 std::vector<G4String> pnames;
1085 for(G4int i = 1; i < (G4int) token.size(); i++)
1086 {
1087 pnames.push_back(token[i]);
1088 }
1089 //
1090 // Attach Filter
1091 mesh->SetFilter(new G4SDParticleFilter(name, pnames));
1092}
1093
1095 G4TokenVec& token)
1096{
1097 G4String& name = token[0];
1098 G4double elow = StoD(token[1]);
1099 G4double ehigh = StoD(token[2]);
1100 G4double unitVal = G4UnitDefinition::GetValueOf(token[3]);
1101 auto filter =
1102 new G4SDParticleWithEnergyFilter(name, elow * unitVal, ehigh * unitVal);
1103 for(G4int i = 4; i < (G4int) token.size(); i++)
1104 {
1105 filter->add(token[i]);
1106 }
1107 mesh->SetFilter(filter);
1108}
1109
1111 G4String& psname,
1112 G4UIcommand* command)
1113{
1114 if(!mesh->FindPrimitiveScorer(psname))
1115 {
1116 return true;
1117 }
1118
1120 ed << "WARNING[" << qTouchCmd->GetCommandPath() << "] : Quantity name, \""
1121 << psname << "\", is already existing.";
1122 command->CommandFailed(ed);
1124 return false;
1125}
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
std::vector< G4String > G4TokenVec
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
virtual void SetUnit(const G4String &unit)
virtual void SetUnit(const G4String &unit)
virtual void SetUnit(const G4String &unit)
virtual void SetUnit(const G4String &unit)
void DivideByArea(G4bool flg=true)
void Weighted(G4bool flg=true)
virtual void SetUnit(const G4String &unit)
void Weighted(G4bool flg=true)
void DivideByArea(G4bool flg=true)
virtual void SetUnit(const G4String &unit)
virtual void SetUnit(const G4String &unit)
void Weighted(G4bool flg=true)
void SetBoundaryFlag(G4bool flg=true)
Definition: G4PSNofStep.hh:57
void Weighted(G4bool flg=true)
virtual void SetUnit(const G4String &unit)
virtual void SetUnit(const G4String &unit)
void Weighted(G4bool flg=true)
void Weighted(G4bool flg=true)
void Weighted(G4bool flg=true)
void Weighted(G4bool flg=true)
virtual void SetUnit(const G4String &unit)
void MultiplyKineticEnergy(G4bool flg=true)
void Weighted(G4bool flg=true)
void DivideByVelocity(G4bool flg=true)
void SetDivCos(G4bool val)
void FParticleWithEnergyCommand(G4VScoringMesh *mesh, G4TokenVec &token)
G4bool CheckMeshPS(G4VScoringMesh *mesh, G4String &psname, G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValues) override
G4String GetCurrentValue(G4UIcommand *) override
G4ScoreQuantityMessenger(G4ScoringManager *SManager)
void FillTokenVec(G4String newValues, G4TokenVec &token)
void FParticleCommand(G4VScoringMesh *mesh, G4TokenVec &token)
G4VScoringMesh * GetCurrentMesh() const
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
const G4String & GetCommandPath() const
Definition: G4UIcommand.hh:137
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:147
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:157
void CommandFailed(G4int errCode, G4ExceptionDescription &ed)
Definition: G4UIcommand.hh:179
G4double StoD(const G4String &s)
G4String DtoS(G4double a)
G4bool StoB(G4String s)
G4int StoI(const G4String &s)
void SetDefaultValue(const char *theDefaultValue)
void SetParameterRange(const char *theRange)
void SetDefaultUnit(const char *theDefaultUnit)
static G4double GetValueOf(const G4String &)
void SetFilter(G4VSDFilter *filter)
void SetNullToCurrentPrimitiveScorer()
G4ThreeVector GetSize() const
MeshShape GetShape() const
G4double GetStartAngle() const
void GetNumberOfSegments(G4int nSegment[3])
G4int GetCopyNumberLevel() const
void SetCurrentPSUnit(const G4String &unit)
void SetCurrentPrimitiveScorer(const G4String &name)
void SetPrimitiveScorer(G4VPrimitiveScorer *ps)
G4String GetCurrentPSUnit()
G4double GetAngleSpan() const
G4bool IsCurrentPrimitiveScorerNull()
G4bool FindPrimitiveScorer(const G4String &psname)
#define DBL_MAX
Definition: templates.hh:62