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