Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes
G4ScoreQuantityMessenger Class Reference

#include <G4ScoreQuantityMessenger.hh>

Inheritance diagram for G4ScoreQuantityMessenger:
G4UImessenger

Public Member Functions

G4bool CommandsShouldBeInMaster () const
 
 G4ScoreQuantityMessenger (G4ScoringManager *SManager)
 
G4String GetCurrentValue (G4UIcommand *)
 
G4bool operator!= (const G4UImessenger &messenger) const
 
G4bool operator== (const G4UImessenger &messenger) const
 
void SetNewValue (G4UIcommand *command, G4String newValues)
 
 ~G4ScoreQuantityMessenger ()
 

Protected Member Functions

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

Protected Attributes

G4UIdirectorybaseDir = nullptr
 
G4String baseDirName = ""
 
G4bool commandsShouldBeInMaster = false
 

Private Member Functions

void FilterCommands ()
 
void QuantityCommands ()
 

Private Attributes

G4UIcmdWithAStringfchargedCmd
 
G4UIdirectoryfilterDir
 
G4UIcommandfkinECmd
 
G4UIcmdWithAStringfneutralCmd
 
G4UIcommandfparticleCmd
 
G4UIcommandfparticleKinECmd
 
G4ScoringManagerfSMan
 
G4UIcommandqCellChgCmd
 
G4UIcommandqCellFluxCmd
 
G4UIcommandqdoseDepCmd
 
G4UIcommandqeDepCmd
 
G4UIcommandqFlatSurfCurrCmd
 
G4UIcommandqFlatSurfFluxCmd
 
G4UIcmdWithoutParameterqGetUnitCmd
 
G4UIcommandqMinKinEAtGeneCmd
 
G4UIcommandqNofCollisionCmd
 
G4UIcommandqnOfSecondaryCmd
 
G4UIcommandqnOfStepCmd
 
G4UIcommandqPassCellCurrCmd
 
G4UIcommandqPassCellFluxCmd
 
G4UIcommandqPassTrackLengthCmd
 
G4UIcommandqPopulationCmd
 
G4UIcmdWithAStringqSetUnitCmd
 
G4UIcommandqStepCheckerCmd
 
G4UIcommandqTerminationCmd
 
G4UIcmdWithAStringqTouchCmd
 
G4UIcommandqTrackCountCmd
 
G4UIcommandqTrackLengthCmd
 
G4UIdirectoryquantityDir
 
G4UIcommandqVolFluxCmd
 

Detailed Description

Definition at line 57 of file G4ScoreQuantityMessenger.hh.

Constructor & Destructor Documentation

◆ G4ScoreQuantityMessenger()

G4ScoreQuantityMessenger::G4ScoreQuantityMessenger ( G4ScoringManager SManager)

◆ ~G4ScoreQuantityMessenger()

G4ScoreQuantityMessenger::~G4ScoreQuantityMessenger ( )

Definition at line 207 of file G4ScoreQuantityMessenger.cc.

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}
G4UIcmdWithoutParameter * qGetUnitCmd

References fchargedCmd, filterDir, fkinECmd, fneutralCmd, fparticleCmd, fparticleKinECmd, qCellChgCmd, qCellFluxCmd, qdoseDepCmd, qeDepCmd, qFlatSurfCurrCmd, qFlatSurfFluxCmd, qGetUnitCmd, qMinKinEAtGeneCmd, qNofCollisionCmd, qnOfSecondaryCmd, qnOfStepCmd, qPassCellCurrCmd, qPassCellFluxCmd, qPassTrackLengthCmd, qPopulationCmd, qSetUnitCmd, qStepCheckerCmd, qTerminationCmd, qTouchCmd, qTrackCountCmd, qTrackLengthCmd, quantityDir, and qVolFluxCmd.

Member Function Documentation

◆ AddUIcommand()

void G4UImessenger::AddUIcommand ( G4UIcommand newCommand)
protectedinherited

Definition at line 149 of file G4UImessenger.cc.

150{
151 G4cerr << "Warning : Old style definition of G4UIcommand <"
152 << newCommand->GetCommandPath() << ">." << G4endl;
153}
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
const G4String & GetCommandPath() const
Definition: G4UIcommand.hh:136

References G4cerr, G4endl, and G4UIcommand::GetCommandPath().

◆ BtoS()

G4String G4UImessenger::BtoS ( G4bool  b)
protectedinherited

Definition at line 98 of file G4UImessenger.cc.

99{
100 G4String vl = "0";
101 if(b)
102 vl = "true";
103 return vl;
104}

◆ CheckMeshPS()

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

Definition at line 860 of file G4ScoreQuantityMessenger.cc.

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
void CommandFailed(G4int errCode, G4ExceptionDescription &ed)
Definition: G4UIcommand.hh:179
void SetNullToCurrentPrimitiveScorer()
G4bool FindPrimitiveScorer(const G4String &psname)

References G4UIcommand::CommandFailed(), G4VScoringMesh::FindPrimitiveScorer(), G4UIcommand::GetCommandPath(), qTouchCmd, and G4VScoringMesh::SetNullToCurrentPrimitiveScorer().

Referenced by SetNewValue().

◆ CommandsShouldBeInMaster()

G4bool G4UImessenger::CommandsShouldBeInMaster ( ) const
inlineinherited

Definition at line 77 of file G4UImessenger.hh.

78 {
80 }
G4bool commandsShouldBeInMaster

References G4UImessenger::commandsShouldBeInMaster.

