source: trunk/source/event/src/G4GeneralParticleSourceMessenger.cc @ 1347

Last change on this file since 1347 was 1347, checked in by garnier, 13 years ago

geant4 tag 9.4

File size: 69.9 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//
28// MODULE:       G4GeneralParticleSourceMessenger.cc
29//
30// Version:      2.0
31// Date:         5/02/04
32// Author:       Fan Lei
33// Organisation: QinetiQ ltd.
34// Customer:     ESA/ESTEC
35//
36///////////////////////////////////////////////////////////////////////////////
37//
38// CHANGE HISTORY
39// --------------
40//
41// Version 2.0, 05/02/2004, Fan Lei, Created.
42//    After changes to version 1.1 as in Geant4 v6.0
43//     - Mutilple particle source definition
44//     - Re-structured commands
45//     - old commands have been retained for backward compatibility, will be
46//       removed in the future.
47//
48///////////////////////////////////////////////////////////////////////////////
49//
50#include "G4Geantino.hh"
51#include "G4ThreeVector.hh"
52#include "G4ParticleTable.hh"
53#include "G4UIdirectory.hh"
54#include "G4UIcmdWithoutParameter.hh"
55#include "G4UIcmdWithAString.hh"
56#include "G4UIcmdWithADoubleAndUnit.hh"
57#include "G4UIcmdWith3Vector.hh"
58#include "G4UIcmdWith3VectorAndUnit.hh"
59#include "G4UIcmdWithAnInteger.hh"
60#include "G4UIcmdWithADouble.hh"
61#include "G4UIcmdWithABool.hh"
62#include "G4ios.hh"
63
64#include "G4Tokenizer.hh"
65#include "G4GeneralParticleSourceMessenger.hh"
66#include "G4SingleParticleSource.hh"
67#include "G4GeneralParticleSource.hh"
68
69///////////////////////////////////////////////////////////////////////////////
70//
71G4GeneralParticleSourceMessenger::G4GeneralParticleSourceMessenger
72  (G4GeneralParticleSource *fPtclGun) 
73    : fGPS(fPtclGun),fShootIon(false)
74{
75  particleTable = G4ParticleTable::GetParticleTable();
76  histtype = "biasx";
77
78  gpsDirectory = new G4UIdirectory("/gps/");
79  gpsDirectory->SetGuidance("General Paricle Source control commands.");
80  //  gpsDirectory->SetGuidance(" The first 9 commands are the same as in G4ParticleGun ");
81
82  // now the commands for mutiple sources
83  sourceDirectory = new G4UIdirectory("/gps/source/");
84  sourceDirectory->SetGuidance("Multiple source control sub-directory");
85
86  addsourceCmd = new G4UIcmdWithADouble("/gps/source/add",this);
87  addsourceCmd->SetGuidance("add a new source defintion to the particle gun with the specified intensity");
88  addsourceCmd->SetParameterName("addsource",false,false);
89  addsourceCmd->SetRange("addsource > 0.");
90
91  listsourceCmd = new G4UIcmdWithoutParameter("/gps/source/list",this);
92  listsourceCmd->SetGuidance("List the defined particle sources");
93
94  clearsourceCmd = new G4UIcmdWithoutParameter("/gps/source/clear",this);
95  clearsourceCmd->SetGuidance("Remove all the defined particle sources");
96
97  getsourceCmd = new G4UIcmdWithoutParameter("/gps/source/show",this);
98  getsourceCmd->SetGuidance("Show the current source index and intensity");
99
100  setsourceCmd = new G4UIcmdWithAnInteger("/gps/source/set",this);
101  setsourceCmd->SetGuidance("set the indexed source as the current one");
102  setsourceCmd->SetGuidance("  so one can change its source definition");
103  setsourceCmd->SetParameterName("setsource",false,false);
104  setsourceCmd->SetRange("setsource >= 0");
105
106  deletesourceCmd = new G4UIcmdWithAnInteger("/gps/source/delete",this);
107  deletesourceCmd->SetGuidance("delete the indexed source from the list");
108  deletesourceCmd->SetParameterName("deletesource",false,false);
109  deletesourceCmd->SetRange("deletesource > 0");
110
111  setintensityCmd = new G4UIcmdWithADouble("/gps/source/intensity",this);
112  setintensityCmd->SetGuidance("reset the current source to the specified intensity");
113  setintensityCmd->SetParameterName("setintensity",false,false);
114  setintensityCmd->SetRange("setintensity > 0."); 
115
116  multiplevertexCmd = new G4UIcmdWithABool("/gps/source/multiplevertex",this);
117  multiplevertexCmd->SetGuidance("true for simulaneous generation mutiple vertex");
118  multiplevertexCmd->SetGuidance("Default is false");
119  multiplevertexCmd->SetParameterName("multiplevertex",true);
120  multiplevertexCmd->SetDefaultValue(false);
121
122  flatsamplingCmd = new G4UIcmdWithABool("/gps/source/flatsampling",this);
123  flatsamplingCmd->SetGuidance("true for appling flat (biased) sampling among the sources");
124  flatsamplingCmd->SetGuidance("Default is false");
125  flatsamplingCmd->SetParameterName("flatsampling",true);
126  flatsamplingCmd->SetDefaultValue(false);
127
128  // below we reproduce commands awailable in G4Particle Gun
129
130  listCmd = new G4UIcmdWithoutParameter("/gps/List",this);
131  listCmd->SetGuidance("List available particles.");
132  listCmd->SetGuidance(" Invoke G4ParticleTable.");
133
134  particleCmd = new G4UIcmdWithAString("/gps/particle",this);
135  particleCmd->SetGuidance("Set particle to be generated.");
136  particleCmd->SetGuidance(" (geantino is default)");
137  particleCmd->SetGuidance(" (ion can be specified for shooting ions)");
138  particleCmd->SetParameterName("particleName",true);
139  particleCmd->SetDefaultValue("geantino");
140  G4String candidateList; 
141  G4int nPtcl = particleTable->entries();
142  for(G4int i=0;i<nPtcl;i++)
143  {
144    candidateList += particleTable->GetParticleName(i);
145    candidateList += " ";
146  }
147  candidateList += "ion ";
148  particleCmd->SetCandidates(candidateList);
149
150  directionCmd = new G4UIcmdWith3Vector("/gps/direction",this);
151  directionCmd->SetGuidance("Set momentum direction.");
152  directionCmd->SetGuidance("Direction needs not to be a unit vector.");
153  directionCmd->SetParameterName("Px","Py","Pz",false,false); 
154  directionCmd->SetRange("Px != 0 || Py != 0 || Pz != 0");
155 
156  energyCmd = new G4UIcmdWithADoubleAndUnit("/gps/energy",this);
157  energyCmd->SetGuidance("Set kinetic energy.");
158  energyCmd->SetParameterName("Energy",false,false);
159  energyCmd->SetDefaultUnit("GeV");
160  //energyCmd->SetUnitCategory("Energy");
161  //energyCmd->SetUnitCandidates("eV keV MeV GeV TeV");
162
163  positionCmd = new G4UIcmdWith3VectorAndUnit("/gps/position",this);
164  positionCmd->SetGuidance("Set starting position of the particle.");
165  positionCmd->SetParameterName("X","Y","Z",false,false);
166  positionCmd->SetDefaultUnit("cm");
167  //positionCmd->SetUnitCategory("Length");
168  //positionCmd->SetUnitCandidates("microm mm cm m km");
169
170  ionCmd = new G4UIcommand("/gps/ion",this);
171  ionCmd->SetGuidance("Set properties of ion to be generated.");
172  ionCmd->SetGuidance("[usage] /gps/ion Z A Q E");
173  ionCmd->SetGuidance("        Z:(int) AtomicNumber");
174  ionCmd->SetGuidance("        A:(int) AtomicMass");
175  ionCmd->SetGuidance("        Q:(int) Charge of Ion (in unit of e)");
176  ionCmd->SetGuidance("        E:(double) Excitation energy (in keV)");
177 
178  G4UIparameter* param;
179  param = new G4UIparameter("Z",'i',false);
180  param->SetDefaultValue("1");
181  ionCmd->SetParameter(param);
182  param = new G4UIparameter("A",'i',false);
183  param->SetDefaultValue("1");
184  ionCmd->SetParameter(param);
185  param = new G4UIparameter("Q",'i',true);
186  param->SetDefaultValue("0");
187  ionCmd->SetParameter(param);
188  param = new G4UIparameter("E",'d',true);
189  param->SetDefaultValue("0.0");
190  ionCmd->SetParameter(param);
191
192
193  timeCmd = new G4UIcmdWithADoubleAndUnit("/gps/time",this);
194  timeCmd->SetGuidance("Set initial time of the particle.");
195  timeCmd->SetParameterName("t0",false,false);
196  timeCmd->SetDefaultUnit("ns");
197  //timeCmd->SetUnitCategory("Time");
198  //timeCmd->SetUnitCandidates("ns ms s");
199 
200  polCmd = new G4UIcmdWith3Vector("/gps/polarization",this);
201  polCmd->SetGuidance("Set polarization.");
202  polCmd->SetParameterName("Px","Py","Pz",false,false); 
203  polCmd->SetRange("Px>=-1.&&Px<=1.&&Py>=-1.&&Py<=1.&&Pz>=-1.&&Pz<=1.");
204
205  numberCmd = new G4UIcmdWithAnInteger("/gps/number",this);
206  numberCmd->SetGuidance("Set number of particles to be generated per vertex.");
207  numberCmd->SetParameterName("N",false,false);
208  numberCmd->SetRange("N>0");
209
210  // verbosity
211  verbosityCmd = new G4UIcmdWithAnInteger("/gps/verbose",this);
212  verbosityCmd->SetGuidance("Set Verbose level for GPS");
213  verbosityCmd->SetGuidance(" 0 : Silent");
214  verbosityCmd->SetGuidance(" 1 : Limited information");
215  verbosityCmd->SetGuidance(" 2 : Detailed information");
216  verbosityCmd->SetParameterName("level",false);
217  verbosityCmd->SetRange("level>=0 && level <=2");
218
219  // now extended commands
220  // Positional ones:
221  positionDirectory = new G4UIdirectory("/gps/pos/");
222  positionDirectory->SetGuidance("Positional commands sub-directory");
223
224  typeCmd1 = new G4UIcmdWithAString("/gps/pos/type",this);
225  typeCmd1->SetGuidance("Sets source distribution type.");
226  typeCmd1->SetGuidance("Either Point, Beam, Plane, Surface or Volume");
227  typeCmd1->SetParameterName("DisType",false,false);
228  typeCmd1->SetDefaultValue("Point");
229  typeCmd1->SetCandidates("Point Beam Plane Surface Volume");
230
231  shapeCmd1 = new G4UIcmdWithAString("/gps/pos/shape",this);
232  shapeCmd1->SetGuidance("Sets source shape for Plan, Surface or Volume type source.");
233  shapeCmd1->SetParameterName("Shape",false,false);
234  shapeCmd1->SetDefaultValue("NULL");
235  shapeCmd1->SetCandidates("Circle Annulus Ellipse Square Rectangle Sphere Ellipsoid Cylinder Para");
236
237  centreCmd1 = new G4UIcmdWith3VectorAndUnit("/gps/pos/centre",this);
238  centreCmd1->SetGuidance("Set centre coordinates of source.");
239  centreCmd1->SetGuidance("   same effect as the /gps/position command");
240  centreCmd1->SetParameterName("X","Y","Z",false,false);
241  centreCmd1->SetDefaultUnit("cm");
242  //  centreCmd1->SetUnitCandidates("micron mm cm m km");
243
244  posrot1Cmd1 = new G4UIcmdWith3Vector("/gps/pos/rot1",this);
245  posrot1Cmd1->SetGuidance("Set the 1st vector defining the rotation matrix'.");
246  posrot1Cmd1->SetGuidance("It does not need to be a unit vector.");
247  posrot1Cmd1->SetParameterName("R1x","R1y","R1z",false,false); 
248  posrot1Cmd1->SetRange("R1x != 0 || R1y != 0 || R1z != 0");
249
250  posrot2Cmd1 = new G4UIcmdWith3Vector("/gps/pos/rot2",this);
251  posrot2Cmd1->SetGuidance("Set the 2nd vector defining the rotation matrix'.");
252  posrot2Cmd1->SetGuidance("It does not need to be a unit vector.");
253  posrot2Cmd1->SetParameterName("R2x","R2y","R2z",false,false); 
254  posrot2Cmd1->SetRange("R2x != 0 || R2y != 0 || R2z != 0");
255
256  halfxCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/halfx",this);
257  halfxCmd1->SetGuidance("Set x half length of source.");
258  halfxCmd1->SetParameterName("Halfx",false,false);
259  halfxCmd1->SetDefaultUnit("cm");
260  //  halfxCmd1->SetUnitCandidates("micron mm cm m km");
261
262  halfyCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/halfy",this);
263  halfyCmd1->SetGuidance("Set y half length of source.");
264  halfyCmd1->SetParameterName("Halfy",false,false);
265  halfyCmd1->SetDefaultUnit("cm");
266  //  halfyCmd1->SetUnitCandidates("micron mm cm m km");
267
268  halfzCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/halfz",this);
269  halfzCmd1->SetGuidance("Set z half length of source.");
270  halfzCmd1->SetParameterName("Halfz",false,false);
271  halfzCmd1->SetDefaultUnit("cm");
272  //  halfzCmd1->SetUnitCandidates("micron mm cm m km");
273
274  radiusCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/radius",this);
275  radiusCmd1->SetGuidance("Set radius of source.");
276  radiusCmd1->SetParameterName("Radius",false,false);
277  radiusCmd1->SetDefaultUnit("cm");
278  //  radiusCmd1->SetUnitCandidates("micron mm cm m km");
279
280  radius0Cmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/inner_radius",this);
281  radius0Cmd1->SetGuidance("Set inner radius of source when required.");
282  radius0Cmd1->SetParameterName("Radius0",false,false);
283  radius0Cmd1->SetDefaultUnit("cm");
284  //  radius0Cmd1->SetUnitCandidates("micron mm cm m km");
285
286  possigmarCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/sigma_r",this);
287  possigmarCmd1->SetGuidance("Set standard deviation in radial of the beam positional profile");
288  possigmarCmd1->SetGuidance(" applicable to Beam type source only");
289  possigmarCmd1->SetParameterName("Sigmar",false,false);
290  possigmarCmd1->SetDefaultUnit("cm");
291  //  possigmarCmd1->SetUnitCandidates("micron mm cm m km");
292
293  possigmaxCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/sigma_x",this);
294  possigmaxCmd1->SetGuidance("Set standard deviation of beam positional profile in x-dir");
295  possigmaxCmd1->SetGuidance(" applicable to Beam type source only");
296  possigmaxCmd1->SetParameterName("Sigmax",false,false);
297  possigmaxCmd1->SetDefaultUnit("cm");
298  //  possigmaxCmd1->SetUnitCandidates("micron mm cm m km");
299
300  possigmayCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/sigma_y",this);
301  possigmayCmd1->SetGuidance("Set standard deviation of beam positional profile in y-dir");
302  possigmayCmd1->SetGuidance(" applicable to Beam type source only");
303  possigmayCmd1->SetParameterName("Sigmay",false,false);
304  possigmayCmd1->SetDefaultUnit("cm");
305  //  possigmayCmd1->SetUnitCandidates("micron mm cm m km");
306
307  paralpCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/paralp",this);
308  paralpCmd1->SetGuidance("Angle from y-axis of y' in Para");
309  paralpCmd1->SetParameterName("paralp",false,false);
310  paralpCmd1->SetDefaultUnit("rad");
311  //  paralpCmd1->SetUnitCandidates("rad deg");
312
313  partheCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/parthe",this);
314  partheCmd1->SetGuidance("Polar angle through centres of z faces");
315  partheCmd1->SetParameterName("parthe",false,false);
316  partheCmd1->SetDefaultUnit("rad");
317  //  partheCmd1->SetUnitCandidates("rad deg");
318
319  parphiCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/parphi",this);
320  parphiCmd1->SetGuidance("Azimuth angle through centres of z faces");
321  parphiCmd1->SetParameterName("parphi",false,false);
322  parphiCmd1->SetDefaultUnit("rad");
323  //  parphiCmd1->SetUnitCandidates("rad deg");
324
325  confineCmd1 = new G4UIcmdWithAString("/gps/pos/confine",this);
326  confineCmd1->SetGuidance("Confine source to volume (NULL to unset).");
327  confineCmd1->SetGuidance("usage: confine VolName");
328  confineCmd1->SetParameterName("VolName",false,false);
329  confineCmd1->SetDefaultValue("NULL");
330
331  // old implementations
332  typeCmd = new G4UIcmdWithAString("/gps/type",this);
333  typeCmd->SetGuidance("Sets source distribution type. (obsolete!)");
334  typeCmd->SetGuidance("Either Point, Beam, Plane, Surface or Volume");
335  typeCmd->SetParameterName("DisType",false,false);
336  typeCmd->SetDefaultValue("Point");
337  typeCmd->SetCandidates("Point Beam Plane Surface Volume");
338
339  shapeCmd = new G4UIcmdWithAString("/gps/shape",this);
340  shapeCmd->SetGuidance("Sets source shape type.(obsolete!)");
341  shapeCmd->SetParameterName("Shape",false,false);
342  shapeCmd->SetDefaultValue("NULL");
343  shapeCmd->SetCandidates("Circle Annulus Ellipse Square Rectangle Sphere Ellipsoid Cylinder Para");
344
345  centreCmd = new G4UIcmdWith3VectorAndUnit("/gps/centre",this);
346  centreCmd->SetGuidance("Set centre coordinates of source.(obsolete!)");
347  centreCmd->SetParameterName("X","Y","Z",false,false);
348  centreCmd->SetDefaultUnit("cm");
349  //  centreCmd->SetUnitCandidates("micron mm cm m km");
350
351  posrot1Cmd = new G4UIcmdWith3Vector("/gps/posrot1",this);
352  posrot1Cmd->SetGuidance("Set rotation matrix of x'.(obsolete!)");
353  posrot1Cmd->SetGuidance("Posrot1 does not need to be a unit vector.");
354  posrot1Cmd->SetParameterName("R1x","R1y","R1z",false,false); 
355  posrot1Cmd->SetRange("R1x != 0 || R1y != 0 || R1z != 0");
356
357  posrot2Cmd = new G4UIcmdWith3Vector("/gps/posrot2",this);
358  posrot2Cmd->SetGuidance("Set rotation matrix of y'.(obsolete!)");
359  posrot2Cmd->SetGuidance("Posrot2 does not need to be a unit vector.");
360  posrot2Cmd->SetParameterName("R2x","R2y","R2z",false,false); 
361  posrot2Cmd->SetRange("R2x != 0 || R2y != 0 || R2z != 0");
362
363  halfxCmd = new G4UIcmdWithADoubleAndUnit("/gps/halfx",this);
364  halfxCmd->SetGuidance("Set x half length of source.(obsolete!)");
365  halfxCmd->SetParameterName("Halfx",false,false);
366  halfxCmd->SetDefaultUnit("cm");
367  //  halfxCmd->SetUnitCandidates("micron mm cm m km");
368
369  halfyCmd = new G4UIcmdWithADoubleAndUnit("/gps/halfy",this);
370  halfyCmd->SetGuidance("Set y half length of source.(obsolete!)");
371  halfyCmd->SetParameterName("Halfy",false,false);
372  halfyCmd->SetDefaultUnit("cm");
373  //  halfyCmd->SetUnitCandidates("micron mm cm m km");
374
375  halfzCmd = new G4UIcmdWithADoubleAndUnit("/gps/halfz",this);
376  halfzCmd->SetGuidance("Set z half length of source.(obsolete!)");
377  halfzCmd->SetParameterName("Halfz",false,false);
378  halfzCmd->SetDefaultUnit("cm");
379  //  halfzCmd->SetUnitCandidates("micron mm cm m km");
380
381  radiusCmd = new G4UIcmdWithADoubleAndUnit("/gps/radius",this);
382  radiusCmd->SetGuidance("Set radius of source.(obsolete!)");
383  radiusCmd->SetParameterName("Radius",false,false);
384  radiusCmd->SetDefaultUnit("cm");
385  //  radiusCmd->SetUnitCandidates("micron mm cm m km");
386
387  radius0Cmd = new G4UIcmdWithADoubleAndUnit("/gps/radius0",this);
388  radius0Cmd->SetGuidance("Set inner radius of source.(obsolete!)");
389  radius0Cmd->SetParameterName("Radius0",false,false);
390  radius0Cmd->SetDefaultUnit("cm");
391  //  radius0Cmd->SetUnitCandidates("micron mm cm m km");
392
393  possigmarCmd = new G4UIcmdWithADoubleAndUnit("/gps/sigmaposr",this);
394  possigmarCmd->SetGuidance("Set standard deviation of beam position in radial(obsolete!)");
395  possigmarCmd->SetParameterName("Sigmar",false,false);
396  possigmarCmd->SetDefaultUnit("cm");
397  // possigmarCmd->SetUnitCandidates("micron mm cm m km");
398
399  possigmaxCmd = new G4UIcmdWithADoubleAndUnit("/gps/sigmaposx",this);
400  possigmaxCmd->SetGuidance("Set standard deviation of beam position in x-dir(obsolete!)");
401  possigmaxCmd->SetParameterName("Sigmax",false,false);
402  possigmaxCmd->SetDefaultUnit("cm");
403  //  possigmaxCmd->SetUnitCandidates("micron mm cm m km");
404
405  possigmayCmd = new G4UIcmdWithADoubleAndUnit("/gps/sigmaposy",this);
406  possigmayCmd->SetGuidance("Set standard deviation of beam position in y-dir(obsolete!)");
407  possigmayCmd->SetParameterName("Sigmay",false,false);
408  possigmayCmd->SetDefaultUnit("cm");
409  //  possigmayCmd->SetUnitCandidates("micron mm cm m km");
410
411  paralpCmd = new G4UIcmdWithADoubleAndUnit("/gps/paralp",this);
412  paralpCmd->SetGuidance("Angle from y-axis of y' in Para(obsolete!)");
413  paralpCmd->SetParameterName("paralp",false,false);
414  paralpCmd->SetDefaultUnit("rad");
415  //  paralpCmd->SetUnitCandidates("rad deg");
416
417  partheCmd = new G4UIcmdWithADoubleAndUnit("/gps/parthe",this);
418  partheCmd->SetGuidance("Polar angle through centres of z faces(obsolete!)");
419  partheCmd->SetParameterName("parthe",false,false);
420  partheCmd->SetDefaultUnit("rad");
421  //  partheCmd->SetUnitCandidates("rad deg");
422
423  parphiCmd = new G4UIcmdWithADoubleAndUnit("/gps/parphi",this);
424  parphiCmd->SetGuidance("Azimuth angle through centres of z faces(obsolete!)");
425  parphiCmd->SetParameterName("parphi",false,false);
426  parphiCmd->SetDefaultUnit("rad");
427  //  parphiCmd->SetUnitCandidates("rad deg");
428
429  confineCmd = new G4UIcmdWithAString("/gps/confine",this);
430  confineCmd->SetGuidance("Confine source to volume (NULL to unset)(obsolete!) .");
431  confineCmd->SetGuidance("usage: confine VolName");
432  confineCmd->SetParameterName("VolName",false,false);
433  confineCmd->SetDefaultValue("NULL");
434
435  // Angular distribution commands
436  angularDirectory = new G4UIdirectory("/gps/ang/");
437  angularDirectory->SetGuidance("Angular commands sub-directory");
438
439  angtypeCmd1 = new G4UIcmdWithAString("/gps/ang/type",this);
440  angtypeCmd1->SetGuidance("Sets angular source distribution type");
441  angtypeCmd1->SetGuidance("Possible variables are: iso, cos, planar, beam1d, beam2d, focused or user");
442  angtypeCmd1->SetParameterName("AngDis",false,false);
443  angtypeCmd1->SetDefaultValue("iso");
444  angtypeCmd1->SetCandidates("iso cos planar beam1d beam2d focused user");
445
446  angrot1Cmd1 = new G4UIcmdWith3Vector("/gps/ang/rot1",this);
447  angrot1Cmd1->SetGuidance("Sets the 1st vector for angular distribution rotation matrix");
448  angrot1Cmd1->SetGuidance("Need not be a unit vector");
449  angrot1Cmd1->SetParameterName("AR1x","AR1y","AR1z",false,false);
450  angrot1Cmd1->SetRange("AR1x != 0 || AR1y != 0 || AR1z != 0");
451
452  angrot2Cmd1 = new G4UIcmdWith3Vector("/gps/ang/rot2",this);
453  angrot2Cmd1->SetGuidance("Sets the 2nd vector for angular distribution rotation matrix");
454  angrot2Cmd1->SetGuidance("Need not be a unit vector");
455  angrot2Cmd1->SetParameterName("AR2x","AR2y","AR2z",false,false);
456  angrot2Cmd1->SetRange("AR2x != 0 || AR2y != 0 || AR2z != 0");
457
458  minthetaCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ang/mintheta",this);
459  minthetaCmd1->SetGuidance("Set minimum theta");
460  minthetaCmd1->SetParameterName("MinTheta",false,false);
461  minthetaCmd1->SetDefaultValue(0.);
462  minthetaCmd1->SetDefaultUnit("rad");
463  //  minthetaCmd1->SetUnitCandidates("rad deg");
464
465  maxthetaCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ang/maxtheta",this);
466  maxthetaCmd1->SetGuidance("Set maximum theta");
467  maxthetaCmd1->SetParameterName("MaxTheta",false,false);
468  maxthetaCmd1->SetDefaultValue(pi);
469  maxthetaCmd1->SetDefaultUnit("rad");
470  //  maxthetaCmd1->SetUnitCandidates("rad deg");
471
472  minphiCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ang/minphi",this);
473  minphiCmd1->SetGuidance("Set minimum phi");
474  minphiCmd1->SetParameterName("MinPhi",false,false);
475  minphiCmd1->SetDefaultUnit("rad");
476  //  minphiCmd1->SetUnitCandidates("rad deg");
477
478  maxphiCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ang/maxphi",this);
479  maxphiCmd1->SetGuidance("Set maximum phi");
480  maxphiCmd1->SetParameterName("MaxPhi",false,false);
481  maxphiCmd1->SetDefaultValue(pi);
482  maxphiCmd1->SetDefaultUnit("rad");
483  //  maxphiCmd1->SetUnitCandidates("rad deg");
484
485  angsigmarCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ang/sigma_r",this);
486  angsigmarCmd1->SetGuidance("Set standard deviation in direction for 1D beam.");
487  angsigmarCmd1->SetParameterName("Sigmara",false,false);
488  angsigmarCmd1->SetDefaultUnit("rad");
489  //  angsigmarCmd1->SetUnitCandidates("rad deg");
490 
491  angsigmaxCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ang/sigma_x",this);
492  angsigmaxCmd1->SetGuidance("Set standard deviation in direction in x-direc. for 2D beam");
493  angsigmaxCmd1->SetParameterName("Sigmaxa",false,false);
494  angsigmaxCmd1->SetDefaultUnit("rad");
495  //  angsigmaxCmd1->SetUnitCandidates("rad deg");
496
497  angsigmayCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ang/sigma_y",this);
498  angsigmayCmd1->SetGuidance("Set standard deviation in direction in y-direc. for 2D beam");
499  angsigmayCmd1->SetParameterName("Sigmaya",false,false);
500  angsigmayCmd1->SetDefaultUnit("rad");
501  //  angsigmayCmd1->SetUnitCandidates("rad deg");
502
503  angfocusCmd = new G4UIcmdWith3VectorAndUnit("/gps/ang/focuspoint",this);
504  angfocusCmd->SetGuidance("Set the focusing point for the beam");
505  angfocusCmd->SetParameterName("x","y","z",false,false);
506  angfocusCmd->SetDefaultUnit("cm");
507  //  angfocusCmd->SetUnitCandidates("micron mm cm m km");
508
509  useuserangaxisCmd1 = new G4UIcmdWithABool("/gps/ang/user_coor",this);
510  useuserangaxisCmd1->SetGuidance("true for using user defined angular co-ordinates");
511  useuserangaxisCmd1->SetGuidance("Default is false");
512  useuserangaxisCmd1->SetParameterName("useuserangaxis",true);
513  useuserangaxisCmd1->SetDefaultValue(false);
514
515  surfnormCmd1 = new G4UIcmdWithABool("/gps/ang/surfnorm",this);
516  surfnormCmd1->SetGuidance("Makes a user-defined distribution with respect to surface normals rather than x,y,z axes.");
517  surfnormCmd1->SetGuidance("Default is false");
518  surfnormCmd1->SetParameterName("surfnorm",true);
519  surfnormCmd1->SetDefaultValue(false);
520
521  // old ones
522  angtypeCmd = new G4UIcmdWithAString("/gps/angtype",this);
523  angtypeCmd->SetGuidance("Sets angular source distribution type (obsolete!)");
524  angtypeCmd->SetGuidance("Possible variables are: iso, cos planar beam1d beam2d or user");
525  angtypeCmd->SetParameterName("AngDis",false,false);
526  angtypeCmd->SetDefaultValue("iso");
527  angtypeCmd->SetCandidates("iso cos planar beam1d beam2d user");
528
529  angrot1Cmd = new G4UIcmdWith3Vector("/gps/angrot1",this);
530  angrot1Cmd->SetGuidance("Sets the x' vector for angular distribution(obsolete!) ");
531  angrot1Cmd->SetGuidance("Need not be a unit vector");
532  angrot1Cmd->SetParameterName("AR1x","AR1y","AR1z",false,false);
533  angrot1Cmd->SetRange("AR1x != 0 || AR1y != 0 || AR1z != 0");
534
535  angrot2Cmd = new G4UIcmdWith3Vector("/gps/angrot2",this);
536  angrot2Cmd->SetGuidance("Sets the y' vector for angular distribution (obsolete!)");
537  angrot2Cmd->SetGuidance("Need not be a unit vector");
538  angrot2Cmd->SetParameterName("AR2x","AR2y","AR2z",false,false);
539  angrot2Cmd->SetRange("AR2x != 0 || AR2y != 0 || AR2z != 0");
540
541  minthetaCmd = new G4UIcmdWithADoubleAndUnit("/gps/mintheta",this);
542  minthetaCmd->SetGuidance("Set minimum theta (obsolete!)");
543  minthetaCmd->SetParameterName("MinTheta",false,false);
544  minthetaCmd->SetDefaultUnit("rad");
545  //  minthetaCmd->SetUnitCandidates("rad deg");
546
547  maxthetaCmd = new G4UIcmdWithADoubleAndUnit("/gps/maxtheta",this);
548  maxthetaCmd->SetGuidance("Set maximum theta (obsolete!)");
549  maxthetaCmd->SetParameterName("MaxTheta",false,false);
550  maxthetaCmd->SetDefaultValue(3.1416);
551  maxthetaCmd->SetDefaultUnit("rad");
552  //  maxthetaCmd->SetUnitCandidates("rad deg");
553
554  minphiCmd = new G4UIcmdWithADoubleAndUnit("/gps/minphi",this);
555  minphiCmd->SetGuidance("Set minimum phi (obsolete!)");
556  minphiCmd->SetParameterName("MinPhi",false,false);
557  minphiCmd->SetDefaultUnit("rad");
558  //  minphiCmd->SetUnitCandidates("rad deg");
559
560  maxphiCmd = new G4UIcmdWithADoubleAndUnit("/gps/maxphi",this);
561  maxphiCmd->SetGuidance("Set maximum phi(obsolete!)");
562  maxphiCmd->SetParameterName("MaxPhi",false,false);
563  maxphiCmd->SetDefaultUnit("rad");
564  //  maxphiCmd->SetUnitCandidates("rad deg");
565
566  angsigmarCmd = new G4UIcmdWithADoubleAndUnit("/gps/sigmaangr",this);
567  angsigmarCmd->SetGuidance("Set standard deviation of beam direction in radial(obsolete!).");
568  angsigmarCmd->SetParameterName("Sigmara",false,false);
569  angsigmarCmd->SetDefaultUnit("rad");
570  //  angsigmarCmd->SetUnitCandidates("rad deg");
571
572  angsigmaxCmd = new G4UIcmdWithADoubleAndUnit("/gps/sigmaangx",this);
573  angsigmaxCmd->SetGuidance("Set standard deviation of beam direction in x-direc(obsolete!).");
574  angsigmaxCmd->SetParameterName("Sigmaxa",false,false);
575  angsigmaxCmd->SetDefaultUnit("rad");
576  //  angsigmaxCmd->SetUnitCandidates("rad deg");
577
578  angsigmayCmd = new G4UIcmdWithADoubleAndUnit("/gps/sigmaangy",this);
579  angsigmayCmd->SetGuidance("Set standard deviation of beam direction in y-direc.(obsolete!)");
580  angsigmayCmd->SetParameterName("Sigmaya",false,false);
581  angsigmayCmd->SetDefaultUnit("rad");
582  //  angsigmayCmd->SetUnitCandidates("rad deg");
583
584  useuserangaxisCmd = new G4UIcmdWithABool("/gps/useuserangaxis",this);
585  useuserangaxisCmd->SetGuidance("true for using user defined angular co-ordinates(obsolete!)");
586  useuserangaxisCmd->SetGuidance("Default is false");
587  useuserangaxisCmd->SetParameterName("useuserangaxis",true);
588  useuserangaxisCmd->SetDefaultValue(false);
589
590  surfnormCmd = new G4UIcmdWithABool("/gps/surfnorm",this);
591  surfnormCmd->SetGuidance("Makes a user-defined distribution with respect to surface normals rather than x,y,z axes (obsolete!).");
592  surfnormCmd->SetGuidance("Default is false");
593  surfnormCmd->SetParameterName("surfnorm",true);
594  surfnormCmd->SetDefaultValue(false);
595
596  // Energy commands
597
598  energyDirectory = new G4UIdirectory("/gps/ene/");
599  energyDirectory->SetGuidance("Spectral commands sub-directory");
600
601  energytypeCmd1 = new G4UIcmdWithAString("/gps/ene/type",this);
602  energytypeCmd1->SetGuidance("Sets energy distribution type");
603  energytypeCmd1->SetParameterName("EnergyDis",false,false);
604  energytypeCmd1->SetDefaultValue("Mono");
605  energytypeCmd1->SetCandidates("Mono Lin Pow Exp Gauss Brem Bbody Cdg User Arb Epn");
606
607  eminCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ene/min",this);
608  eminCmd1->SetGuidance("Sets minimum energy");
609  eminCmd1->SetParameterName("emin",false,false);
610  eminCmd1->SetDefaultUnit("keV");
611  //  eminCmd1->SetUnitCandidates("eV keV MeV GeV TeV PeV");
612
613  emaxCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ene/max",this);
614  emaxCmd1->SetGuidance("Sets maximum energy");
615  emaxCmd1->SetParameterName("emax",false,false);
616  emaxCmd1->SetDefaultUnit("keV");
617  //  emaxCmd1->SetUnitCandidates("eV keV MeV GeV TeV PeV");
618
619  monoenergyCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ene/mono",this);
620  monoenergyCmd1->SetGuidance("Sets a monocromatic energy (same as  gps/energy)");
621  monoenergyCmd1->SetParameterName("monoenergy",false,false);
622  monoenergyCmd1->SetDefaultUnit("keV");
623  //  monoenergyCmd1->SetUnitCandidates("eV keV MeV GeV TeV PeV");
624
625  engsigmaCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ene/sigma",this);
626  engsigmaCmd1->SetGuidance("Sets the standard deviation for Gaussian energy dist.");
627  engsigmaCmd1->SetParameterName("Sigmae",false,false);
628  engsigmaCmd1->SetDefaultUnit("keV");
629  //  engsigmaCmd1->SetUnitCandidates("eV keV MeV GeV TeV PeV");
630
631  alphaCmd1 = new G4UIcmdWithADouble("/gps/ene/alpha",this);
632  alphaCmd1->SetGuidance("Sets Alpha (index) for power-law energy dist.");
633  alphaCmd1->SetParameterName("alpha",false,false);
634 
635  tempCmd1 = new G4UIcmdWithADouble("/gps/ene/temp",this);
636  tempCmd1->SetGuidance("Sets the temperature for Brem and BBody distributions (in Kelvin)");
637  tempCmd1->SetParameterName("temp",false,false);
638
639  ezeroCmd1 = new G4UIcmdWithADouble("/gps/ene/ezero",this);
640  ezeroCmd1->SetGuidance("Sets E_0 for exponential distribution (in MeV)");
641  ezeroCmd1->SetParameterName("ezero",false,false);
642
643  gradientCmd1 = new G4UIcmdWithADouble("/gps/ene/gradient",this);
644  gradientCmd1->SetGuidance("Sets the gradient for Lin distribution (in 1/MeV)");
645  gradientCmd1->SetParameterName("gradient",false,false);
646
647  interceptCmd1 = new G4UIcmdWithADouble("/gps/ene/intercept",this);
648  interceptCmd1->SetGuidance("Sets the intercept for Lin distributions (in MeV)");
649  interceptCmd1->SetParameterName("intercept",false,false);
650
651  arbeintCmd1 = new G4UIcmdWithADouble("/gps/ene/biasAlpha",this);
652  arbeintCmd1->SetGuidance("Set the power-law index for the energy sampling distri. )");
653  arbeintCmd1->SetParameterName("arbeint",false,false);
654
655  calculateCmd1 = new G4UIcmdWithoutParameter("/gps/ene/calculate",this);
656  calculateCmd1->SetGuidance("Calculates the distributions for Cdg and BBody");
657
658  energyspecCmd1 = new G4UIcmdWithABool("/gps/ene/emspec",this);
659  energyspecCmd1->SetGuidance("True for energy and false for momentum spectra");
660  energyspecCmd1->SetParameterName("energyspec",true);
661  energyspecCmd1->SetDefaultValue(true);
662
663  diffspecCmd1 = new G4UIcmdWithABool("/gps/ene/diffspec",this);
664  diffspecCmd1->SetGuidance("True for differential and flase for integral spectra");
665  diffspecCmd1->SetParameterName("diffspec",true);
666  diffspecCmd1->SetDefaultValue(true);
667
668  //old ones
669  energytypeCmd = new G4UIcmdWithAString("/gps/energytype",this);
670  energytypeCmd->SetGuidance("Sets energy distribution type (obsolete!)");
671  energytypeCmd->SetParameterName("EnergyDis",false,false);
672  energytypeCmd->SetDefaultValue("Mono");
673  energytypeCmd->SetCandidates("Mono Lin Pow Exp Gauss Brem Bbody Cdg User Arb Epn");
674
675  eminCmd = new G4UIcmdWithADoubleAndUnit("/gps/emin",this);
676  eminCmd->SetGuidance("Sets Emin (obsolete!)");
677  eminCmd->SetParameterName("emin",false,false);
678  eminCmd->SetDefaultUnit("keV");
679  //  eminCmd->SetUnitCandidates("eV keV MeV GeV TeV PeV");
680
681  emaxCmd = new G4UIcmdWithADoubleAndUnit("/gps/emax",this);
682  emaxCmd->SetGuidance("Sets Emax (obsolete!)");
683  emaxCmd->SetParameterName("emax",false,false);
684  emaxCmd->SetDefaultUnit("keV");
685  //  emaxCmd->SetUnitCandidates("eV keV MeV GeV TeV PeV");
686
687  monoenergyCmd = new G4UIcmdWithADoubleAndUnit("/gps/monoenergy",this);
688  monoenergyCmd->SetGuidance("Sets Monoenergy (obsolete, use gps/energy instead!)");
689  monoenergyCmd->SetParameterName("monoenergy",false,false);
690  monoenergyCmd->SetDefaultUnit("keV");
691  //  monoenergyCmd->SetUnitCandidates("eV keV MeV GeV TeV PeV");
692
693  engsigmaCmd = new G4UIcmdWithADoubleAndUnit("/gps/sigmae",this);
694  engsigmaCmd->SetGuidance("Sets the standard deviation for Gaussian energy dist.(obsolete!)");
695  engsigmaCmd->SetParameterName("Sigmae",false,false);
696  engsigmaCmd->SetDefaultUnit("keV");
697  //  engsigmaCmd->SetUnitCandidates("eV keV MeV GeV TeV PeV");
698
699  alphaCmd = new G4UIcmdWithADouble("/gps/alpha",this);
700  alphaCmd->SetGuidance("Sets Alpha (index) for power-law energy dist(obsolete!).");
701  alphaCmd->SetParameterName("alpha",false,false);
702 
703  tempCmd = new G4UIcmdWithADouble("/gps/temp",this);
704  tempCmd->SetGuidance("Sets the temperature for Brem and BBody (in Kelvin)(obsolete!)");
705  tempCmd->SetParameterName("temp",false,false);
706
707  ezeroCmd = new G4UIcmdWithADouble("/gps/ezero",this);
708  ezeroCmd->SetGuidance("Sets ezero exponential distributions (in MeV)(obsolete!)");
709  ezeroCmd->SetParameterName("ezero",false,false);
710
711  gradientCmd = new G4UIcmdWithADouble("/gps/gradient",this);
712  gradientCmd->SetGuidance("Sets the gradient for Lin distributions (in 1/MeV)(obsolete!)");
713  gradientCmd->SetParameterName("gradient",false,false);
714
715  interceptCmd = new G4UIcmdWithADouble("/gps/intercept",this);
716  interceptCmd->SetGuidance("Sets the intercept for Lin distributions (in MeV)(obsolete!)");
717  interceptCmd->SetParameterName("intercept",false,false);
718
719  calculateCmd = new G4UIcmdWithoutParameter("/gps/calculate",this);
720  calculateCmd->SetGuidance("Calculates distributions for Cdg and BBody(obsolete!)");
721
722  energyspecCmd = new G4UIcmdWithABool("/gps/energyspec",this);
723  energyspecCmd->SetGuidance("True for energy and false for momentum spectra(obsolete!)");
724  energyspecCmd->SetParameterName("energyspec",true);
725  energyspecCmd->SetDefaultValue(true);
726
727  diffspecCmd = new G4UIcmdWithABool("/gps/diffspec",this);
728  diffspecCmd->SetGuidance("True for differential and flase for integral spectra(obsolete!)");
729  diffspecCmd->SetParameterName("diffspec",true);
730  diffspecCmd->SetDefaultValue(true);
731
732  // Biasing + histograms in general
733  histDirectory = new G4UIdirectory("/gps/hist/");
734  histDirectory->SetGuidance("Histogram, biasing commands sub-directory");
735
736  histnameCmd1 = new G4UIcmdWithAString("/gps/hist/type",this);
737  histnameCmd1->SetGuidance("Sets histogram type");
738  histnameCmd1->SetParameterName("HistType",false,false);
739  histnameCmd1->SetDefaultValue("biasx");
740  histnameCmd1->SetCandidates("biasx biasy biasz biast biasp biase biaspt biaspp theta phi energy arb epn");
741
742  resethistCmd1 = new G4UIcmdWithAString("/gps/hist/reset",this);
743  resethistCmd1->SetGuidance("Reset (clean) the histogram ");
744  resethistCmd1->SetParameterName("HistType",false,false);
745  resethistCmd1->SetDefaultValue("energy");
746  resethistCmd1->SetCandidates("biasx biasy biasz biast biasp biase biaspt biaspp theta phi energy arb epn");
747
748  histpointCmd1 = new G4UIcmdWith3Vector("/gps/hist/point",this);
749  histpointCmd1->SetGuidance("Allows user to define a histogram");
750  histpointCmd1->SetGuidance("Enter: Ehi Weight");
751  histpointCmd1->SetParameterName("Ehi","Weight","Junk",true,true);
752  histpointCmd1->SetRange("Ehi >= 0. && Weight >= 0.");
753
754  histfileCmd1 = new G4UIcmdWithAString("/gps/hist/file",this);
755  histfileCmd1->SetGuidance("import the arb energy hist in an ASCII file");
756  histfileCmd1->SetParameterName("HistFile",false,false);
757
758  arbintCmd1 = new G4UIcmdWithAString("/gps/hist/inter",this);
759  arbintCmd1->SetGuidance("Sets the interpolation method for arbitrary distribution.");
760  arbintCmd1->SetParameterName("int",false,false);
761  arbintCmd1->SetDefaultValue("Lin");
762  arbintCmd1->SetCandidates("Lin Log Exp Spline");
763
764  // old ones
765  histnameCmd = new G4UIcmdWithAString("/gps/histname",this);
766  histnameCmd->SetGuidance("Sets histogram type (obsolete!)");
767  histnameCmd->SetParameterName("HistType",false,false);
768  histnameCmd->SetDefaultValue("biasx");
769  histnameCmd->SetCandidates("biasx biasy biasz biast biasp biase biaspt biaspp theta phi energy arb epn");
770
771  // re-set the histograms
772  resethistCmd = new G4UIcmdWithAString("/gps/resethist",this);
773  resethistCmd->SetGuidance("Re-Set the histogram (obsolete!)");
774  resethistCmd->SetParameterName("HistType",false,false);
775  resethistCmd->SetDefaultValue("energy");
776  resethistCmd->SetCandidates("biasx biasy biasz biast biasp biase biaspt biaspp theta phi energy arb epn");
777
778  histpointCmd = new G4UIcmdWith3Vector("/gps/histpoint",this);
779  histpointCmd->SetGuidance("Allows user to define a histogram (obsolete!)");
780  histpointCmd->SetGuidance("Enter: Ehi Weight");
781  histpointCmd->SetParameterName("Ehi","Weight","Junk",false,false);
782  histpointCmd->SetRange("Ehi >= 0. && Weight >= 0.");
783
784  arbintCmd = new G4UIcmdWithAString("/gps/arbint",this);
785  arbintCmd->SetGuidance("Sets Arbitrary Interpolation type.(obsolete!) ");
786  arbintCmd->SetParameterName("int",false,false);
787  arbintCmd->SetDefaultValue("NULL");
788  arbintCmd->SetCandidates("Lin Log Exp Spline");
789
790}
791
792G4GeneralParticleSourceMessenger::~G4GeneralParticleSourceMessenger()
793{
794  delete positionDirectory;
795  delete typeCmd;
796  delete shapeCmd;
797  delete centreCmd;
798  delete posrot1Cmd;
799  delete posrot2Cmd;
800  delete halfxCmd;
801  delete halfyCmd;
802  delete halfzCmd;
803  delete radiusCmd;
804  delete radius0Cmd;
805  delete possigmarCmd;
806  delete possigmaxCmd;
807  delete possigmayCmd;
808  delete paralpCmd;
809  delete partheCmd;
810  delete parphiCmd;
811  delete confineCmd;
812  delete typeCmd1;
813  delete shapeCmd1;
814  delete centreCmd1;
815  delete posrot1Cmd1;
816  delete posrot2Cmd1;
817  delete halfxCmd1;
818  delete halfyCmd1;
819  delete halfzCmd1;
820  delete radiusCmd1;
821  delete radius0Cmd1;
822  delete possigmarCmd1;
823  delete possigmaxCmd1;
824  delete possigmayCmd1;
825  delete paralpCmd1;
826  delete partheCmd1;
827  delete parphiCmd1;
828  delete confineCmd1;
829
830  delete angularDirectory;
831  delete angtypeCmd;
832  delete angrot1Cmd;
833  delete angrot2Cmd;
834  delete minthetaCmd;
835  delete maxthetaCmd;
836  delete minphiCmd;
837  delete maxphiCmd;
838  delete angsigmarCmd;
839  delete angsigmaxCmd;
840  delete angsigmayCmd;
841  delete useuserangaxisCmd;
842  delete surfnormCmd;
843  delete angtypeCmd1;
844  delete angrot1Cmd1;
845  delete angrot2Cmd1;
846  delete minthetaCmd1;
847  delete maxthetaCmd1;
848  delete minphiCmd1;
849  delete maxphiCmd1;
850  delete angsigmarCmd1;
851  delete angsigmaxCmd1;
852  delete angfocusCmd;
853  delete useuserangaxisCmd1;
854  delete surfnormCmd1;
855
856  delete energyDirectory;
857  delete energytypeCmd;
858  delete eminCmd;
859  delete emaxCmd;
860  delete monoenergyCmd;
861  delete engsigmaCmd;
862  delete alphaCmd;
863  delete tempCmd;
864  delete ezeroCmd;
865  delete gradientCmd;
866  delete interceptCmd;
867  delete calculateCmd;
868  delete energyspecCmd;
869  delete diffspecCmd;
870  delete energytypeCmd1;
871  delete eminCmd1;
872  delete emaxCmd1;
873  delete monoenergyCmd1;
874  delete engsigmaCmd1;
875  delete alphaCmd1;
876  delete tempCmd1;
877  delete ezeroCmd1;
878  delete gradientCmd1;
879  delete interceptCmd1;
880  delete arbeintCmd1;
881  delete calculateCmd1;
882  delete energyspecCmd1;
883  delete diffspecCmd1;
884
885  delete histDirectory;
886  delete histnameCmd;
887  delete resethistCmd;
888  delete histpointCmd;
889  delete arbintCmd;
890  delete histnameCmd1;
891  delete resethistCmd1;
892  delete histpointCmd1;
893  delete histfileCmd1;
894  delete arbintCmd1;
895
896  delete verbosityCmd;
897  delete ionCmd;
898  delete particleCmd;
899  delete timeCmd;
900  delete polCmd;
901  delete numberCmd;
902  delete positionCmd;
903  delete directionCmd;
904  delete energyCmd;
905  delete listCmd;
906
907  delete sourceDirectory;
908  delete addsourceCmd;
909  delete listsourceCmd;
910  delete clearsourceCmd;
911  delete getsourceCmd;
912  delete setsourceCmd;
913  delete setintensityCmd;
914  delete deletesourceCmd;
915  delete multiplevertexCmd;
916  delete flatsamplingCmd;
917
918  delete gpsDirectory;
919 
920}
921
922void G4GeneralParticleSourceMessenger::SetNewValue(G4UIcommand *command, G4String newValues)
923{
924  if(command == typeCmd)
925    {
926      fParticleGun->GetPosDist()->SetPosDisType(newValues);
927      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
928             << " The command is obsolete and will be removed soon." << G4endl
929             << " Please try to use the new structured commands!" << G4endl;
930    }
931  else if(command == shapeCmd)
932    {
933      fParticleGun->GetPosDist()->SetPosDisShape(newValues);
934      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
935             << " The command is obsolete and will be removed soon." << G4endl
936             << " Please try to use the new structured commands!" << G4endl;
937    }
938  else if(command == centreCmd)
939    {
940      fParticleGun->GetPosDist()->SetCentreCoords(centreCmd->GetNew3VectorValue(newValues));
941      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
942             << " The command is obsolete and will be removed soon." << G4endl
943             << " Please try to use the new structured commands!" << G4endl;
944    }
945  else if(command == posrot1Cmd)
946    {
947      fParticleGun->GetPosDist()->SetPosRot1(posrot1Cmd->GetNew3VectorValue(newValues));
948      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
949             << " The command is obsolete and will be removed soon." << G4endl
950             << " Please try to use the new structured commands!" << G4endl;
951    }
952  else if(command == posrot2Cmd)
953    {
954      fParticleGun->GetPosDist()->SetPosRot2(posrot2Cmd->GetNew3VectorValue(newValues));
955      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
956             << " The command is obsolete and will be removed soon." << G4endl
957             << " Please try to use the new structured commands!" << G4endl;
958    }
959  else if(command == halfxCmd)
960    {
961      fParticleGun->GetPosDist()->SetHalfX(halfxCmd->GetNewDoubleValue(newValues));
962      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
963             << " The command is obsolete and will be removed soon." << G4endl
964             << " Please try to use the new structured commands!" << G4endl;
965    }
966  else if(command == halfyCmd)
967    {
968      fParticleGun->GetPosDist()->SetHalfY(halfyCmd->GetNewDoubleValue(newValues));
969      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
970             << " The command is obsolete and will be removed soon." << G4endl
971             << " Please try to use the new structured commands!" << G4endl;
972    }
973  else if(command == halfzCmd)
974    {
975      fParticleGun->GetPosDist()->SetHalfZ(halfzCmd->GetNewDoubleValue(newValues));
976      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
977             << " The command is obsolete and will be removed soon." << G4endl
978             << " Please try to use the new structured commands!" << G4endl;
979    }
980  else if(command == radiusCmd)
981    {
982      fParticleGun->GetPosDist()->SetRadius(radiusCmd->GetNewDoubleValue(newValues));
983      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
984             << " The command is obsolete and will be removed soon." << G4endl
985             << " Please try to use the new structured commands!" << G4endl;
986    }
987  else if(command == radius0Cmd)
988    {
989      fParticleGun->GetPosDist()->SetRadius0(radius0Cmd->GetNewDoubleValue(newValues));
990      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
991             << " The command is obsolete and will be removed soon." << G4endl
992             << " Please try to use the new structured commands!" << G4endl;
993    }
994  else if(command == possigmarCmd)
995    {
996      fParticleGun->GetPosDist()->SetBeamSigmaInR(possigmarCmd->GetNewDoubleValue(newValues));
997      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
998             << " The command is obsolete and will be removed soon." << G4endl
999             << " Please try to use the new structured commands!" << G4endl;
1000    }
1001  else if(command == possigmaxCmd)
1002    {
1003      fParticleGun->GetPosDist()->SetBeamSigmaInX(possigmaxCmd->GetNewDoubleValue(newValues));
1004      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1005             << " The command is obsolete and will be removed soon." << G4endl
1006             << " Please try to use the new structured commands!" << G4endl;
1007    }
1008  else if(command == possigmayCmd)
1009    {
1010      fParticleGun->GetPosDist()->SetBeamSigmaInY(possigmayCmd->GetNewDoubleValue(newValues));
1011      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1012             << " The command is obsolete and will be removed soon." << G4endl
1013             << " Please try to use the new structured commands!" << G4endl;
1014    }
1015  else if(command == paralpCmd)
1016    {
1017      fParticleGun->GetPosDist()->SetParAlpha(paralpCmd->GetNewDoubleValue(newValues));
1018      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1019             << " The command is obsolete and will be removed soon." << G4endl
1020             << " Please try to use the new structured commands!" << G4endl;
1021    }
1022  else if(command == partheCmd)
1023    {
1024      fParticleGun->GetPosDist()->SetParTheta(partheCmd->GetNewDoubleValue(newValues));
1025      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1026             << " The command is obsolete and will be removed soon." << G4endl
1027             << " Please try to use the new structured commands!" << G4endl;
1028    }
1029  else if(command == parphiCmd)
1030    {
1031      fParticleGun->GetPosDist()->SetParPhi(parphiCmd->GetNewDoubleValue(newValues));
1032      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1033             << " The command is obsolete and will be removed soon." << G4endl
1034             << " Please try to use the new structured commands!" << G4endl;
1035    }
1036  else if(command == confineCmd)
1037    {
1038      fParticleGun->GetPosDist()->ConfineSourceToVolume(newValues);
1039      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1040             << " The command is obsolete and will be removed soon." << G4endl
1041             << " Please try to use the new structured commands!" << G4endl;
1042    }
1043  else if(command == angtypeCmd)
1044    {
1045      fParticleGun->GetAngDist()->SetAngDistType(newValues);
1046      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1047             << " The command is obsolete and will be removed soon." << G4endl
1048             << " Please try to use the new structured commands!" << G4endl;
1049    }
1050  else if(command == angrot1Cmd)
1051    {
1052      G4String a = "angref1";
1053      fParticleGun->GetAngDist()->DefineAngRefAxes(a,angrot1Cmd->GetNew3VectorValue(newValues));
1054      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1055             << " The command is obsolete and will be removed soon." << G4endl
1056             << " Please try to use the new structured commands!" << G4endl;
1057    }
1058  else if(command == angrot2Cmd)
1059    {
1060      G4String a = "angref2";
1061      fParticleGun->GetAngDist()->DefineAngRefAxes(a,angrot2Cmd->GetNew3VectorValue(newValues));
1062      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1063             << " The command is obsolete and will be removed soon." << G4endl
1064             << " Please try to use the new structured commands!" << G4endl;
1065    }
1066  else if(command == minthetaCmd)
1067    {
1068      fParticleGun->GetAngDist()->SetMinTheta(minthetaCmd->GetNewDoubleValue(newValues));
1069      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1070             << " The command is obsolete and will be removed soon." << G4endl
1071             << " Please try to use the new structured commands!" << G4endl;
1072    }
1073  else if(command == minphiCmd)
1074    {
1075      fParticleGun->GetAngDist()->SetMinPhi(minphiCmd->GetNewDoubleValue(newValues));
1076      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1077             << " The command is obsolete and will be removed soon." << G4endl
1078             << " Please try to use the new structured commands!" << G4endl;
1079    }
1080  else if(command == maxthetaCmd)
1081    {
1082      fParticleGun->GetAngDist()->SetMaxTheta(maxthetaCmd->GetNewDoubleValue(newValues));
1083      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1084             << " The command is obsolete and will be removed soon." << G4endl
1085             << " Please try to use the new structured commands!" << G4endl;
1086    }
1087  else if(command == maxphiCmd)
1088    {
1089      fParticleGun->GetAngDist()->SetMaxPhi(maxphiCmd->GetNewDoubleValue(newValues));
1090      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1091             << " The command is obsolete and will be removed soon." << G4endl
1092             << " Please try to use the new structured commands!" << G4endl;
1093    }
1094  else if(command == angsigmarCmd)
1095    {
1096      fParticleGun->GetAngDist()->SetBeamSigmaInAngR(angsigmarCmd->GetNewDoubleValue(newValues));
1097      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1098             << " The command is obsolete and will be removed soon." << G4endl
1099             << " Please try to use the new structured commands!" << G4endl;
1100    }
1101  else if(command == angsigmaxCmd)
1102    {
1103      fParticleGun->GetAngDist()->SetBeamSigmaInAngX(angsigmaxCmd->GetNewDoubleValue(newValues));
1104      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1105             << " The command is obsolete and will be removed soon." << G4endl
1106             << " Please try to use the new structured commands!" << G4endl;
1107    }
1108  else if(command == angsigmayCmd)
1109    {
1110      fParticleGun->GetAngDist()->SetBeamSigmaInAngY(angsigmayCmd->GetNewDoubleValue(newValues));
1111      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1112             << " The command is obsolete and will be removed soon." << G4endl
1113             << " Please try to use the new structured commands!" << G4endl;
1114    }
1115  else if(command == useuserangaxisCmd)
1116    {
1117      fParticleGun->GetAngDist()->SetUseUserAngAxis(useuserangaxisCmd->GetNewBoolValue(newValues));
1118      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1119             << " The command is obsolete and will be removed soon." << G4endl
1120             << " Please try to use the new structured commands!" << G4endl;
1121    }
1122  else if(command == surfnormCmd)
1123    {
1124      fParticleGun->GetAngDist()->SetUserWRTSurface(surfnormCmd->GetNewBoolValue(newValues));
1125      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1126             << " The command is obsolete and will be removed soon." << G4endl
1127             << " Please try to use the new structured commands!" << G4endl;
1128    }
1129  else if(command == energytypeCmd)
1130    {
1131      fParticleGun->GetEneDist()->SetEnergyDisType(newValues);
1132      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1133             << " The command is obsolete and will be removed soon." << G4endl
1134             << " Please try to use the new structured commands!" << G4endl;
1135    }
1136  else if(command == eminCmd)
1137    {
1138      fParticleGun->GetEneDist()->SetEmin(eminCmd->GetNewDoubleValue(newValues));
1139      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1140             << " The command is obsolete and will be removed soon." << G4endl
1141             << " Please try to use the new structured commands!" << G4endl;
1142    }
1143  else if(command == emaxCmd)
1144    {
1145      fParticleGun->GetEneDist()->SetEmax(emaxCmd->GetNewDoubleValue(newValues));
1146      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1147             << " The command is obsolete and will be removed soon." << G4endl
1148             << " Please try to use the new structured commands!" << G4endl;
1149    }
1150  else if(command == monoenergyCmd)
1151    {
1152      fParticleGun->GetEneDist()->SetMonoEnergy(monoenergyCmd->GetNewDoubleValue(newValues));
1153      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1154             << " The command is obsolete and will be removed soon." << G4endl
1155             << " Please try to use the new structured commands!" << G4endl;
1156    }
1157  else if(command == engsigmaCmd)
1158    {
1159      fParticleGun->GetEneDist()->SetBeamSigmaInE(engsigmaCmd->GetNewDoubleValue(newValues));
1160      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1161             << " The command is obsolete and will be removed soon." << G4endl
1162             << " Please try to use the new structured commands!" << G4endl;
1163    }
1164  else if(command == alphaCmd)
1165    {
1166      fParticleGun->GetEneDist()->SetAlpha(alphaCmd->GetNewDoubleValue(newValues));
1167      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1168             << " The command is obsolete and will be removed soon." << G4endl
1169             << " Please try to use the new structured commands!" << G4endl;
1170    }
1171  else if(command == tempCmd)
1172    {
1173      fParticleGun->GetEneDist()->SetTemp(tempCmd->GetNewDoubleValue(newValues));
1174      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1175             << " The command is obsolete and will be removed soon." << G4endl
1176             << " Please try to use the new structured commands!" << G4endl;
1177    }
1178  else if(command == ezeroCmd)
1179    {
1180      fParticleGun->GetEneDist()->SetEzero(ezeroCmd->GetNewDoubleValue(newValues));
1181      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1182             << " The command is obsolete and will be removed soon." << G4endl
1183             << " Please try to use the new structured commands!" << G4endl;
1184    }
1185  else if(command == gradientCmd)
1186    {
1187      fParticleGun->GetEneDist()->SetGradient(gradientCmd->GetNewDoubleValue(newValues));
1188      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1189             << " The command is obsolete and will be removed soon." << G4endl
1190             << " Please try to use the new structured commands!" << G4endl;
1191    }
1192  else if(command == interceptCmd)
1193    {
1194      fParticleGun->GetEneDist()->SetInterCept(interceptCmd->GetNewDoubleValue(newValues));
1195      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1196             << " The command is obsolete and will be removed soon." << G4endl
1197             << " Please try to use the new structured commands!" << G4endl;
1198    }
1199  else if(command == calculateCmd)
1200    {
1201      fParticleGun->GetEneDist()->Calculate();
1202      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1203             << " The command is obsolete and will be removed soon." << G4endl
1204             << " Please try to use the new structured commands!" << G4endl;
1205    }
1206  else if(command == energyspecCmd)
1207    {
1208      fParticleGun->GetEneDist()->InputEnergySpectra(energyspecCmd->GetNewBoolValue(newValues));
1209      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1210             << " The command is obsolete and will be removed soon." << G4endl
1211             << " Please try to use the new structured commands!" << G4endl;
1212    }
1213  else if(command == diffspecCmd)
1214    {
1215      fParticleGun->GetEneDist()->InputDifferentialSpectra(diffspecCmd->GetNewBoolValue(newValues));
1216      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1217             << " The command is obsolete and will be removed soon." << G4endl
1218             << " Please try to use the new structured commands!" << G4endl;
1219    }
1220  else if(command == histnameCmd)
1221    {
1222      histtype = newValues;
1223      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1224             << " The command is obsolete and will be removed soon." << G4endl
1225             << " Please try to use the new structured commands!" << G4endl;
1226    }
1227  else if(command == histpointCmd)
1228    {
1229      if(histtype == "biasx")
1230        fParticleGun->GetBiasRndm()->SetXBias(histpointCmd->GetNew3VectorValue(newValues));
1231      if(histtype == "biasy")
1232        fParticleGun->GetBiasRndm()->SetYBias(histpointCmd->GetNew3VectorValue(newValues));
1233      if(histtype == "biasz")
1234        fParticleGun->GetBiasRndm()->SetZBias(histpointCmd->GetNew3VectorValue(newValues));
1235      if(histtype == "biast")
1236        fParticleGun->GetBiasRndm()->SetThetaBias(histpointCmd->GetNew3VectorValue(newValues));
1237      if(histtype == "biasp")
1238        fParticleGun->GetBiasRndm()->SetPhiBias(histpointCmd->GetNew3VectorValue(newValues));
1239      if(histtype == "biase")
1240        fParticleGun->GetBiasRndm()->SetEnergyBias(histpointCmd->GetNew3VectorValue(newValues));
1241      if(histtype == "theta")
1242        fParticleGun->GetAngDist()->UserDefAngTheta(histpointCmd->GetNew3VectorValue(newValues));
1243      if(histtype == "phi")
1244        fParticleGun->GetAngDist()->UserDefAngPhi(histpointCmd->GetNew3VectorValue(newValues));
1245      if(histtype == "energy")
1246        fParticleGun->GetEneDist()->UserEnergyHisto(histpointCmd->GetNew3VectorValue(newValues));
1247      if(histtype == "arb")
1248        fParticleGun->GetEneDist()->ArbEnergyHisto(histpointCmd->GetNew3VectorValue(newValues));
1249      if(histtype == "epn")
1250        fParticleGun->GetEneDist()->EpnEnergyHisto(histpointCmd->GetNew3VectorValue(newValues));
1251      G4cout << " G4GeneralParticleSourceMessenger - Warning: The command is obsolete and will be removed soon. Please try to use the new structured commands!" << G4endl;
1252    }
1253  else if(command == resethistCmd)
1254    {
1255      if(newValues == "theta" || newValues == "phi") {
1256        fParticleGun->GetAngDist()->ReSetHist(newValues);
1257      } else if (newValues == "energy" || newValues == "arb" || newValues == "epn") {
1258        fParticleGun->GetEneDist()->ReSetHist(newValues);
1259      } else {
1260        fParticleGun->GetBiasRndm()->ReSetHist(newValues);
1261      }
1262      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1263             << " The command is obsolete and will be removed soon." << G4endl
1264             << " Please try to use the new structured commands!" << G4endl;
1265    }
1266  else if(command == arbintCmd)
1267    {
1268      fParticleGun->GetEneDist()->ArbInterpolate(newValues);
1269      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1270             << " The command is obsolete and will be removed soon." << G4endl
1271             << " Please try to use the new structured commands!" << G4endl;
1272    }
1273  else if( command==directionCmd )
1274    { 
1275      fParticleGun->GetAngDist()->SetAngDistType("planar");
1276      fParticleGun->GetAngDist()->SetParticleMomentumDirection(directionCmd->GetNew3VectorValue(newValues));
1277    }
1278  else if( command==energyCmd )
1279    {   
1280      fParticleGun->GetEneDist()->SetEnergyDisType("Mono");
1281      fParticleGun->GetEneDist()->SetMonoEnergy(energyCmd->GetNewDoubleValue(newValues));
1282    }
1283  else if( command==positionCmd )
1284    { 
1285      fParticleGun->GetPosDist()->SetPosDisType("Point");   
1286      fParticleGun->GetPosDist()->SetCentreCoords(positionCmd->GetNew3VectorValue(newValues));
1287    }
1288  else if(command == verbosityCmd)
1289    {
1290      fParticleGun->SetVerbosity(verbosityCmd->GetNewIntValue(newValues));
1291    }
1292  else if( command==particleCmd )
1293    {
1294      if (newValues =="ion") {
1295        fShootIon = true;
1296      } else {
1297        fShootIon = false;
1298        G4ParticleDefinition* pd = particleTable->FindParticle(newValues);
1299        if(pd != NULL)
1300          { fParticleGun->SetParticleDefinition( pd ); }
1301      }
1302    }
1303  else if( command==timeCmd )
1304    { fParticleGun->SetParticleTime(timeCmd->GetNewDoubleValue(newValues)); }
1305  else if( command==polCmd )
1306    { fParticleGun->SetParticlePolarization(polCmd->GetNew3VectorValue(newValues)); }
1307  else if( command==numberCmd )
1308    { fParticleGun->SetNumberOfParticles(numberCmd->GetNewIntValue(newValues)); }
1309  else if( command==ionCmd )
1310    { IonCommand(newValues); }
1311  else if( command==listCmd ){ 
1312    particleTable->DumpTable(); 
1313  }
1314  else if( command==addsourceCmd )
1315    { 
1316      fGPS->AddaSource(addsourceCmd->GetNewDoubleValue(newValues));
1317    }
1318  else if( command==listsourceCmd )
1319    { 
1320      fGPS->ListSource();
1321    }
1322  else if( command==clearsourceCmd )
1323    { 
1324      fGPS->ClearAll();
1325    }
1326  else if( command==getsourceCmd )
1327    { 
1328      G4cout << " Current source index:" << fGPS->GetCurrentSourceIndex() 
1329             << " ; Intensity:" << fGPS->GetCurrentSourceIntensity() << G4endl;
1330    }
1331  else if( command==setsourceCmd )
1332    { 
1333      fGPS->SetCurrentSourceto(setsourceCmd->GetNewIntValue(newValues));
1334    }
1335  else if( command==setintensityCmd )
1336    { 
1337      fGPS->SetCurrentSourceIntensity(setintensityCmd->GetNewDoubleValue(newValues));
1338    }
1339  else if( command==deletesourceCmd )
1340    { 
1341      fGPS->DeleteaSource(deletesourceCmd->GetNewIntValue(newValues));
1342    }
1343  else if(command == multiplevertexCmd)
1344    {
1345      fGPS->SetMultipleVertex(multiplevertexCmd->GetNewBoolValue(newValues));
1346    }
1347  else if(command == flatsamplingCmd)
1348    {
1349      fGPS->SetFlatSampling(flatsamplingCmd->GetNewBoolValue(newValues));
1350    }
1351  //
1352  // new implementations
1353  //
1354  //
1355  else if(command == typeCmd1)
1356    {
1357      fParticleGun->GetPosDist()->SetPosDisType(newValues);
1358    }
1359  else if(command == shapeCmd1)
1360    {
1361      fParticleGun->GetPosDist()->SetPosDisShape(newValues);
1362    }
1363  else if(command == centreCmd1)
1364    {
1365      fParticleGun->GetPosDist()->SetCentreCoords(centreCmd1->GetNew3VectorValue(newValues));
1366    }
1367  else if(command == posrot1Cmd1)
1368    {
1369      fParticleGun->GetPosDist()->SetPosRot1(posrot1Cmd1->GetNew3VectorValue(newValues));
1370    }
1371  else if(command == posrot2Cmd1)
1372    {
1373      fParticleGun->GetPosDist()->SetPosRot2(posrot2Cmd1->GetNew3VectorValue(newValues));
1374    }
1375  else if(command == halfxCmd1)
1376    {
1377      fParticleGun->GetPosDist()->SetHalfX(halfxCmd1->GetNewDoubleValue(newValues));
1378    }
1379  else if(command == halfyCmd1)
1380    {
1381      fParticleGun->GetPosDist()->SetHalfY(halfyCmd1->GetNewDoubleValue(newValues));
1382    }
1383  else if(command == halfzCmd1)
1384    {
1385      fParticleGun->GetPosDist()->SetHalfZ(halfzCmd1->GetNewDoubleValue(newValues));
1386    }
1387  else if(command == radiusCmd1)
1388    {
1389      fParticleGun->GetPosDist()->SetRadius(radiusCmd1->GetNewDoubleValue(newValues));
1390    }
1391  else if(command == radius0Cmd1)
1392    {
1393      fParticleGun->GetPosDist()->SetRadius0(radius0Cmd1->GetNewDoubleValue(newValues));
1394    }
1395  else if(command == possigmarCmd1)
1396    {
1397      fParticleGun->GetPosDist()->SetBeamSigmaInR(possigmarCmd1->GetNewDoubleValue(newValues));
1398    }
1399  else if(command == possigmaxCmd1)
1400    {
1401      fParticleGun->GetPosDist()->SetBeamSigmaInX(possigmaxCmd1->GetNewDoubleValue(newValues));
1402    }
1403  else if(command == possigmayCmd1)
1404    {
1405      fParticleGun->GetPosDist()->SetBeamSigmaInY(possigmayCmd1->GetNewDoubleValue(newValues));
1406    }
1407  else if(command == paralpCmd1)
1408    {
1409      fParticleGun->GetPosDist()->SetParAlpha(paralpCmd1->GetNewDoubleValue(newValues));
1410    }
1411  else if(command == partheCmd1)
1412    {
1413      fParticleGun->GetPosDist()->SetParTheta(partheCmd1->GetNewDoubleValue(newValues));
1414    }
1415  else if(command == parphiCmd1)
1416    {
1417      fParticleGun->GetPosDist()->SetParPhi(parphiCmd1->GetNewDoubleValue(newValues));
1418    }
1419  else if(command == confineCmd1)
1420    {
1421      fParticleGun->GetPosDist()->ConfineSourceToVolume(newValues);
1422    }
1423  else if(command == angtypeCmd1)
1424    {
1425      fParticleGun->GetAngDist()->SetAngDistType(newValues);
1426    }
1427  else if(command == angrot1Cmd1)
1428    {
1429      G4String a = "angref1";
1430      fParticleGun->GetAngDist()->DefineAngRefAxes(a,angrot1Cmd1->GetNew3VectorValue(newValues));
1431    }
1432  else if(command == angrot2Cmd1)
1433    {
1434      G4String a = "angref2";
1435      fParticleGun->GetAngDist()->DefineAngRefAxes(a,angrot2Cmd1->GetNew3VectorValue(newValues));
1436    }
1437  else if(command == minthetaCmd1)
1438    {
1439      fParticleGun->GetAngDist()->SetMinTheta(minthetaCmd1->GetNewDoubleValue(newValues));
1440    }
1441  else if(command == minphiCmd1)
1442    {
1443      fParticleGun->GetAngDist()->SetMinPhi(minphiCmd1->GetNewDoubleValue(newValues));
1444    }
1445  else if(command == maxthetaCmd1)
1446    {
1447      fParticleGun->GetAngDist()->SetMaxTheta(maxthetaCmd1->GetNewDoubleValue(newValues));
1448    }
1449  else if(command == maxphiCmd1)
1450    {
1451      fParticleGun->GetAngDist()->SetMaxPhi(maxphiCmd1->GetNewDoubleValue(newValues));
1452    }
1453  else if(command == angsigmarCmd1)
1454    {
1455      fParticleGun->GetAngDist()->SetBeamSigmaInAngR(angsigmarCmd1->GetNewDoubleValue(newValues));
1456    }
1457  else if(command == angsigmaxCmd1)
1458    {
1459      fParticleGun->GetAngDist()->SetBeamSigmaInAngX(angsigmaxCmd1->GetNewDoubleValue(newValues));
1460    }
1461  else if(command == angsigmayCmd1)
1462    {
1463      fParticleGun->GetAngDist()->SetBeamSigmaInAngY(angsigmayCmd1->GetNewDoubleValue(newValues));
1464    }
1465  else if(command == angfocusCmd)
1466    {
1467      fParticleGun->GetAngDist()->SetFocusPoint(angfocusCmd->GetNew3VectorValue(newValues));
1468    }
1469  else if(command == useuserangaxisCmd1)
1470    {
1471      fParticleGun->GetAngDist()->SetUseUserAngAxis(useuserangaxisCmd1->GetNewBoolValue(newValues));
1472    }
1473  else if(command == surfnormCmd1)
1474    {
1475      fParticleGun->GetAngDist()->SetUserWRTSurface(surfnormCmd1->GetNewBoolValue(newValues));
1476    }
1477  else if(command == energytypeCmd1)
1478    {
1479      fParticleGun->GetEneDist()->SetEnergyDisType(newValues);
1480    }
1481  else if(command == eminCmd1)
1482    {
1483      fParticleGun->GetEneDist()->SetEmin(eminCmd1->GetNewDoubleValue(newValues));
1484    }
1485  else if(command == emaxCmd1)
1486    {
1487      fParticleGun->GetEneDist()->SetEmax(emaxCmd1->GetNewDoubleValue(newValues));
1488    }
1489  else if(command == monoenergyCmd1)
1490    {
1491      fParticleGun->GetEneDist()->SetMonoEnergy(monoenergyCmd1->GetNewDoubleValue(newValues));
1492    }
1493  else if(command == engsigmaCmd1)
1494    {
1495      fParticleGun->GetEneDist()->SetBeamSigmaInE(engsigmaCmd1->GetNewDoubleValue(newValues));
1496    }
1497  else if(command == alphaCmd1)
1498    {
1499      fParticleGun->GetEneDist()->SetAlpha(alphaCmd1->GetNewDoubleValue(newValues));
1500    }
1501  else if(command == tempCmd1)
1502    {
1503      fParticleGun->GetEneDist()->SetTemp(tempCmd1->GetNewDoubleValue(newValues));
1504    }
1505  else if(command == ezeroCmd1)
1506    {
1507      fParticleGun->GetEneDist()->SetEzero(ezeroCmd1->GetNewDoubleValue(newValues));
1508    }
1509  else if(command == gradientCmd1)
1510    {
1511      fParticleGun->GetEneDist()->SetGradient(gradientCmd1->GetNewDoubleValue(newValues));
1512    }
1513  else if(command == interceptCmd1)
1514    {
1515      fParticleGun->GetEneDist()->SetInterCept(interceptCmd1->GetNewDoubleValue(newValues));
1516    }
1517  else if(command == arbeintCmd1)
1518    {
1519      fParticleGun->GetEneDist()->SetBiasAlpha(arbeintCmd1->GetNewDoubleValue(newValues));
1520    }
1521  else if(command == calculateCmd1)
1522    {
1523      fParticleGun->GetEneDist()->Calculate();
1524    }
1525  else if(command == energyspecCmd1)
1526    {
1527      fParticleGun->GetEneDist()->InputEnergySpectra(energyspecCmd1->GetNewBoolValue(newValues));
1528    }
1529  else if(command == diffspecCmd1)
1530    {
1531      fParticleGun->GetEneDist()->InputDifferentialSpectra(diffspecCmd1->GetNewBoolValue(newValues));
1532    }
1533  else if(command == histnameCmd1)
1534    {
1535      histtype = newValues;
1536    }
1537  else if(command == histfileCmd1)
1538    {
1539      histtype = "arb";
1540      fParticleGun->GetEneDist()->ArbEnergyHistoFile(newValues);
1541    }
1542  else if(command == histpointCmd1)
1543    {
1544      if(histtype == "biasx")
1545        fParticleGun->GetBiasRndm()->SetXBias(histpointCmd1->GetNew3VectorValue(newValues));
1546      if(histtype == "biasy")
1547        fParticleGun->GetBiasRndm()->SetYBias(histpointCmd1->GetNew3VectorValue(newValues));
1548      if(histtype == "biasz")
1549        fParticleGun->GetBiasRndm()->SetZBias(histpointCmd1->GetNew3VectorValue(newValues));
1550      if(histtype == "biast")
1551        fParticleGun->GetBiasRndm()->SetThetaBias(histpointCmd1->GetNew3VectorValue(newValues));
1552      if(histtype == "biasp")
1553        fParticleGun->GetBiasRndm()->SetPhiBias(histpointCmd1->GetNew3VectorValue(newValues));
1554      if(histtype == "biaspt")
1555        fParticleGun->GetBiasRndm()->SetPosThetaBias(histpointCmd1->GetNew3VectorValue(newValues));
1556      if(histtype == "biaspp")
1557        fParticleGun->GetBiasRndm()->SetPosPhiBias(histpointCmd1->GetNew3VectorValue(newValues));
1558      if(histtype == "biase")
1559        fParticleGun->GetBiasRndm()->SetEnergyBias(histpointCmd1->GetNew3VectorValue(newValues));
1560      if(histtype == "theta")
1561        fParticleGun->GetAngDist()->UserDefAngTheta(histpointCmd1->GetNew3VectorValue(newValues));
1562      if(histtype == "phi")
1563        fParticleGun->GetAngDist()->UserDefAngPhi(histpointCmd1->GetNew3VectorValue(newValues));
1564      if(histtype == "energy")
1565        fParticleGun->GetEneDist()->UserEnergyHisto(histpointCmd1->GetNew3VectorValue(newValues));
1566      if(histtype == "arb")
1567        fParticleGun->GetEneDist()->ArbEnergyHisto(histpointCmd1->GetNew3VectorValue(newValues));
1568      if(histtype == "epn")
1569        fParticleGun->GetEneDist()->EpnEnergyHisto(histpointCmd1->GetNew3VectorValue(newValues));
1570    }
1571  else if(command == resethistCmd1)
1572    {
1573      if(newValues == "theta" || newValues == "phi") {
1574        fParticleGun->GetAngDist()->ReSetHist(newValues);
1575      } else if (newValues == "energy" || newValues == "arb" || newValues == "epn") {
1576        fParticleGun->GetEneDist()->ReSetHist(newValues);
1577      } else {
1578        fParticleGun->GetBiasRndm()->ReSetHist(newValues);
1579      }
1580    }
1581  else if(command == arbintCmd1)
1582    {
1583      fParticleGun->GetEneDist()->ArbInterpolate(newValues);
1584    }
1585  else 
1586    {
1587      G4cout << "Error entering command" << G4endl;
1588    }
1589}
1590
1591G4String G4GeneralParticleSourceMessenger::GetCurrentValue(G4UIcommand *)
1592{
1593  G4String cv;
1594 
1595  //  if( command==directionCmd )
1596  //  { cv = directionCmd->ConvertToString(fParticleGun->GetParticleMomentumDirection()); }
1597  //  else if( command==energyCmd )
1598  //  { cv = energyCmd->ConvertToString(fParticleGun->GetParticleEnergy(),"GeV"); }
1599  //  else if( command==positionCmd )
1600  //  { cv = positionCmd->ConvertToString(fParticleGun->GetParticlePosition(),"cm"); }
1601  //  else if( command==timeCmd )
1602  //  { cv = timeCmd->ConvertToString(fParticleGun->GetParticleTime(),"ns"); }
1603  //  else if( command==polCmd )
1604  //  { cv = polCmd->ConvertToString(fParticleGun->GetParticlePolarization()); }
1605  //  else if( command==numberCmd )
1606  //  { cv = numberCmd->ConvertToString(fParticleGun->GetNumberOfParticles()); }
1607 
1608  cv = "Not implemented yet";
1609
1610  return cv;
1611}
1612
1613void G4GeneralParticleSourceMessenger::IonCommand(G4String newValues)
1614{
1615  fShootIon = true;
1616
1617  if (fShootIon)
1618  {
1619    G4Tokenizer next( newValues );
1620    // check argument
1621    fAtomicNumber = StoI(next());
1622    fAtomicMass = StoI(next());
1623    G4String sQ = next();
1624    if (sQ.isNull())
1625    {
1626        fIonCharge = fAtomicNumber;
1627    }
1628    else
1629    {
1630        fIonCharge = StoI(sQ);
1631        sQ = next();
1632        if (sQ.isNull())
1633      {
1634          fIonExciteEnergy = 0.0;
1635      }
1636      else
1637      {
1638          fIonExciteEnergy = StoD(sQ) * keV;
1639      }
1640    }
1641    G4ParticleDefinition* ion;
1642    ion =  particleTable->GetIon( fAtomicNumber, fAtomicMass, fIonExciteEnergy);
1643    if (ion==0)
1644    {
1645      G4cout << "Ion with Z=" << fAtomicNumber;
1646      G4cout << " A=" << fAtomicMass << "is not be defined" << G4endl;   
1647    }
1648    else
1649    {
1650      fParticleGun->SetParticleDefinition(ion);
1651      fParticleGun->SetParticleCharge(fIonCharge*eplus);
1652    }
1653  }
1654  else
1655  {
1656    G4cout << "Set /gps/particle to ion before using /gps/ion command";
1657    G4cout << G4endl; 
1658  }
1659}
Note: See TracBrowser for help on using the repository browser.