Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4ParticleGunMessenger Class Reference

#include <G4ParticleGunMessenger.hh>

Inheritance diagram for G4ParticleGunMessenger:
G4UImessenger

Public Member Functions

 G4ParticleGunMessenger (G4ParticleGun *fPtclGun)
 
 ~G4ParticleGunMessenger ()
 
void SetNewValue (G4UIcommand *command, G4String newValues)
 
G4String GetCurrentValue (G4UIcommand *command)
 
- Public Member Functions inherited from G4UImessenger
 G4UImessenger ()
 
 G4UImessenger (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
virtual ~G4UImessenger ()
 
G4bool operator== (const G4UImessenger &messenger) const
 

Additional Inherited Members

- Protected Member Functions inherited from G4UImessenger
G4String ItoS (G4int i)
 
G4String DtoS (G4double a)
 
G4String BtoS (G4bool b)
 
G4int StoI (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)
 
- Protected Attributes inherited from G4UImessenger
G4UIdirectorybaseDir
 
G4String baseDirName
 

Detailed Description

Definition at line 54 of file G4ParticleGunMessenger.hh.

Constructor & Destructor Documentation

G4ParticleGunMessenger::G4ParticleGunMessenger ( G4ParticleGun fPtclGun)

Definition at line 47 of file G4ParticleGunMessenger.cc.

References python.hepunit::cm, G4Geantino::Geantino(), G4ParticleDefinition::GetDecayTable(), G4ParticleTable::GetIterator(), G4ParticleDefinition::GetParticleName(), G4ParticleTable::GetParticleTable(), python.hepunit::GeV, G4ParticleDefinition::IsShortLived(), ns, G4ParticleTableIterator< K, V >::reset(), G4UIcmdWithAString::SetCandidates(), G4UIcmdWithADoubleAndUnit::SetDefaultUnit(), G4UIcmdWith3VectorAndUnit::SetDefaultUnit(), G4UIcmdWithAString::SetDefaultValue(), G4UIparameter::SetDefaultValue(), G4UIcommand::SetGuidance(), G4UIcommand::SetParameter(), G4UIcmdWithAString::SetParameterName(), G4UIcmdWithAnInteger::SetParameterName(), G4UIcmdWith3Vector::SetParameterName(), G4UIcmdWithADoubleAndUnit::SetParameterName(), G4UIcmdWith3VectorAndUnit::SetParameterName(), G4ParticleGun::SetParticleDefinition(), G4ParticleGun::SetParticleEnergy(), G4ParticleGun::SetParticleMomentumDirection(), G4VPrimaryGenerator::SetParticlePosition(), G4VPrimaryGenerator::SetParticleTime(), G4UIcommand::SetRange(), and G4ParticleTableIterator< K, V >::value().

48  :fParticleGun(fPtclGun),fShootIon(false),
49  fAtomicNumber(0),fAtomicMass(0),fIonCharge(0),fIonExciteEnergy(0.0),fIonEnergyLevel(0)
50 {
51  particleTable = G4ParticleTable::GetParticleTable();
52 
53  gunDirectory = new G4UIdirectory("/gun/");
54  gunDirectory->SetGuidance("Particle Gun control commands.");
55 
56  listCmd = new G4UIcmdWithoutParameter("/gun/List",this);
57  listCmd->SetGuidance("List available particles.");
58  listCmd->SetGuidance(" Invoke G4ParticleTable.");
59 
60  particleCmd = new G4UIcmdWithAString("/gun/particle",this);
61  particleCmd->SetGuidance("Set particle to be generated.");
62  particleCmd->SetGuidance(" (geantino is default)");
63  particleCmd->SetGuidance(" (ion can be specified for shooting ions)");
64  particleCmd->SetParameterName("particleName",true);
65  particleCmd->SetDefaultValue("geantino");
66  G4String candidateList;
67  G4ParticleTable::G4PTblDicIterator* itr = particleTable->GetIterator();
68  itr->reset();
69  while( (*itr)() )
70  {
71  G4ParticleDefinition* pd = itr->value();
72  if( !(pd->IsShortLived()) || pd->GetDecayTable() )
73  {
74  candidateList += pd->GetParticleName();
75  candidateList += " ";
76  }
77  }
78  candidateList += "ion ";
79  particleCmd->SetCandidates(candidateList);
80 
81  directionCmd = new G4UIcmdWith3Vector("/gun/direction",this);
82  directionCmd->SetGuidance("Set momentum direction.");
83  directionCmd->SetGuidance("Direction needs not to be a unit vector.");
84  directionCmd->SetParameterName("ex","ey","ez",true,true);
85  directionCmd->SetRange("ex != 0 || ey != 0 || ez != 0");
86 
87  energyCmd = new G4UIcmdWithADoubleAndUnit("/gun/energy",this);
88  energyCmd->SetGuidance("Set kinetic energy.");
89  energyCmd->SetParameterName("Energy",true,true);
90  energyCmd->SetDefaultUnit("GeV");
91  //energyCmd->SetUnitCategory("Energy");
92  //energyCmd->SetUnitCandidates("eV keV MeV GeV TeV");
93 
94  momCmd = new G4UIcmdWith3VectorAndUnit("/gun/momentum",this);
95  momCmd->SetGuidance("Set momentum. This command is equivalent to two commands /gun/direction and /gun/momentumAmp");
96  momCmd->SetParameterName("px","py","pz",true,true);
97  momCmd->SetRange("px != 0 || py != 0 || pz != 0");
98  momCmd->SetDefaultUnit("GeV");
99 
100  momAmpCmd = new G4UIcmdWithADoubleAndUnit("/gun/momentumAmp",this);
101  momAmpCmd->SetGuidance("Set absolute value of momentum.");
102  momAmpCmd->SetGuidance("Direction should be set by /gun/direction command.");
103  momAmpCmd->SetGuidance("This command should be used alternatively with /gun/energy.");
104  momAmpCmd->SetParameterName("Momentum",true,true);
105  momAmpCmd->SetDefaultUnit("GeV");
106 
107  positionCmd = new G4UIcmdWith3VectorAndUnit("/gun/position",this);
108  positionCmd->SetGuidance("Set starting position of the particle.");
109  positionCmd->SetParameterName("X","Y","Z",true,true);
110  positionCmd->SetDefaultUnit("cm");
111  //positionCmd->SetUnitCategory("Length");
112  //positionCmd->SetUnitCandidates("microm mm cm m km");
113 
114  timeCmd = new G4UIcmdWithADoubleAndUnit("/gun/time",this);
115  timeCmd->SetGuidance("Set initial time of the particle.");
116  timeCmd->SetParameterName("t0",true,true);
117  timeCmd->SetDefaultUnit("ns");
118  //timeCmd->SetUnitCategory("Time");
119  //timeCmd->SetUnitCandidates("ns ms s");
120 
121  polCmd = new G4UIcmdWith3Vector("/gun/polarization",this);
122  polCmd->SetGuidance("Set polarization.");
123  polCmd->SetParameterName("Px","Py","Pz",true,true);
124  polCmd->SetRange("Px>=-1.&&Px<=1.&&Py>=-1.&&Py<=1.&&Pz>=-1.&&Pz<=1.");
125 
126  numberCmd = new G4UIcmdWithAnInteger("/gun/number",this);
127  numberCmd->SetGuidance("Set number of particles to be generated.");
128  numberCmd->SetParameterName("N",true,true);
129  numberCmd->SetRange("N>0");
130 
131  ionCmd = new G4UIcommand("/gun/ion",this);
132  ionCmd->SetGuidance("Set properties of ion to be generated.");
133  ionCmd->SetGuidance("[usage] /gun/ion Z A Q");
134  ionCmd->SetGuidance(" Z:(int) AtomicNumber");
135  ionCmd->SetGuidance(" A:(int) AtomicMass");
136  ionCmd->SetGuidance(" Q:(int) Charge of Ion (in unit of e)");
137  ionCmd->SetGuidance(" E:(double) Excitation energy (in keV)");
138 
139  G4UIparameter* param;
140  param = new G4UIparameter("Z",'i',false);
141  ionCmd->SetParameter(param);
142  param = new G4UIparameter("A",'i',false);
143  ionCmd->SetParameter(param);
144  param = new G4UIparameter("Q",'i',true);
145  param->SetDefaultValue(-1);
146  ionCmd->SetParameter(param);
147  param = new G4UIparameter("E",'d',true);
148  param->SetDefaultValue(0.0);
149  ionCmd->SetParameter(param);
150 
151  ionLvlCmd = new G4UIcommand("/gun/ionL",this);
152  ionLvlCmd->SetGuidance("Set properties of ion to be generated.");
153  ionLvlCmd->SetGuidance("[usage] /gun/ion Z A Q I");
154  ionLvlCmd->SetGuidance(" Z:(int) AtomicNumber");
155  ionLvlCmd->SetGuidance(" A:(int) AtomicMass");
156  ionLvlCmd->SetGuidance(" Q:(int) Charge of Ion (in unit of e)");
157  ionLvlCmd->SetGuidance(" I:(int) Level number of metastable state (0 = ground)");
158 
159  G4UIparameter* paraml;
160  paraml = new G4UIparameter("Z",'i',false);
161  ionLvlCmd->SetParameter(paraml);
162  paraml = new G4UIparameter("A",'i',false);
163  ionLvlCmd->SetParameter(paraml);
164  paraml = new G4UIparameter("Q",'i',true);
165  paraml->SetDefaultValue(-1);
166  ionLvlCmd->SetParameter(paraml);
167  paraml = new G4UIparameter("I",'i',true);
168  paraml->SetDefaultValue("0");
169  ionLvlCmd->SetParameter(paraml);
170 
171  // set initial value to G4ParticleGun
172  fParticleGun->SetParticleDefinition( G4Geantino::Geantino() );
173  fParticleGun->SetParticleMomentumDirection( G4ThreeVector(1.0,0.0,0.0) );
174  fParticleGun->SetParticleEnergy( 1.0*GeV );
175  fParticleGun->SetParticlePosition(G4ThreeVector(0.0*cm, 0.0*cm, 0.0*cm));
176  fParticleGun->SetParticleTime( 0.0*ns );
177 }
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:152
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
CLHEP::Hep3Vector G4ThreeVector
void SetDefaultUnit(const char *defUnit)
void SetParameterName(const char *theNameX, const char *theNameY, const char *theNameZ, G4bool omittable, G4bool currentAsDefault=false)
void SetDefaultValue(const char *theDefaultValue)
void SetParticleMomentumDirection(G4ParticleMomentum aMomentumDirection)
void SetParameterName(const char *theNameX, const char *theNameY, const char *theNameZ, G4bool omittable, G4bool currentAsDefault=false)
const G4String & GetParticleName() const
G4DecayTable * GetDecayTable() const
void SetParticlePosition(G4ThreeVector aPosition)
void reset(G4bool ifSkipIon=true)
void SetRange(const char *rs)
Definition: G4UIcommand.hh:125
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:161
void SetParticleEnergy(G4double aKineticEnergy)
static G4ParticleTable * GetParticleTable()
void SetDefaultValue(const char *defVal)
void SetDefaultUnit(const char *defUnit)
static G4Geantino * Geantino()
Definition: G4Geantino.cc:87
void SetCandidates(const char *candidateList)
G4PTblDicIterator * GetIterator() const
#define ns
Definition: xmlparse.cc:597
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetParticleDefinition(G4ParticleDefinition *aParticleDefinition)
void SetParticleTime(G4double aTime)
G4ParticleGunMessenger::~G4ParticleGunMessenger ( )

Definition at line 179 of file G4ParticleGunMessenger.cc.

180 {
181  delete listCmd;
182  delete particleCmd;
183  delete directionCmd;
184  delete energyCmd;
185  delete momCmd;
186  delete momAmpCmd;
187  delete positionCmd;
188  delete timeCmd;
189  delete polCmd;
190  delete numberCmd;
191  delete ionCmd;
192  delete ionLvlCmd;
193  delete gunDirectory;
194 }

Member Function Documentation

G4String G4ParticleGunMessenger::GetCurrentValue ( G4UIcommand command)
virtual

Reimplemented from G4UImessenger.

Definition at line 232 of file G4ParticleGunMessenger.cc.

References G4UIcommand::ConvertToString(), G4cerr, G4endl, G4ParticleGun::GetNumberOfParticles(), G4ParticleGun::GetParticleDefinition(), G4ParticleGun::GetParticleEnergy(), G4ParticleGun::GetParticleMomentum(), G4ParticleGun::GetParticleMomentumDirection(), G4ParticleDefinition::GetParticleName(), G4ParticleGun::GetParticlePolarization(), G4VPrimaryGenerator::GetParticlePosition(), G4VPrimaryGenerator::GetParticleTime(), and G4UImessenger::ItoS().

233 {
234  G4String cv;
235 
236  if( command==directionCmd )
237  { cv = directionCmd->ConvertToString(fParticleGun->GetParticleMomentumDirection()); }
238  else if( command==particleCmd )
239  { cv = fParticleGun->GetParticleDefinition()->GetParticleName(); }
240  else if( command==energyCmd )
241  {
242  G4double ene = fParticleGun->GetParticleEnergy();
243  if(ene == 0.)
244  { G4cerr << " G4ParticleGun: was defined in terms of momentum." << G4endl; }
245  else
246  { cv = energyCmd->ConvertToString(ene,"GeV"); }
247  }
248  else if( command==momCmd || command==momAmpCmd )
249  {
250  G4double mom = fParticleGun->GetParticleMomentum();
251  if(mom == 0.)
252  { G4cerr << " G4ParticleGun: was defined in terms of kinetic energy." << G4endl; }
253  else
254  {
255  if( command==momCmd )
256  { cv = momCmd->ConvertToString(mom*(fParticleGun->GetParticleMomentumDirection()),"GeV"); }
257  else
258  { cv = momAmpCmd->ConvertToString(mom,"GeV"); }
259  }
260  }
261  else if( command==positionCmd )
262  { cv = positionCmd->ConvertToString(fParticleGun->GetParticlePosition(),"cm"); }
263  else if( command==timeCmd )
264  { cv = timeCmd->ConvertToString(fParticleGun->GetParticleTime(),"ns"); }
265  else if( command==polCmd )
266  { cv = polCmd->ConvertToString(fParticleGun->GetParticlePolarization()); }
267  else if( command==numberCmd )
268  { cv = numberCmd->ConvertToString(fParticleGun->GetNumberOfParticles()); }
269  else if( command==ionCmd )
270  {
271  if (fShootIon) {
272  cv = ItoS(fAtomicNumber) + " " + ItoS(fAtomicMass) + " ";
273  cv += ItoS(fIonCharge);
274  } else {
275  cv = "";
276  }
277  }
278  return cv;
279 }
G4double GetParticleMomentum() const
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:357
G4ThreeVector GetParticlePosition()
const G4String & GetParticleName() const
G4ParticleMomentum GetParticleMomentumDirection() const
G4ThreeVector GetParticlePolarization() const
G4String ItoS(G4int i)
G4int GetNumberOfParticles() const
G4ParticleDefinition * GetParticleDefinition() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double GetParticleEnergy() const
G4GLOB_DLL std::ostream G4cerr
void G4ParticleGunMessenger::SetNewValue ( G4UIcommand command,
G4String  newValues 
)
virtual

Reimplemented from G4UImessenger.

Definition at line 196 of file G4ParticleGunMessenger.cc.

References G4ParticleTable::DumpTable(), G4ParticleTable::FindParticle(), G4UIcmdWith3Vector::GetNew3VectorValue(), G4UIcmdWith3VectorAndUnit::GetNew3VectorValue(), G4UIcmdWithADoubleAndUnit::GetNewDoubleValue(), G4UIcmdWithAnInteger::GetNewIntValue(), G4ParticleGun::SetNumberOfParticles(), G4ParticleGun::SetParticleDefinition(), G4ParticleGun::SetParticleEnergy(), G4ParticleGun::SetParticleMomentum(), G4ParticleGun::SetParticleMomentumDirection(), G4ParticleGun::SetParticlePolarization(), G4VPrimaryGenerator::SetParticlePosition(), and G4VPrimaryGenerator::SetParticleTime().

197 {
198  if (command==listCmd) {
199  particleTable->DumpTable();
200  } else if (command==particleCmd) {
201  if (newValues =="ion") {
202  fShootIon = true;
203  } else {
204  fShootIon = false;
205  G4ParticleDefinition* pd = particleTable->FindParticle(newValues);
206  if(pd != 0)
207  { fParticleGun->SetParticleDefinition( pd ); }
208  }
209 
210  } else if( command==directionCmd )
211  { fParticleGun->SetParticleMomentumDirection(directionCmd->GetNew3VectorValue(newValues)); }
212  else if( command==energyCmd )
213  { fParticleGun->SetParticleEnergy(energyCmd->GetNewDoubleValue(newValues)); }
214  else if( command==momCmd )
215  { fParticleGun->SetParticleMomentum(momCmd->GetNew3VectorValue(newValues)); }
216  else if( command==momAmpCmd )
217  { fParticleGun->SetParticleMomentum(momAmpCmd->GetNewDoubleValue(newValues)); }
218  else if( command==positionCmd )
219  { fParticleGun->SetParticlePosition(positionCmd->GetNew3VectorValue(newValues)); }
220  else if( command==timeCmd )
221  { fParticleGun->SetParticleTime(timeCmd->GetNewDoubleValue(newValues)); }
222  else if( command==polCmd )
223  { fParticleGun->SetParticlePolarization(polCmd->GetNew3VectorValue(newValues)); }
224  else if( command==numberCmd )
225  { fParticleGun->SetNumberOfParticles(numberCmd->GetNewIntValue(newValues)); }
226  else if( command==ionCmd )
227  { IonCommand(newValues); }
228  else if( command==ionLvlCmd )
229  { IonLevelCommand(newValues); }
230 }
void SetParticleMomentum(G4double aMomentum)
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
void DumpTable(const G4String &particle_name="ALL")
static G4int GetNewIntValue(const char *paramString)
static G4ThreeVector GetNew3VectorValue(const char *paramString)
void SetParticleMomentumDirection(G4ParticleMomentum aMomentumDirection)
static G4double GetNewDoubleValue(const char *paramString)
void SetParticlePolarization(G4ThreeVector aVal)
void SetParticlePosition(G4ThreeVector aPosition)
static G4ThreeVector GetNew3VectorValue(const char *paramString)
void SetNumberOfParticles(G4int i)
void SetParticleEnergy(G4double aKineticEnergy)
void SetParticleDefinition(G4ParticleDefinition *aParticleDefinition)
void SetParticleTime(G4double aTime)

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