Referenced by G4UIcommand::G4UIcommandCommonConstructorCode().

◆ CreateCommand()

template<typename T >
T * G4UImessenger::CreateCommand ( const G4String cname,
const G4String dsc 
)
protectedinherited

Definition at line 110 of file G4UImessenger.hh.

111{
112 G4String path;
113 if(cname[0] != '/')
114 {
115 path = baseDirName + cname;
116 if(path[0] != '/')
117 path = "/" + path;
118 }
119
120 T* command = new T(path.c_str(), this);
121 command->SetGuidance(dsc.c_str());
122
123 return command;
124}
G4String baseDirName

References G4UImessenger::baseDirName.

◆ CreateDirectory()

void G4UImessenger::CreateDirectory ( const G4String path,
const G4String dsc,
G4bool  commandsToBeBroadcasted = true 
)
protectedinherited

Definition at line 156 of file G4UImessenger.cc.

158{
160
161 G4String fullpath = path;
162 if(fullpath.back() != '/')
163 fullpath.append("/");
164
165 G4UIcommandTree* tree = ui->GetTree()->FindCommandTree(fullpath.c_str());
166 if(tree != nullptr)
167 {
168 baseDirName = tree->GetPathName();
169 }
170 else
171 {
172 baseDir = new G4UIdirectory(fullpath.c_str(), commandsToBeBroadcasted);
173 baseDirName = fullpath;
174 baseDir->SetGuidance(dsc.c_str());
175 }
176}
const G4String & GetPathName() const
G4UIcommandTree * FindCommandTree(const char *commandPath)
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:156
G4UIcommandTree * GetTree() const
Definition: G4UImanager.hh:186
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77
G4UIdirectory * baseDir

References G4UImessenger::baseDir, G4UImessenger::baseDirName, G4UIcommandTree::FindCommandTree(), G4UIcommandTree::GetPathName(), G4UImanager::GetTree(), G4UImanager::GetUIpointer(), and G4UIcommand::SetGuidance().

Referenced by G4MoleculeShootMessenger::G4MoleculeShootMessenger(), and G4UImessenger::G4UImessenger().

◆ DtoS()

G4String G4UImessenger::DtoS ( G4double  a)
protectedinherited

Definition at line 90 of file G4UImessenger.cc.

91{
92 std::ostringstream os;
93 os << a;
94 return G4String(os.str());
95}

Referenced by FilterCommands(), and G4UIcontrolMessenger::SetNewValue().

◆ FillTokenVec()

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

Definition at line 815 of file G4ScoreQuantityMessenger.cc.

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}

References anonymous_namespace{G4MTcoutDestination.cc}::empty.

Referenced by SetNewValue().

◆ FilterCommands()

void G4ScoreQuantityMessenger::FilterCommands ( )
private

Definition at line 124 of file G4ScoreQuantityMessenger.cc.

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}
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:146
G4String DtoS(G4double a)
void SetDefaultValue(const char *theDefaultValue)
void SetDefaultUnit(const char *theDefaultUnit)
#define DBL_MAX
Definition: templates.hh:62

References DBL_MAX, G4UImessenger::DtoS(), fchargedCmd, filterDir, fkinECmd, fneutralCmd, fparticleCmd, fparticleKinECmd, G4UIparameter::SetDefaultUnit(), G4UIparameter::SetDefaultValue(), G4UIcommand::SetGuidance(), G4UIcommand::SetParameter(), and G4UIcmdWithAString::SetParameterName().

Referenced by G4ScoreQuantityMessenger().

◆ FParticleCommand()

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

Definition at line 826 of file G4ScoreQuantityMessenger.cc.

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}
int G4int
Definition: G4Types.hh:85
void SetFilter(G4VSDFilter *filter)
const char * name(G4int ptype)

References G4InuclParticleNames::name(), and G4VScoringMesh::SetFilter().

Referenced by SetNewValue().

◆ FParticleWithEnergyCommand()

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

Definition at line 844 of file G4ScoreQuantityMessenger.cc.

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}
double G4double
Definition: G4Types.hh:83
void add(const G4String &particleName)
G4double StoD(G4String s)
static G4double GetValueOf(const G4String &)

References G4SDParticleWithEnergyFilter::add(), G4UnitDefinition::GetValueOf(), G4InuclParticleNames::name(), G4VScoringMesh::SetFilter(), and G4UImessenger::StoD().

Referenced by SetNewValue().

◆ GetCurrentValue()

G4String G4ScoreQuantityMessenger::GetCurrentValue ( G4UIcommand )
virtual

Reimplemented from G4UImessenger.

Definition at line 808 of file G4ScoreQuantityMessenger.cc.

809{
810 G4String val;
811
812 return val;
813}

◆ ItoS()

G4String G4UImessenger::ItoS ( G4int  i)
protectedinherited

Definition at line 82 of file G4UImessenger.cc.

83{
84 std::ostringstream os;
85 os << i;
86 return G4String(os.str());
87}

Referenced by G4GenericMessenger::DeclareMethod(), and G4ParticleGunMessenger::GetCurrentValue().

◆ operator!=()

G4bool G4UImessenger::operator!= ( const G4UImessenger messenger) const
inherited

Definition at line 76 of file G4UImessenger.cc.

77{
78 return this != &messenger;
79}

◆ operator==()

G4bool G4UImessenger::operator== ( const G4UImessenger messenger) const
inherited

Definition at line 70 of file G4UImessenger.cc.

