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

#include <G4ScoreQuantityMessenger.hh>

+ Inheritance diagram for G4ScoreQuantityMessenger:

Public Member Functions

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

Protected Member Functions

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

Additional Inherited Members

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

Detailed Description

Definition at line 58 of file G4ScoreQuantityMessenger.hh.

Constructor & Destructor Documentation

◆ G4ScoreQuantityMessenger()

G4ScoreQuantityMessenger::G4ScoreQuantityMessenger ( G4ScoringManager SManager)

Definition at line 117 of file G4ScoreQuantityMessenger.cc.

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

◆ ~G4ScoreQuantityMessenger()

G4ScoreQuantityMessenger::~G4ScoreQuantityMessenger ( )

Definition at line 200 of file G4ScoreQuantityMessenger.cc.

201{
202 delete quantityDir;
203 delete qTouchCmd;
204 delete qGetUnitCmd;
205 delete qSetUnitCmd;
206
207 //
208 delete qCellChgCmd;
209 delete qCellFluxCmd;
210 delete qPassCellFluxCmd;
211 delete qeDepCmd;
212 delete qdoseDepCmd;
213 delete qnOfStepCmd;
214 delete qnOfSecondaryCmd;
215 //
216 delete qTrackLengthCmd;
217 delete qPassCellCurrCmd;
218 delete qPassTrackLengthCmd;
219 delete qFlatSurfCurrCmd;
220 delete qFlatSurfFluxCmd;
221 delete qVolFluxCmd;
222// delete qSphereSurfCurrCmd;
223// delete qSphereSurfFluxCmd;
224// delete qCylSurfCurrCmd;
225// delete qCylSurfFluxCmd;
226 delete qNofCollisionCmd;
227 delete qPopulationCmd;
228 delete qTrackCountCmd;
229 delete qTerminationCmd;
230 delete qMinKinEAtGeneCmd;
231 //
232 delete qStepCheckerCmd;
233 //
234 delete filterDir;
235 delete fchargedCmd;
236 delete fneutralCmd;
237 delete fkinECmd;
238 delete fparticleCmd;
239 delete fparticleKinECmd;
240}

Member Function Documentation

◆ CheckMeshPS()

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

Definition at line 660 of file G4ScoreQuantityMessenger.cc.

661 {
662 if(!mesh->FindPrimitiveScorer(psname)) {
663 return true;
664 } else {
666 ed << "WARNING[" << qTouchCmd->GetCommandPath()
667 << "] : Quantity name, \"" << psname << "\", is already existing.";
668 command->CommandFailed(ed);
670 return false;
671 }
672}
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
const G4String & GetCommandPath() const
Definition: G4UIcommand.hh:136
void CommandFailed(G4int errCode, G4ExceptionDescription &ed)
Definition: G4UIcommand.hh:179
void SetNullToCurrentPrimitiveScorer()
G4bool FindPrimitiveScorer(const G4String &psname)

Referenced by SetNewValue().

◆ FillTokenVec()

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

Definition at line 622 of file G4ScoreQuantityMessenger.cc.

622 {
623
624 G4Tokenizer next(newValues);
625 G4String val;
626 while ( !(val = next()).isNull() ) { // Loop checking 12.18.2015 M.Asai
627 token.push_back(val);
628 }
629}

Referenced by SetNewValue().

◆ FParticleCommand()

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

Definition at line 632 of file G4ScoreQuantityMessenger.cc.

632 {
633 //
634 // Filter name
635 G4String name = token[0];
636 //
637 // particle list
638 std::vector<G4String> pnames;
639 for ( G4int i = 1; i<(G4int)token.size(); i++){
640 pnames.push_back(token[i]);
641 }
642 //
643 // Attach Filter
644 mesh->SetFilter(new G4SDParticleFilter(name,pnames));
645}
int G4int
Definition: G4Types.hh:85
void SetFilter(G4VSDFilter *filter)
const char * name(G4int ptype)

Referenced by SetNewValue().

◆ FParticleWithEnergyCommand()

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

Definition at line 647 of file G4ScoreQuantityMessenger.cc.

