Geant4-11
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{
122}
123
125{
126 G4UIparameter* param;
127
128 //
129 // Filter commands
130 filterDir = new G4UIdirectory("/score/filter/");
131 filterDir->SetGuidance(" Scoring filter commands.");
132 //
133 fchargedCmd = new G4UIcmdWithAString("/score/filter/charged", this);
134 fchargedCmd->SetGuidance("Charged particle filter.");
135 fchargedCmd->SetParameterName("fname", false);
136 //
137 fneutralCmd = new G4UIcmdWithAString("/score/filter/neutral", this);
138 fneutralCmd->SetGuidance("Neutral particle filter.");
139 fneutralCmd->SetParameterName("fname", false);
140 //
141 fkinECmd = new G4UIcommand("/score/filter/kineticEnergy", this);
142 fkinECmd->SetGuidance("Kinetic energy filter.");
144 "[usage] /score/filter/kineticEnergy fname Elow Ehigh unit");
145 fkinECmd->SetGuidance(" fname :(String) Filter Name ");
146 fkinECmd->SetGuidance(" Elow :(Double) Lower edge of kinetic energy");
147 fkinECmd->SetGuidance(" Ehigh :(Double) Higher edge of kinetic energy");
148 fkinECmd->SetGuidance(" unit :(String) unit of given kinetic energy");
149 param = new G4UIparameter("fname", 's', false);
150 fkinECmd->SetParameter(param);
151 param = new G4UIparameter("elow", 'd', true);
152 param->SetDefaultValue("0.0");
153 fkinECmd->SetParameter(param);
154 param = new G4UIparameter("ehigh", 'd', true);
155 fkinECmd->SetParameter(param);
156 G4String smax = DtoS(DBL_MAX);
157 param->SetDefaultValue(smax);
158 param = new G4UIparameter("unit", 's', true);
159 param->SetDefaultUnit("keV");
160 fkinECmd->SetParameter(param);
161 //
162 fparticleCmd = new G4UIcommand("/score/filter/particle", this);
163 fparticleCmd->SetGuidance("Particle filter.");
164 fparticleCmd->SetGuidance("[usage] /score/filter/particle fname p0 .. pn");
165 fparticleCmd->SetGuidance(" fname :(String) Filter Name ");
166 fparticleCmd->SetGuidance(" p0 .. pn :(String) particle names");
167 param = new G4UIparameter("fname", 's', false);
169 param = new G4UIparameter("particlelist", 's', false);
170 param->SetDefaultValue("");
172 //
173 //
174 //
176 new G4UIcommand("/score/filter/particleWithKineticEnergy", this);
177 fparticleKinECmd->SetGuidance("Particle with kinetic energy filter.");
179 "[usage] /score/filter/particleWithKineticEnergy fname Elow Ehigh unit p0 "
180 ".. pn");
181 fparticleKinECmd->SetGuidance(" fname :(String) Filter Name ");
183 " Elow :(Double) Lower edge of kinetic energy");
185 " Ehigh :(Double) Higher edge of kinetic energy");
187 " unit :(String) unit of given kinetic energy");
188 fparticleKinECmd->SetGuidance(" p0 .. pn :(String) particle names");
189 param = new G4UIparameter("fname", 's', false);
191 param = new G4UIparameter("elow", 'd', false);
192 param->SetDefaultValue("0.0");
194 param = new G4UIparameter("ehigh", 'd', true);
195 param->SetDefaultValue(smax);
197 param = new G4UIparameter("unit", 's', true);
198 param->SetDefaultUnit("keV");
200 param = new G4UIparameter("particlelist", 's', false);
201 param->SetDefaultValue("");
203 //
204 //
205}
206
208{
209 delete quantityDir;
210 delete qTouchCmd;
211 delete qGetUnitCmd;
212 delete qSetUnitCmd;
213
214 //
215 delete qCellChgCmd;
216 delete qCellFluxCmd;
217 delete qPassCellFluxCmd;
218 delete qeDepCmd;
219 delete qdoseDepCmd;
220 delete qnOfStepCmd;
221 delete qnOfSecondaryCmd;
222 //
223 delete qTrackLengthCmd;
224 delete qPassCellCurrCmd;
225 delete qPassTrackLengthCmd;
226 delete qFlatSurfCurrCmd;
227 delete qFlatSurfFluxCmd;
228 delete qVolFluxCmd;
229 // delete qSphereSurfCurrCmd;
230 // delete qSphereSurfFluxCmd;
231 // delete qCylSurfCurrCmd;
232 // delete qCylSurfFluxCmd;
233 delete qNofCollisionCmd;
234 delete qPopulationCmd;
235 delete qTrackCountCmd;
236 delete qTerminationCmd;
237 delete qMinKinEAtGeneCmd;
238 //
239 delete qStepCheckerCmd;
240 //
241 delete filterDir;
242 delete fchargedCmd;
243 delete fneutralCmd;
244 delete fkinECmd;
245 delete fparticleCmd;
246 delete fparticleKinECmd;
247}
248
250 G4String newVal)
251{
252 using MeshShape = G4VScoringMesh::MeshShape;
253
255
256 //
257 // Get Current Mesh
258 //
260 if(!mesh)
261 {
262 ed << "ERROR : No mesh is currently open. Open/create a mesh first. "
263 "Command ignored.";
264 command->CommandFailed(ed);
265 return;
266 }
267 // Mesh type
268 auto shape = mesh->GetShape();
269 // Tokens
270 G4TokenVec token;
271 FillTokenVec(newVal, token);
272 //
273 // Commands for Current Mesh
274 if(command == qTouchCmd)
275 {
276 mesh->SetCurrentPrimitiveScorer(newVal);
277 }
278 else if(command == qGetUnitCmd)
279 {
280 G4cout << "Unit: " << mesh->GetCurrentPSUnit() << G4endl;
281 }
282 else if(command == qSetUnitCmd)
283 {
284 mesh->SetCurrentPSUnit(newVal);
285 }
286 else if(command == qCellChgCmd)
287 {
288 if(CheckMeshPS(mesh, token[0], command))
289 {
290 G4PSCellCharge* ps = nullptr;
291 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
292 {
293 ps = new G4PSCellCharge(token[0], mesh->GetCopyNumberLevel());
294 }
295 else
296 {
297 ps = new G4PSCellCharge3D(token[0]);
298 }
299 ps->SetUnit(token[1]);
300 mesh->SetPrimitiveScorer(ps);
301 }
302 }
303 else if(command == qCellFluxCmd)
304 {
305 if(CheckMeshPS(mesh, token[0], command))
306 {
307 G4PSCellFlux* ps = nullptr;
308 if(shape == MeshShape::box)
309 {
310 ps = new G4PSCellFlux3D(token[0]);
311 }
312 else if(shape == MeshShape::cylinder)
313 {
315 new G4PSCellFluxForCylinder3D(token[0]);
316 pps->SetCylinderSize(mesh->GetSize(),mesh->GetStartAngle(),mesh->GetAngleSpan());
317 G4int nSeg[3];
318 mesh->GetNumberOfSegments(nSeg);
319 pps->SetNumberOfSegments(nSeg);
320 ps = pps;
321 }
322 else if(shape == MeshShape::realWorldLogVol)
323 {
324 ed << "Cell flux for real world volume is not yet supported. Command "
325 "ignored.";
326 command->CommandFailed(ed);
327 return;
328 }
329 else if(shape == MeshShape::probe)
330 {
331 ps = new G4PSCellFlux(token[0]);
332 }
333 ps->SetUnit(token[1]);
334 mesh->SetPrimitiveScorer(ps);
335 }
336 }
337 else if(command == qPassCellFluxCmd)
338 {
339 if(CheckMeshPS(mesh, token[0], command))
340 {
341 G4PSPassageCellFlux* ps = nullptr;
342 if(shape == MeshShape::box)
343 {
344 ps = new G4PSPassageCellFlux3D(token[0]);
345 }
346 else if(shape == MeshShape::cylinder)
347 {
350 pps->SetCylinderSize(mesh->GetSize(),mesh->GetStartAngle(),mesh->GetAngleSpan());
351 G4int nSeg[3];
352 mesh->GetNumberOfSegments(nSeg);
353 pps->SetNumberOfSegments(nSeg);
354 ps = pps;
355 }
356 else if(shape == MeshShape::realWorldLogVol)
357 {
358 ed << "Passing cell flux for real world volume is not yet supported. "
359 "Command ignored.";
360 command->CommandFailed(ed);
361 return;
362 }
363 else if(shape == MeshShape::probe)
364 {
365 ps = new G4PSPassageCellFlux(token[0]);
366 }
367 ps->SetUnit(token[1]);
368 mesh->SetPrimitiveScorer(ps);
369 }
370 }
371 else if(command == qeDepCmd)
372 {
373 if(CheckMeshPS(mesh, token[0], command))
374 {
375 G4PSEnergyDeposit* ps = nullptr;
376 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
377 {
378 ps = new G4PSEnergyDeposit(token[0], mesh->GetCopyNumberLevel());
379 }
380 else
381 {
382 ps = new G4PSEnergyDeposit3D(token[0]);
383 }
384 ps->SetUnit(token[1]);
385 mesh->SetPrimitiveScorer(ps);
386 }
387 }
388 else if(command == qdoseDepCmd)
389 {
390 if(CheckMeshPS(mesh, token[0], command))
391 {
392 G4PSDoseDeposit* ps = nullptr;
393 if(shape == MeshShape::box)
394 {
395 ps = new G4PSDoseDeposit3D(token[0]);
396 }
397 else if(shape == MeshShape::cylinder)
398 {
400 new G4PSDoseDepositForCylinder3D(token[0]);
401 pps->SetUnit(token[1]);
402 pps->SetCylinderSize(mesh->GetSize(),mesh->GetStartAngle(),mesh->GetAngleSpan());
403 G4int nSeg[3];
404 mesh->GetNumberOfSegments(nSeg);
405 pps->SetNumberOfSegments(nSeg);
406 ps = pps;
407 }
408 else if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
409 {
410 ps = new G4PSDoseDeposit(token[0], mesh->GetCopyNumberLevel());
411 }
412 ps->SetUnit(token[1]);
413 mesh->SetPrimitiveScorer(ps);
414 }
415 }
416 else if(command == qnOfStepCmd)
417 {
418 if(CheckMeshPS(mesh, token[0], command))
419 {
420 G4PSNofStep* ps = nullptr;
421 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
422 {
423 ps = new G4PSNofStep(token[0], mesh->GetCopyNumberLevel());
424 }
425 else
426 {
427 ps = new G4PSNofStep3D(token[0]);
428 }
429 ps->SetBoundaryFlag(StoB(token[1]));
430 mesh->SetPrimitiveScorer(ps);
431 }
432 }
433 else if(command == qnOfSecondaryCmd)
434 {
435 if(CheckMeshPS(mesh, token[0], command))
436 {
437 G4PSNofSecondary* ps = nullptr;
438 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
439 {
440 ps = new G4PSNofSecondary(token[0], mesh->GetCopyNumberLevel());
441 }
442 else
443 {
444 ps = new G4PSNofSecondary3D(token[0]);
445 }
446 mesh->SetPrimitiveScorer(ps);
447 }
448 }
449 else if(command == qTrackLengthCmd)
450 {
451 if(CheckMeshPS(mesh, token[0], command))
452 {
453 G4PSTrackLength* ps = nullptr;
454 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
455 {
456 ps = new G4PSTrackLength(token[0], mesh->GetCopyNumberLevel());
457 }
458 else
459 {
460 ps = new G4PSTrackLength3D(token[0]);
461 }
462 ps->Weighted(StoB(token[1]));
463 ps->MultiplyKineticEnergy(StoB(token[2]));
464 ps->DivideByVelocity(StoB(token[3]));
465 ps->SetUnit(token[4]);
466 mesh->SetPrimitiveScorer(ps);
467 }
468 }
469 else if(command == qPassCellCurrCmd)
470 {
471 if(CheckMeshPS(mesh, token[0], command))
472 {
473 G4PSPassageCellCurrent* ps = nullptr;
474 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
475 {
476 ps = new G4PSPassageCellCurrent(token[0], mesh->GetCopyNumberLevel());
477 }
478 else
479 {
480 ps = new G4PSPassageCellCurrent3D(token[0]);
481 }
482 ps->Weighted(StoB(token[1]));
483 mesh->SetPrimitiveScorer(ps);
484 }
485 }
486 else if(command == qPassTrackLengthCmd)
487 {
488 if(CheckMeshPS(mesh, token[0], command))
489 {
490 G4PSPassageTrackLength* ps = nullptr;
491 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
492 {
493 ps = new G4PSPassageTrackLength(token[0], mesh->GetCopyNumberLevel());
494 }
495 else
496 {
497 ps = new G4PSPassageTrackLength3D(token[0]);
498 }
499 ps->Weighted(StoB(token[1]));
500 ps->SetUnit(token[2]);
501 mesh->SetPrimitiveScorer(ps);
502 }
503 }
504 else if(command == qFlatSurfCurrCmd)
505 {
506 if(CheckMeshPS(mesh, token[0], command))
507 {
508 G4PSFlatSurfaceCurrent* ps = nullptr;
509 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
510 {
511 ps = new G4PSFlatSurfaceCurrent(token[0], StoI(token[1]),
512 mesh->GetCopyNumberLevel());
513 }
514 else
515 {
516 ps = new G4PSFlatSurfaceCurrent3D(token[0], StoI(token[1]));
517 }
518 ps->Weighted(StoB(token[2]));
519 ps->DivideByArea(StoB(token[3]));
520 if(StoB(token[3]))
521 {
522 ps->SetUnit(token[4]);
523 }
524 else
525 {
526 ps->SetUnit("");
527 }
528 mesh->SetPrimitiveScorer(ps);
529 }
530 }
531 else if(command == qFlatSurfFluxCmd)
532 {
533 if(CheckMeshPS(mesh, token[0], command))
534 {
535 G4PSFlatSurfaceFlux* ps = nullptr;
536 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
537 {
538 ps = new G4PSFlatSurfaceFlux(token[0], StoI(token[1]),
539 mesh->GetCopyNumberLevel());
540 }
541 else
542 {
543 ps = new G4PSFlatSurfaceFlux3D(token[0], StoI(token[1]));
544 }
545 ps->Weighted(StoB(token[2]));
546 ps->DivideByArea(StoB(token[3]));
547 if(StoB(token[3]))
548 {
549 ps->SetUnit(token[4]);
550 }
551 else
552 {
553 ps->SetUnit("");
554 }
555 mesh->SetPrimitiveScorer(ps);
556 }
557 }
558 else if(command == qVolFluxCmd)
559 {
560 if(CheckMeshPS(mesh, token[0], command))
561 {
562 G4PSVolumeFlux* ps = nullptr;
563 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
564 {
565 ps = new G4PSVolumeFlux(token[0], StoI(token[2]),
566 mesh->GetCopyNumberLevel());
567 }
568 else
569 {
570 ps = new G4PSVolumeFlux3D(token[0], StoI(token[2]));
571 }
572 ps->SetDivCos(StoI(token[1]));
573 mesh->SetPrimitiveScorer(ps);
574 }
575
576 // } else if(command== qSphereSurfCurrCmd){
577 // if( CheckMeshPS(mesh, token[0],command )) {
578 // G4PSSphereSurfaceCurrent3D* ps =
579 // new G4PSSphereSurfaceCurrent3D(token[0],StoI(token[1]));
580 // ps->Weighted(StoB(token[2]));
581 // ps->DivideByArea(StoB(token[3]));
582 // if ( StoB(token[3]) ){
583 // ps->SetUnit(token[4]);
584 // }else{
585 // ps->SetUnit("");
586 // }
587 // mesh->SetPrimitiveScorer(ps);
588 // }
589 // } else if(command== qSphereSurfFluxCmd){
590 // if( CheckMeshPS(mesh,token[0],command)) {
591 // G4PSSphereSurfaceFlux3D* ps = new
592 // G4PSSphereSurfaceFlux3D(token[0], StoI(token[1]));
593 // ps->Weighted(StoB(token[2]));
594 // ps->DivideByArea(StoB(token[3]));
595 // if ( StoB(token[3]) ){
596 // ps->SetUnit(token[4]);
597 // }else{
598 // ps->SetUnit("");
599 // }
600 // mesh->SetPrimitiveScorer(ps);
601 // }
602 // } else if(command== qCylSurfCurrCmd){
603 // if( CheckMeshPS(mesh, token[0],command ) ) {
604 // G4PSCylinderSurfaceCurrent3D* ps =
605 // new G4PSCylinderSurfaceCurrent3D(token[0],StoI(token[1]));
606 // ps->Weighted(StoB(token[2]));
607 // ps->DivideByArea(StoB(token[3]));
608 // if ( StoB(token[3]) ){
609 // ps->SetUnit(token[4]);
610 // }else{
611 // ps->SetUnit("");
612 // }
613 // ps->SetUnit(token[4]);
614 // mesh->SetPrimitiveScorer(ps);
615 // }
616 // } else if(command== qCylSurfFluxCmd){
617 // if( CheckMeshPS(mesh, token[0],command ) {
618 // G4PSCylinerSurfaceFlux3D* ps =new
619 // G4PSCylinderSurfaceFlux3D(token[0], StoI(token[1]));
620 // ps->Weighted(StoB(token[2]));
621 // ps->DivideByArea(StoB(token[3]));
622 // if ( StoB(token[3]) ){
623 // ps->SetUnit(token[4]);
624 // }else{
625 // ps->SetUnit("");
626 // }
627 // mesh->SetPrimitiveScorer(ps);
628 // }
629 }
630 else if(command == qNofCollisionCmd)
631 {
632 if(CheckMeshPS(mesh, token[0], command))
633 {
634 G4PSNofCollision* ps = nullptr;
635 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
636 {
637 ps = new G4PSNofCollision3D(token[0], mesh->GetCopyNumberLevel());
638 }
639 else
640 {
641 ps = new G4PSNofCollision3D(token[0]);
642 }
643 ps->Weighted(StoB(token[1]));
644 mesh->SetPrimitiveScorer(ps);
645 }
646 }
647 else if(command == qPopulationCmd)
648 {
649 if(CheckMeshPS(mesh, token[0], command))
650 {
651 G4PSPopulation* ps = nullptr;
652 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
653 {
654 ps = new G4PSPopulation(token[0], mesh->GetCopyNumberLevel());
655 }
656 else
657 {
658 ps = new G4PSPopulation3D(token[0]);
659 }
660 ps->Weighted(StoB(token[1]));
661 mesh->SetPrimitiveScorer(ps);
662 }
663 }
664 else if(command == qTrackCountCmd)
665 {
666 if(CheckMeshPS(mesh, token[0], command))
667 {
668 G4PSTrackCounter* ps = nullptr;
669 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
670 {
671 ps = new G4PSTrackCounter(token[0], StoI(token[1]),
672 mesh->GetCopyNumberLevel());
673 }
674 else
675 {
676 ps = new G4PSTrackCounter3D(token[0], StoI(token[1]));
677 }
678 ps->Weighted(StoB(token[2]));
679 mesh->SetPrimitiveScorer(ps);
680 }
681 }
682 else if(command == qTerminationCmd)
683 {
684 if(CheckMeshPS(mesh, token[0], command))
685 {
686 G4PSTermination* ps = nullptr;
687 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
688 {
689 ps = new G4PSTermination(token[0], mesh->GetCopyNumberLevel());
690 }
691 else
692 {
693 ps = new G4PSTermination3D(token[0]);
694 }
695 ps->Weighted(StoB(token[1]));
696 mesh->SetPrimitiveScorer(ps);
697 }
698 }
699 else if(command == qMinKinEAtGeneCmd)
700 {
701 if(CheckMeshPS(mesh, token[0], command))
702 {
703 G4PSMinKinEAtGeneration* ps = nullptr;
704 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
705 {
706 ps = new G4PSMinKinEAtGeneration(token[0], mesh->GetCopyNumberLevel());
707 }
708 else
709 {
710 ps = new G4PSMinKinEAtGeneration3D(token[0]);
711 }
712 ps->SetUnit(token[1]);
713 mesh->SetPrimitiveScorer(ps);
714 }
715 }
716 else if(command == qStepCheckerCmd)
717 {
718 if(CheckMeshPS(mesh, token[0], command))
719 {
720 G4PSStepChecker* ps = nullptr;
721 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
722 {
723 ps = new G4PSStepChecker(token[0], mesh->GetCopyNumberLevel());
724 }
725 else
726 {
727 ps = new G4PSStepChecker3D(token[0]);
728 }
729 mesh->SetPrimitiveScorer(ps);
730 }
731
732 //
733 // Filters
734 //
735 }
736 else if(command == fchargedCmd)
737 {
739 {
740 mesh->SetFilter(new G4SDChargedFilter(token[0]));
741 }
742 else
743 {
744 ed << "WARNING[" << fchargedCmd->GetCommandPath()
745 << "] : Current quantity is not set. Set or touch a quantity first.";
746 command->CommandFailed(ed);
747 }
748 }
749 else if(command == fneutralCmd)
750 {
752 {
753 mesh->SetFilter(new G4SDNeutralFilter(token[0]));
754 }
755 else
756 {
757 ed << "WARNING[" << fneutralCmd->GetCommandPath()
758 << "] : Current quantity is not set. Set or touch a quantity first.";
759 command->CommandFailed(ed);
760 }
761 }
762 else if(command == fkinECmd)
763 {
765 {
766 G4String& name = token[0];
767 G4double elow = StoD(token[1]);
768 G4double ehigh = StoD(token[2]);
769 G4double unitVal = G4UnitDefinition::GetValueOf(token[3]);
770 mesh->SetFilter(
771 new G4SDKineticEnergyFilter(name, elow * unitVal, ehigh * unitVal));
772 }
773 else
774 {
775 ed << "WARNING[" << fkinECmd->GetCommandPath()
776 << "] : Current quantity is not set. Set or touch a quantity first.";
777 command->CommandFailed(ed);
778 }
779 }
780 else if(command == fparticleKinECmd)
781 {
783 {
784 FParticleWithEnergyCommand(mesh, token);
785 }
786 else
787 {
788 ed << "WARNING[" << fparticleKinECmd->GetCommandPath()
789 << "] : Current quantity is not set. Set or touch a quantity first.";
790 command->CommandFailed(ed);
791 }
792 }
793 else if(command == fparticleCmd)
794 {
796 {
797 FParticleCommand(mesh, token);
798 }
799 else
800 {
801 ed << "WARNING[" << fparticleCmd->GetCommandPath()
802 << "] : Current quantity is not set. Set or touch a quantity first.";
803 command->CommandFailed(ed);
804 }
805 }
806}
807
809{
810 G4String val;
811
812 return val;
813}
814
816 G4TokenVec& token)
817{
818 G4Tokenizer next(newValues);
819 G4String val;
820 while(!(val = next()).empty())
821 { // Loop checking 12.18.2015 M.Asai
822 token.push_back(val);
823 }
824}
825
827 G4TokenVec& token)
828{
829 //
830 // Filter name
831 G4String name = token[0];
832 //
833 // particle list
834 std::vector<G4String> pnames;
835 for(G4int i = 1; i < (G4int) token.size(); i++)
836 {
837 pnames.push_back(token[i]);
838 }
839 //
840 // Attach Filter
841 mesh->SetFilter(new G4SDParticleFilter(name, pnames));
842}
843
845 G4TokenVec& token)
846{
847 G4String& name = token[0];
848 G4double elow = StoD(token[1]);
849 G4double ehigh = StoD(token[2]);
850 G4double unitVal = G4UnitDefinition::GetValueOf(token[3]);
852 new G4SDParticleWithEnergyFilter(name, elow * unitVal, ehigh * unitVal);
853 for(G4int i = 4; i < (G4int) token.size(); i++)
854 {
855 filter->add(token[i]);
856 }
857 mesh->SetFilter(filter);
858}
859
861 G4String& psname,
862 G4UIcommand* command)
863{
864 if(!mesh->FindPrimitiveScorer(psname))
865 {
866 return true;
867 }
868 else
869 {
871 ed << "WARNING[" << qTouchCmd->GetCommandPath() << "] : Quantity name, \""
872 << psname << "\", is already existing.";
873 command->CommandFailed(ed);
875 return false;
876 }
877}
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
static constexpr double ps
Definition: G4SIunits.hh:157
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
void SetCylinderSize(G4ThreeVector cylSize, G4double startAng, G4double angSpan)
void SetCylinderSize(G4ThreeVector cylSize, G4double startAng, G4double angSpan)
virtual void SetUnit(const G4String &unit)
void SetCylinderSize(G4ThreeVector cylSize, G4double startAng, G4double angSpan)
void add(const G4String &particleName)
void FParticleWithEnergyCommand(G4VScoringMesh *mesh, G4TokenVec &token)
G4bool CheckMeshPS(G4VScoringMesh *mesh, G4String &psname, G4UIcommand *command)
G4String GetCurrentValue(G4UIcommand *)
void SetNewValue(G4UIcommand *command, G4String newValues)
G4ScoreQuantityMessenger(G4ScoringManager *SManager)
void FillTokenVec(G4String newValues, G4TokenVec &token)
void FParticleCommand(G4VScoringMesh *mesh, G4TokenVec &token)
G4UIcmdWithoutParameter * qGetUnitCmd
G4VScoringMesh * GetCurrentMesh() const
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
const G4String & GetCommandPath() const
Definition: G4UIcommand.hh:136
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:146
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:156
void CommandFailed(G4int errCode, G4ExceptionDescription &ed)
Definition: G4UIcommand.hh:179
G4int StoI(G4String s)
G4String DtoS(G4double a)
G4double StoD(G4String s)
G4bool StoB(G4String s)
void SetDefaultValue(const char *theDefaultValue)
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)
const char * name(G4int ptype)
#define DBL_MAX
Definition: templates.hh:62