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