71{
72 return this == &messenger;
73}

◆ QuantityCommands()

void G4ScoreQuantityMessenger::QuantityCommands ( )
private

Definition at line 85 of file G4ScoreQuantityMessengerQCmd.cc.

86{
87 G4UIparameter* param;
88
89 //
90 // Quantity commands
91 quantityDir = new G4UIdirectory("/score/quantity/");
92 quantityDir->SetGuidance("Scoring quantity of the mesh.");
93 //
94 qTouchCmd = new G4UIcmdWithAString("/score/quantity/touch", this);
96 "Assign previously defined quantity to the current quantity.");
97 qTouchCmd->SetParameterName("qname", false);
98 //
99 qGetUnitCmd = new G4UIcmdWithoutParameter("/score/quantity/get/unit", this);
100 qGetUnitCmd->SetGuidance("Print output unit of the current quantity.");
101 //
102 qSetUnitCmd = new G4UIcmdWithAString("/score/quantity/set/unit", this);
103 qSetUnitCmd->SetGuidance("Set output unit of the current quantity.");
104 qSetUnitCmd->SetParameterName("unit", false);
105
106 // Primitive Scorers
107 qeDepCmd = new G4UIcommand("/score/quantity/energyDeposit", this);
108 qeDepCmd->SetGuidance("Energy deposit scorer.");
109 qeDepCmd->SetGuidance("[usage] /score/quantity/energyDeposit qname unit");
110 qeDepCmd->SetGuidance(" qname :(String) scorer name");
111 qeDepCmd->SetGuidance(" unit :(String) unit");
112 param = new G4UIparameter("qname", 's', false);
113 qeDepCmd->SetParameter(param);
114 param = new G4UIparameter("unit", 's', true);
115 param->SetDefaultUnit("MeV");
116 qeDepCmd->SetParameter(param);
117 //
118 qCellChgCmd = new G4UIcommand("/score/quantity/cellCharge", this);
119 qCellChgCmd->SetGuidance("Cell charge scorer.");
120 qCellChgCmd->SetGuidance("[usage] /score/quantity/cellCharge qname unit");
121 qCellChgCmd->SetGuidance(" qname :(String) scorer name");
122 qCellChgCmd->SetGuidance(" unit :(String) unit");
123 param = new G4UIparameter("qname", 's', false);
125 param = new G4UIparameter("unit", 's', true);
126 param->SetDefaultUnit("e+");
128 //
129 qCellFluxCmd = new G4UIcommand("/score/quantity/cellFlux", this);
130 qCellFluxCmd->SetGuidance("Cell flux scorer.");
131 qCellFluxCmd->SetGuidance("[usage] /score/quantity/cellFlux qname unit");
132 qCellFluxCmd->SetGuidance(" qname :(String) scorer name");
133 qCellFluxCmd->SetGuidance(" unit :(String) unit");
134 param = new G4UIparameter("qname", 's', false);
136 param = new G4UIparameter("unit", 's', true);
137 param->SetDefaultValue("percm2");
139 //
140 qPassCellFluxCmd = new G4UIcommand("/score/quantity/passageCellFlux", this);
141 qPassCellFluxCmd->SetGuidance("Passage cell flux scorer");
143 "[usage] /score/quantity/passageCellFlux qname unit");
144 qPassCellFluxCmd->SetGuidance(" qname :(String) scorer name");
145 qPassCellFluxCmd->SetGuidance(" unit :(String) unit");
146 param = new G4UIparameter("qname", 's', false);
148 param = new G4UIparameter("unit", 's', true);
149 param->SetDefaultValue("percm2");
151 //
152 qdoseDepCmd = new G4UIcommand("/score/quantity/doseDeposit", this);
153 qdoseDepCmd->SetGuidance("Dose deposit scorer.");
154 qdoseDepCmd->SetGuidance("[usage] /score/quantity/doseDeposit qname unit");
155 qdoseDepCmd->SetGuidance(" qname :(String) scorer name");
156 qdoseDepCmd->SetGuidance(" unit :(String) unit");
157 param = new G4UIparameter("qname", 's', false);
159 param = new G4UIparameter("unit", 's', true);
160 param->SetDefaultUnit("Gy");
162 //
163 qnOfStepCmd = new G4UIcommand("/score/quantity/nOfStep", this);
164 qnOfStepCmd->SetGuidance("Number of step scorer.");
165 qnOfStepCmd->SetGuidance("[usage] /score/quantity/nOfStep qname");
166 qnOfStepCmd->SetGuidance("[usage] /score/quantity/nOfStep qname bflag");
167 qnOfStepCmd->SetGuidance(" qname :(String) scorer name");
168 qnOfStepCmd->SetGuidance(" bflag :(Bool) Skip zero step ");
169 qnOfStepCmd->SetGuidance(" at geometry boundary if true");
170 param = new G4UIparameter("qname", 's', false);
172 param = new G4UIparameter("bflag", 'b', true);
173 param->SetDefaultValue("false");
175 //
176 qnOfSecondaryCmd = new G4UIcommand("/score/quantity/nOfSecondary", this);
177 qnOfSecondaryCmd->SetGuidance("Number of secondary scorer.");
178 qnOfSecondaryCmd->SetGuidance("[usage] /score/quantity/nOfSecondary qname");
179 qnOfSecondaryCmd->SetGuidance(" qname :(String) scorer name");
180 param = new G4UIparameter("qname", 's', false);
182 //
183 qTrackLengthCmd = new G4UIcommand("/score/quantity/trackLength", this);
184 qTrackLengthCmd->SetGuidance("Track length scorer.");
186 "[usage] /score/quantity/trackLength qname wflag kflag vflag unit");
187 qTrackLengthCmd->SetGuidance(" qname :(String) scorer name");
188 qTrackLengthCmd->SetGuidance(" wflag :(Bool) weighted");
189 qTrackLengthCmd->SetGuidance(" kflag :(Bool) multiply kinetic energy");
190 qTrackLengthCmd->SetGuidance(" vflag :(Bool) divide by velocity");
191 qTrackLengthCmd->SetGuidance(" unit :(String) unit");
192 param = new G4UIparameter("qname", 's', false);
194 param = new G4UIparameter("wflag", 'b', true);
195 param->SetDefaultValue("false");
197 param = new G4UIparameter("kflag", 'b', true);
198 param->SetDefaultValue("false");
200 param = new G4UIparameter("vflag", 'b', true);
201 param->SetDefaultValue("false");
203 param = new G4UIparameter("unit", 's', true);
204 param->SetDefaultValue("mm");
206 //
208 new G4UIcommand("/score/quantity/passageCellCurrent", this);
209 qPassCellCurrCmd->SetGuidance("Passage cell current scorer.");
211 "[usage] /score/quantity/passageCellCurrent qname wflag");
212 qPassCellCurrCmd->SetGuidance(" qname :(String) scorer name");
213 qPassCellCurrCmd->SetGuidance(" wflag :(Bool) weighted");
214 param = new G4UIparameter("qname", 's', false);
216 param = new G4UIparameter("wflag", 'b', true);
217 param->SetDefaultValue("true");
219 //
221 new G4UIcommand("/score/quantity/passageTrackLength", this);
222 qPassTrackLengthCmd->SetGuidance("Passage track length scorer.");
224 "[usage] /score/quantity/passageTrackLength qname wflag unit");
225 qPassTrackLengthCmd->SetGuidance(" qname :(String) scorer name");
226 qPassTrackLengthCmd->SetGuidance(" wflag :(Bool) weighted");
227 qPassTrackLengthCmd->SetGuidance(" unit :(Bool) unit");
228 param = new G4UIparameter("qname", 's', false);
230 param = new G4UIparameter("wflag", 'b', true);
231 param->SetDefaultValue("true");
233 param = new G4UIparameter("unit", 's', true);
234 param->SetDefaultUnit("mm");
236 //
238 new G4UIcommand("/score/quantity/flatSurfaceCurrent", this);
239 qFlatSurfCurrCmd->SetGuidance("Flat surface current Scorer.");
241 "[usage] /score/quantity/flatSurfaceCurrent qname dflag wflag aflag unit");
242 qFlatSurfCurrCmd->SetGuidance(" qname :(String) scorer name");
243 qFlatSurfCurrCmd->SetGuidance(" dflag :(Int) direction flag");
244 qFlatSurfCurrCmd->SetGuidance(" : 0 = Both In and Out");
245 qFlatSurfCurrCmd->SetGuidance(" : 1 = In only");
246 qFlatSurfCurrCmd->SetGuidance(" : 2 = Out only");
247 qFlatSurfCurrCmd->SetGuidance(" wflag :(Bool) weighted");
248 qFlatSurfCurrCmd->SetGuidance(" aflag :(Bool) divide by area");
249 qFlatSurfCurrCmd->SetGuidance(" unit :(String) unit");
250 param = new G4UIparameter("qname", 's', false);
252 param = new G4UIparameter("dflag", 'i', true);
253 param->SetDefaultValue("0");
255 param = new G4UIparameter("wflag", 'b', true);
256 param->SetDefaultValue("true");
258 param = new G4UIparameter("aflag", 'b', true);
259 param->SetDefaultValue("true");
261 param = new G4UIparameter("unit", 's', true);
262 param->SetDefaultValue("percm2");
264 //
265 qFlatSurfFluxCmd = new G4UIcommand("/score/quantity/flatSurfaceFlux", this);
266 qFlatSurfFluxCmd->SetGuidance("Flat surface flux scorer.");
268 "[usage] /score/quantity/flatSurfaceFlux qname dflag unit");
269 qFlatSurfFluxCmd->SetGuidance(" qname :(String) scorer name");
270 qFlatSurfFluxCmd->SetGuidance(" dflag :(Int) direction flag");
271 qFlatSurfFluxCmd->SetGuidance(" : 0 = Both In and Out");
272 qFlatSurfFluxCmd->SetGuidance(" : 1 = In only");
273 qFlatSurfFluxCmd->SetGuidance(" : 2 = Out only");
274 qFlatSurfFluxCmd->SetGuidance(" wflag :(Bool) weighted");
275 qFlatSurfFluxCmd->SetGuidance(" aflag :(Bool) divide by area");
276 qFlatSurfFluxCmd->SetGuidance(" unit :(String) unit");
277 param = new G4UIparameter("qname", 's', false);
279 param = new G4UIparameter("dflag", 'i', true);
280 param->SetDefaultValue("0");
282 param = new G4UIparameter("wflag", 'b', true);
283 param->SetDefaultValue("true");
285 param = new G4UIparameter("aflag", 'b', true);
286 param->SetDefaultValue("true");
288 param = new G4UIparameter("unit", 's', true);
289 param->SetDefaultValue("percm2");
291 //
292
293 qVolFluxCmd = new G4UIcommand("/score/quantity/volumeFlux", this);
294 qVolFluxCmd->SetGuidance("Volume flux scorer.");
296 "This scorer scores the number of particles getting into the volume "
297 "without normalized by the surface area.");
299 "[usage] /score/quantity/volumeFlux qname divcos dflag");
300 qVolFluxCmd->SetGuidance(" qname :(String) scorer name");
301 qVolFluxCmd->SetGuidance(" divcos :(Bool) divide by cos(theta), where theta "
302 "is the incident angle (default : false)");
304 " dflag :(Int) direction, 1 : inward (default), 2 : outward");
305 param = new G4UIparameter("qname", 's', false);
307 param = new G4UIparameter("divcos", 'b', true);
308 param->SetDefaultValue(false);
310 param = new G4UIparameter("dflag", 'i', true);
311 param->SetParameterRange("dflag>=1 && dflag<=2");
312 param->SetDefaultValue(1);
314
315 // qSphereSurfCurrCmd = new
316 // G4UIcommand("/score/quantity/sphereSurfaceCurrent",this);
317 // qSphereSurfCurrCmd->SetGuidance("Sphere surface current Scorer.");
318 // qSphereSurfCurrCmd->
319 // SetGuidance("[usage] /score/quantity/sphereSurfaceCurrent qname dflag
320 // wflag aflag unit");
321 // qSphereSurfCurrCmd->SetGuidance(" qname :(String) scorer name");
322 // qSphereSurfCurrCmd->SetGuidance(" dflag :(Int) direction flag");
323 // qSphereSurfCurrCmd->SetGuidance(" : 0 = Both In and Out");
324 // qSphereSurfCurrCmd->SetGuidance(" : 1 = In only");
325 // qSphereSurfCurrCmd->SetGuidance(" : 2 = Out only");
326 // qSphereSurfCurrCmd->SetGuidance(" wflag :(Bool) Weighted");
327 // qSphereSurfCurrCmd->SetGuidance(" aflag :(Bool) DivideByArea");
328 // qSphereSurfCurrCmd->SetGuidance(" unit :(String) unit");
329 // param = new G4UIparameter("qname",'s',false);
330 // qSphereSurfCurrCmd->SetParameter(param);
331 // param = new G4UIparameter("dflag",'i',true);
332 // param->SetDefaultValue("0");
333 // qSphereSurfCurrCmd->SetParameter(param);
334 // param = new G4UIparameter("wflag",'b',true);
335 // param->SetDefaultValue("true");
336 // qSphereSurfCurrCmd->SetParameter(param);
337 // param = new G4UIparameter("aflag",'b',true);
338 // param->SetDefaultValue("true");
339 // qSphereSurfCurrCmd->SetParameter(param);
340 // param = new G4UIparameter("unit",'s',true);
341 // param->SetDefaultValue("percm2");
342 // qSphereSurfCurrCmd->SetParameter(param);
343
344 //
345 // qSphereSurfFluxCmd = new
346 // G4UIcommand("/score/quantity/sphereSurfaceFlux",this);
347 // qSphereSurfFluxCmd->SetGuidance("Sphere surface Flux Scorer.");
348 // qSphereSurfFluxCmd->
349 // SetGuidance("[usage] /score/quantity/sphereSurfaceFlux qname dflag unit");
350 // qSphereSurfFluxCmd->SetGuidance(" qname :(String) scorer name");
351 // qSphereSurfFluxCmd->SetGuidance(" dflag :(Int) direction flag");
352 // qSphereSurfFluxCmd->SetGuidance(" : 0 = Both In and Out");
353 // qSphereSurfFluxCmd->SetGuidance(" : 1 = In only");
354 // qSphereSurfFluxCmd->SetGuidance(" : 2 = Out only");
355 // qSphereSurfFluxCmd->SetGuidance(" wflag :(Bool) weighted");
356 // qSphereSurfFluxCmd->SetGuidance(" aflag :(Bool) divide by area");
357 // qSphereSurfFluxCmd->SetGuidance(" unit :(String) unit");
358 // param = new G4UIparameter("qname",'s',false);
359 // qSphereSurfFluxCmd->SetParameter(param);
360 // param = new G4UIparameter("dflag",'i',true);
361 // param->SetDefaultValue("0");
362 // qSphereSurfFluxCmd->SetParameter(param);
363 // param = new G4UIparameter("wflag",'b',true);
364 // param->SetDefaultValue("true");
365 // qSphereSurfFluxCmd->SetParameter(param);
366 // param = new G4UIparameter("aflag",'b',true);
367 // param->SetDefaultValue("true");
368 // qSphereSurfFluxCmd->SetParameter(param);
369 // param = new G4UIparameter("unit",'s',true);
370 // param->SetDefaultValue("percm2");
371 // qSphereSurfFluxCmd->SetParameter(param);
372
373 //
374 // qCylSurfCurrCmd = new
375 // G4UIcommand("/score/quantity/cylinderSurfaceCurrent",this);
376 // qCylSurfCurrCmd->SetGuidance("Cylinder surface current Scorer.");
377 // qCylSurfCurrCmd->
378 // SetGuidance("[usage] /score/quantity/cylinderSurfaceCurrent qname dflag
379 // wflag aflag unit"); qCylSurfCurrCmd->SetGuidance(" qname :(String)
380 // scorer name"); qCylSurfCurrCmd->SetGuidance(" dflag :(Int) direction
381 // flag"); qCylSurfCurrCmd->SetGuidance(" : 0 = Both In and Out");
382 // qCylSurfCurrCmd->SetGuidance(" : 1 = In only");
383 // qCylSurfCurrCmd->SetGuidance(" : 2 = Out only");
384 // qCylSurfCurrCmd->SetGuidance(" wflag :(Bool) Weighted");
385 // qCylSurfCurrCmd->SetGuidance(" aflag :(Bool) DivideByArea");
386 // qCylSurfCurrCmd->SetGuidance(" unit :(String) unit");
387 // param = new G4UIparameter("qname",'s',false);
388 // qCylSurfCurrCmd->SetParameter(param);
389 // param = new G4UIparameter("dflag",'i',true);
390 // param->SetDefaultValue("0");
391 // qCylSurfCurrCmd->SetParameter(param);
392 // param = new G4UIparameter("wflag",'b',true);
393 // param->SetDefaultValue("true");
394 // qCylSurfCurrCmd->SetParameter(param);
395 // param = new G4UIparameter("aflag",'b',true);
396 // param->SetDefaultValue("true");
397 // qCylSurfCurrCmd->SetParameter(param);
398 // param = new G4UIparameter("unit",'s',true);
399 // param->SetDefaultValue("percm2");
400 // qCylSurfCurrCmd->SetParameter(param);
401 //
402 // qCylSurfFluxCmd = new
403 // G4UIcommand("/score/quantity/cylinderSurfaceFlux",this);
404 // qCylSurfFluxCmd->SetGuidance("Cylinder surface Flux Scorer.");
405 // qCylSurfFluxCmd->
406 // SetGuidance("[usage] /score/quantity/cylinderSurfaceFlux qname dflag
407 // unit"); qCylSurfFluxCmd->SetGuidance(" qname :(String) scorer name");
408 // qCylSurfFluxCmd->SetGuidance(" dflag :(Int) direction flag");
409 // qCylSurfFluxCmd->SetGuidance(" : 0 = Both In and Out");
410 // qCylSurfFluxCmd->SetGuidance(" : 1 = In only");
411 // qCylSurfFluxCmd->SetGuidance(" : 2 = Out only");
412 // qCylSurfFluxCmd->SetGuidance(" wflag :(Bool) weighted");
413 // qCylSurfFluxCmd->SetGuidance(" aflag :(Bool) divide by area");
414 // qCylSurfFluxCmd->SetGuidance(" unit :(String) unit");
415 // param = new G4UIparameter("qname",'s',false);
416 // qCylSurfFluxCmd->SetParameter(param);
417 // param = new G4UIparameter("dflag",'i',true);
418 // param->SetDefaultValue("0");
419 // qCylSurfFluxCmd->SetParameter(param);
420 // param = new G4UIparameter("wflag",'b',true);
421 // param->SetDefaultValue("true");
422 // qCylSurfFluxCmd->SetParameter(param);
423 // param = new G4UIparameter("aflag",'b',true);
424 // param->SetDefaultValue("true");
425 // qCylSurfFluxCmd->SetParameter(param);
426 // param = new G4UIparameter("unit",'s',true);
427 // param->SetDefaultValue("percm2");
428 // qCylSurfFluxCmd->SetParameter(param);
429 //
430 //
431 qNofCollisionCmd = new G4UIcommand("/score/quantity/nOfCollision", this);
432 qNofCollisionCmd->SetGuidance("Number of collision scorer.");
434 "[usage] /score/quantity/nOfCollision qname wflag");
435 qNofCollisionCmd->SetGuidance(" qname :(String) scorer name");
436 param = new G4UIparameter("qname", 's', false);
438 param = new G4UIparameter("wflag", 'b', true);
439 param->SetDefaultValue("false");
441 //
442 qPopulationCmd = new G4UIcommand("/score/quantity/population", this);
443 qPopulationCmd->SetGuidance("Population scorer.");
444 qPopulationCmd->SetGuidance("[usage] /score/quantity/population qname wflag");
445 qPopulationCmd->SetGuidance(" qname :(String) scorer name");
446 qPopulationCmd->SetGuidance(" wflag :(Bool) weighted");
447 param = new G4UIparameter("qname", 's', false);
449 param = new G4UIparameter("wflag", 'b', true);
450 param->SetDefaultValue("false");
452
453 //
454 qTrackCountCmd = new G4UIcommand("/score/quantity/nOfTrack", this);
455 qTrackCountCmd->SetGuidance("Number of track scorer.");
457 "[usage] /score/quantity/nOfTrack qname dflag wflag");
458 qTrackCountCmd->SetGuidance(" qname :(String) scorer name");
459 qTrackCountCmd->SetGuidance(" dflag :(Int) direction");
460 qTrackCountCmd->SetGuidance(" : 0 = Both In and Out");
461 qTrackCountCmd->SetGuidance(" : 1 = In only");
462 qTrackCountCmd->SetGuidance(" : 2 = Out only");
463 qTrackCountCmd->SetGuidance(" wflag :(Bool) weighted");
464 param = new G4UIparameter("qname", 's', false);
466 param = new G4UIparameter("dflag", 'i', true);
467 param->SetDefaultValue("0");
469 param = new G4UIparameter("wflag", 'b', true);
470 param->SetDefaultValue("false");
472
473 //
474 qTerminationCmd = new G4UIcommand("/score/quantity/nOfTerminatedTrack", this);
475 qTerminationCmd->SetGuidance("Number of terminated tracks scorer.");
477 "[usage] /score/quantity/nOfTerminatedTrack qname wflag");
478 qTerminationCmd->SetGuidance(" qname :(String) scorer name");
479 qTerminationCmd->SetGuidance(" wflag :(Bool) weighted");
480 param = new G4UIparameter("qname", 's', false);
482 param = new G4UIparameter("wflag", 'b', true);
483 param->SetDefaultValue("false");
485
486 //
488 new G4UIcommand("/score/quantity/minKinEAtGeneration", this);
489 qMinKinEAtGeneCmd->SetGuidance("Min Kinetic Energy at Generation");
491 "[usage] /score/quantity/minKinEAtGeneration qname unit");
492 qMinKinEAtGeneCmd->SetGuidance(" qname :(String) scorer name");
493 qMinKinEAtGeneCmd->SetGuidance(" unit :(String) unit name");
494 param = new G4UIparameter("qname", 's', false);
496 param = new G4UIparameter("unit", 's', true);
497 param->SetDefaultUnit("MeV");
499 //
500 qStepCheckerCmd = new G4UIcommand("/score/quantity/stepChecker", this);
501 qStepCheckerCmd->SetGuidance("Display a comment when this PS is invoked");
502 qStepCheckerCmd->SetGuidance("[usage] /score/quantity/stepChecker qname");
503 qStepCheckerCmd->SetGuidance(" qname :(String) scorer name");
504 param = new G4UIparameter("qname", 's', false);
506}
void SetParameterRange(const char *theRange)

