source: trunk/source/digits_hits/utils/src/G4ScoreQuantityMessenger.cc @ 1340

Last change on this file since 1340 was 1340, checked in by garnier, 14 years ago

update ti head

File size: 16.8 KB
Line 
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,v 1.10 2010/11/03 08:28:42 taso Exp $
28// GEANT4 tag $Name:  $
29//
30// ---------------------------------------------------------------------
31// Modifications
32// 08-Oct-2010 T.Aso remove unit of G4PSPassageCellCurrent.
33// ---------------------------------------------------------------------
34
35#include "G4ScoreQuantityMessenger.hh"
36#include "G4ScoringManager.hh"
37#include "G4VScoringMesh.hh"
38
39#include "G4PSCellCharge3D.hh"
40#include "G4PSCellFlux3D.hh"
41#include "G4PSPassageCellFlux3D.hh"
42#include "G4PSEnergyDeposit3D.hh"
43#include "G4PSDoseDeposit3D.hh"
44#include "G4PSNofStep3D.hh"
45#include "G4PSNofSecondary3D.hh"
46//
47#include "G4PSTrackLength3D.hh"
48#include "G4PSPassageCellCurrent3D.hh"
49#include "G4PSPassageTrackLength3D.hh"
50#include "G4PSFlatSurfaceCurrent3D.hh"
51#include "G4PSFlatSurfaceFlux3D.hh"
52#include "G4PSSphereSurfaceCurrent3D.hh"
53#include "G4PSSphereSurfaceFlux3D.hh"
54#include "G4PSCylinderSurfaceCurrent3D.hh"
55#include "G4PSCylinderSurfaceFlux3D.hh"
56#include "G4PSNofCollision3D.hh"
57#include "G4PSPopulation3D.hh"
58#include "G4PSTrackCounter3D.hh"
59#include "G4PSTermination3D.hh"
60#include "G4PSMinKinEAtGeneration3D.hh"
61
62#include "G4SDChargedFilter.hh"
63#include "G4SDNeutralFilter.hh"
64#include "G4SDKineticEnergyFilter.hh"
65#include "G4SDParticleFilter.hh"
66#include "G4SDParticleWithEnergyFilter.hh"
67
68#include "G4UIdirectory.hh"
69#include "G4UIcmdWithoutParameter.hh"
70#include "G4UIcmdWithAnInteger.hh"
71#include "G4UIcmdWithAString.hh"
72#include "G4UIcmdWithABool.hh"
73#include "G4UIcmdWithADoubleAndUnit.hh"
74#include "G4UIcmdWith3VectorAndUnit.hh"
75#include "G4UIcommand.hh"
76#include "G4Tokenizer.hh"
77#include "G4UnitsTable.hh"
78
79G4ScoreQuantityMessenger::G4ScoreQuantityMessenger(G4ScoringManager* SManager)
80:fSMan(SManager)
81{
82  QuantityCommands();
83  FilterCommands();
84}
85
86void G4ScoreQuantityMessenger::FilterCommands()
87{
88  G4UIparameter* param;
89
90  //
91  // Filter commands
92  filterDir = new G4UIdirectory("/score/filter/");
93  filterDir->SetGuidance("  Scoring filter commands.");
94  //
95  fchargedCmd = new G4UIcmdWithAString("/score/filter/charged",this);
96  fchargedCmd->SetGuidance("Charged particle filter.");
97  fchargedCmd->SetParameterName("fname",false);
98  //
99  fneutralCmd = new G4UIcmdWithAString("/score/filter/neutral",this);
100  fneutralCmd->SetGuidance("Neutral particle filter.");
101  fneutralCmd->SetParameterName("fname",false);
102  //
103  fkinECmd = new G4UIcommand("/score/filter/kineticEnergy",this);
104  fkinECmd->SetGuidance("Kinetic energy filter.");
105  fkinECmd->SetGuidance("[usage] /score/filter/kineticEnergy fname Elow Ehigh unit");
106  fkinECmd->SetGuidance("  fname     :(String) Filter Name ");
107  fkinECmd->SetGuidance("  Elow      :(Double) Lower edge of kinetic energy");
108  fkinECmd->SetGuidance("  Ehigh     :(Double) Higher edge of kinetic energy");
109  fkinECmd->SetGuidance("  unit      :(String) unit of given kinetic energy");
110  param = new G4UIparameter("fname",'s',false);
111  fkinECmd->SetParameter(param);
112  param = new G4UIparameter("elow",'d',true);
113  param->SetDefaultValue("0.0");
114  fkinECmd->SetParameter(param);
115  param = new G4UIparameter("ehigh",'d',false);
116  fkinECmd->SetParameter(param);
117  G4String smax = DtoS(DBL_MAX);
118  param->SetDefaultValue(smax);
119  param = new G4UIparameter("unit",'s',false);
120  param->SetDefaultValue("keV");
121  fkinECmd->SetParameter(param);
122  //
123  fparticleCmd = new G4UIcommand("/score/filter/particle",this);
124  fparticleCmd->SetGuidance("Particle filter.");
125  fparticleCmd->SetGuidance("[usage] /score/filter/particle fname p0 .. pn");
126  fparticleCmd->SetGuidance("  fname     :(String) Filter Name ");
127  fparticleCmd->SetGuidance("  p0 .. pn  :(String) particle names");
128  param = new G4UIparameter("fname",'s',false);
129  fparticleCmd->SetParameter(param);
130  param = new G4UIparameter("particlelist",'s',false);
131  param->SetDefaultValue("");
132  fparticleCmd->SetParameter(param);
133  //
134  //
135  //
136  fparticleKinECmd = new G4UIcommand("/score/filter/particleWithKineticEnergy",this);
137  fparticleKinECmd->SetGuidance("Particle with kinetic energy filter.");
138  fparticleKinECmd->SetGuidance("[usage] /score/filter/particleWithKineticEnergy fname Elow Ehigh unit p0 .. pn");
139  fparticleKinECmd->SetGuidance("  fname     :(String) Filter Name ");
140  fparticleKinECmd->SetGuidance("  Elow      :(Double) Lower edge of kinetic energy");
141  fparticleKinECmd->SetGuidance("  Ehigh     :(Double) Higher edge of kinetic energy");
142  fparticleKinECmd->SetGuidance("  unit      :(String) unit of given kinetic energy");
143  fparticleKinECmd->SetGuidance("  p0 .. pn  :(String) particle names");
144  param = new G4UIparameter("fname",'s',false);
145  fparticleKinECmd->SetParameter(param);
146  param = new G4UIparameter("elow",'d',false);
147  param->SetDefaultValue("0.0");
148  fparticleKinECmd->SetParameter(param);
149  param = new G4UIparameter("ehigh",'d',true);
150  param->SetDefaultValue(smax);
151  fparticleKinECmd->SetParameter(param);
152  param = new G4UIparameter("unit",'s',true);
153  param->SetDefaultValue("keV");
154  fparticleKinECmd->SetParameter(param);
155  param = new G4UIparameter("particlelist",'s',false);
156  param->SetDefaultValue("");
157  fparticleKinECmd->SetParameter(param);
158  //
159  //
160}
161
162G4ScoreQuantityMessenger::~G4ScoreQuantityMessenger()
163{
164    delete         quantityDir;
165    delete         qTouchCmd;
166    delete         qGetUnitCmd;
167    delete         qSetUnitCmd;
168
169    //
170    delete    qCellChgCmd;
171    delete    qCellFluxCmd;
172    delete    qPassCellFluxCmd;
173    delete    qeDepCmd;
174    delete    qdoseDepCmd;
175    delete    qnOfStepCmd;
176    delete    qnOfSecondaryCmd;
177    //
178    delete          qTrackLengthCmd;
179    delete          qPassCellCurrCmd;
180    delete          qPassTrackLengthCmd;
181    delete          qFlatSurfCurrCmd;
182    delete          qFlatSurfFluxCmd;
183//    delete          qSphereSurfCurrCmd;
184//    delete          qSphereSurfFluxCmd;
185//    delete          qCylSurfCurrCmd;
186//    delete          qCylSurfFluxCmd;
187    delete          qNofCollisionCmd;
188    delete          qPopulationCmd;
189    delete          qTrackCountCmd;
190    delete          qTerminationCmd;
191    delete          qMinKinEAtGeneCmd;
192    //
193    delete   filterDir;
194    delete   fchargedCmd;
195    delete   fneutralCmd;
196    delete   fkinECmd;
197    delete   fparticleCmd;
198    delete   fparticleKinECmd;
199}
200
201void G4ScoreQuantityMessenger::SetNewValue(G4UIcommand * command,G4String newVal)
202{
203      //
204      // Get Current Mesh
205      //
206      G4VScoringMesh* mesh = fSMan->GetCurrentMesh();
207      if(!mesh)
208      {
209        G4cerr << "ERROR : No mesh is currently open. Open/create a mesh first. Command ignored." << G4endl;
210        return;
211      }
212      // Tokens
213      G4TokenVec token;
214      FillTokenVec(newVal,token);
215      //
216      // Commands for Current Mesh
217      if(command==qTouchCmd) {
218              mesh->SetCurrentPrimitiveScorer(newVal);
219      } else if(command == qGetUnitCmd ){
220          G4cout << "Unit:  "<< mesh->GetCurrentPSUnit() <<G4endl;
221      } else if(command == qSetUnitCmd ){
222          mesh->SetCurrentPSUnit(newVal);
223      } else if(command== qCellChgCmd) {
224          if ( CheckMeshPS(mesh,token[0]) ){
225              G4PSCellCharge3D* ps = new G4PSCellCharge3D(token[0]);
226              ps->SetUnit(token[1]);
227              mesh->SetPrimitiveScorer(ps);
228          }
229      } else if(command== qCellFluxCmd) {
230          if ( CheckMeshPS(mesh,token[0]) ){
231              G4PSCellFlux3D* ps = new G4PSCellFlux3D(token[0]);
232              ps->SetUnit(token[1]);
233              mesh->SetPrimitiveScorer(ps);
234          }
235      } else if(command== qPassCellFluxCmd) {
236          if ( CheckMeshPS(mesh,token[0]) ){
237              G4PSPassageCellFlux3D* ps = new G4PSPassageCellFlux3D(token[0]);
238              ps->SetUnit(token[1]);
239              mesh->SetPrimitiveScorer(ps);
240          }
241      } else if(command==qeDepCmd) {
242          if ( CheckMeshPS(mesh,token[0]) ){
243              G4PSEnergyDeposit3D* ps =new G4PSEnergyDeposit3D(token[0]);
244              ps->SetUnit(token[1]);
245              mesh->SetPrimitiveScorer(ps);
246          }
247      } else if(command== qdoseDepCmd) {
248          if ( CheckMeshPS(mesh,token[0]) ){
249              G4PSDoseDeposit3D* ps = new G4PSDoseDeposit3D(token[0]);
250              ps->SetUnit(token[1]);
251              mesh->SetPrimitiveScorer(ps);
252          }
253      } else if(command== qnOfStepCmd) {
254          if ( CheckMeshPS(mesh,token[0]) ){
255              G4PSNofStep3D* ps = new G4PSNofStep3D(token[0]);
256              mesh->SetPrimitiveScorer(ps);
257          }
258      } else if(command== qnOfSecondaryCmd) {
259          if ( CheckMeshPS(mesh,token[0]) ){
260              G4PSNofSecondary3D* ps =new G4PSNofSecondary3D(token[0]);
261              mesh->SetPrimitiveScorer(ps);
262          }
263      } else if(command== qTrackLengthCmd) {
264          if ( CheckMeshPS(mesh,token[0]) ){
265              G4PSTrackLength3D* ps = new G4PSTrackLength3D(token[0]);
266              ps->Weighted(StoB(token[1]));
267              ps->MultiplyKineticEnergy(StoB(token[2]));
268              ps->DivideByVelocity(StoB(token[3]));
269              ps->SetUnit(token[4]);
270              mesh->SetPrimitiveScorer(ps);
271          }
272      } else if(command== qPassCellCurrCmd){
273          if( CheckMeshPS(mesh,token[0]) ) {
274              G4PSPassageCellCurrent* ps = new G4PSPassageCellCurrent3D(token[0]);
275              ps->Weighted(StoB(token[1]));
276              //ps->SetUnit(token[2]);
277              mesh->SetPrimitiveScorer(ps);
278          }
279      } else if(command== qPassTrackLengthCmd){
280          if( CheckMeshPS(mesh,token[0]) ) {
281              G4PSPassageTrackLength* ps = new G4PSPassageTrackLength3D(token[0]);
282              ps->Weighted(StoB(token[1]));
283              mesh->SetPrimitiveScorer(ps);
284          }
285      } else if(command== qFlatSurfCurrCmd){
286          if( CheckMeshPS(mesh,token[0])) {
287              G4PSFlatSurfaceCurrent3D* ps = 
288                new G4PSFlatSurfaceCurrent3D(token[0],StoI(token[1]));
289              ps->Weighted(StoB(token[2]));
290              ps->DivideByArea(StoB(token[3]));
291              ps->SetUnit(token[4]);
292              mesh->SetPrimitiveScorer(ps);
293          }
294      } else if(command== qFlatSurfFluxCmd){
295          if( CheckMeshPS(mesh, token[0] )) {
296              G4PSFlatSurfaceFlux3D* ps = new G4PSFlatSurfaceFlux3D(token[0],StoI(token[1]));
297              ps->SetUnit(token[2]);
298              mesh->SetPrimitiveScorer(ps);
299          }
300//    } else if(command== qSphereSurfCurrCmd){
301//        if( CheckMeshPS(mesh, token[0] )) {
302//            G4PSSphereSurfaceCurrent3D* ps =
303//              new G4PSSphereSurfaceCurrent3D(token[0],StoI(token[1]));
304//            ps->Weighted(StoB(token[2]));
305//            ps->DivideByArea(StoB(token[3]));
306//            ps->SetUnit(token[4]);
307//            mesh->SetPrimitiveScorer(ps);
308//        }
309//   } else if(command== qSphereSurfFluxCmd){
310//        if( CheckMeshPS(mesh,token[0])) {
311//            G4PSSphereSurfaceFlux3D* ps = new G4PSSphereSurfaceFlux3D(token[0], StoI(token[1]));
312//            ps->SetUnit(token[2]);       
313//            mesh->SetPrimitiveScorer(ps);
314//        }
315//   } else if(command== qCylSurfCurrCmd){
316//        if( CheckMeshPS(mesh, token[0] ) ) {
317//            G4PSCylinderSurfaceCurrent3D* ps =
318//              new G4PSCylinderSurfaceCurrent3D(token[0],StoI(token[1]));
319//            ps->Weighted(StoB(token[2]));
320//            ps->DivideByArea(StoB(token[3]));
321//            ps->SetUnit(token[4]);
322//            mesh->SetPrimitiveScorer(ps);
323//        }
324//   } else if(command== qCylSurfFluxCmd){
325//        if( CheckMeshPS(mesh, token[0] ) {
326//            G4PSCylinerSurfaceFlux3D* ps =new G4PSCylinderSurfaceFlux3D(token[0], StoI(token[1]));
327//            ps->SetUnit(token[2]);
328//            mesh->SetPrimitiveScorer(ps);
329//        }
330      } else if(command== qNofCollisionCmd){
331          if( CheckMeshPS(mesh,token[0])) {
332              G4PSNofCollision3D* ps =new G4PSNofCollision3D(token[0]); 
333              ps->Weighted(StoB(token[1]));
334              mesh->SetPrimitiveScorer(ps);
335          }
336      } else if(command== qPopulationCmd){
337          if( CheckMeshPS(mesh,token[0]) ) {
338              G4PSPopulation3D* ps =new G4PSPopulation3D(token[0]); 
339              ps->Weighted(StoB(token[1]));
340              mesh->SetPrimitiveScorer(ps);
341          }
342      } else if(command== qTrackCountCmd){
343          if( CheckMeshPS(mesh,token[0])) {
344              G4PSTrackCounter3D* ps =new G4PSTrackCounter3D(token[0],StoI(token[1])); 
345              mesh->SetPrimitiveScorer(ps);
346          }
347      } else if(command== qTerminationCmd){
348          if( CheckMeshPS(mesh,token[0])) {
349              G4PSTermination3D* ps =new G4PSTermination3D(token[0]); 
350              ps->Weighted(StoB(token[1]));
351              mesh->SetPrimitiveScorer(ps);
352          }
353
354      } else if(command== qMinKinEAtGeneCmd){
355          if( CheckMeshPS(mesh,token[0]) ){
356              G4PSMinKinEAtGeneration3D* ps =new G4PSMinKinEAtGeneration3D(token[0]); 
357              ps->SetUnit(token[1]);
358              mesh->SetPrimitiveScorer(ps);
359          }
360
361            //
362            // Filters
363            //
364      }else if(command== fchargedCmd){
365            if(!mesh->IsCurrentPrimitiveScorerNull()) {
366              mesh->SetFilter(new G4SDChargedFilter(token[0])); 
367            } else {
368              G4cout << "WARNING[" << fchargedCmd->GetCommandPath()
369                     << "] : Current quantity is not set. Set or touch a quantity first." << G4endl;
370            }
371      }else if(command== fneutralCmd){
372            if(!mesh->IsCurrentPrimitiveScorerNull()) {
373              mesh->SetFilter(new G4SDNeutralFilter(token[0])); 
374            } else {
375              G4cout << "WARNING[" << fneutralCmd->GetCommandPath()
376                     << "] : Current quantity is not set. Set or touch a quantity first." << G4endl;
377            }
378      }else if(command== fkinECmd){
379            if(!mesh->IsCurrentPrimitiveScorerNull()) {
380              G4String& name = token[0];
381              G4double elow  = StoD(token[1]);
382              G4double ehigh = StoD(token[2]);
383              mesh->SetFilter(new G4SDKineticEnergyFilter(name,elow,ehigh));
384            } else {
385              G4cout << "WARNING[" << fkinECmd->GetCommandPath()
386                     << "] : Current quantity is not set. Set or touch a quantity first." << G4endl;
387            }
388      }else if(command== fparticleKinECmd){
389            if(!mesh->IsCurrentPrimitiveScorerNull()) {
390              FParticleWithEnergyCommand(mesh,token); 
391            } else {
392              G4cout << "WARNING[" << fparticleKinECmd->GetCommandPath()
393                     << "] : Current quantity is not set. Set or touch a quantity first." << G4endl;
394            }
395      } else if(command==fparticleCmd) {
396            if(!mesh->IsCurrentPrimitiveScorerNull()) {
397              FParticleCommand(mesh,token);
398            } else {
399              G4cout << "WARNING[" << fparticleCmd->GetCommandPath()
400                     << "] : Current quantity is not set. Set or touch a quantity first." << G4endl;
401            }
402      }
403}
404
405G4String G4ScoreQuantityMessenger::GetCurrentValue(G4UIcommand * /*command*/)
406{
407  G4String val;
408
409  return val;
410}
411
412void G4ScoreQuantityMessenger::FillTokenVec(G4String newValues, G4TokenVec& token){
413
414    G4Tokenizer next(newValues);
415    G4String val;
416    while ( !(val = next()).isNull() ) {
417        token.push_back(val);
418    }
419}
420
421
422void G4ScoreQuantityMessenger::FParticleCommand(G4VScoringMesh* mesh, G4TokenVec& token){
423    //
424    // Filter name
425    G4String name = token[0];
426    //
427    // particle list
428    std::vector<G4String> pnames;
429    for ( G4int i = 1; i<(G4int)token.size(); i++){
430        pnames.push_back(token[i]);
431    }
432    //
433    // Attach Filter
434    mesh->SetFilter(new G4SDParticleFilter(name,pnames));
435}   
436
437void G4ScoreQuantityMessenger::FParticleWithEnergyCommand(G4VScoringMesh* mesh,G4TokenVec& token){
438    G4String& name = token[0];
439    G4double  elow = StoD(token[1]);
440    G4double  ehigh= StoD(token[2]);
441    G4double  unitVal = G4UnitDefinition::GetValueOf(token[3]);
442    G4SDParticleWithEnergyFilter* filter = 
443        new G4SDParticleWithEnergyFilter(name,elow*unitVal,ehigh*unitVal);
444    for ( G4int i = 4; i < (G4int)token.size(); i++){
445        filter->add(token[i]);
446    }
447    mesh->SetFilter(filter);
448}
449
450G4bool G4ScoreQuantityMessenger::CheckMeshPS(G4VScoringMesh* mesh, G4String& psname){
451    if(!mesh->FindPrimitiveScorer(psname)) {
452        return true;
453    } else {
454        G4cout << "WARNING[" << qTouchCmd->GetCommandPath()
455               << "] : Quantity name, \"" << psname << "\", is already existing." << G4endl;
456        mesh->SetNullToCurrentPrimitiveScorer();
457        return false;
458    }
459}
Note: See TracBrowser for help on using the repository browser.