647 {
648 G4String& name = token[0];
649 G4double elow = StoD(token[1]);
650 G4double ehigh= StoD(token[2]);
651 G4double unitVal = G4UnitDefinition::GetValueOf(token[3]);
653 new G4SDParticleWithEnergyFilter(name,elow*unitVal,ehigh*unitVal);
654 for ( G4int i = 4; i < (G4int)token.size(); i++){
655 filter->add(token[i]);
656 }
657 mesh->SetFilter(filter);
658}
double G4double
Definition: G4Types.hh:83
void add(const G4String &particleName)
G4double StoD(G4String s)
static G4double GetValueOf(const G4String &)

Referenced by SetNewValue().

◆ GetCurrentValue()

G4String G4ScoreQuantityMessenger::GetCurrentValue ( G4UIcommand )
virtual

Reimplemented from G4UImessenger.

Definition at line 615 of file G4ScoreQuantityMessenger.cc.

616{
617 G4String val;
618
619 return val;
620}

◆ SetNewValue()

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

Reimplemented from G4UImessenger.

Definition at line 242 of file G4ScoreQuantityMessenger.cc.

243{
244 using MeshShape = G4VScoringMesh::MeshShape;
245
247
248 //
249 // Get Current Mesh
250 //
251 G4VScoringMesh* mesh = fSMan->GetCurrentMesh();
252 if(!mesh)
253 {
254 ed << "ERROR : No mesh is currently open. Open/create a mesh first. Command ignored.";
255 command->CommandFailed(ed);
256 return;
257 }
258 // Mesh type
259 auto shape = mesh->GetShape();
260 // Tokens
261 G4TokenVec token;
262 FillTokenVec(newVal,token);
263 //
264 // Commands for Current Mesh
265 if(command==qTouchCmd) {
266 mesh->SetCurrentPrimitiveScorer(newVal);
267 } else if(command == qGetUnitCmd ){
268 G4cout << "Unit: "<< mesh->GetCurrentPSUnit() <<G4endl;
269 } else if(command == qSetUnitCmd ){
270 mesh->SetCurrentPSUnit(newVal);
271 } else if(command== qCellChgCmd) {
272 if ( CheckMeshPS(mesh,token[0],command) ){
273 G4PSCellCharge* ps = nullptr;
274 if(shape==MeshShape::realWorldLogVol || shape==MeshShape::probe)
275 { ps = new G4PSCellCharge(token[0],mesh->GetCopyNumberLevel()); }
276 else
277 { ps = new G4PSCellCharge3D(token[0]); }
278 ps->SetUnit(token[1]);
279 mesh->SetPrimitiveScorer(ps);
280 }
281 } else if(command== qCellFluxCmd) {
282 if ( CheckMeshPS(mesh,token[0],command) ){
283 G4PSCellFlux* ps = nullptr;
284 if( shape==MeshShape::box ) {
285 ps = new G4PSCellFlux3D(token[0]);
286 } else if( shape==MeshShape::cylinder ) {
288 G4ThreeVector msize = mesh->GetSize(); // gevin in R Z N/A
289 pps->SetCylinderSize(msize[0],msize[1]); // given in dr dz
290 G4int nSeg[3];
291 mesh->GetNumberOfSegments(nSeg);
292 pps->SetNumberOfSegments(nSeg);
293 ps = pps;
294 } else if(shape==MeshShape::realWorldLogVol) {
295 ed<<"Cell flux for real world volume is not yet supported. Command ignored.";
296 command->CommandFailed(ed);
297 return;
298 } else if(shape==MeshShape::probe) {
299 ps = new G4PSCellFlux(token[0]);
300 }
301 ps->SetUnit(token[1]);
302 mesh->SetPrimitiveScorer(ps);
303 }
304 } else if(command== qPassCellFluxCmd) {
305 if ( CheckMeshPS(mesh,token[0],command) ){
306 G4PSPassageCellFlux* ps = nullptr;
307 if( shape==MeshShape::box ) {
308 ps = new G4PSPassageCellFlux3D(token[0]);
309 } else if( shape==MeshShape::cylinder ) {
311 G4ThreeVector msize = mesh->GetSize(); // gevin in R Z N/A
312 pps->SetCylinderSize(msize[0],msize[1]); // given in dr dz
313 G4int nSeg[3];
314 mesh->GetNumberOfSegments(nSeg);
315 pps->SetNumberOfSegments(nSeg);
316 ps = pps;
317 } else if(shape==MeshShape::realWorldLogVol) {
318 ed<<"Passing cell flux for real world volume is not yet supported. Command ignored.";
319 command->CommandFailed(ed);
320 return;
321 } else if(shape==MeshShape::probe) {
322 ps = new G4PSPassageCellFlux(token[0]);
323 }
324 ps->SetUnit(token[1]);
325 mesh->SetPrimitiveScorer(ps);
326 }
327 } else if(command==qeDepCmd) {
328 if ( CheckMeshPS(mesh,token[0],command) ){
329 G4PSEnergyDeposit* ps = nullptr;
330 if(shape==MeshShape::realWorldLogVol || shape==MeshShape::probe)
331 { ps =new G4PSEnergyDeposit(token[0],mesh->GetCopyNumberLevel()); }
332 else
333 { ps =new G4PSEnergyDeposit3D(token[0]); }
334 ps->SetUnit(token[1]);
335 mesh->SetPrimitiveScorer(ps);
336 }
337 } else if(command== qdoseDepCmd) {
338 if ( CheckMeshPS(mesh,token[0],command) ){
339 G4PSDoseDeposit* ps = nullptr;
340 if( shape==MeshShape::box ) {
341 ps = new G4PSDoseDeposit3D(token[0]);
342 } else if( shape==MeshShape::cylinder ) {
344 pps->SetUnit(token[1]);
345 G4ThreeVector msize = mesh->GetSize(); // gevin in R Z N/A
346 pps->SetCylinderSize(msize[0],msize[1]); // given in dr dz
347 G4int nSeg[3];
348 mesh->GetNumberOfSegments(nSeg);
349 pps->SetNumberOfSegments(nSeg);
350 ps = pps;
351 } else if(shape==MeshShape::realWorldLogVol || shape==MeshShape::probe) {
352 ps = new G4PSDoseDeposit(token[0],mesh->GetCopyNumberLevel());
353 }
354 ps->SetUnit(token[1]);
355 mesh->SetPrimitiveScorer(ps);
356 }
357 } else if(command== qnOfStepCmd) {
358 if ( CheckMeshPS(mesh,token[0],command) ){
359 G4PSNofStep* ps = nullptr;
360 if(shape==MeshShape::realWorldLogVol || shape==MeshShape::probe)
361 { ps = new G4PSNofStep(token[0],mesh->GetCopyNumberLevel()); }
362 else
363 { ps = new G4PSNofStep3D(token[0]); }
364 ps->SetBoundaryFlag(StoB(token[1]));
365 mesh->SetPrimitiveScorer(ps);
366 }
367 } else if(command== qnOfSecondaryCmd) {
368 if ( CheckMeshPS(mesh,token[0],command) ){
369 G4PSNofSecondary* ps = nullptr;
370 if(shape==MeshShape::realWorldLogVol || shape==MeshShape::probe)
371 { ps = new G4PSNofSecondary(token[0],mesh->GetCopyNumberLevel()); }
372 else
373 { ps = new G4PSNofSecondary3D(token[0]); }
374 mesh->SetPrimitiveScorer(ps);
375 }
376 } else if(command== qTrackLengthCmd) {
377 if ( CheckMeshPS(mesh,token[0],command) ){
378 G4PSTrackLength* ps = nullptr;
379 if(shape==MeshShape::realWorldLogVol || shape==MeshShape::probe)
380 { ps = new G4PSTrackLength(token[0],mesh->GetCopyNumberLevel()); }
381 else
382 { ps = new G4PSTrackLength3D(token[0]); }
383 ps->Weighted(StoB(token[1]));
384 ps->MultiplyKineticEnergy(StoB(token[2]));
385 ps->DivideByVelocity(StoB(token[3]));
386 ps->SetUnit(token[4]);
387 mesh->SetPrimitiveScorer(ps);
388 }
389 } else if(command== qPassCellCurrCmd){
390 if( CheckMeshPS(mesh,token[0],command) ) {
391 G4PSPassageCellCurrent* ps = nullptr;
392 if(shape==MeshShape::realWorldLogVol || shape==MeshShape::probe)
393 { ps = new G4PSPassageCellCurrent(token[0],mesh->GetCopyNumberLevel()); }
394 else
395 { ps = new G4PSPassageCellCurrent3D(token[0]); }
396 ps->Weighted(StoB(token[1]));
397 mesh->SetPrimitiveScorer(ps);
398 }
399 } else if(command== qPassTrackLengthCmd){
400 if( CheckMeshPS(mesh,token[0],command) ) {
401 G4PSPassageTrackLength* ps = nullptr;
402 if(shape==MeshShape::realWorldLogVol || shape==MeshShape::probe)
403 { ps = new G4PSPassageTrackLength(token[0],mesh->GetCopyNumberLevel()); }
404 else
405 { ps = new G4PSPassageTrackLength3D(token[0]); }
406 ps->Weighted(StoB(token[1]));
407 ps->SetUnit(token[2]);
408 mesh->SetPrimitiveScorer(ps);
409 }
410 } else if(command== qFlatSurfCurrCmd){
411 if( CheckMeshPS(mesh,token[0],command)) {
412 G4PSFlatSurfaceCurrent* ps = nullptr;
413 if(shape==MeshShape::realWorldLogVol || shape==MeshShape::probe)
414 { ps = new G4PSFlatSurfaceCurrent(token[0],StoI(token[1]),mesh->GetCopyNumberLevel()); }
415 else
416 { ps = new G4PSFlatSurfaceCurrent3D(token[0],StoI(token[1])); }
417 ps->Weighted(StoB(token[2]));
418 ps->DivideByArea(StoB(token[3]));
419 if ( StoB(token[3]) ){
420 ps->SetUnit(token[4]);
421 }else{
422 ps->SetUnit("");
423 }
424 mesh->SetPrimitiveScorer(ps);
425 }
426 } else if(command== qFlatSurfFluxCmd){
427 if( CheckMeshPS(mesh, token[0],command )) {
428 G4PSFlatSurfaceFlux* ps = nullptr;
429 if(shape==MeshShape::realWorldLogVol || shape==MeshShape::probe)
430 { ps = new G4PSFlatSurfaceFlux(token[0],StoI(token[1]),mesh->GetCopyNumberLevel()); }
431 else
432 { ps = new G4PSFlatSurfaceFlux3D(token[0],StoI(token[1])); }
433 ps->Weighted(StoB(token[2]));
434 ps->DivideByArea(StoB(token[3]));
435 if ( StoB(token[3]) ){
436 ps->SetUnit(token[4]);
437 }else{
438 ps->SetUnit("");
439 }
440 mesh->SetPrimitiveScorer(ps);
441 }
442 } else if(command== qVolFluxCmd) {
443 if( CheckMeshPS(mesh, token[0],command )) {
444 G4PSVolumeFlux* ps = nullptr;
445 if(shape==MeshShape::realWorldLogVol || shape==MeshShape::probe)
446 { ps = new G4PSVolumeFlux(token[0],StoI(token[2]),mesh->GetCopyNumberLevel()); }
447 else
448 { ps = new G4PSVolumeFlux3D(token[0],StoI(token[2])); }
449 ps->SetDivCos(StoI(token[1]));
450 mesh->SetPrimitiveScorer(ps);
451 }
452
453// } else if(command== qSphereSurfCurrCmd){
454// if( CheckMeshPS(mesh, token[0],command )) {
455// G4PSSphereSurfaceCurrent3D* ps =
456// new G4PSSphereSurfaceCurrent3D(token[0],StoI(token[1]));
457// ps->Weighted(StoB(token[2]));
458// ps->DivideByArea(StoB(token[3]));
459// if ( StoB(token[3]) ){
460// ps->SetUnit(token[4]);
461// }else{
462// ps->SetUnit("");
463// }
464// mesh->SetPrimitiveScorer(ps);
465// }
466// } else if(command== qSphereSurfFluxCmd){
467// if( CheckMeshPS(mesh,token[0],command)) {
468// G4PSSphereSurfaceFlux3D* ps = new G4PSSphereSurfaceFlux3D(token[0], StoI(token[1]));
469// ps->Weighted(StoB(token[2]));
470// ps->DivideByArea(StoB(token[3]));
471// if ( StoB(token[3]) ){
472// ps->SetUnit(token[4]);
473// }else{
474// ps->SetUnit("");
475// }
476// mesh->SetPrimitiveScorer(ps);
477// }
478// } else if(command== qCylSurfCurrCmd){
479// if( CheckMeshPS(mesh, token[0],command ) ) {
480// G4PSCylinderSurfaceCurrent3D* ps =
481// new G4PSCylinderSurfaceCurrent3D(token[0],StoI(token[1]));
482// ps->Weighted(StoB(token[2]));
483// ps->DivideByArea(StoB(token[3]));
484// if ( StoB(token[3]) ){
485// ps->SetUnit(token[4]);
486// }else{
487// ps->SetUnit("");
488// }
489// ps->SetUnit(token[4]);
490// mesh->SetPrimitiveScorer(ps);
491// }
492// } else if(command== qCylSurfFluxCmd){
493// if( CheckMeshPS(mesh, token[0],command ) {
494// G4PSCylinerSurfaceFlux3D* ps =new G4PSCylinderSurfaceFlux3D(token[0], StoI(token[1]));
495// ps->Weighted(StoB(token[2]));
496// ps->DivideByArea(StoB(token[3]));
497// if ( StoB(token[3]) ){
498// ps->SetUnit(token[4]);
499// }else{
500// ps->SetUnit("");
501// }
502// mesh->SetPrimitiveScorer(ps);
503// }
504 } else if(command== qNofCollisionCmd){
505 if( CheckMeshPS(mesh,token[0],command)) {
506 G4PSNofCollision* ps = nullptr;
507 if(shape==MeshShape::realWorldLogVol || shape==MeshShape::probe)
508 { ps = new G4PSNofCollision3D(token[0],mesh->GetCopyNumberLevel()); }
509 else
510 { ps = new G4PSNofCollision3D(token[0]); }
511 ps->Weighted(StoB(token[1]));
512 mesh->SetPrimitiveScorer(ps);
513 }
514 } else if(command== qPopulationCmd){
515 if( CheckMeshPS(mesh,token[0],command) ) {
516 G4PSPopulation* ps = nullptr;
517 if(shape==MeshShape::realWorldLogVol || shape==MeshShape::probe)
518 { ps = new G4PSPopulation(token[0],mesh->GetCopyNumberLevel()); }
519 else
520 { ps = new G4PSPopulation3D(token[0]); }
521 ps->Weighted(StoB(token[1]));
522 mesh->SetPrimitiveScorer(ps);
523 }
524 } else if(command== qTrackCountCmd){
525 if( CheckMeshPS(mesh,token[0],command)) {
526 G4PSTrackCounter* ps = nullptr;
527 if(shape==MeshShape::realWorldLogVol || shape==MeshShape::probe)
528 { ps = new G4PSTrackCounter(token[0],StoI(token[1]),mesh->GetCopyNumberLevel()); }
529 else
530 { ps = new G4PSTrackCounter3D(token[0],StoI(token[1])); }
531 ps->Weighted(StoB(token[2]));
532 mesh->SetPrimitiveScorer(ps);
533 }
534 } else if(command== qTerminationCmd){
535 if( CheckMeshPS(mesh,token[0],command)) {
536 G4PSTermination* ps = nullptr;
537 if(shape==MeshShape::realWorldLogVol || shape==MeshShape::probe)
538 { ps = new G4PSTermination(token[0],mesh->GetCopyNumberLevel()); }
539 else
540 { ps = new G4PSTermination3D(token[0]); }
541 ps->Weighted(StoB(token[1]));
542 mesh->SetPrimitiveScorer(ps);
543 }
544
545 } else if(command== qMinKinEAtGeneCmd){
546 if( CheckMeshPS(mesh,token[0],command) ){
547 G4PSMinKinEAtGeneration* ps = nullptr;
548 if(shape==MeshShape::realWorldLogVol || shape==MeshShape::probe)
549 { ps = new G4PSMinKinEAtGeneration(token[0],mesh->GetCopyNumberLevel()); }
550 else
551 { ps = new G4PSMinKinEAtGeneration3D(token[0]); }
552 ps->SetUnit(token[1]);
553 mesh->SetPrimitiveScorer(ps);
554 }
555 } else if(command== qStepCheckerCmd){
556 if( CheckMeshPS(mesh,token[0],command) ){
557 G4PSStepChecker* ps = nullptr;
558 if(shape==MeshShape::realWorldLogVol || shape==MeshShape::probe)
559 { ps = new G4PSStepChecker(token[0],mesh->GetCopyNumberLevel()); }
560 else
561 { ps = new G4PSStepChecker3D(token[0]); }
562 mesh->SetPrimitiveScorer(ps);
563 }
564
565 //
566 // Filters
567 //
568 }else if(command== fchargedCmd){
569 if(!mesh->IsCurrentPrimitiveScorerNull()) {
570 mesh->SetFilter(new G4SDChargedFilter(token[0]));
571 } else {
572 ed << "WARNING[" << fchargedCmd->GetCommandPath()
573 << "] : Current quantity is not set. Set or touch a quantity first.";
574 command->CommandFailed(ed);
575 }
576 }else if(command== fneutralCmd){
577 if(!mesh->IsCurrentPrimitiveScorerNull()) {
578 mesh->SetFilter(new G4SDNeutralFilter(token[0]));
579 } else {
580 ed << "WARNING[" << fneutralCmd->GetCommandPath()
581 << "] : Current quantity is not set. Set or touch a quantity first.";
582 command->CommandFailed(ed);
583 }
584 }else if(command== fkinECmd){
585 if(!mesh->IsCurrentPrimitiveScorerNull()) {
586 G4String& name = token[0];
587 G4double elow = StoD(token[1]);
588 G4double ehigh = StoD(token[2]);
589 G4double unitVal = G4UnitDefinition::GetValueOf(token[3]);
590 mesh->SetFilter(new G4SDKineticEnergyFilter(name,elow*unitVal,ehigh*unitVal));
591 } else {
592 ed << "WARNING[" << fkinECmd->GetCommandPath()
593 << "] : Current quantity is not set. Set or touch a quantity first.";
594 command->CommandFailed(ed);
595 }
596 }else if(command== fparticleKinECmd){
597 if(!mesh->IsCurrentPrimitiveScorerNull()) {
598 FParticleWithEnergyCommand(mesh,token);
599 } else {
600 ed << "WARNING[" << fparticleKinECmd->GetCommandPath()
601 << "] : Current quantity is not set. Set or touch a quantity first.";
602 command->CommandFailed(ed);
603 }
604 } else if(command==fparticleCmd) {
605 if(!mesh->IsCurrentPrimitiveScorerNull()) {
606 FParticleCommand(mesh,token);
607 } else {
608 ed << "WARNING[" << fparticleCmd->GetCommandPath()
609 << "] : Current quantity is not set. Set or touch a quantity first.";
610 command->CommandFailed(ed);
611 }
612 }
613}
std::vector< G4String > G4TokenVec
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
virtual void SetUnit(const G4String &unit)
void SetCylinderSize(G4double dr, G4double dz)
virtual void SetUnit(const G4String &unit)
void SetCylinderSize(G4double dr, G4double dz)
virtual void SetUnit(const G4String &unit)
virtual void SetUnit(const G4String &unit)
void DivideByArea(G4bool flg=true)
void Weighted(G4bool flg=true)
virtual void SetUnit(const G4String &unit)
void Weighted(G4bool flg=true)
void DivideByArea(G4bool flg=true)
virtual void SetUnit(const G4String &unit)
virtual void SetUnit(const G4String &unit)
void Weighted(G4bool flg=true)
void SetBoundaryFlag(G4bool flg=true)
Definition: G4PSNofStep.hh:71
void Weighted(G4bool flg=true)
void SetCylinderSize(G4double dr, G4double dz)
virtual void SetUnit(const G4String &unit)
virtual void SetUnit(const G4String &unit)
void Weighted(G4bool flg=true)
void Weighted(G4bool flg=true)
void Weighted(G4bool flg=true)
void Weighted(G4bool flg=true)
virtual void SetUnit(const G4String &unit)
void MultiplyKineticEnergy(G4bool flg=true)
void Weighted(G4bool flg=true)
void DivideByVelocity(G4bool flg=true)
void SetDivCos(G4bool val)
void FParticleWithEnergyCommand(G4VScoringMesh *mesh, G4TokenVec &token)
G4bool CheckMeshPS(G4VScoringMesh *mesh, 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
void GetNumberOfSegments(G4int nSegment[3])
G4int GetCopyNumberLevel() const
void SetCurrentPSUnit(const G4String &unit)
void SetCurrentPrimitiveScorer(const G4String &name)
void SetPrimitiveScorer(G4VPrimitiveScorer *ps)
G4String GetCurrentPSUnit()
G4bool IsCurrentPrimitiveScorerNull()

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