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