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

Last change on this file since 1245 was 1196, checked in by garnier, 15 years ago

update CVS release candidate geant4.9.3.01

File size: 69.0 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",true,true); 
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",true,true);
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",true,true);
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",true,true);
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",true,true); 
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true); 
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",true,true); 
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true); 
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",true,true); 
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
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",true,true);
638
639  ezeroCmd1 = new G4UIcmdWithADouble("/gps/ene/ezero",this);
640  ezeroCmd1->SetGuidance("Sets E_0 for exponential distribution (in MeV)");
641  ezeroCmd1->SetParameterName("ezero",true,true);
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",true,true);
646
647  interceptCmd1 = new G4UIcmdWithADouble("/gps/ene/intercept",this);
648  interceptCmd1->SetGuidance("Sets the intercept for Lin distributions (in MeV)");
649  interceptCmd1->SetParameterName("intercept",true,true);
650
651  calculateCmd1 = new G4UIcmdWithoutParameter("/gps/ene/calculate",this);
652  calculateCmd1->SetGuidance("Calculates the distributions for Cdg and BBody");
653
654  energyspecCmd1 = new G4UIcmdWithABool("/gps/ene/emspec",this);
655  energyspecCmd1->SetGuidance("True for energy and false for momentum spectra");
656  energyspecCmd1->SetParameterName("energyspec",true);
657  energyspecCmd1->SetDefaultValue(true);
658
659  diffspecCmd1 = new G4UIcmdWithABool("/gps/ene/diffspec",this);
660  diffspecCmd1->SetGuidance("True for differential and flase for integral spectra");
661  diffspecCmd1->SetParameterName("diffspec",true);
662  diffspecCmd1->SetDefaultValue(true);
663
664  //old ones
665  energytypeCmd = new G4UIcmdWithAString("/gps/energytype",this);
666  energytypeCmd->SetGuidance("Sets energy distribution type (obsolete!)");
667  energytypeCmd->SetParameterName("EnergyDis",true,true);
668  energytypeCmd->SetDefaultValue("Mono");
669  energytypeCmd->SetCandidates("Mono Lin Pow Exp Gauss Brem Bbody Cdg User Arb Epn");
670
671  eminCmd = new G4UIcmdWithADoubleAndUnit("/gps/emin",this);
672  eminCmd->SetGuidance("Sets Emin (obsolete!)");
673  eminCmd->SetParameterName("emin",true,true);
674  eminCmd->SetDefaultUnit("keV");
675  //  eminCmd->SetUnitCandidates("eV keV MeV GeV TeV PeV");
676
677  emaxCmd = new G4UIcmdWithADoubleAndUnit("/gps/emax",this);
678  emaxCmd->SetGuidance("Sets Emax (obsolete!)");
679  emaxCmd->SetParameterName("emax",true,true);
680  emaxCmd->SetDefaultUnit("keV");
681  //  emaxCmd->SetUnitCandidates("eV keV MeV GeV TeV PeV");
682
683  monoenergyCmd = new G4UIcmdWithADoubleAndUnit("/gps/monoenergy",this);
684  monoenergyCmd->SetGuidance("Sets Monoenergy (obsolete, use gps/energy instead!)");
685  monoenergyCmd->SetParameterName("monoenergy",true,true);
686  monoenergyCmd->SetDefaultUnit("keV");
687  //  monoenergyCmd->SetUnitCandidates("eV keV MeV GeV TeV PeV");
688
689  engsigmaCmd = new G4UIcmdWithADoubleAndUnit("/gps/sigmae",this);
690  engsigmaCmd->SetGuidance("Sets the standard deviation for Gaussian energy dist.(obsolete!)");
691  engsigmaCmd->SetParameterName("Sigmae",true,true);
692  engsigmaCmd->SetDefaultUnit("keV");
693  //  engsigmaCmd->SetUnitCandidates("eV keV MeV GeV TeV PeV");
694
695  alphaCmd = new G4UIcmdWithADouble("/gps/alpha",this);
696  alphaCmd->SetGuidance("Sets Alpha (index) for power-law energy dist(obsolete!).");
697  alphaCmd->SetParameterName("alpha",true,true);
698 
699  tempCmd = new G4UIcmdWithADouble("/gps/temp",this);
700  tempCmd->SetGuidance("Sets the temperature for Brem and BBody (in Kelvin)(obsolete!)");
701  tempCmd->SetParameterName("temp",true,true);
702
703  ezeroCmd = new G4UIcmdWithADouble("/gps/ezero",this);
704  ezeroCmd->SetGuidance("Sets ezero exponential distributions (in MeV)(obsolete!)");
705  ezeroCmd->SetParameterName("ezero",true,true);
706
707  gradientCmd = new G4UIcmdWithADouble("/gps/gradient",this);
708  gradientCmd->SetGuidance("Sets the gradient for Lin distributions (in 1/MeV)(obsolete!)");
709  gradientCmd->SetParameterName("gradient",true,true);
710
711  interceptCmd = new G4UIcmdWithADouble("/gps/intercept",this);
712  interceptCmd->SetGuidance("Sets the intercept for Lin distributions (in MeV)(obsolete!)");
713  interceptCmd->SetParameterName("intercept",true,true);
714
715  calculateCmd = new G4UIcmdWithoutParameter("/gps/calculate",this);
716  calculateCmd->SetGuidance("Calculates distributions for Cdg and BBody(obsolete!)");
717
718  energyspecCmd = new G4UIcmdWithABool("/gps/energyspec",this);
719  energyspecCmd->SetGuidance("True for energy and false for momentum spectra(obsolete!)");
720  energyspecCmd->SetParameterName("energyspec",true);
721  energyspecCmd->SetDefaultValue(true);
722
723  diffspecCmd = new G4UIcmdWithABool("/gps/diffspec",this);
724  diffspecCmd->SetGuidance("True for differential and flase for integral spectra(obsolete!)");
725  diffspecCmd->SetParameterName("diffspec",true);
726  diffspecCmd->SetDefaultValue(true);
727
728  // Biasing + histograms in general
729  histDirectory = new G4UIdirectory("/gps/hist/");
730  histDirectory->SetGuidance("Histogram, biasing commands sub-directory");
731
732  histnameCmd1 = new G4UIcmdWithAString("/gps/hist/type",this);
733  histnameCmd1->SetGuidance("Sets histogram type");
734  histnameCmd1->SetParameterName("HistType",true,true);
735  histnameCmd1->SetDefaultValue("biasx");
736  histnameCmd1->SetCandidates("biasx biasy biasz biast biasp biase biaspt biaspp theta phi energy arb epn");
737
738  resethistCmd1 = new G4UIcmdWithAString("/gps/hist/reset",this);
739  resethistCmd1->SetGuidance("Reset (clean) the histogram ");
740  resethistCmd1->SetParameterName("HistType",true,true);
741  resethistCmd1->SetDefaultValue("energy");
742  resethistCmd1->SetCandidates("biasx biasy biasz biast biasp biase biaspt biaspp theta phi energy arb epn");
743
744  histpointCmd1 = new G4UIcmdWith3Vector("/gps/hist/point",this);
745  histpointCmd1->SetGuidance("Allows user to define a histogram");
746  histpointCmd1->SetGuidance("Enter: Ehi Weight");
747  histpointCmd1->SetParameterName("Ehi","Weight","Junk",true,true);
748  histpointCmd1->SetRange("Ehi >= 0. && Weight >= 0.");
749
750  arbintCmd1 = new G4UIcmdWithAString("/gps/hist/inter",this);
751  arbintCmd1->SetGuidance("Sets the interpolation method for arbitrary distribution.");
752  arbintCmd1->SetParameterName("int",true,true);
753  arbintCmd1->SetDefaultValue("Lin");
754  arbintCmd1->SetCandidates("Lin Log Exp Spline");
755
756  // old ones
757  histnameCmd = new G4UIcmdWithAString("/gps/histname",this);
758  histnameCmd->SetGuidance("Sets histogram type (obsolete!)");
759  histnameCmd->SetParameterName("HistType",true,true);
760  histnameCmd->SetDefaultValue("biasx");
761  histnameCmd->SetCandidates("biasx biasy biasz biast biasp biase biaspt biaspp theta phi energy arb epn");
762
763  // re-set the histograms
764  resethistCmd = new G4UIcmdWithAString("/gps/resethist",this);
765  resethistCmd->SetGuidance("Re-Set the histogram (obsolete!)");
766  resethistCmd->SetParameterName("HistType",true,true);
767  resethistCmd->SetDefaultValue("energy");
768  resethistCmd->SetCandidates("biasx biasy biasz biast biasp biase biaspt biaspp theta phi energy arb epn");
769
770  histpointCmd = new G4UIcmdWith3Vector("/gps/histpoint",this);
771  histpointCmd->SetGuidance("Allows user to define a histogram (obsolete!)");
772  histpointCmd->SetGuidance("Enter: Ehi Weight");
773  histpointCmd->SetParameterName("Ehi","Weight","Junk",true,true);
774  histpointCmd->SetRange("Ehi >= 0. && Weight >= 0.");
775
776  arbintCmd = new G4UIcmdWithAString("/gps/arbint",this);
777  arbintCmd->SetGuidance("Sets Arbitrary Interpolation type.(obsolete!) ");
778  arbintCmd->SetParameterName("int",true,true);
779  arbintCmd->SetDefaultValue("NULL");
780  arbintCmd->SetCandidates("Lin Log Exp Spline");
781
782}
783
784G4GeneralParticleSourceMessenger::~G4GeneralParticleSourceMessenger()
785{
786  delete positionDirectory;
787  delete typeCmd;
788  delete shapeCmd;
789  delete centreCmd;
790  delete posrot1Cmd;
791  delete posrot2Cmd;
792  delete halfxCmd;
793  delete halfyCmd;
794  delete halfzCmd;
795  delete radiusCmd;
796  delete radius0Cmd;
797  delete possigmarCmd;
798  delete possigmaxCmd;
799  delete possigmayCmd;
800  delete paralpCmd;
801  delete partheCmd;
802  delete parphiCmd;
803  delete confineCmd;
804  delete typeCmd1;
805  delete shapeCmd1;
806  delete centreCmd1;
807  delete posrot1Cmd1;
808  delete posrot2Cmd1;
809  delete halfxCmd1;
810  delete halfyCmd1;
811  delete halfzCmd1;
812  delete radiusCmd1;
813  delete radius0Cmd1;
814  delete possigmarCmd1;
815  delete possigmaxCmd1;
816  delete possigmayCmd1;
817  delete paralpCmd1;
818  delete partheCmd1;
819  delete parphiCmd1;
820  delete confineCmd1;
821
822  delete angularDirectory;
823  delete angtypeCmd;
824  delete angrot1Cmd;
825  delete angrot2Cmd;
826  delete minthetaCmd;
827  delete maxthetaCmd;
828  delete minphiCmd;
829  delete maxphiCmd;
830  delete angsigmarCmd;
831  delete angsigmaxCmd;
832  delete angsigmayCmd;
833  delete useuserangaxisCmd;
834  delete surfnormCmd;
835  delete angtypeCmd1;
836  delete angrot1Cmd1;
837  delete angrot2Cmd1;
838  delete minthetaCmd1;
839  delete maxthetaCmd1;
840  delete minphiCmd1;
841  delete maxphiCmd1;
842  delete angsigmarCmd1;
843  delete angsigmaxCmd1;
844  delete angfocusCmd;
845  delete useuserangaxisCmd1;
846  delete surfnormCmd1;
847
848  delete energyDirectory;
849  delete energytypeCmd;
850  delete eminCmd;
851  delete emaxCmd;
852  delete monoenergyCmd;
853  delete engsigmaCmd;
854  delete alphaCmd;
855  delete tempCmd;
856  delete ezeroCmd;
857  delete gradientCmd;
858  delete interceptCmd;
859  delete calculateCmd;
860  delete energyspecCmd;
861  delete diffspecCmd;
862  delete energytypeCmd1;
863  delete eminCmd1;
864  delete emaxCmd1;
865  delete monoenergyCmd1;
866  delete engsigmaCmd1;
867  delete alphaCmd1;
868  delete tempCmd1;
869  delete ezeroCmd1;
870  delete gradientCmd1;
871  delete interceptCmd1;
872  delete calculateCmd1;
873  delete energyspecCmd1;
874  delete diffspecCmd1;
875
876  delete histDirectory;
877  delete histnameCmd;
878  delete resethistCmd;
879  delete histpointCmd;
880  delete arbintCmd;
881  delete histnameCmd1;
882  delete resethistCmd1;
883  delete histpointCmd1;
884  delete arbintCmd1;
885
886  delete verbosityCmd;
887  delete ionCmd;
888  delete particleCmd;
889  delete timeCmd;
890  delete polCmd;
891  delete numberCmd;
892  delete positionCmd;
893  delete directionCmd;
894  delete energyCmd;
895  delete listCmd;
896
897  delete sourceDirectory;
898  delete addsourceCmd;
899  delete listsourceCmd;
900  delete clearsourceCmd;
901  delete getsourceCmd;
902  delete setsourceCmd;
903  delete setintensityCmd;
904  delete deletesourceCmd;
905  delete multiplevertexCmd;
906  delete flatsamplingCmd;
907
908  delete gpsDirectory;
909 
910}
911
912void G4GeneralParticleSourceMessenger::SetNewValue(G4UIcommand *command, G4String newValues)
913{
914  if(command == typeCmd)
915    {
916      fParticleGun->GetPosDist()->SetPosDisType(newValues);
917      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
918             << " The command is obsolete and will be removed soon." << G4endl
919             << " Please try to use the new structured commands!" << G4endl;
920    }
921  else if(command == shapeCmd)
922    {
923      fParticleGun->GetPosDist()->SetPosDisShape(newValues);
924      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
925             << " The command is obsolete and will be removed soon." << G4endl
926             << " Please try to use the new structured commands!" << G4endl;
927    }
928  else if(command == centreCmd)
929    {
930      fParticleGun->GetPosDist()->SetCentreCoords(centreCmd->GetNew3VectorValue(newValues));
931      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
932             << " The command is obsolete and will be removed soon." << G4endl
933             << " Please try to use the new structured commands!" << G4endl;
934    }
935  else if(command == posrot1Cmd)
936    {
937      fParticleGun->GetPosDist()->SetPosRot1(posrot1Cmd->GetNew3VectorValue(newValues));
938      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
939             << " The command is obsolete and will be removed soon." << G4endl
940             << " Please try to use the new structured commands!" << G4endl;
941    }
942  else if(command == posrot2Cmd)
943    {
944      fParticleGun->GetPosDist()->SetPosRot2(posrot2Cmd->GetNew3VectorValue(newValues));
945      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
946             << " The command is obsolete and will be removed soon." << G4endl
947             << " Please try to use the new structured commands!" << G4endl;
948    }
949  else if(command == halfxCmd)
950    {
951      fParticleGun->GetPosDist()->SetHalfX(halfxCmd->GetNewDoubleValue(newValues));
952      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
953             << " The command is obsolete and will be removed soon." << G4endl
954             << " Please try to use the new structured commands!" << G4endl;
955    }
956  else if(command == halfyCmd)
957    {
958      fParticleGun->GetPosDist()->SetHalfY(halfyCmd->GetNewDoubleValue(newValues));
959      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
960             << " The command is obsolete and will be removed soon." << G4endl
961             << " Please try to use the new structured commands!" << G4endl;
962    }
963  else if(command == halfzCmd)
964    {
965      fParticleGun->GetPosDist()->SetHalfZ(halfzCmd->GetNewDoubleValue(newValues));
966      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
967             << " The command is obsolete and will be removed soon." << G4endl
968             << " Please try to use the new structured commands!" << G4endl;
969    }
970  else if(command == radiusCmd)
971    {
972      fParticleGun->GetPosDist()->SetRadius(radiusCmd->GetNewDoubleValue(newValues));
973      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
974             << " The command is obsolete and will be removed soon." << G4endl
975             << " Please try to use the new structured commands!" << G4endl;
976    }
977  else if(command == radius0Cmd)
978    {
979      fParticleGun->GetPosDist()->SetRadius0(radius0Cmd->GetNewDoubleValue(newValues));
980      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
981             << " The command is obsolete and will be removed soon." << G4endl
982             << " Please try to use the new structured commands!" << G4endl;
983    }
984  else if(command == possigmarCmd)
985    {
986      fParticleGun->GetPosDist()->SetBeamSigmaInR(possigmarCmd->GetNewDoubleValue(newValues));
987      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
988             << " The command is obsolete and will be removed soon." << G4endl
989             << " Please try to use the new structured commands!" << G4endl;
990    }
991  else if(command == possigmaxCmd)
992    {
993      fParticleGun->GetPosDist()->SetBeamSigmaInX(possigmaxCmd->GetNewDoubleValue(newValues));
994      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
995             << " The command is obsolete and will be removed soon." << G4endl
996             << " Please try to use the new structured commands!" << G4endl;
997    }
998  else if(command == possigmayCmd)
999    {
1000      fParticleGun->GetPosDist()->SetBeamSigmaInY(possigmayCmd->GetNewDoubleValue(newValues));
1001      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1002             << " The command is obsolete and will be removed soon." << G4endl
1003             << " Please try to use the new structured commands!" << G4endl;
1004    }
1005  else if(command == paralpCmd)
1006    {
1007      fParticleGun->GetPosDist()->SetParAlpha(paralpCmd->GetNewDoubleValue(newValues));
1008      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1009             << " The command is obsolete and will be removed soon." << G4endl
1010             << " Please try to use the new structured commands!" << G4endl;
1011    }
1012  else if(command == partheCmd)
1013    {
1014      fParticleGun->GetPosDist()->SetParTheta(partheCmd->GetNewDoubleValue(newValues));
1015      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1016             << " The command is obsolete and will be removed soon." << G4endl
1017             << " Please try to use the new structured commands!" << G4endl;
1018    }
1019  else if(command == parphiCmd)
1020    {
1021      fParticleGun->GetPosDist()->SetParPhi(parphiCmd->GetNewDoubleValue(newValues));
1022      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1023             << " The command is obsolete and will be removed soon." << G4endl
1024             << " Please try to use the new structured commands!" << G4endl;
1025    }
1026  else if(command == confineCmd)
1027    {
1028      fParticleGun->GetPosDist()->ConfineSourceToVolume(newValues);
1029      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1030             << " The command is obsolete and will be removed soon." << G4endl
1031             << " Please try to use the new structured commands!" << G4endl;
1032    }
1033  else if(command == angtypeCmd)
1034    {
1035      fParticleGun->GetAngDist()->SetAngDistType(newValues);
1036      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1037             << " The command is obsolete and will be removed soon." << G4endl
1038             << " Please try to use the new structured commands!" << G4endl;
1039    }
1040  else if(command == angrot1Cmd)
1041    {
1042      G4String a = "angref1";
1043      fParticleGun->GetAngDist()->DefineAngRefAxes(a,angrot1Cmd->GetNew3VectorValue(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 == angrot2Cmd)
1049    {
1050      G4String a = "angref2";
1051      fParticleGun->GetAngDist()->DefineAngRefAxes(a,angrot2Cmd->GetNew3VectorValue(newValues));
1052      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1053             << " The command is obsolete and will be removed soon." << G4endl
1054             << " Please try to use the new structured commands!" << G4endl;
1055    }
1056  else if(command == minthetaCmd)
1057    {
1058      fParticleGun->GetAngDist()->SetMinTheta(minthetaCmd->GetNewDoubleValue(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 == minphiCmd)
1064    {
1065      fParticleGun->GetAngDist()->SetMinPhi(minphiCmd->GetNewDoubleValue(newValues));
1066      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1067             << " The command is obsolete and will be removed soon." << G4endl
1068             << " Please try to use the new structured commands!" << G4endl;
1069    }
1070  else if(command == maxthetaCmd)
1071    {
1072      fParticleGun->GetAngDist()->SetMaxTheta(maxthetaCmd->GetNewDoubleValue(newValues));
1073      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1074             << " The command is obsolete and will be removed soon." << G4endl
1075             << " Please try to use the new structured commands!" << G4endl;
1076    }
1077  else if(command == maxphiCmd)
1078    {
1079      fParticleGun->GetAngDist()->SetMaxPhi(maxphiCmd->GetNewDoubleValue(newValues));
1080      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1081             << " The command is obsolete and will be removed soon." << G4endl
1082             << " Please try to use the new structured commands!" << G4endl;
1083    }
1084  else if(command == angsigmarCmd)
1085    {
1086      fParticleGun->GetAngDist()->SetBeamSigmaInAngR(angsigmarCmd->GetNewDoubleValue(newValues));
1087      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1088             << " The command is obsolete and will be removed soon." << G4endl
1089             << " Please try to use the new structured commands!" << G4endl;
1090    }
1091  else if(command == angsigmaxCmd)
1092    {
1093      fParticleGun->GetAngDist()->SetBeamSigmaInAngX(angsigmaxCmd->GetNewDoubleValue(newValues));
1094      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1095             << " The command is obsolete and will be removed soon." << G4endl
1096             << " Please try to use the new structured commands!" << G4endl;
1097    }
1098  else if(command == angsigmayCmd)
1099    {
1100      fParticleGun->GetAngDist()->SetBeamSigmaInAngY(angsigmayCmd->GetNewDoubleValue(newValues));
1101      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1102             << " The command is obsolete and will be removed soon." << G4endl
1103             << " Please try to use the new structured commands!" << G4endl;
1104    }
1105  else if(command == useuserangaxisCmd)
1106    {
1107      fParticleGun->GetAngDist()->SetUseUserAngAxis(useuserangaxisCmd->GetNewBoolValue(newValues));
1108      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1109             << " The command is obsolete and will be removed soon." << G4endl
1110             << " Please try to use the new structured commands!" << G4endl;
1111    }
1112  else if(command == surfnormCmd)
1113    {
1114      fParticleGun->GetAngDist()->SetUserWRTSurface(surfnormCmd->GetNewBoolValue(newValues));
1115      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1116             << " The command is obsolete and will be removed soon." << G4endl
1117             << " Please try to use the new structured commands!" << G4endl;
1118    }
1119  else if(command == energytypeCmd)
1120    {
1121      fParticleGun->GetEneDist()->SetEnergyDisType(newValues);
1122      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1123             << " The command is obsolete and will be removed soon." << G4endl
1124             << " Please try to use the new structured commands!" << G4endl;
1125    }
1126  else if(command == eminCmd)
1127    {
1128      fParticleGun->GetEneDist()->SetEmin(eminCmd->GetNewDoubleValue(newValues));
1129      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1130             << " The command is obsolete and will be removed soon." << G4endl
1131             << " Please try to use the new structured commands!" << G4endl;
1132    }
1133  else if(command == emaxCmd)
1134    {
1135      fParticleGun->GetEneDist()->SetEmax(emaxCmd->GetNewDoubleValue(newValues));
1136      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1137             << " The command is obsolete and will be removed soon." << G4endl
1138             << " Please try to use the new structured commands!" << G4endl;
1139    }
1140  else if(command == monoenergyCmd)
1141    {
1142      fParticleGun->GetEneDist()->SetMonoEnergy(monoenergyCmd->GetNewDoubleValue(newValues));
1143      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1144             << " The command is obsolete and will be removed soon." << G4endl
1145             << " Please try to use the new structured commands!" << G4endl;
1146    }
1147  else if(command == engsigmaCmd)
1148    {
1149      fParticleGun->GetEneDist()->SetBeamSigmaInE(engsigmaCmd->GetNewDoubleValue(newValues));
1150      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1151             << " The command is obsolete and will be removed soon." << G4endl
1152             << " Please try to use the new structured commands!" << G4endl;
1153    }
1154  else if(command == alphaCmd)
1155    {
1156      fParticleGun->GetEneDist()->SetAlpha(alphaCmd->GetNewDoubleValue(newValues));
1157      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1158             << " The command is obsolete and will be removed soon." << G4endl
1159             << " Please try to use the new structured commands!" << G4endl;
1160    }
1161  else if(command == tempCmd)
1162    {
1163      fParticleGun->GetEneDist()->SetTemp(tempCmd->GetNewDoubleValue(newValues));
1164      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1165             << " The command is obsolete and will be removed soon." << G4endl
1166             << " Please try to use the new structured commands!" << G4endl;
1167    }
1168  else if(command == ezeroCmd)
1169    {
1170      fParticleGun->GetEneDist()->SetEzero(ezeroCmd->GetNewDoubleValue(newValues));
1171      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1172             << " The command is obsolete and will be removed soon." << G4endl
1173             << " Please try to use the new structured commands!" << G4endl;
1174    }
1175  else if(command == gradientCmd)
1176    {
1177      fParticleGun->GetEneDist()->SetGradient(gradientCmd->GetNewDoubleValue(newValues));
1178      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1179             << " The command is obsolete and will be removed soon." << G4endl
1180             << " Please try to use the new structured commands!" << G4endl;
1181    }
1182  else if(command == interceptCmd)
1183    {
1184      fParticleGun->GetEneDist()->SetInterCept(interceptCmd->GetNewDoubleValue(newValues));
1185      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1186             << " The command is obsolete and will be removed soon." << G4endl
1187             << " Please try to use the new structured commands!" << G4endl;
1188    }
1189  else if(command == calculateCmd)
1190    {
1191      fParticleGun->GetEneDist()->Calculate();
1192      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1193             << " The command is obsolete and will be removed soon." << G4endl
1194             << " Please try to use the new structured commands!" << G4endl;
1195    }
1196  else if(command == energyspecCmd)
1197    {
1198      fParticleGun->GetEneDist()->InputEnergySpectra(energyspecCmd->GetNewBoolValue(newValues));
1199      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1200             << " The command is obsolete and will be removed soon." << G4endl
1201             << " Please try to use the new structured commands!" << G4endl;
1202    }
1203  else if(command == diffspecCmd)
1204    {
1205      fParticleGun->GetEneDist()->InputDifferentialSpectra(diffspecCmd->GetNewBoolValue(newValues));
1206      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1207             << " The command is obsolete and will be removed soon." << G4endl
1208             << " Please try to use the new structured commands!" << G4endl;
1209    }
1210  else if(command == histnameCmd)
1211    {
1212      histtype = newValues;
1213      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1214             << " The command is obsolete and will be removed soon." << G4endl
1215             << " Please try to use the new structured commands!" << G4endl;
1216    }
1217  else if(command == histpointCmd)
1218    {
1219      if(histtype == "biasx")
1220        fParticleGun->GetBiasRndm()->SetXBias(histpointCmd->GetNew3VectorValue(newValues));
1221      if(histtype == "biasy")
1222        fParticleGun->GetBiasRndm()->SetYBias(histpointCmd->GetNew3VectorValue(newValues));
1223      if(histtype == "biasz")
1224        fParticleGun->GetBiasRndm()->SetZBias(histpointCmd->GetNew3VectorValue(newValues));
1225      if(histtype == "biast")
1226        fParticleGun->GetBiasRndm()->SetThetaBias(histpointCmd->GetNew3VectorValue(newValues));
1227      if(histtype == "biasp")
1228        fParticleGun->GetBiasRndm()->SetPhiBias(histpointCmd->GetNew3VectorValue(newValues));
1229      if(histtype == "biase")
1230        fParticleGun->GetBiasRndm()->SetEnergyBias(histpointCmd->GetNew3VectorValue(newValues));
1231      if(histtype == "theta")
1232        fParticleGun->GetAngDist()->UserDefAngTheta(histpointCmd->GetNew3VectorValue(newValues));
1233      if(histtype == "phi")
1234        fParticleGun->GetAngDist()->UserDefAngPhi(histpointCmd->GetNew3VectorValue(newValues));
1235      if(histtype == "energy")
1236        fParticleGun->GetEneDist()->UserEnergyHisto(histpointCmd->GetNew3VectorValue(newValues));
1237      if(histtype == "arb")
1238        fParticleGun->GetEneDist()->ArbEnergyHisto(histpointCmd->GetNew3VectorValue(newValues));
1239      if(histtype == "epn")
1240        fParticleGun->GetEneDist()->EpnEnergyHisto(histpointCmd->GetNew3VectorValue(newValues));
1241      G4cout << " G4GeneralParticleSourceMessenger - Warning: The command is obsolete and will be removed soon. Please try to use the new structured commands!" << G4endl;
1242    }
1243  else if(command == resethistCmd)
1244    {
1245      if(newValues == "theta" || newValues == "phi") {
1246        fParticleGun->GetAngDist()->ReSetHist(newValues);
1247      } else if (newValues == "energy" || newValues == "arb" || newValues == "epn") {
1248        fParticleGun->GetEneDist()->ReSetHist(newValues);
1249      } else {
1250        fParticleGun->GetBiasRndm()->ReSetHist(newValues);
1251      }
1252      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1253             << " The command is obsolete and will be removed soon." << G4endl
1254             << " Please try to use the new structured commands!" << G4endl;
1255    }
1256  else if(command == arbintCmd)
1257    {
1258      fParticleGun->GetEneDist()->ArbInterpolate(newValues);
1259      G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1260             << " The command is obsolete and will be removed soon." << G4endl
1261             << " Please try to use the new structured commands!" << G4endl;
1262    }
1263  else if( command==directionCmd )
1264    { 
1265      fParticleGun->GetAngDist()->SetAngDistType("planar");
1266      fParticleGun->GetAngDist()->SetParticleMomentumDirection(directionCmd->GetNew3VectorValue(newValues));
1267    }
1268  else if( command==energyCmd )
1269    {   
1270      fParticleGun->GetEneDist()->SetEnergyDisType("Mono");
1271      fParticleGun->GetEneDist()->SetMonoEnergy(energyCmd->GetNewDoubleValue(newValues));
1272    }
1273  else if( command==positionCmd )
1274    { 
1275      fParticleGun->GetPosDist()->SetPosDisType("Point");   
1276      fParticleGun->GetPosDist()->SetCentreCoords(positionCmd->GetNew3VectorValue(newValues));
1277    }
1278  else if(command == verbosityCmd)
1279    {
1280      fParticleGun->SetVerbosity(verbosityCmd->GetNewIntValue(newValues));
1281    }
1282  else if( command==particleCmd )
1283    {
1284      if (newValues =="ion") {
1285        fShootIon = true;
1286      } else {
1287        fShootIon = false;
1288        G4ParticleDefinition* pd = particleTable->FindParticle(newValues);
1289        if(pd != NULL)
1290          { fParticleGun->SetParticleDefinition( pd ); }
1291      }
1292    }
1293  else if( command==timeCmd )
1294    { fParticleGun->SetParticleTime(timeCmd->GetNewDoubleValue(newValues)); }
1295  else if( command==polCmd )
1296    { fParticleGun->SetParticlePolarization(polCmd->GetNew3VectorValue(newValues)); }
1297  else if( command==numberCmd )
1298    { fParticleGun->SetNumberOfParticles(numberCmd->GetNewIntValue(newValues)); }
1299  else if( command==ionCmd )
1300    { IonCommand(newValues); }
1301  else if( command==listCmd ){ 
1302    particleTable->DumpTable(); 
1303  }
1304  else if( command==addsourceCmd )
1305    { 
1306      fGPS->AddaSource(addsourceCmd->GetNewDoubleValue(newValues));
1307    }
1308  else if( command==listsourceCmd )
1309    { 
1310      fGPS->ListSource();
1311    }
1312  else if( command==clearsourceCmd )
1313    { 
1314      fGPS->ClearAll();
1315    }
1316  else if( command==getsourceCmd )
1317    { 
1318      G4cout << " Current source index:" << fGPS->GetCurrentSourceIndex() 
1319             << " ; Intensity:" << fGPS->GetCurrentSourceIntensity() << G4endl;
1320    }
1321  else if( command==setsourceCmd )
1322    { 
1323      fGPS->SetCurrentSourceto(setsourceCmd->GetNewIntValue(newValues));
1324    }
1325  else if( command==setintensityCmd )
1326    { 
1327      fGPS->SetCurrentSourceIntensity(setintensityCmd->GetNewDoubleValue(newValues));
1328    }
1329  else if( command==deletesourceCmd )
1330    { 
1331      fGPS->DeleteaSource(deletesourceCmd->GetNewIntValue(newValues));
1332    }
1333  else if(command == multiplevertexCmd)
1334    {
1335      fGPS->SetMultipleVertex(multiplevertexCmd->GetNewBoolValue(newValues));
1336    }
1337  else if(command == flatsamplingCmd)
1338    {
1339      fGPS->SetFlatSampling(flatsamplingCmd->GetNewBoolValue(newValues));
1340    }
1341  //
1342  // new implementations
1343  //
1344  //
1345  else if(command == typeCmd1)
1346    {
1347      fParticleGun->GetPosDist()->SetPosDisType(newValues);
1348    }
1349  else if(command == shapeCmd1)
1350    {
1351      fParticleGun->GetPosDist()->SetPosDisShape(newValues);
1352    }
1353  else if(command == centreCmd1)
1354    {
1355      fParticleGun->GetPosDist()->SetCentreCoords(centreCmd1->GetNew3VectorValue(newValues));
1356    }
1357  else if(command == posrot1Cmd1)
1358    {
1359      fParticleGun->GetPosDist()->SetPosRot1(posrot1Cmd1->GetNew3VectorValue(newValues));
1360    }
1361  else if(command == posrot2Cmd1)
1362    {
1363      fParticleGun->GetPosDist()->SetPosRot2(posrot2Cmd1->GetNew3VectorValue(newValues));
1364    }
1365  else if(command == halfxCmd1)
1366    {
1367      fParticleGun->GetPosDist()->SetHalfX(halfxCmd1->GetNewDoubleValue(newValues));
1368    }
1369  else if(command == halfyCmd1)
1370    {
1371      fParticleGun->GetPosDist()->SetHalfY(halfyCmd1->GetNewDoubleValue(newValues));
1372    }
1373  else if(command == halfzCmd1)
1374    {
1375      fParticleGun->GetPosDist()->SetHalfZ(halfzCmd1->GetNewDoubleValue(newValues));
1376    }
1377  else if(command == radiusCmd1)
1378    {
1379      fParticleGun->GetPosDist()->SetRadius(radiusCmd1->GetNewDoubleValue(newValues));
1380    }
1381  else if(command == radius0Cmd1)
1382    {
1383      fParticleGun->GetPosDist()->SetRadius0(radius0Cmd1->GetNewDoubleValue(newValues));
1384    }
1385  else if(command == possigmarCmd1)
1386    {
1387      fParticleGun->GetPosDist()->SetBeamSigmaInR(possigmarCmd1->GetNewDoubleValue(newValues));
1388    }
1389  else if(command == possigmaxCmd1)
1390    {
1391      fParticleGun->GetPosDist()->SetBeamSigmaInX(possigmaxCmd1->GetNewDoubleValue(newValues));
1392    }
1393  else if(command == possigmayCmd1)
1394    {
1395      fParticleGun->GetPosDist()->SetBeamSigmaInY(possigmayCmd1->GetNewDoubleValue(newValues));
1396    }
1397  else if(command == paralpCmd1)
1398    {
1399      fParticleGun->GetPosDist()->SetParAlpha(paralpCmd1->GetNewDoubleValue(newValues));
1400    }
1401  else if(command == partheCmd1)
1402    {
1403      fParticleGun->GetPosDist()->SetParTheta(partheCmd1->GetNewDoubleValue(newValues));
1404    }
1405  else if(command == parphiCmd1)
1406    {
1407      fParticleGun->GetPosDist()->SetParPhi(parphiCmd1->GetNewDoubleValue(newValues));
1408    }
1409  else if(command == confineCmd1)
1410    {
1411      fParticleGun->GetPosDist()->ConfineSourceToVolume(newValues);
1412    }
1413  else if(command == angtypeCmd1)
1414    {
1415      fParticleGun->GetAngDist()->SetAngDistType(newValues);
1416    }
1417  else if(command == angrot1Cmd1)
1418    {
1419      G4String a = "angref1";
1420      fParticleGun->GetAngDist()->DefineAngRefAxes(a,angrot1Cmd1->GetNew3VectorValue(newValues));
1421    }
1422  else if(command == angrot2Cmd1)
1423    {
1424      G4String a = "angref2";
1425      fParticleGun->GetAngDist()->DefineAngRefAxes(a,angrot2Cmd1->GetNew3VectorValue(newValues));
1426    }
1427  else if(command == minthetaCmd1)
1428    {
1429      fParticleGun->GetAngDist()->SetMinTheta(minthetaCmd1->GetNewDoubleValue(newValues));
1430    }
1431  else if(command == minphiCmd1)
1432    {
1433      fParticleGun->GetAngDist()->SetMinPhi(minphiCmd1->GetNewDoubleValue(newValues));
1434    }
1435  else if(command == maxthetaCmd1)
1436    {
1437      fParticleGun->GetAngDist()->SetMaxTheta(maxthetaCmd1->GetNewDoubleValue(newValues));
1438    }
1439  else if(command == maxphiCmd1)
1440    {
1441      fParticleGun->GetAngDist()->SetMaxPhi(maxphiCmd1->GetNewDoubleValue(newValues));
1442    }
1443  else if(command == angsigmarCmd1)
1444    {
1445      fParticleGun->GetAngDist()->SetBeamSigmaInAngR(angsigmarCmd1->GetNewDoubleValue(newValues));
1446    }
1447  else if(command == angsigmaxCmd1)
1448    {
1449      fParticleGun->GetAngDist()->SetBeamSigmaInAngX(angsigmaxCmd1->GetNewDoubleValue(newValues));
1450    }
1451  else if(command == angsigmayCmd1)
1452    {
1453      fParticleGun->GetAngDist()->SetBeamSigmaInAngY(angsigmayCmd1->GetNewDoubleValue(newValues));
1454    }
1455  else if(command == angfocusCmd)
1456    {
1457      fParticleGun->GetAngDist()->SetFocusPoint(angfocusCmd->GetNew3VectorValue(newValues));
1458    }
1459  else if(command == useuserangaxisCmd1)
1460    {
1461      fParticleGun->GetAngDist()->SetUseUserAngAxis(useuserangaxisCmd1->GetNewBoolValue(newValues));
1462    }
1463  else if(command == surfnormCmd1)
1464    {
1465      fParticleGun->GetAngDist()->SetUserWRTSurface(surfnormCmd1->GetNewBoolValue(newValues));
1466    }
1467  else if(command == energytypeCmd1)
1468    {
1469      fParticleGun->GetEneDist()->SetEnergyDisType(newValues);
1470    }
1471  else if(command == eminCmd1)
1472    {
1473      fParticleGun->GetEneDist()->SetEmin(eminCmd1->GetNewDoubleValue(newValues));
1474    }
1475  else if(command == emaxCmd1)
1476    {
1477      fParticleGun->GetEneDist()->SetEmax(emaxCmd1->GetNewDoubleValue(newValues));
1478    }
1479  else if(command == monoenergyCmd1)
1480    {
1481      fParticleGun->GetEneDist()->SetMonoEnergy(monoenergyCmd1->GetNewDoubleValue(newValues));
1482    }
1483  else if(command == engsigmaCmd1)
1484    {
1485      fParticleGun->GetEneDist()->SetBeamSigmaInE(engsigmaCmd1->GetNewDoubleValue(newValues));
1486    }
1487  else if(command == alphaCmd1)
1488    {
1489      fParticleGun->GetEneDist()->SetAlpha(alphaCmd1->GetNewDoubleValue(newValues));
1490    }
1491  else if(command == tempCmd1)
1492    {
1493      fParticleGun->GetEneDist()->SetTemp(tempCmd1->GetNewDoubleValue(newValues));
1494    }
1495  else if(command == ezeroCmd1)
1496    {
1497      fParticleGun->GetEneDist()->SetEzero(ezeroCmd1->GetNewDoubleValue(newValues));
1498    }
1499  else if(command == gradientCmd1)
1500    {
1501      fParticleGun->GetEneDist()->SetGradient(gradientCmd1->GetNewDoubleValue(newValues));
1502    }
1503  else if(command == interceptCmd1)
1504    {
1505      fParticleGun->GetEneDist()->SetInterCept(interceptCmd1->GetNewDoubleValue(newValues));
1506    }
1507  else if(command == calculateCmd1)
1508    {
1509      fParticleGun->GetEneDist()->Calculate();
1510    }
1511  else if(command == energyspecCmd1)
1512    {
1513      fParticleGun->GetEneDist()->InputEnergySpectra(energyspecCmd1->GetNewBoolValue(newValues));
1514    }
1515  else if(command == diffspecCmd1)
1516    {
1517      fParticleGun->GetEneDist()->InputDifferentialSpectra(diffspecCmd1->GetNewBoolValue(newValues));
1518    }
1519  else if(command == histnameCmd1)
1520    {
1521      histtype = newValues;
1522    }
1523  else if(command == histpointCmd1)
1524    {
1525      if(histtype == "biasx")
1526        fParticleGun->GetBiasRndm()->SetXBias(histpointCmd1->GetNew3VectorValue(newValues));
1527      if(histtype == "biasy")
1528        fParticleGun->GetBiasRndm()->SetYBias(histpointCmd1->GetNew3VectorValue(newValues));
1529      if(histtype == "biasz")
1530        fParticleGun->GetBiasRndm()->SetZBias(histpointCmd1->GetNew3VectorValue(newValues));
1531      if(histtype == "biast")
1532        fParticleGun->GetBiasRndm()->SetThetaBias(histpointCmd1->GetNew3VectorValue(newValues));
1533      if(histtype == "biasp")
1534        fParticleGun->GetBiasRndm()->SetPhiBias(histpointCmd1->GetNew3VectorValue(newValues));
1535      if(histtype == "biaspt")
1536        fParticleGun->GetBiasRndm()->SetPosThetaBias(histpointCmd1->GetNew3VectorValue(newValues));
1537      if(histtype == "biaspp")
1538        fParticleGun->GetBiasRndm()->SetPosPhiBias(histpointCmd1->GetNew3VectorValue(newValues));
1539      if(histtype == "biase")
1540        fParticleGun->GetBiasRndm()->SetEnergyBias(histpointCmd1->GetNew3VectorValue(newValues));
1541      if(histtype == "theta")
1542        fParticleGun->GetAngDist()->UserDefAngTheta(histpointCmd1->GetNew3VectorValue(newValues));
1543      if(histtype == "phi")
1544        fParticleGun->GetAngDist()->UserDefAngPhi(histpointCmd1->GetNew3VectorValue(newValues));
1545      if(histtype == "energy")
1546        fParticleGun->GetEneDist()->UserEnergyHisto(histpointCmd1->GetNew3VectorValue(newValues));
1547      if(histtype == "arb")
1548        fParticleGun->GetEneDist()->ArbEnergyHisto(histpointCmd1->GetNew3VectorValue(newValues));
1549      if(histtype == "epn")
1550        fParticleGun->GetEneDist()->EpnEnergyHisto(histpointCmd1->GetNew3VectorValue(newValues));
1551    }
1552  else if(command == resethistCmd1)
1553    {
1554      if(newValues == "theta" || newValues == "phi") {
1555        fParticleGun->GetAngDist()->ReSetHist(newValues);
1556      } else if (newValues == "energy" || newValues == "arb" || newValues == "epn") {
1557        fParticleGun->GetEneDist()->ReSetHist(newValues);
1558      } else {
1559        fParticleGun->GetBiasRndm()->ReSetHist(newValues);
1560      }
1561    }
1562  else if(command == arbintCmd1)
1563    {
1564      fParticleGun->GetEneDist()->ArbInterpolate(newValues);
1565    }
1566  else 
1567    {
1568      G4cout << "Error entering command" << G4endl;
1569    }
1570}
1571
1572G4String G4GeneralParticleSourceMessenger::GetCurrentValue(G4UIcommand *)
1573{
1574  G4String cv;
1575 
1576  //  if( command==directionCmd )
1577  //  { cv = directionCmd->ConvertToString(fParticleGun->GetParticleMomentumDirection()); }
1578  //  else if( command==energyCmd )
1579  //  { cv = energyCmd->ConvertToString(fParticleGun->GetParticleEnergy(),"GeV"); }
1580  //  else if( command==positionCmd )
1581  //  { cv = positionCmd->ConvertToString(fParticleGun->GetParticlePosition(),"cm"); }
1582  //  else if( command==timeCmd )
1583  //  { cv = timeCmd->ConvertToString(fParticleGun->GetParticleTime(),"ns"); }
1584  //  else if( command==polCmd )
1585  //  { cv = polCmd->ConvertToString(fParticleGun->GetParticlePolarization()); }
1586  //  else if( command==numberCmd )
1587  //  { cv = numberCmd->ConvertToString(fParticleGun->GetNumberOfParticles()); }
1588 
1589  cv = "Not implemented yet";
1590
1591  return cv;
1592}
1593
1594void G4GeneralParticleSourceMessenger::IonCommand(G4String newValues)
1595{
1596  fShootIon = true;
1597
1598  if (fShootIon)
1599  {
1600    G4Tokenizer next( newValues );
1601    // check argument
1602    fAtomicNumber = StoI(next());
1603    fAtomicMass = StoI(next());
1604    G4String sQ = next();
1605    if (sQ.isNull())
1606    {
1607        fIonCharge = fAtomicNumber;
1608    }
1609    else
1610    {
1611        fIonCharge = StoI(sQ);
1612        sQ = next();
1613        if (sQ.isNull())
1614      {
1615          fIonExciteEnergy = 0.0;
1616      }
1617      else
1618      {
1619          fIonExciteEnergy = StoD(sQ) * keV;
1620      }
1621    }
1622    G4ParticleDefinition* ion;
1623    ion =  particleTable->GetIon( fAtomicNumber, fAtomicMass, fIonExciteEnergy);
1624    if (ion==0)
1625    {
1626      G4cout << "Ion with Z=" << fAtomicNumber;
1627      G4cout << " A=" << fAtomicMass << "is not be defined" << G4endl;   
1628    }
1629    else
1630    {
1631      fParticleGun->SetParticleDefinition(ion);
1632      fParticleGun->SetParticleCharge(fIonCharge*eplus);
1633    }
1634  }
1635  else
1636  {
1637    G4cout << "Set /gps/particle to ion before using /gps/ion command";
1638    G4cout << G4endl; 
1639  }
1640}
Note: See TracBrowser for help on using the repository browser.