Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4ScoreQuantityMessenger.cc
이 파일의 문서화 페이지로 가기
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 108220 2018-01-19 15:06:38Z 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 {
96 }
97 
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);
158  param = new G4UIparameter("elow",'d',false);
159  param->SetDefaultValue("0.0");
161  param = new G4UIparameter("ehigh",'d',true);
162  param->SetDefaultValue(smax);
164  param = new G4UIparameter("unit",'s',true);
165  param->SetDefaultValue("keV");
167  param = new G4UIparameter("particlelist",'s',false);
168  param->SetDefaultValue("");
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 {
218  //
219  // Get Current Mesh
220  //
222  if(!mesh)
223  {
224  ed << "ERROR : No mesh is currently open. Open/create a mesh first. Command ignored.";
225  command->CommandFailed(ed);
226  return;
227  }
228  // Tokens
229  G4TokenVec token;
230  FillTokenVec(newVal,token);
231  //
232  // Commands for Current Mesh
233  if(command==qTouchCmd) {
234  mesh->SetCurrentPrimitiveScorer(newVal);
235  } else if(command == qGetUnitCmd ){
236  G4cout << "Unit: "<< mesh->GetCurrentPSUnit() <<G4endl;
237  } else if(command == qSetUnitCmd ){
238  mesh->SetCurrentPSUnit(newVal);
239  } else if(command== qCellChgCmd) {
240  if ( CheckMeshPS(mesh,token[0],command) ){
241  G4PSCellCharge3D* ps = new G4PSCellCharge3D(token[0]);
242  ps->SetUnit(token[1]);
243  mesh->SetPrimitiveScorer(ps);
244  }
245  } else if(command== qCellFluxCmd) {
246  if ( CheckMeshPS(mesh,token[0],command) ){
247  if( mesh->GetShape()==boxMesh ) {
248  G4PSCellFlux3D* ps = new G4PSCellFlux3D(token[0]);
249  ps->SetUnit(token[1]);
250  mesh->SetPrimitiveScorer(ps);
251  } else if( mesh->GetShape()==cylinderMesh ) {
253  new G4PSCellFluxForCylinder3D(token[0]);
254  ps->SetUnit(token[1]);
255  G4ThreeVector msize = mesh->GetSize(); // gevin in R Z N/A
256  ps->SetCylinderSize(msize[0],msize[1]); // given in dr dz
257  G4int nSeg[3];
258  mesh->GetNumberOfSegments(nSeg);
259  ps->SetNumberOfSegments(nSeg);
260  mesh->SetPrimitiveScorer(ps);
261  }
262  }
263  } else if(command== qPassCellFluxCmd) {
264  if ( CheckMeshPS(mesh,token[0],command) ){
265  if( mesh->GetShape()==boxMesh ) {
267  ps->SetUnit(token[1]);
268  mesh->SetPrimitiveScorer(ps);
269  } else if( mesh->GetShape()==cylinderMesh ) {
271  new G4PSPassageCellFluxForCylinder3D(token[0]);
272  ps->SetUnit(token[1]);
273  G4ThreeVector msize = mesh->GetSize(); // gevin in R Z N/A
274  ps->SetCylinderSize(msize[0],msize[1]); // given in dr dz
275  G4int nSeg[3];
276  mesh->GetNumberOfSegments(nSeg);
277  ps->SetNumberOfSegments(nSeg);
278  mesh->SetPrimitiveScorer(ps);
279  }
280  }
281  } else if(command==qeDepCmd) {
282  if ( CheckMeshPS(mesh,token[0],command) ){
284  ps->SetUnit(token[1]);
285  mesh->SetPrimitiveScorer(ps);
286  }
287  } else if(command== qdoseDepCmd) {
288  if ( CheckMeshPS(mesh,token[0],command) ){
289  if( mesh->GetShape()==boxMesh ) {
290  G4PSDoseDeposit3D* ps = new G4PSDoseDeposit3D(token[0]);
291  ps->SetUnit(token[1]);
292  mesh->SetPrimitiveScorer(ps);
293  } else if( mesh->GetShape()==cylinderMesh ) {
295  new G4PSDoseDepositForCylinder3D(token[0]);
296  ps->SetUnit(token[1]);
297  G4ThreeVector msize = mesh->GetSize(); // gevin in R Z N/A
298  ps->SetCylinderSize(msize[0],msize[1]); // given in dr dz
299  G4int nSeg[3];
300  mesh->GetNumberOfSegments(nSeg);
301  ps->SetNumberOfSegments(nSeg);
302  mesh->SetPrimitiveScorer(ps);
303  }
304  }
305  } else if(command== qnOfStepCmd) {
306  if ( CheckMeshPS(mesh,token[0],command) ){
307  G4PSNofStep3D* ps = new G4PSNofStep3D(token[0]);
308  ps->SetBoundaryFlag(StoB(token[1]));
309  mesh->SetPrimitiveScorer(ps);
310  }
311  } else if(command== qnOfSecondaryCmd) {
312  if ( CheckMeshPS(mesh,token[0],command) ){
313  G4PSNofSecondary3D* ps =new G4PSNofSecondary3D(token[0]);
314  mesh->SetPrimitiveScorer(ps);
315  }
316  } else if(command== qTrackLengthCmd) {
317  if ( CheckMeshPS(mesh,token[0],command) ){
318  G4PSTrackLength3D* ps = new G4PSTrackLength3D(token[0]);
319  ps->Weighted(StoB(token[1]));
320  ps->MultiplyKineticEnergy(StoB(token[2]));
321  ps->DivideByVelocity(StoB(token[3]));
322  ps->SetUnit(token[4]);
323  mesh->SetPrimitiveScorer(ps);
324  }
325  } else if(command== qPassCellCurrCmd){
326  if( CheckMeshPS(mesh,token[0],command) ) {
328  ps->Weighted(StoB(token[1]));
329  mesh->SetPrimitiveScorer(ps);
330  }
331  } else if(command== qPassTrackLengthCmd){
332  if( CheckMeshPS(mesh,token[0],command) ) {
334  ps->Weighted(StoB(token[1]));
335  ps->SetUnit(token[2]);
336  mesh->SetPrimitiveScorer(ps);
337  }
338  } else if(command== qFlatSurfCurrCmd){
339  if( CheckMeshPS(mesh,token[0],command)) {
341  new G4PSFlatSurfaceCurrent3D(token[0],StoI(token[1]));
342  ps->Weighted(StoB(token[2]));
343  ps->DivideByArea(StoB(token[3]));
344  if ( StoB(token[3]) ){
345  ps->SetUnit(token[4]);
346  }else{
347  ps->SetUnit("");
348  }
349  mesh->SetPrimitiveScorer(ps);
350  }
351  } else if(command== qFlatSurfFluxCmd){
352  if( CheckMeshPS(mesh, token[0],command )) {
353  G4PSFlatSurfaceFlux3D* ps = new G4PSFlatSurfaceFlux3D(token[0],StoI(token[1]));
354  ps->Weighted(StoB(token[2]));
355  ps->DivideByArea(StoB(token[3]));
356  if ( StoB(token[3]) ){
357  ps->SetUnit(token[4]);
358  }else{
359  ps->SetUnit("");
360  }
361  mesh->SetPrimitiveScorer(ps);
362  }
363 // } else if(command== qSphereSurfCurrCmd){
364 // if( CheckMeshPS(mesh, token[0],command )) {
365 // G4PSSphereSurfaceCurrent3D* ps =
366 // new G4PSSphereSurfaceCurrent3D(token[0],StoI(token[1]));
367 // ps->Weighted(StoB(token[2]));
368 // ps->DivideByArea(StoB(token[3]));
369 // if ( StoB(token[3]) ){
370 // ps->SetUnit(token[4]);
371 // }else{
372 // ps->SetUnit("");
373 // }
374 // mesh->SetPrimitiveScorer(ps);
375 // }
376 // } else if(command== qSphereSurfFluxCmd){
377 // if( CheckMeshPS(mesh,token[0],command)) {
378 // G4PSSphereSurfaceFlux3D* ps = new G4PSSphereSurfaceFlux3D(token[0], StoI(token[1]));
379 // ps->Weighted(StoB(token[2]));
380 // ps->DivideByArea(StoB(token[3]));
381 // if ( StoB(token[3]) ){
382 // ps->SetUnit(token[4]);
383 // }else{
384 // ps->SetUnit("");
385 // }
386 // mesh->SetPrimitiveScorer(ps);
387 // }
388 // } else if(command== qCylSurfCurrCmd){
389 // if( CheckMeshPS(mesh, token[0],command ) ) {
390 // G4PSCylinderSurfaceCurrent3D* ps =
391 // new G4PSCylinderSurfaceCurrent3D(token[0],StoI(token[1]));
392 // ps->Weighted(StoB(token[2]));
393 // ps->DivideByArea(StoB(token[3]));
394 // if ( StoB(token[3]) ){
395 // ps->SetUnit(token[4]);
396 // }else{
397 // ps->SetUnit("");
398 // }
399 // ps->SetUnit(token[4]);
400 // mesh->SetPrimitiveScorer(ps);
401 // }
402 // } else if(command== qCylSurfFluxCmd){
403 // if( CheckMeshPS(mesh, token[0],command ) {
404 // G4PSCylinerSurfaceFlux3D* ps =new G4PSCylinderSurfaceFlux3D(token[0], StoI(token[1]));
405 // ps->Weighted(StoB(token[2]));
406 // ps->DivideByArea(StoB(token[3]));
407 // if ( StoB(token[3]) ){
408 // ps->SetUnit(token[4]);
409 // }else{
410 // ps->SetUnit("");
411 // }
412 // mesh->SetPrimitiveScorer(ps);
413 // }
414  } else if(command== qNofCollisionCmd){
415  if( CheckMeshPS(mesh,token[0],command)) {
416  G4PSNofCollision3D* ps =new G4PSNofCollision3D(token[0]);
417  ps->Weighted(StoB(token[1]));
418  mesh->SetPrimitiveScorer(ps);
419  }
420  } else if(command== qPopulationCmd){
421  if( CheckMeshPS(mesh,token[0],command) ) {
422  G4PSPopulation3D* ps =new G4PSPopulation3D(token[0]);
423  ps->Weighted(StoB(token[1]));
424  mesh->SetPrimitiveScorer(ps);
425  }
426  } else if(command== qTrackCountCmd){
427  if( CheckMeshPS(mesh,token[0],command)) {
428  G4PSTrackCounter3D* ps =new G4PSTrackCounter3D(token[0],StoI(token[1]));
429  ps->Weighted(StoB(token[2]));
430  mesh->SetPrimitiveScorer(ps);
431  }
432  } else if(command== qTerminationCmd){
433  if( CheckMeshPS(mesh,token[0],command)) {
434  G4PSTermination3D* ps =new G4PSTermination3D(token[0]);
435  ps->Weighted(StoB(token[1]));
436  mesh->SetPrimitiveScorer(ps);
437  }
438 
439  } else if(command== qMinKinEAtGeneCmd){
440  if( CheckMeshPS(mesh,token[0],command) ){
442  ps->SetUnit(token[1]);
443  mesh->SetPrimitiveScorer(ps);
444  }
445  } else if(command== qStepCheckerCmd){
446  if( CheckMeshPS(mesh,token[0],command) ){
447  G4PSStepChecker3D* ps =new G4PSStepChecker3D(token[0]);
448  mesh->SetPrimitiveScorer(ps);
449  }
450 
451  //
452  // Filters
453  //
454  }else if(command== fchargedCmd){
455  if(!mesh->IsCurrentPrimitiveScorerNull()) {
456  mesh->SetFilter(new G4SDChargedFilter(token[0]));
457  } else {
458  ed << "WARNING[" << fchargedCmd->GetCommandPath()
459  << "] : Current quantity is not set. Set or touch a quantity first.";
460  command->CommandFailed(ed);
461  }
462  }else if(command== fneutralCmd){
463  if(!mesh->IsCurrentPrimitiveScorerNull()) {
464  mesh->SetFilter(new G4SDNeutralFilter(token[0]));
465  } else {
466  ed << "WARNING[" << fneutralCmd->GetCommandPath()
467  << "] : Current quantity is not set. Set or touch a quantity first.";
468  command->CommandFailed(ed);
469  }
470  }else if(command== fkinECmd){
471  if(!mesh->IsCurrentPrimitiveScorerNull()) {
472  G4String& name = token[0];
473  G4double elow = StoD(token[1]);
474  G4double ehigh = StoD(token[2]);
475  G4double unitVal = G4UnitDefinition::GetValueOf(token[3]);
476  mesh->SetFilter(new G4SDKineticEnergyFilter(name,elow*unitVal,ehigh*unitVal));
477  } else {
478  ed << "WARNING[" << fkinECmd->GetCommandPath()
479  << "] : Current quantity is not set. Set or touch a quantity first.";
480  command->CommandFailed(ed);
481  }
482  }else if(command== fparticleKinECmd){
483  if(!mesh->IsCurrentPrimitiveScorerNull()) {
484  FParticleWithEnergyCommand(mesh,token);
485  } else {
486  ed << "WARNING[" << fparticleKinECmd->GetCommandPath()
487  << "] : Current quantity is not set. Set or touch a quantity first.";
488  command->CommandFailed(ed);
489  }
490  } else if(command==fparticleCmd) {
491  if(!mesh->IsCurrentPrimitiveScorerNull()) {
492  FParticleCommand(mesh,token);
493  } else {
494  ed << "WARNING[" << fparticleCmd->GetCommandPath()
495  << "] : Current quantity is not set. Set or touch a quantity first.";
496  command->CommandFailed(ed);
497  }
498  }
499 }
500 
502 {
503  G4String val;
504 
505  return val;
506 }
507 
509 
510  G4Tokenizer next(newValues);
511  G4String val;
512  while ( !(val = next()).isNull() ) { // Loop checking 12.18.2015 M.Asai
513  token.push_back(val);
514  }
515 }
516 
517 
519  //
520  // Filter name
521  G4String name = token[0];
522  //
523  // particle list
524  std::vector<G4String> pnames;
525  for ( G4int i = 1; i<(G4int)token.size(); i++){
526  pnames.push_back(token[i]);
527  }
528  //
529  // Attach Filter
530  mesh->SetFilter(new G4SDParticleFilter(name,pnames));
531 }
532 
534  G4String& name = token[0];
535  G4double elow = StoD(token[1]);
536  G4double ehigh= StoD(token[2]);
537  G4double unitVal = G4UnitDefinition::GetValueOf(token[3]);
539  new G4SDParticleWithEnergyFilter(name,elow*unitVal,ehigh*unitVal);
540  for ( G4int i = 4; i < (G4int)token.size(); i++){
541  filter->add(token[i]);
542  }
543  mesh->SetFilter(filter);
544 }
545 
547  G4UIcommand* command){
548  if(!mesh->FindPrimitiveScorer(psname)) {
549  return true;
550  } else {
552  ed << "WARNING[" << qTouchCmd->GetCommandPath()
553  << "] : Quantity name, \"" << psname << "\", is already existing.";
554  command->CommandFailed(ed);
556  return false;
557  }
558 }
G4ThreeVector GetSize() const
void Weighted(G4bool flg=true)
G4int StoI(G4String s)
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:161
const XML_Char * name
Definition: expat.h:151
void FParticleCommand(G4VScoringMesh *mesh, G4TokenVec &token)
G4String GetCurrentPSUnit()
void SetBoundaryFlag(G4bool flg=true)
Definition: G4PSNofStep.hh:70
void CommandFailed(G4int errCode, G4ExceptionDescription &ed)
Definition: G4UIcommand.hh:202
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:45
void Weighted(G4bool flg=true)
virtual void SetUnit(const G4String &unit)
G4UIcmdWithoutParameter * qGetUnitCmd
#define G4endl
Definition: G4ios.hh:61
void Weighted(G4bool flg=true)
std::vector< G4String > G4TokenVec
void SetCylinderSize(G4double dr, G4double dz)
void add(const G4String &particleName)
virtual void SetUnit(const G4String &unit)
void DivideByVelocity(G4bool flg=true)
MeshShape GetShape() const
virtual void SetUnit(const G4String &unit)
virtual void SetUnit(const G4String &unit)
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:152
G4ScoreQuantityMessenger(G4ScoringManager *SManager)
G4VScoringMesh * GetCurrentMesh() const
G4String GetCurrentValue(G4UIcommand *)
void SetDefaultValue(const char *theDefaultValue)
static constexpr double ps
Definition: G4SIunits.hh:172
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
void Weighted(G4bool flg=true)
void GetNumberOfSegments(G4int nSegment[3])
virtual void SetUnit(const G4String &unit)
virtual void SetUnit(const G4String &unit)
virtual void SetUnit(const G4String &unit)
static G4double GetValueOf(const G4String &)
void SetCylinderSize(G4double dr, G4double dz)
void SetPrimitiveScorer(G4VPrimitiveScorer *ps)
virtual void SetUnit(const G4String &unit)
virtual void SetUnit(const G4String &unit)
void Weighted(G4bool flg=true)
const G4String & GetCommandPath() const
Definition: G4UIcommand.hh:139
void SetNullToCurrentPrimitiveScorer()
void Weighted(G4bool flg=true)
void SetFilter(G4VSDFilter *filter)
int G4int
Definition: G4Types.hh:78
void SetCurrentPSUnit(const G4String &unit)
void FParticleWithEnergyCommand(G4VScoringMesh *mesh, G4TokenVec &token)
void FillTokenVec(G4String newValues, G4TokenVec &token)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void Weighted(G4bool flg=true)
G4bool CheckMeshPS(G4VScoringMesh *mesh, G4String &psname, G4UIcommand *command)
void SetCylinderSize(G4double dr, G4double dz)
G4double StoD(G4String s)
G4bool FindPrimitiveScorer(const G4String &psname)
G4GLOB_DLL std::ostream G4cout
const G4int smax
G4bool StoB(G4String s)
void SetNewValue(G4UIcommand *command, G4String newValues)
void DivideByArea(G4bool flg=true)
void MultiplyKineticEnergy(G4bool flg=true)
void Weighted(G4bool flg=true)
virtual void SetUnit(const G4String &unit)
void SetCurrentPrimitiveScorer(const G4String &name)
void Weighted(G4bool flg=true)
#define DBL_MAX
Definition: templates.hh:83
void DivideByArea(G4bool flg=true)
G4String DtoS(G4double a)
G4bool IsCurrentPrimitiveScorerNull()