Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ScoreQuantityMessenger Class Reference

#include <G4ScoreQuantityMessenger.hh>

+ Inheritance diagram for G4ScoreQuantityMessenger:

Public Member Functions

 G4ScoreQuantityMessenger (G4ScoringManager *SManager)
 
 ~G4ScoreQuantityMessenger () override
 
void SetNewValue (G4UIcommand *command, G4String newValues) override
 
G4String GetCurrentValue (G4UIcommand *) override
 
- Public Member Functions inherited from G4UImessenger
 G4UImessenger ()=default
 
 G4UImessenger (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
virtual ~G4UImessenger ()
 
G4bool CommandsShouldBeInMaster () const
 

Protected Member Functions

void FillTokenVec (const G4String &newValues, G4TokenVec &token)
 
void FParticleCommand (G4VScoringMesh *mesh, G4TokenVec &token)
 
void FParticleWithEnergyCommand (G4VScoringMesh *mesh, G4TokenVec &token)
 
G4bool CheckMeshPS (G4VScoringMesh *mesh, const G4String &psname, G4UIcommand *command)
 
- Protected Member Functions inherited from G4UImessenger
G4String ItoS (G4int i)
 
G4String LtoS (G4long l)
 
G4String DtoS (G4double a)
 
G4String BtoS (G4bool b)
 
G4int StoI (const G4String &s)
 
G4long StoL (const G4String &s)
 
G4double StoD (const G4String &s)
 
G4bool StoB (const G4String &s)
 
void AddUIcommand (G4UIcommand *newCommand)
 
void CreateDirectory (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
template<typename T>
T * CreateCommand (const G4String &cname, const G4String &dsc)
 

Additional Inherited Members

- Protected Attributes inherited from G4UImessenger
G4UIdirectorybaseDir = nullptr
 
G4String baseDirName = ""
 
G4bool commandsShouldBeInMaster = false
 

Detailed Description

Definition at line 52 of file G4ScoreQuantityMessenger.hh.

Constructor & Destructor Documentation

◆ G4ScoreQuantityMessenger()

G4ScoreQuantityMessenger::G4ScoreQuantityMessenger ( G4ScoringManager * SManager)

Definition at line 115 of file G4ScoreQuantityMessenger.cc.

116 : fSMan(SManager)
117{
118 QuantityCommands();
119 FilterCommands();
120}

◆ ~G4ScoreQuantityMessenger()

G4ScoreQuantityMessenger::~G4ScoreQuantityMessenger ( )
override

Definition at line 512 of file G4ScoreQuantityMessenger.cc.

513{
514 delete quantityDir;
515 delete qTouchCmd;
516 delete qGetUnitCmd;
517 delete qSetUnitCmd;
518
519 //
520 delete qCellChgCmd;
521 delete qCellFluxCmd;
522 delete qPassCellFluxCmd;
523 delete qeDepCmd;
524 delete qdoseDepCmd;
525 delete qnOfStepCmd;
526 delete qnOfSecondaryCmd;
527 //
528 delete qTrackLengthCmd;
529 delete qPassCellCurrCmd;
530 delete qPassTrackLengthCmd;
531 delete qFlatSurfCurrCmd;
532 delete qFlatSurfFluxCmd;
533 delete qVolFluxCmd;
534
535 delete qNofCollisionCmd;
536 delete qPopulationCmd;
537 delete qTrackCountCmd;
538 delete qTerminationCmd;
539 delete qMinKinEAtGeneCmd;
540 //
541 delete qStepCheckerCmd;
542 //
543 delete filterDir;
544 delete fchargedCmd;
545 delete fneutralCmd;
546 delete fkinECmd;
547 delete fparticleCmd;
548 delete fparticleKinECmd;
549}

Member Function Documentation

◆ CheckMeshPS()

G4bool G4ScoreQuantityMessenger::CheckMeshPS ( G4VScoringMesh * mesh,
const G4String & psname,
G4UIcommand * command )
protected

Definition at line 1108 of file G4ScoreQuantityMessenger.cc.

1111{
1112 if(!mesh->FindPrimitiveScorer(psname))
1113 {
1114 return true;
1115 }
1116
1118 ed << "WARNING[" << qTouchCmd->GetCommandPath() << "] : Quantity name, \""
1119 << psname << "\", is already existing.";
1120 command->CommandFailed(ed);
1122 return false;
1123}
std::ostringstream G4ExceptionDescription
void CommandFailed(G4int errCode, G4ExceptionDescription &ed)
void SetNullToCurrentPrimitiveScorer()
G4bool FindPrimitiveScorer(const G4String &psname)

Referenced by SetNewValue().

◆ FillTokenVec()

void G4ScoreQuantityMessenger::FillTokenVec ( const G4String & newValues,
G4TokenVec & token )
protected

Definition at line 1063 of file G4ScoreQuantityMessenger.cc.

1065{
1066 G4Tokenizer next(newValues);
1067 G4String val;
1068 while(!(val = next()).empty())
1069 { // Loop checking 12.18.2015 M.Asai
1070 token.push_back(val);
1071 }
1072}

Referenced by SetNewValue().

◆ FParticleCommand()

void G4ScoreQuantityMessenger::FParticleCommand ( G4VScoringMesh * mesh,
G4TokenVec & token )
protected

Definition at line 1074 of file G4ScoreQuantityMessenger.cc.

1076{
1077 //
1078 // Filter name
1079 G4String name = token[0];
1080 //
1081 // particle list
1082 std::vector<G4String> pnames;
1083 for(G4int i = 1; i < (G4int) token.size(); i++)
1084 {
1085 pnames.push_back(token[i]);
1086 }
1087 //
1088 // Attach Filter
1089 mesh->SetFilter(new G4SDParticleFilter(name, pnames));
1090}
int G4int
Definition G4Types.hh:85
void SetFilter(G4VSDFilter *filter)
const char * name(G4int ptype)

Referenced by SetNewValue().

◆ FParticleWithEnergyCommand()

void G4ScoreQuantityMessenger::FParticleWithEnergyCommand ( G4VScoringMesh * mesh,
G4TokenVec & token )
protected

Definition at line 1092 of file G4ScoreQuantityMessenger.cc.

1094{
1095 G4String& name = token[0];
1096 G4double elow = StoD(token[1]);
1097 G4double ehigh = StoD(token[2]);
1098 G4double unitVal = G4UnitDefinition::GetValueOf(token[3]);
1099 auto filter =
1100 new G4SDParticleWithEnergyFilter(name, elow * unitVal, ehigh * unitVal);
1101 for(G4int i = 4; i < (G4int) token.size(); i++)
1102 {
1103 filter->add(token[i]);
1104 }
1105 mesh->SetFilter(filter);
1106}
double G4double
Definition G4Types.hh:83
G4double StoD(const G4String &s)
static G4double GetValueOf(const G4String &)

Referenced by SetNewValue().

◆ GetCurrentValue()

G4String G4ScoreQuantityMessenger::GetCurrentValue ( G4UIcommand * )
overridevirtual

Reimplemented from G4UImessenger.

Definition at line 1056 of file G4ScoreQuantityMessenger.cc.

1057{
1058 G4String val;
1059
1060 return val;
1061}

◆ SetNewValue()

void G4ScoreQuantityMessenger::SetNewValue ( G4UIcommand * command,
G4String newValues )
overridevirtual

Reimplemented from G4UImessenger.

Definition at line 551 of file G4ScoreQuantityMessenger.cc.

553{
554 using MeshShape = G4VScoringMesh::MeshShape;
555
557
558 //
559 // Get Current Mesh
560 //
561 G4VScoringMesh* mesh = fSMan->GetCurrentMesh();
562 if(mesh == nullptr)
563 {
564 ed << "ERROR : No mesh is currently open. Open/create a mesh first. "
565 "Command ignored.";
566 command->CommandFailed(ed);
567 return;
568 }
569 // Mesh type
570 auto shape = mesh->GetShape();
571 // Tokens
572 G4TokenVec token;
573 FillTokenVec(newVal, token);
574 //
575 // Commands for Current Mesh
576 if(command == qTouchCmd)
577 {
578 mesh->SetCurrentPrimitiveScorer(newVal);
579 }
580 else if(command == qGetUnitCmd)
581 {
582 G4cout << "Unit: " << mesh->GetCurrentPSUnit() << G4endl;
583 }
584 else if(command == qSetUnitCmd)
585 {
586 mesh->SetCurrentPSUnit(newVal);
587 }
588 else if(command == qCellChgCmd)
589 {
590 if(CheckMeshPS(mesh, token[0], command))
591 {
592 G4PSCellCharge* ps = nullptr;
593 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
594 {
595 ps = new G4PSCellCharge(token[0], mesh->GetCopyNumberLevel());
596 }
597 else
598 {
599 ps = new G4PSCellCharge3D(token[0]);
600 }
601 ps->SetUnit(token[1]);
602 mesh->SetPrimitiveScorer(ps);
603 }
604 }
605 else if(command == qCellFluxCmd)
606 {
607 if(CheckMeshPS(mesh, token[0], command))
608 {
609 G4PSCellFlux* ps = nullptr;
610 if(shape == MeshShape::box)
611 {
612 ps = new G4PSCellFlux3D(token[0]);
613 }
614 else if(shape == MeshShape::cylinder)
615 {
616 auto pps =
617 new G4PSCellFluxForCylinder3D(token[0]);
618 pps->SetCylinderSize(mesh->GetSize(),mesh->GetStartAngle(),mesh->GetAngleSpan());
619 G4int nSeg[3];
620 mesh->GetNumberOfSegments(nSeg);
621 pps->SetNumberOfSegments(nSeg);
622 ps = pps;
623 }
624 else if(shape == MeshShape::realWorldLogVol)
625 {
626 ed << "Cell flux for real world volume is not yet supported. Command "
627 "ignored.";
628 command->CommandFailed(ed);
629 return;
630 }
631 else if(shape == MeshShape::probe)
632 {
633 ps = new G4PSCellFlux(token[0]);
634 }
635 ps->SetUnit(token[1]);
636 mesh->SetPrimitiveScorer(ps);
637 }
638 }
639 else if(command == qPassCellFluxCmd)
640 {
641 if(CheckMeshPS(mesh, token[0], command))
642 {
643 G4PSPassageCellFlux* ps = nullptr;
644 if(shape == MeshShape::box)
645 {
646 ps = new G4PSPassageCellFlux3D(token[0]);
647 }
648 else if(shape == MeshShape::cylinder)
649 {
650 auto pps =
651 new G4PSPassageCellFluxForCylinder3D(token[0]);
652 pps->SetCylinderSize(mesh->GetSize(),mesh->GetStartAngle(),mesh->GetAngleSpan());
653 G4int nSeg[3];
654 mesh->GetNumberOfSegments(nSeg);
655 pps->SetNumberOfSegments(nSeg);
656 ps = pps;
657 }
658 else if(shape == MeshShape::realWorldLogVol)
659 {
660 ed << "Passing cell flux for real world volume is not yet supported. "
661 "Command ignored.";
662 command->CommandFailed(ed);
663 return;
664 }
665 else if(shape == MeshShape::probe)
666 {
667 ps = new G4PSPassageCellFlux(token[0]);
668 }
669 ps->SetUnit(token[1]);
670 mesh->SetPrimitiveScorer(ps);
671 }
672 }
673 else if(command == qeDepCmd)
674 {
675 if(CheckMeshPS(mesh, token[0], command))
676 {
677 G4PSEnergyDeposit* ps = nullptr;
678 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
679 {
680 ps = new G4PSEnergyDeposit(token[0], mesh->GetCopyNumberLevel());
681 }
682 else
683 {
684 ps = new G4PSEnergyDeposit3D(token[0]);
685 }
686 ps->SetUnit(token[1]);
687 mesh->SetPrimitiveScorer(ps);
688 }
689 }
690 else if(command == qdoseDepCmd)
691 {
692 if(CheckMeshPS(mesh, token[0], command))
693 {
694 G4PSDoseDeposit* ps = nullptr;
695 if(shape == MeshShape::box)
696 {
697 ps = new G4PSDoseDeposit3D(token[0]);
698 }
699 else if(shape == MeshShape::cylinder)
700 {
701 auto pps =
702 new G4PSDoseDepositForCylinder3D(token[0]);
703 pps->SetUnit(token[1]);
704 pps->SetCylinderSize(mesh->GetSize(),mesh->GetStartAngle(),mesh->GetAngleSpan());
705 G4int nSeg[3];
706 mesh->GetNumberOfSegments(nSeg);
707 pps->SetNumberOfSegments(nSeg);
708 ps = pps;
709 }
710 else if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
711 {
712 ps = new G4PSDoseDeposit(token[0], mesh->GetCopyNumberLevel());
713 }
714 ps->SetUnit(token[1]);
715 mesh->SetPrimitiveScorer(ps);
716 }
717 }
718 else if(command == qnOfStepCmd)
719 {
720 if(CheckMeshPS(mesh, token[0], command))
721 {
722 G4PSNofStep* ps = nullptr;
723 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
724 {
725 ps = new G4PSNofStep(token[0], mesh->GetCopyNumberLevel());
726 }
727 else
728 {
729 ps = new G4PSNofStep3D(token[0]);
730 }
731 ps->SetBoundaryFlag(StoB(token[1]));
732 mesh->SetPrimitiveScorer(ps);
733 }
734 }
735 else if(command == qnOfSecondaryCmd)
736 {
737 if(CheckMeshPS(mesh, token[0], command))
738 {
739 G4PSNofSecondary* ps = nullptr;
740 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
741 {
742 ps = new G4PSNofSecondary(token[0], mesh->GetCopyNumberLevel());
743 }
744 else
745 {
746 ps = new G4PSNofSecondary3D(token[0]);
747 }
748 mesh->SetPrimitiveScorer(ps);
749 }
750 }
751 else if(command == qTrackLengthCmd)
752 {
753 if(CheckMeshPS(mesh, token[0], command))
754 {
755 G4PSTrackLength* ps = nullptr;
756 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
757 {
758 ps = new G4PSTrackLength(token[0], mesh->GetCopyNumberLevel());
759 }
760 else
761 {
762 ps = new G4PSTrackLength3D(token[0]);
763 }
764 ps->Weighted(StoB(token[1]));
765 ps->MultiplyKineticEnergy(StoB(token[2]));
766 ps->DivideByVelocity(StoB(token[3]));
767 ps->SetUnit(token[4]);
768 mesh->SetPrimitiveScorer(ps);
769 }
770 }
771 else if(command == qPassCellCurrCmd)
772 {
773 if(CheckMeshPS(mesh, token[0], command))
774 {
775 G4PSPassageCellCurrent* ps = nullptr;
776 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
777 {
778 ps = new G4PSPassageCellCurrent(token[0], mesh->GetCopyNumberLevel());
779 }
780 else
781 {
782 ps = new G4PSPassageCellCurrent3D(token[0]);
783 }
784 ps->Weighted(StoB(token[1]));
785 mesh->SetPrimitiveScorer(ps);
786 }
787 }
788 else if(command == qPassTrackLengthCmd)
789 {
790 if(CheckMeshPS(mesh, token[0], command))
791 {
792 G4PSPassageTrackLength* ps = nullptr;
793 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
794 {
795 ps = new G4PSPassageTrackLength(token[0], mesh->GetCopyNumberLevel());
796 }
797 else
798 {
799 ps = new G4PSPassageTrackLength3D(token[0]);
800 }
801 ps->Weighted(StoB(token[1]));
802 ps->SetUnit(token[2]);
803 mesh->SetPrimitiveScorer(ps);
804 }
805 }
806 else if(command == qFlatSurfCurrCmd)
807 {
808 if(CheckMeshPS(mesh, token[0], command))
809 {
810 G4PSFlatSurfaceCurrent* ps = nullptr;
811 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
812 {
813 ps = new G4PSFlatSurfaceCurrent(token[0], StoI(token[1]),
814 mesh->GetCopyNumberLevel());
815 }
816 else
817 {
818 ps = new G4PSFlatSurfaceCurrent3D(token[0], StoI(token[1]));
819 }
820 ps->Weighted(StoB(token[2]));
821 ps->DivideByArea(StoB(token[3]));
822 if(StoB(token[3]))
823 {
824 ps->SetUnit(token[4]);
825 }
826 else
827 {
828 ps->SetUnit("");
829 }
830 mesh->SetPrimitiveScorer(ps);
831 }
832 }
833 else if(command == qFlatSurfFluxCmd)
834 {
835 if(CheckMeshPS(mesh, token[0], command))
836 {
837 G4PSFlatSurfaceFlux* ps = nullptr;
838 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
839 {
840 ps = new G4PSFlatSurfaceFlux(token[0], StoI(token[1]),
841 mesh->GetCopyNumberLevel());
842 }
843 else
844 {
845 ps = new G4PSFlatSurfaceFlux3D(token[0], StoI(token[1]));
846 }
847 ps->Weighted(StoB(token[2]));
848 ps->DivideByArea(StoB(token[3]));
849 if(StoB(token[3]))
850 {
851 ps->SetUnit(token[4]);
852 }
853 else
854 {
855 ps->SetUnit("");
856 }
857 mesh->SetPrimitiveScorer(ps);
858 }
859 }
860 else if(command == qVolFluxCmd)
861 {
862 if(CheckMeshPS(mesh, token[0], command))
863 {
864 G4PSVolumeFlux* ps = nullptr;
865 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
866 {
867 ps = new G4PSVolumeFlux(token[0], StoI(token[2]),
868 mesh->GetCopyNumberLevel());
869 }
870 else
871 {
872 ps = new G4PSVolumeFlux3D(token[0], StoI(token[2]));
873 }
874 ps->SetDivCos(StoI(token[1]) != 0);
875 mesh->SetPrimitiveScorer(ps);
876 }
877 }
878 else if(command == qNofCollisionCmd)
879 {
880 if(CheckMeshPS(mesh, token[0], command))
881 {
882 G4PSNofCollision* ps = nullptr;
883 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
884 {
885 ps = new G4PSNofCollision3D(token[0], mesh->GetCopyNumberLevel());
886 }
887 else
888 {
889 ps = new G4PSNofCollision3D(token[0]);
890 }
891 ps->Weighted(StoB(token[1]));
892 mesh->SetPrimitiveScorer(ps);
893 }
894 }
895 else if(command == qPopulationCmd)
896 {
897 if(CheckMeshPS(mesh, token[0], command))
898 {
899 G4PSPopulation* ps = nullptr;
900 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
901 {
902 ps = new G4PSPopulation(token[0], mesh->GetCopyNumberLevel());
903 }
904 else
905 {
906 ps = new G4PSPopulation3D(token[0]);
907 }
908 ps->Weighted(StoB(token[1]));
909 mesh->SetPrimitiveScorer(ps);
910 }
911 }
912 else if(command == qTrackCountCmd)
913 {
914 if(CheckMeshPS(mesh, token[0], command))
915 {
916 G4PSTrackCounter* ps = nullptr;
917 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
918 {
919 ps = new G4PSTrackCounter(token[0], StoI(token[1]),
920 mesh->GetCopyNumberLevel());
921 }
922 else
923 {
924 ps = new G4PSTrackCounter3D(token[0], StoI(token[1]));
925 }
926 ps->Weighted(StoB(token[2]));
927 mesh->SetPrimitiveScorer(ps);
928 }
929 }
930 else if(command == qTerminationCmd)
931 {
932 if(CheckMeshPS(mesh, token[0], command))
933 {
934 G4PSTermination* ps = nullptr;
935 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
936 {
937 ps = new G4PSTermination(token[0], mesh->GetCopyNumberLevel());
938 }
939 else
940 {
941 ps = new G4PSTermination3D(token[0]);
942 }
943 ps->Weighted(StoB(token[1]));
944 mesh->SetPrimitiveScorer(ps);
945 }
946 }
947 else if(command == qMinKinEAtGeneCmd)
948 {
949 if(CheckMeshPS(mesh, token[0], command))
950 {
951 G4PSMinKinEAtGeneration* ps = nullptr;
952 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
953 {
954 ps = new G4PSMinKinEAtGeneration(token[0], mesh->GetCopyNumberLevel());
955 }
956 else
957 {
958 ps = new G4PSMinKinEAtGeneration3D(token[0]);
959 }
960 ps->SetUnit(token[1]);
961 mesh->SetPrimitiveScorer(ps);
962 }
963 }
964 else if(command == qStepCheckerCmd)
965 {
966 if(CheckMeshPS(mesh, token[0], command))
967 {
968 G4PSStepChecker* ps = nullptr;
969 if(shape == MeshShape::realWorldLogVol || shape == MeshShape::probe)
970 {
971 ps = new G4PSStepChecker(token[0], mesh->GetCopyNumberLevel());
972 }
973 else
974 {
975 ps = new G4PSStepChecker3D(token[0]);
976 }
977 mesh->SetPrimitiveScorer(ps);
978 }
979
980 //
981 // Filters
982 //
983 }
984 else if(command == fchargedCmd)
985 {
987 {
988 mesh->SetFilter(new G4SDChargedFilter(token[0]));
989 }
990 else
991 {
992 ed << "WARNING[" << fchargedCmd->GetCommandPath()
993 << "] : Current quantity is not set. Set or touch a quantity first.";
994 command->CommandFailed(ed);
995 }
996 }
997 else if(command == fneutralCmd)
998 {
1000 {
1001 mesh->SetFilter(new G4SDNeutralFilter(token[0]));
1002 }
1003 else
1004 {
1005 ed << "WARNING[" << fneutralCmd->GetCommandPath()
1006 << "] : Current quantity is not set. Set or touch a quantity first.";
1007 command->CommandFailed(ed);
1008 }
1009 }
1010 else if(command == fkinECmd)
1011 {
1012 if(!mesh->IsCurrentPrimitiveScorerNull())
1013 {
1014 G4String& name = token[0];
1015 G4double elow = StoD(token[1]);
1016 G4double ehigh = StoD(token[2]);
1017 G4double unitVal = G4UnitDefinition::GetValueOf(token[3]);
1018 mesh->SetFilter(
1019 new G4SDKineticEnergyFilter(name, elow * unitVal, ehigh * unitVal));
1020 }
1021 else
1022 {
1023 ed << "WARNING[" << fkinECmd->GetCommandPath()
1024 << "] : Current quantity is not set. Set or touch a quantity first.";
1025 command->CommandFailed(ed);
1026 }
1027 }
1028 else if(command == fparticleKinECmd)
1029 {
1030 if(!mesh->IsCurrentPrimitiveScorerNull())
1031 {
1032 FParticleWithEnergyCommand(mesh, token);
1033 }
1034 else
1035 {
1036 ed << "WARNING[" << fparticleKinECmd->GetCommandPath()
1037 << "] : Current quantity is not set. Set or touch a quantity first.";
1038 command->CommandFailed(ed);
1039 }
1040 }
1041 else if(command == fparticleCmd)
1042 {
1043 if(!mesh->IsCurrentPrimitiveScorerNull())
1044 {
1045 FParticleCommand(mesh, token);
1046 }
1047 else
1048 {
1049 ed << "WARNING[" << fparticleCmd->GetCommandPath()
1050 << "] : Current quantity is not set. Set or touch a quantity first.";
1051 command->CommandFailed(ed);
1052 }
1053 }
1054}
std::vector< G4String > G4TokenVec
#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 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)
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, const G4String &psname, G4UIcommand *command)
void FillTokenVec(const G4String &newValues, G4TokenVec &token)
void FParticleCommand(G4VScoringMesh *mesh, G4TokenVec &token)
G4bool StoB(const G4String &s)
G4int StoI(const G4String &s)
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()

The documentation for this class was generated from the following files: