Geant4 11.2.2
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 (G4String newValues, G4TokenVec &token)
 
void FParticleCommand (G4VScoringMesh *mesh, G4TokenVec &token)
 
void FParticleWithEnergyCommand (G4VScoringMesh *mesh, G4TokenVec &token)
 
G4bool CheckMeshPS (G4VScoringMesh *mesh, 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 (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 57 of file G4ScoreQuantityMessenger.hh.

Constructor & Destructor Documentation

◆ G4ScoreQuantityMessenger()

G4ScoreQuantityMessenger::G4ScoreQuantityMessenger ( G4ScoringManager * SManager)

Definition at line 117 of file G4ScoreQuantityMessenger.cc.

118 : fSMan(SManager)
119{
120 QuantityCommands();
121 FilterCommands();
122}

◆ ~G4ScoreQuantityMessenger()

G4ScoreQuantityMessenger::~G4ScoreQuantityMessenger ( )
override

Definition at line 514 of file G4ScoreQuantityMessenger.cc.

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}

Member Function Documentation

◆ CheckMeshPS()

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

Definition at line 1110 of file G4ScoreQuantityMessenger.cc.

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
const G4String & GetCommandPath() const
void CommandFailed(G4int errCode, G4ExceptionDescription &ed)
void SetNullToCurrentPrimitiveScorer()
G4bool FindPrimitiveScorer(const G4String &psname)

Referenced by SetNewValue().

◆ FillTokenVec()

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

Definition at line 1065 of file G4ScoreQuantityMessenger.cc.

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}

Referenced by SetNewValue().

◆ FParticleCommand()

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

Definition at line 1076 of file G4ScoreQuantityMessenger.cc.

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}
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 1094 of file G4ScoreQuantityMessenger.cc.

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}
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 1058 of file G4ScoreQuantityMessenger.cc.

1059{
1060 G4String val;
1061
1062 return val;
1063}

◆ SetNewValue()

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

Reimplemented from G4UImessenger.

Definition at line 553 of file G4ScoreQuantityMessenger.cc.

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}
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 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, G4String &psname, G4UIcommand *command)
void FillTokenVec(G4String newValues, G4TokenVec &token)
void FParticleCommand(G4VScoringMesh *mesh, G4TokenVec &token)
G4VScoringMesh * GetCurrentMesh() const
G4bool StoB(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: