Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ScoreQuantityMessenger.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4ScoreQuantityMessenger.cc 73819 2013-09-12 15:52:52Z gcosmo $
28 //
29 // ---------------------------------------------------------------------
30 // Modifications
31 // 08-Oct-2010 T.Aso remove unit of G4PSPassageCellCurrent.
32 // 24-Mar-2011 T.Aso Add StepChecker for debugging.
33 // 24-Mar-2011 T.Aso Size and segmentation for replicated cylinder.
34 // 01-Jun-2012 T.Aso Support weighted/dividedByArea options
35 // in flatCurrent and flatFulx commands.
36 // 27-Mar-2013 T.Aso Unit option in the kineticEnergy filter was
37 // supported.
38 //
39 // ---------------------------------------------------------------------
40 
42 #include "G4ScoringManager.hh"
43 #include "G4VScoringMesh.hh"
44 
45 #include "G4PSCellCharge3D.hh"
46 #include "G4PSCellFlux3D.hh"
48 #include "G4PSPassageCellFlux3D.hh"
50 #include "G4PSEnergyDeposit3D.hh"
51 #include "G4PSDoseDeposit3D.hh"
53 #include "G4PSNofStep3D.hh"
54 #include "G4PSNofSecondary3D.hh"
55 //
56 #include "G4PSTrackLength3D.hh"
60 #include "G4PSFlatSurfaceFlux3D.hh"
65 #include "G4PSNofCollision3D.hh"
66 #include "G4PSPopulation3D.hh"
67 #include "G4PSTrackCounter3D.hh"
68 #include "G4PSTermination3D.hh"
70 //
71 // For debug purpose
72 #include "G4PSStepChecker3D.hh"
73 
74 #include "G4SDChargedFilter.hh"
75 #include "G4SDNeutralFilter.hh"
77 #include "G4SDParticleFilter.hh"
79 
80 #include "G4UIdirectory.hh"
82 #include "G4UIcmdWithAnInteger.hh"
83 #include "G4UIcmdWithAString.hh"
84 #include "G4UIcmdWithABool.hh"
87 #include "G4UIcommand.hh"
88 #include "G4Tokenizer.hh"
89 #include "G4UnitsTable.hh"
90 
92 :fSMan(SManager)
93 {
94  QuantityCommands();
95  FilterCommands();
96 }
97 
98 void G4ScoreQuantityMessenger::FilterCommands()
99 {
100  G4UIparameter* param;
101 
102  //
103  // Filter commands
104  filterDir = new G4UIdirectory("/score/filter/");
105  filterDir->SetGuidance(" Scoring filter commands.");
106  //
107  fchargedCmd = new G4UIcmdWithAString("/score/filter/charged",this);
108  fchargedCmd->SetGuidance("Charged particle filter.");
109  fchargedCmd->SetParameterName("fname",false);
110  //
111  fneutralCmd = new G4UIcmdWithAString("/score/filter/neutral",this);
112  fneutralCmd->SetGuidance("Neutral particle filter.");
113  fneutralCmd->SetParameterName("fname",false);
114  //
115  fkinECmd = new G4UIcommand("/score/filter/kineticEnergy",this);
116  fkinECmd->SetGuidance("Kinetic energy filter.");
117  fkinECmd->SetGuidance("[usage] /score/filter/kineticEnergy fname Elow Ehigh unit");
118  fkinECmd->SetGuidance(" fname :(String) Filter Name ");
119  fkinECmd->SetGuidance(" Elow :(Double) Lower edge of kinetic energy");
120  fkinECmd->SetGuidance(" Ehigh :(Double) Higher edge of kinetic energy");
121  fkinECmd->SetGuidance(" unit :(String) unit of given kinetic energy");
122  param = new G4UIparameter("fname",'s',false);
123  fkinECmd->SetParameter(param);
124  param = new G4UIparameter("elow",'d',true);
125  param->SetDefaultValue("0.0");
126  fkinECmd->SetParameter(param);
127  param = new G4UIparameter("ehigh",'d',true);
128  fkinECmd->SetParameter(param);
130  param->SetDefaultValue(smax);
131  param = new G4UIparameter("unit",'s',true);
132  param->SetDefaultValue("keV");
133  fkinECmd->SetParameter(param);
134  //
135  fparticleCmd = new G4UIcommand("/score/filter/particle",this);
136  fparticleCmd->SetGuidance("Particle filter.");
137  fparticleCmd->SetGuidance("[usage] /score/filter/particle fname p0 .. pn");
138  fparticleCmd->SetGuidance(" fname :(String) Filter Name ");
139  fparticleCmd->SetGuidance(" p0 .. pn :(String) particle names");
140  param = new G4UIparameter("fname",'s',false);
141  fparticleCmd->SetParameter(param);
142  param = new G4UIparameter("particlelist",'s',false);
143  param->SetDefaultValue("");
144  fparticleCmd->SetParameter(param);
145  //
146  //
147  //
148  fparticleKinECmd = new G4UIcommand("/score/filter/particleWithKineticEnergy",this);
149  fparticleKinECmd->SetGuidance("Particle with kinetic energy filter.");
150  fparticleKinECmd->SetGuidance("[usage] /score/filter/particleWithKineticEnergy fname Elow Ehigh unit p0 .. pn");
151  fparticleKinECmd->SetGuidance(" fname :(String) Filter Name ");
152  fparticleKinECmd->SetGuidance(" Elow :(Double) Lower edge of kinetic energy");
153  fparticleKinECmd->SetGuidance(" Ehigh :(Double) Higher edge of kinetic energy");
154  fparticleKinECmd->SetGuidance(" unit :(String) unit of given kinetic energy");
155  fparticleKinECmd->SetGuidance(" p0 .. pn :(String) particle names");
156  param = new G4UIparameter("fname",'s',false);
157  fparticleKinECmd->SetParameter(param);
158  param = new G4UIparameter("elow",'d',false);
159  param->SetDefaultValue("0.0");
160  fparticleKinECmd->SetParameter(param);
161  param = new G4UIparameter("ehigh",'d',true);
162  param->SetDefaultValue(smax);
163  fparticleKinECmd->SetParameter(param);
164  param = new G4UIparameter("unit",'s',true);
165  param->SetDefaultValue("keV");
166  fparticleKinECmd->SetParameter(param);
167  param = new G4UIparameter("particlelist",'s',false);
168  param->SetDefaultValue("");
169  fparticleKinECmd->SetParameter(param);
170  //
171  //
172 }
173 
175 {
176  delete quantityDir;
177  delete qTouchCmd;
178  delete qGetUnitCmd;
179  delete qSetUnitCmd;
180 
181  //
182  delete qCellChgCmd;
183  delete qCellFluxCmd;
184  delete qPassCellFluxCmd;
185  delete qeDepCmd;
186  delete qdoseDepCmd;
187  delete qnOfStepCmd;
188  delete qnOfSecondaryCmd;
189  //
190  delete qTrackLengthCmd;
191  delete qPassCellCurrCmd;
192  delete qPassTrackLengthCmd;
193  delete qFlatSurfCurrCmd;
194  delete qFlatSurfFluxCmd;
195 // delete qSphereSurfCurrCmd;
196 // delete qSphereSurfFluxCmd;
197 // delete qCylSurfCurrCmd;
198 // delete qCylSurfFluxCmd;
199  delete qNofCollisionCmd;
200  delete qPopulationCmd;
201  delete qTrackCountCmd;
202  delete qTerminationCmd;
203  delete qMinKinEAtGeneCmd;
204  //
205  delete qStepCheckerCmd;
206  //
207  delete filterDir;
208  delete fchargedCmd;
209  delete fneutralCmd;
210  delete fkinECmd;
211  delete fparticleCmd;
212  delete fparticleKinECmd;
213 }
214 
216 {
217  //
218  // Get Current Mesh
219  //
220  G4VScoringMesh* mesh = fSMan->GetCurrentMesh();
221  if(!mesh)
222  {
223  G4cerr << "ERROR : No mesh is currently open. Open/create a mesh first. Command ignored." << G4endl;
224  return;
225  }
226  // Tokens
227  G4TokenVec token;
228  FillTokenVec(newVal,token);
229  //
230  // Commands for Current Mesh
231  if(command==qTouchCmd) {
232  mesh->SetCurrentPrimitiveScorer(newVal);
233  } else if(command == qGetUnitCmd ){
234  G4cout << "Unit: "<< mesh->GetCurrentPSUnit() <<G4endl;
235  } else if(command == qSetUnitCmd ){
236  mesh->SetCurrentPSUnit(newVal);
237  } else if(command== qCellChgCmd) {
238  if ( CheckMeshPS(mesh,token[0]) ){
239  G4PSCellCharge3D* ps = new G4PSCellCharge3D(token[0]);
240  ps->SetUnit(token[1]);
241  mesh->SetPrimitiveScorer(ps);
242  }
243  } else if(command== qCellFluxCmd) {
244  if ( CheckMeshPS(mesh,token[0]) ){
245  if( mesh->GetShape()==boxMesh ) {
246  G4PSCellFlux3D* ps = new G4PSCellFlux3D(token[0]);
247  ps->SetUnit(token[1]);
248  mesh->SetPrimitiveScorer(ps);
249  } else if( mesh->GetShape()==cylinderMesh ) {
251  new G4PSCellFluxForCylinder3D(token[0]);
252  ps->SetUnit(token[1]);
253  G4ThreeVector msize = mesh->GetSize(); // gevin in R Z N/A
254  ps->SetCylinderSize(msize[0],msize[1]); // given in dr dz
255  G4int nSeg[3];
256  mesh->GetNumberOfSegments(nSeg);
257  ps->SetNumberOfSegments(nSeg);
258  mesh->SetPrimitiveScorer(ps);
259  }
260  }
261  } else if(command== qPassCellFluxCmd) {
262  if ( CheckMeshPS(mesh,token[0]) ){
263  if( mesh->GetShape()==boxMesh ) {
264  G4PSPassageCellFlux3D* ps = new G4PSPassageCellFlux3D(token[0]);
265  ps->SetUnit(token[1]);
266  mesh->SetPrimitiveScorer(ps);
267  } else if( mesh->GetShape()==cylinderMesh ) {
269  new G4PSPassageCellFluxForCylinder3D(token[0]);
270  ps->SetUnit(token[1]);
271  G4ThreeVector msize = mesh->GetSize(); // gevin in R Z N/A
272  ps->SetCylinderSize(msize[0],msize[1]); // given in dr dz
273  G4int nSeg[3];
274  mesh->GetNumberOfSegments(nSeg);
275  ps->SetNumberOfSegments(nSeg);
276  mesh->SetPrimitiveScorer(ps);
277  }
278  }
279  } else if(command==qeDepCmd) {
280  if ( CheckMeshPS(mesh,token[0]) ){
281  G4PSEnergyDeposit3D* ps =new G4PSEnergyDeposit3D(token[0]);
282  ps->SetUnit(token[1]);
283  mesh->SetPrimitiveScorer(ps);
284  }
285  } else if(command== qdoseDepCmd) {
286  if ( CheckMeshPS(mesh,token[0]) ){
287  if( mesh->GetShape()==boxMesh ) {
288  G4PSDoseDeposit3D* ps = new G4PSDoseDeposit3D(token[0]);
289  ps->SetUnit(token[1]);
290  mesh->SetPrimitiveScorer(ps);
291  } else if( mesh->GetShape()==cylinderMesh ) {
293  new G4PSDoseDepositForCylinder3D(token[0]);
294  ps->SetUnit(token[1]);
295  G4ThreeVector msize = mesh->GetSize(); // gevin in R Z N/A
296  ps->SetCylinderSize(msize[0],msize[1]); // given in dr dz
297  G4int nSeg[3];
298  mesh->GetNumberOfSegments(nSeg);
299  ps->SetNumberOfSegments(nSeg);
300  mesh->SetPrimitiveScorer(ps);
301  }
302  }
303  } else if(command== qnOfStepCmd) {
304  if ( CheckMeshPS(mesh,token[0]) ){
305  G4PSNofStep3D* ps = new G4PSNofStep3D(token[0]);
306  ps->SetBoundaryFlag(StoB(token[1]));
307  mesh->SetPrimitiveScorer(ps);
308  }
309  } else if(command== qnOfSecondaryCmd) {
310  if ( CheckMeshPS(mesh,token[0]) ){
311  G4PSNofSecondary3D* ps =new G4PSNofSecondary3D(token[0]);
312  mesh->SetPrimitiveScorer(ps);
313  }
314  } else if(command== qTrackLengthCmd) {
315  if ( CheckMeshPS(mesh,token[0]) ){
316  G4PSTrackLength3D* ps = new G4PSTrackLength3D(token[0]);
317  ps->Weighted(StoB(token[1]));
318  ps->MultiplyKineticEnergy(StoB(token[2]));
319  ps->DivideByVelocity(StoB(token[3]));
320  ps->SetUnit(token[4]);
321  mesh->SetPrimitiveScorer(ps);
322  }
323  } else if(command== qPassCellCurrCmd){
324  if( CheckMeshPS(mesh,token[0]) ) {
326  ps->Weighted(StoB(token[1]));
327  mesh->SetPrimitiveScorer(ps);
328  }
329  } else if(command== qPassTrackLengthCmd){
330  if( CheckMeshPS(mesh,token[0]) ) {
332  ps->Weighted(StoB(token[1]));
333  ps->SetUnit(token[2]);
334  mesh->SetPrimitiveScorer(ps);
335  }
336  } else if(command== qFlatSurfCurrCmd){
337  if( CheckMeshPS(mesh,token[0])) {
339  new G4PSFlatSurfaceCurrent3D(token[0],StoI(token[1]));
340  ps->Weighted(StoB(token[2]));
341  ps->DivideByArea(StoB(token[3]));
342  if ( StoB(token[3]) ){
343  ps->SetUnit(token[4]);
344  }else{
345  ps->SetUnit("");
346  }
347  mesh->SetPrimitiveScorer(ps);
348  }
349  } else if(command== qFlatSurfFluxCmd){
350  if( CheckMeshPS(mesh, token[0] )) {
351  G4PSFlatSurfaceFlux3D* ps = new G4PSFlatSurfaceFlux3D(token[0],StoI(token[1]));
352  ps->Weighted(StoB(token[2]));
353  ps->DivideByArea(StoB(token[3]));
354  if ( StoB(token[3]) ){
355  ps->SetUnit(token[4]);
356  }else{
357  ps->SetUnit("");
358  }
359  mesh->SetPrimitiveScorer(ps);
360  }
361 // } else if(command== qSphereSurfCurrCmd){
362 // if( CheckMeshPS(mesh, token[0] )) {
363 // G4PSSphereSurfaceCurrent3D* ps =
364 // new G4PSSphereSurfaceCurrent3D(token[0],StoI(token[1]));
365 // ps->Weighted(StoB(token[2]));
366 // ps->DivideByArea(StoB(token[3]));
367 // if ( StoB(token[3]) ){
368 // ps->SetUnit(token[4]);
369 // }else{
370 // ps->SetUnit("");
371 // }
372 // mesh->SetPrimitiveScorer(ps);
373 // }
374 // } else if(command== qSphereSurfFluxCmd){
375 // if( CheckMeshPS(mesh,token[0])) {
376 // G4PSSphereSurfaceFlux3D* ps = new G4PSSphereSurfaceFlux3D(token[0], StoI(token[1]));
377 // ps->Weighted(StoB(token[2]));
378 // ps->DivideByArea(StoB(token[3]));
379 // if ( StoB(token[3]) ){
380 // ps->SetUnit(token[4]);
381 // }else{
382 // ps->SetUnit("");
383 // }
384 // mesh->SetPrimitiveScorer(ps);
385 // }
386 // } else if(command== qCylSurfCurrCmd){
387 // if( CheckMeshPS(mesh, token[0] ) ) {
388 // G4PSCylinderSurfaceCurrent3D* ps =
389 // new G4PSCylinderSurfaceCurrent3D(token[0],StoI(token[1]));
390 // ps->Weighted(StoB(token[2]));
391 // ps->DivideByArea(StoB(token[3]));
392 // if ( StoB(token[3]) ){
393 // ps->SetUnit(token[4]);
394 // }else{
395 // ps->SetUnit("");
396 // }
397 // ps->SetUnit(token[4]);
398 // mesh->SetPrimitiveScorer(ps);
399 // }
400 // } else if(command== qCylSurfFluxCmd){
401 // if( CheckMeshPS(mesh, token[0] ) {
402 // G4PSCylinerSurfaceFlux3D* ps =new G4PSCylinderSurfaceFlux3D(token[0], StoI(token[1]));
403 // ps->Weighted(StoB(token[2]));
404 // ps->DivideByArea(StoB(token[3]));
405 // if ( StoB(token[3]) ){
406 // ps->SetUnit(token[4]);
407 // }else{
408 // ps->SetUnit("");
409 // }
410 // mesh->SetPrimitiveScorer(ps);
411 // }
412  } else if(command== qNofCollisionCmd){
413  if( CheckMeshPS(mesh,token[0])) {
414  G4PSNofCollision3D* ps =new G4PSNofCollision3D(token[0]);
415  ps->Weighted(StoB(token[1]));
416  mesh->SetPrimitiveScorer(ps);
417  }
418  } else if(command== qPopulationCmd){
419  if( CheckMeshPS(mesh,token[0]) ) {
420  G4PSPopulation3D* ps =new G4PSPopulation3D(token[0]);
421  ps->Weighted(StoB(token[1]));
422  mesh->SetPrimitiveScorer(ps);
423  }
424  } else if(command== qTrackCountCmd){
425  if( CheckMeshPS(mesh,token[0])) {
426  G4PSTrackCounter3D* ps =new G4PSTrackCounter3D(token[0],StoI(token[1]));
427  ps->Weighted(StoB(token[2]));
428  mesh->SetPrimitiveScorer(ps);
429  }
430  } else if(command== qTerminationCmd){
431  if( CheckMeshPS(mesh,token[0])) {
432  G4PSTermination3D* ps =new G4PSTermination3D(token[0]);
433  ps->Weighted(StoB(token[1]));
434  mesh->SetPrimitiveScorer(ps);
435  }
436 
437  } else if(command== qMinKinEAtGeneCmd){
438  if( CheckMeshPS(mesh,token[0]) ){
440  ps->SetUnit(token[1]);
441  mesh->SetPrimitiveScorer(ps);
442  }
443  } else if(command== qStepCheckerCmd){
444  if( CheckMeshPS(mesh,token[0]) ){
445  G4PSStepChecker3D* ps =new G4PSStepChecker3D(token[0]);
446  mesh->SetPrimitiveScorer(ps);
447  }
448 
449  //
450  // Filters
451  //
452  }else if(command== fchargedCmd){
453  if(!mesh->IsCurrentPrimitiveScorerNull()) {
454  mesh->SetFilter(new G4SDChargedFilter(token[0]));
455  } else {
456  G4cout << "WARNING[" << fchargedCmd->GetCommandPath()
457  << "] : Current quantity is not set. Set or touch a quantity first." << G4endl;
458  }
459  }else if(command== fneutralCmd){
460  if(!mesh->IsCurrentPrimitiveScorerNull()) {
461  mesh->SetFilter(new G4SDNeutralFilter(token[0]));
462  } else {
463  G4cout << "WARNING[" << fneutralCmd->GetCommandPath()
464  << "] : Current quantity is not set. Set or touch a quantity first." << G4endl;
465  }
466  }else if(command== fkinECmd){
467  if(!mesh->IsCurrentPrimitiveScorerNull()) {
468  G4String& name = token[0];
469  G4double elow = StoD(token[1]);
470  G4double ehigh = StoD(token[2]);
471  G4double unitVal = G4UnitDefinition::GetValueOf(token[3]);
472  mesh->SetFilter(new G4SDKineticEnergyFilter(name,elow*unitVal,ehigh*unitVal));
473  } else {
474  G4cout << "WARNING[" << fkinECmd->GetCommandPath()
475  << "] : Current quantity is not set. Set or touch a quantity first." << G4endl;
476  }
477  }else if(command== fparticleKinECmd){
478  if(!mesh->IsCurrentPrimitiveScorerNull()) {
479  FParticleWithEnergyCommand(mesh,token);
480  } else {
481  G4cout << "WARNING[" << fparticleKinECmd->GetCommandPath()
482  << "] : Current quantity is not set. Set or touch a quantity first." << G4endl;
483  }
484  } else if(command==fparticleCmd) {
485  if(!mesh->IsCurrentPrimitiveScorerNull()) {
486  FParticleCommand(mesh,token);
487  } else {
488  G4cout << "WARNING[" << fparticleCmd->GetCommandPath()
489  << "] : Current quantity is not set. Set or touch a quantity first." << G4endl;
490  }
491  }
492 }
493 
495 {
496  G4String val;
497 
498  return val;
499 }
500 
502 
503  G4Tokenizer next(newValues);
504  G4String val;
505  while ( !(val = next()).isNull() ) {
506  token.push_back(val);
507  }
508 }
509 
510 
512  //
513  // Filter name
514  G4String name = token[0];
515  //
516  // particle list
517  std::vector<G4String> pnames;
518  for ( G4int i = 1; i<(G4int)token.size(); i++){
519  pnames.push_back(token[i]);
520  }
521  //
522  // Attach Filter
523  mesh->SetFilter(new G4SDParticleFilter(name,pnames));
524 }
525 
527  G4String& name = token[0];
528  G4double elow = StoD(token[1]);
529  G4double ehigh= StoD(token[2]);
530  G4double unitVal = G4UnitDefinition::GetValueOf(token[3]);
532  new G4SDParticleWithEnergyFilter(name,elow*unitVal,ehigh*unitVal);
533  for ( G4int i = 4; i < (G4int)token.size(); i++){
534  filter->add(token[i]);
535  }
536  mesh->SetFilter(filter);
537 }
538 
540  if(!mesh->FindPrimitiveScorer(psname)) {
541  return true;
542  } else {
543  G4cout << "WARNING[" << qTouchCmd->GetCommandPath()
544  << "] : Quantity name, \"" << psname << "\", is already existing." << G4endl;
546  return false;
547  }
548 }
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void GetNumberOfSegments(G4int nSegment[3])
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:152
pid_t filter
Definition: tracer.cxx:30
void FParticleCommand(G4VScoringMesh *mesh, G4TokenVec &token)
void SetCylinderSize(G4double dr, G4double dz)
std::vector< G4String > G4TokenVec
G4bool FindPrimitiveScorer(const G4String &psname)
void SetPrimitiveScorer(G4VPrimitiveScorer *ps)
void SetCurrentPrimitiveScorer(const G4String &name)
virtual void SetUnit(const G4String &unit)
void SetDefaultValue(const char *theDefaultValue)
void FillTokenVec(G4String newValues, G4TokenVec &token)
const XML_Char * name
void SetCylinderSize(G4double dr, G4double dz)
void SetBoundaryFlag(G4bool flg=true)
Definition: G4PSNofStep.hh:70
void SetFilter(G4VSDFilter *filter)
G4ThreeVector GetSize() const
int G4int
Definition: G4Types.hh:78
void DivideByArea(G4bool flg=true)
virtual void SetUnit(const G4String &unit)
const G4int smax
void Weighted(G4bool flg=true)
void Weighted(G4bool flg=true)
G4bool StoB(G4String s)
void SetNullToCurrentPrimitiveScorer()
G4bool IsCurrentPrimitiveScorerNull()
virtual void SetUnit(const G4String &unit)
void SetCylinderSize(G4double dr, G4double dz)
static G4double GetValueOf(const G4String &)
void Weighted(G4bool flg=true)
G4GLOB_DLL std::ostream G4cout
void Weighted(G4bool flg=true)
MeshShape GetShape() const
bool G4bool
Definition: G4Types.hh:79
G4ScoreQuantityMessenger(G4ScoringManager *SManager)
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:161
virtual void SetUnit(const G4String &unit)
G4bool CheckMeshPS(G4VScoringMesh *mesh, G4String &psname)
G4String GetCurrentValue(G4UIcommand *)
const G4String & GetCommandPath() const
Definition: G4UIcommand.hh:139
G4String DtoS(G4double a)
void SetCurrentPSUnit(const G4String &unit)
virtual void SetUnit(const G4String &unit)
G4int StoI(G4String s)
void Weighted(G4bool flg=true)
G4VScoringMesh * GetCurrentMesh() const
void DivideByArea(G4bool flg=true)
G4double StoD(G4String s)
void SetNewValue(G4UIcommand *command, G4String newValues)
void FParticleWithEnergyCommand(G4VScoringMesh *mesh, G4TokenVec &token)
void MultiplyKineticEnergy(G4bool flg=true)
void Weighted(G4bool flg=true)
void Weighted(G4bool flg=true)
#define G4endl
Definition: G4ios.hh:61
G4String GetCurrentPSUnit()
virtual void SetUnit(const G4String &unit)
virtual void SetUnit(const G4String &unit)
double G4double
Definition: G4Types.hh:76
virtual void SetUnit(const G4String &unit)
virtual void SetUnit(const G4String &unit)
void Weighted(G4bool flg=true)
#define DBL_MAX
Definition: templates.hh:83
void Weighted(G4bool flg=true)
void add(const G4String &particleName)
void DivideByVelocity(G4bool flg=true)
virtual void SetUnit(const G4String &unit)
G4GLOB_DLL std::ostream G4cerr