References qCellChgCmd, qCellFluxCmd, qdoseDepCmd, qeDepCmd, qFlatSurfCurrCmd, qFlatSurfFluxCmd, qGetUnitCmd, qMinKinEAtGeneCmd, qNofCollisionCmd, qnOfSecondaryCmd, qnOfStepCmd, qPassCellCurrCmd, qPassCellFluxCmd, qPassTrackLengthCmd, qPopulationCmd, qSetUnitCmd, qStepCheckerCmd, qTerminationCmd, qTouchCmd, qTrackCountCmd, qTrackLengthCmd, quantityDir, qVolFluxCmd, G4UIparameter::SetDefaultUnit(), G4UIparameter::SetDefaultValue(), G4UIcommand::SetGuidance(), G4UIcommand::SetParameter(), G4UIcmdWithAString::SetParameterName(), and G4UIparameter::SetParameterRange().

Referenced by G4ScoreQuantityMessenger().

◆ SetNewValue()

void G4ScoreQuantityMessenger::SetNewValue ( G4UIcommand command,
G4String  newValues 
)
virtual

Reimplemented from G4UImessenger.

Definition at line 249 of file G4ScoreQuantityMessenger.cc.

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}
static constexpr double ps
Definition: G4SIunits.hh:157
std::vector< G4String > G4TokenVec
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 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
G4int StoI(G4String s)
G4bool StoB(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()

References CheckMeshPS(), G4UIcommand::CommandFailed(), fchargedCmd, FillTokenVec(), fkinECmd, fneutralCmd, fparticleCmd, FParticleCommand(), fparticleKinECmd, FParticleWithEnergyCommand(), fSMan, G4cout, G4endl, G4VScoringMesh::GetAngleSpan(), G4UIcommand::GetCommandPath(), G4VScoringMesh::GetCopyNumberLevel(), G4ScoringManager::GetCurrentMesh(), G4VScoringMesh::GetCurrentPSUnit(), G4VScoringMesh::GetNumberOfSegments(), G4VScoringMesh::GetShape(), G4VScoringMesh::GetSize(), G4VScoringMesh::GetStartAngle(), G4UnitDefinition::GetValueOf(), G4VScoringMesh::IsCurrentPrimitiveScorerNull(), G4InuclParticleNames::name(), ps, qCellChgCmd, qCellFluxCmd, qdoseDepCmd, qeDepCmd, qFlatSurfCurrCmd, qFlatSurfFluxCmd, qGetUnitCmd, qMinKinEAtGeneCmd, qNofCollisionCmd, qnOfSecondaryCmd, qnOfStepCmd, qPassCellCurrCmd, qPassCellFluxCmd, qPassTrackLengthCmd, qPopulationCmd, qSetUnitCmd, qStepCheckerCmd, qTerminationCmd, qTouchCmd, qTrackCountCmd, qTrackLengthCmd, qVolFluxCmd, G4VScoringMesh::SetCurrentPrimitiveScorer(), G4VScoringMesh::SetCurrentPSUnit(), G4PSCellFluxForCylinder3D::SetCylinderSize(), G4PSDoseDepositForCylinder3D::SetCylinderSize(), G4PSPassageCellFluxForCylinder3D::SetCylinderSize(), G4VScoringMesh::SetFilter(), G4PSCellFluxForCylinder3D::SetNumberOfSegments(), G4PSDoseDepositForCylinder3D::SetNumberOfSegments(), G4PSPassageCellFluxForCylinder3D::SetNumberOfSegments(), G4VScoringMesh::SetPrimitiveScorer(), G4PSDoseDeposit::SetUnit(), G4UImessenger::StoB(), G4UImessenger::StoD(), and G4UImessenger::StoI().

◆ StoB()

G4bool G4UImessenger::StoB ( G4String  s)
protectedinherited

Definition at line 137 of file G4UImessenger.cc.

138{
140 G4bool vl = false;
141 if(v == "Y" || v == "YES" || v == "1" || v == "T" || v == "TRUE")
142 {
143 vl = true;
144 }
145 return vl;
146}
bool G4bool
Definition: G4Types.hh:86
G4String to_upper_copy(G4String str)
Return uppercase copy of string.

References G4StrUtil::to_upper_copy().

Referenced by G4LocalThreadCoutMessenger::SetNewValue(), G4CascadeParamMessenger::SetNewValue(), SetNewValue(), and G4ScoringMessenger::SetNewValue().

◆ StoD()

G4double G4UImessenger::StoD ( G4String  s)
protectedinherited

◆ StoI()

G4int G4UImessenger::StoI ( G4String  s)
protectedinherited

◆ StoL()

G4long G4UImessenger::StoL ( G4String  s)
protectedinherited

Definition at line 117 of file G4UImessenger.cc.

118{
119 G4long vl;
120 const char* t = str;
121 std::istringstream is(t);
122 is >> vl;
123 return vl;
124}
long G4long
Definition: G4Types.hh:87

Referenced by G4RunMessenger::SetNewValue().

Field Documentation

◆ baseDir

G4UIdirectory* G4UImessenger::baseDir = nullptr
protectedinherited

◆ baseDirName

G4String G4UImessenger::baseDirName = ""
protectedinherited

◆ commandsShouldBeInMaster

G4bool G4UImessenger::commandsShouldBeInMaster = false
protectedinherited

◆ fchargedCmd

G4UIcmdWithAString* G4ScoreQuantityMessenger::fchargedCmd
private

◆ filterDir

G4UIdirectory* G4ScoreQuantityMessenger::filterDir
private

Definition at line 118 of file G4ScoreQuantityMessenger.hh.

Referenced by FilterCommands(), and ~G4ScoreQuantityMessenger().

◆ fkinECmd

G4UIcommand* G4ScoreQuantityMessenger::fkinECmd
private

◆ fneutralCmd

G4UIcmdWithAString* G4ScoreQuantityMessenger::fneutralCmd
private

◆ fparticleCmd

G4UIcommand* G4ScoreQuantityMessenger::fparticleCmd
private

◆ fparticleKinECmd

G4UIcommand* G4ScoreQuantityMessenger::fparticleKinECmd
private

◆ fSMan

G4ScoringManager* G4ScoreQuantityMessenger::fSMan
private

Definition at line 82 of file G4ScoreQuantityMessenger.hh.

Referenced by SetNewValue().

◆ qCellChgCmd

G4UIcommand* G4ScoreQuantityMessenger::qCellChgCmd
private

◆ qCellFluxCmd

G4UIcommand* G4ScoreQuantityMessenger::qCellFluxCmd
private

◆ qdoseDepCmd

G4UIcommand* G4ScoreQuantityMessenger::qdoseDepCmd
private

◆ qeDepCmd

G4UIcommand* G4ScoreQuantityMessenger::qeDepCmd
private

◆ qFlatSurfCurrCmd

G4UIcommand* G4ScoreQuantityMessenger::qFlatSurfCurrCmd
private

◆ qFlatSurfFluxCmd

G4UIcommand* G4ScoreQuantityMessenger::qFlatSurfFluxCmd
private

◆ qGetUnitCmd

G4UIcmdWithoutParameter* G4ScoreQuantityMessenger::qGetUnitCmd
private

◆ qMinKinEAtGeneCmd

G4UIcommand* G4ScoreQuantityMessenger::qMinKinEAtGeneCmd
private

◆ qNofCollisionCmd

G4UIcommand* G4ScoreQuantityMessenger::qNofCollisionCmd
private

◆ qnOfSecondaryCmd

G4UIcommand* G4ScoreQuantityMessenger::qnOfSecondaryCmd
private

◆ qnOfStepCmd

G4UIcommand* G4ScoreQuantityMessenger::qnOfStepCmd
private

◆ qPassCellCurrCmd

G4UIcommand* G4ScoreQuantityMessenger::qPassCellCurrCmd
private

◆ qPassCellFluxCmd

G4UIcommand* G4ScoreQuantityMessenger::qPassCellFluxCmd
private

◆ qPassTrackLengthCmd

G4UIcommand* G4ScoreQuantityMessenger::qPassTrackLengthCmd
private

◆ qPopulationCmd

G4UIcommand* G4ScoreQuantityMessenger::qPopulationCmd
private

◆ qSetUnitCmd

G4UIcmdWithAString* G4ScoreQuantityMessenger::qSetUnitCmd
private

◆ qStepCheckerCmd

G4UIcommand* G4ScoreQuantityMessenger::qStepCheckerCmd
private

◆ qTerminationCmd

G4UIcommand* G4ScoreQuantityMessenger::qTerminationCmd
private

◆ qTouchCmd

G4UIcmdWithAString* G4ScoreQuantityMessenger::qTouchCmd
private

◆ qTrackCountCmd

G4UIcommand* G4ScoreQuantityMessenger::qTrackCountCmd
private

◆ qTrackLengthCmd

G4UIcommand* G4ScoreQuantityMessenger::qTrackLengthCmd
private

◆ quantityDir

G4UIdirectory* G4ScoreQuantityMessenger::quantityDir
private

Definition at line 85 of file G4ScoreQuantityMessenger.hh.

Referenced by QuantityCommands(), and ~G4ScoreQuantityMessenger().

◆ qVolFluxCmd

G4UIcommand* G4ScoreQuantityMessenger::qVolFluxCmd
private

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