Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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",false,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",false,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",false,false);
430 minphiCmd1->SetDefaultUnit("rad");
431 // minphiCmd1->SetUnitCandidates("rad deg");
432
433 maxphiCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ang/maxphi",this);
434 maxphiCmd1->SetGuidance("Set maximum phi");
435 maxphiCmd1->SetParameterName("MaxPhi",false,false);
436 maxphiCmd1->SetDefaultValue(pi);
437 maxphiCmd1->SetDefaultUnit("rad");
438 // maxphiCmd1->SetUnitCandidates("rad deg");
439
440 angsigmarCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ang/sigma_r",this);
441 angsigmarCmd1->SetGuidance("Set standard deviation in direction for 1D beam.");
442 angsigmarCmd1->SetParameterName("Sigmara",false,false);
443 angsigmarCmd1->SetDefaultUnit("rad");
444 // angsigmarCmd1->SetUnitCandidates("rad deg");
445
446 angsigmaxCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ang/sigma_x",this);
447 angsigmaxCmd1->SetGuidance("Set standard deviation in direction in x-direc. for 2D beam");
448 angsigmaxCmd1->SetParameterName("Sigmaxa",false,false);
449 angsigmaxCmd1->SetDefaultUnit("rad");
450 // angsigmaxCmd1->SetUnitCandidates("rad deg");
451
452 angsigmayCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ang/sigma_y",this);
453 angsigmayCmd1->SetGuidance("Set standard deviation in direction in y-direc. for 2D beam");
454 angsigmayCmd1->SetParameterName("Sigmaya",false,false);
455 angsigmayCmd1->SetDefaultUnit("rad");
456 // angsigmayCmd1->SetUnitCandidates("rad deg");
457
458 angfocusCmd = new G4UIcmdWith3VectorAndUnit("/gps/ang/focuspoint",this);
459 angfocusCmd->SetGuidance("Set the focusing point for the beam");
460 angfocusCmd->SetParameterName("x","y","z",false,false);
461 angfocusCmd->SetDefaultUnit("cm");
462 // angfocusCmd->SetUnitCandidates("micron mm cm m km");
463
464 useuserangaxisCmd1 = new G4UIcmdWithABool("/gps/ang/user_coor",this);
465 useuserangaxisCmd1->SetGuidance("True for using user defined angular co-ordinates");
466 useuserangaxisCmd1->SetGuidance(" Default is false");
467 useuserangaxisCmd1->SetParameterName("useuserangaxis",true);
468 useuserangaxisCmd1->SetDefaultValue(false);
469
470 surfnormCmd1 = new G4UIcmdWithABool("/gps/ang/surfnorm",this);
471 surfnormCmd1->SetGuidance("Makes a user-defined distribution with respect to surface normals rather than x,y,z axes.");
472 surfnormCmd1->SetGuidance(" Default is false");
473 surfnormCmd1->SetParameterName("surfnorm",true);
474 surfnormCmd1->SetDefaultValue(false);
475
476 // Energy commands
477 //
478 energyDirectory = new G4UIdirectory("/gps/ene/");
479 energyDirectory->SetGuidance("Spectral commands sub-directory");
480
481 energytypeCmd1 = new G4UIcmdWithAString("/gps/ene/type",this);
482 energytypeCmd1->SetGuidance("Sets energy distribution type");
483 energytypeCmd1->SetParameterName("EnergyDis",false,false);
484 energytypeCmd1->SetDefaultValue("Mono");
485 energytypeCmd1->SetCandidates("Mono Lin Pow Exp CPow Gauss Brem Bbody Cdg User Arb Epn LW");
486
487 eminCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ene/min",this);
488 eminCmd1->SetGuidance("Sets minimum energy");
489 eminCmd1->SetParameterName("emin",false,false);
490 eminCmd1->SetDefaultUnit("keV");
491 // eminCmd1->SetUnitCandidates("eV keV MeV GeV TeV PeV");
492
493 emaxCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ene/max",this);
494 emaxCmd1->SetGuidance("Sets maximum energy");
495 emaxCmd1->SetParameterName("emax",false,false);
496 emaxCmd1->SetDefaultUnit("keV");
497 // emaxCmd1->SetUnitCandidates("eV keV MeV GeV TeV PeV");
498
499 monoenergyCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ene/mono",this);
500 monoenergyCmd1->SetGuidance("Sets a monocromatic energy (same as gps/energy)");
501 monoenergyCmd1->SetParameterName("monoenergy",false,false);
502 monoenergyCmd1->SetDefaultUnit("keV");
503 // monoenergyCmd1->SetUnitCandidates("eV keV MeV GeV TeV PeV");
504
505 engsigmaCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ene/sigma",this);
506 engsigmaCmd1->SetGuidance("Sets the standard deviation for Gaussian energy dist.");
507 engsigmaCmd1->SetParameterName("Sigmae",false,false);
508 engsigmaCmd1->SetDefaultUnit("keV");
509 // engsigmaCmd1->SetUnitCandidates("eV keV MeV GeV TeV PeV");
510
511 alphaCmd1 = new G4UIcmdWithADouble("/gps/ene/alpha",this);
512 alphaCmd1->SetGuidance("Sets Alpha (index) for power-law energy dist.");
513 alphaCmd1->SetParameterName("alpha",false,false);
514
515 tempCmd1 = new G4UIcmdWithADouble("/gps/ene/temp",this);
516 tempCmd1->SetGuidance("Sets the temperature for Brem and BBody distributions (in Kelvin)");
517 tempCmd1->SetParameterName("temp",false,false);
518
519 ezeroCmd1 = new G4UIcmdWithADouble("/gps/ene/ezero",this);
520 ezeroCmd1->SetGuidance("Sets E_0 for exponential distribution (in MeV)");
521 ezeroCmd1->SetParameterName("ezero",false,false);
522
523 gradientCmd1 = new G4UIcmdWithADouble("/gps/ene/gradient",this);
524 gradientCmd1->SetGuidance("Sets the gradient for Lin distribution (in 1/MeV)");
525 gradientCmd1->SetParameterName("gradient",false,false);
526
527 interceptCmd1 = new G4UIcmdWithADouble("/gps/ene/intercept",this);
528 interceptCmd1->SetGuidance("Sets the intercept for Lin distributions (in MeV)");
529 interceptCmd1->SetParameterName("intercept",false,false);
530
531 arbeintCmd1 = new G4UIcmdWithADouble("/gps/ene/biasAlpha",this);
532 arbeintCmd1->SetGuidance("Sets the power-law index for the energy sampling distri. )");
533 arbeintCmd1->SetParameterName("arbeint",false,false);
534
535 calculateCmd1 = new G4UIcmdWithoutParameter("/gps/ene/calculate",this);
536 calculateCmd1->SetGuidance("Calculates the distributions for Cdg and BBody");
537
538 energyspecCmd1 = new G4UIcmdWithABool("/gps/ene/emspec",this);
539 energyspecCmd1->SetGuidance("True for energy and false for momentum spectra");
540 energyspecCmd1->SetParameterName("energyspec",true);
541 energyspecCmd1->SetDefaultValue(true);
542
543 diffspecCmd1 = new G4UIcmdWithABool("/gps/ene/diffspec",this);
544 diffspecCmd1->SetGuidance("True for differential and flase for integral spectra");
545 diffspecCmd1->SetParameterName("diffspec",true);
546 diffspecCmd1->SetDefaultValue(true);
547
548 applyEnergyWeightCmd1 = new G4UIcmdWithABool("/gps/ene/applyEneWeight",this);
549 applyEnergyWeightCmd1->SetGuidance("Apply energy weight.");
550 applyEnergyWeightCmd1->SetGuidance("- Instead of using the Arb type histogram for sampling the energy spectrum,");
551 applyEnergyWeightCmd1->SetGuidance(" energy is sampled by Linear distribution, and the Arb type histogram is");
552 applyEnergyWeightCmd1->SetGuidance(" used for the weight of the generated particle.");
553 applyEnergyWeightCmd1->SetGuidance("- \"/gps/ene/type LW\" automatically applies this command.");
554 applyEnergyWeightCmd1->SetGuidance("- If this command has to be explicitly used, \"/gps/ene/type\" distribution mush be Lin.");
555 applyEnergyWeightCmd1->SetParameterName("flag",true);
556 applyEnergyWeightCmd1->SetDefaultValue(true);
557
558 // Biasing + histograms in general
559 //
560 histDirectory = new G4UIdirectory("/gps/hist/");
561 histDirectory->SetGuidance("Histogram, biasing commands sub-directory");
562
563 histnameCmd1 = new G4UIcmdWithAString("/gps/hist/type",this);
564 histnameCmd1->SetGuidance("Sets histogram type");
565 histnameCmd1->SetParameterName("HistType",false,false);
566 histnameCmd1->SetDefaultValue("biasx");
567 histnameCmd1->SetCandidates("biasx biasy biasz biast biasp biase biaspt biaspp theta phi energy arb epn");
568
569 resethistCmd1 = new G4UIcmdWithAString("/gps/hist/reset",this);
570 resethistCmd1->SetGuidance("Reset (clean) the histogram ");
571 resethistCmd1->SetParameterName("HistType",false,false);
572 resethistCmd1->SetDefaultValue("energy");
573 resethistCmd1->SetCandidates("biasx biasy biasz biast biasp biase biaspt biaspp theta phi energy arb epn");
574
575 histpointCmd1 = new G4UIcmdWith3Vector("/gps/hist/point",this);
576 histpointCmd1->SetGuidance("Allows user to define a histogram");
577 histpointCmd1->SetGuidance(" Enter: Ehi Weight");
578 histpointCmd1->SetParameterName("Ehi","Weight","Junk",true,true);
579 histpointCmd1->SetRange("Ehi >= 0. && Weight >= 0.");
580
581 histfileCmd1 = new G4UIcmdWithAString("/gps/hist/file",this);
582 histfileCmd1->SetGuidance("Imports the arb energy hist in an ASCII file");
583 histfileCmd1->SetParameterName("HistFile",false,false);
584
585 arbintCmd1 = new G4UIcmdWithAString("/gps/hist/inter",this);
586 arbintCmd1->SetGuidance("Sets the interpolation method for arbitrary distribution.");
587 arbintCmd1->SetGuidance("Spline interpolation may not be applicable for some distributions.");
588 arbintCmd1->SetParameterName("int",false,false);
589 arbintCmd1->SetDefaultValue("Lin");
590 arbintCmd1->SetCandidates("Lin Log Exp Spline");
591}
592
593G4GeneralParticleSourceMessenger::~G4GeneralParticleSourceMessenger()
594{
595 delete positionDirectory;
596 delete typeCmd1;
597 delete shapeCmd1;
598 delete centreCmd1;
599 delete posrot1Cmd1;
600 delete posrot2Cmd1;
601 delete halfxCmd1;
602 delete halfyCmd1;
603 delete halfzCmd1;
604 delete radiusCmd1;
605 delete radius0Cmd1;
606 delete possigmarCmd1;
607 delete possigmaxCmd1;
608 delete possigmayCmd1;
609 delete paralpCmd1;
610 delete partheCmd1;
611 delete parphiCmd1;
612 delete confineCmd1;
613
614 delete angularDirectory;
615 delete angtypeCmd1;
616 delete angrot1Cmd1;
617 delete angrot2Cmd1;
618 delete minthetaCmd1;
619 delete maxthetaCmd1;
620 delete minphiCmd1;
621 delete maxphiCmd1;
622 delete angsigmarCmd1;
623 delete angsigmaxCmd1;
624 delete angsigmayCmd1;
625 delete angfocusCmd;
626 delete useuserangaxisCmd1;
627 delete surfnormCmd1;
628
629 delete energyDirectory;
630 delete energytypeCmd1;
631 delete eminCmd1;
632 delete emaxCmd1;
633 delete monoenergyCmd1;
634 delete engsigmaCmd1;
635 delete alphaCmd1;
636 delete tempCmd1;
637 delete ezeroCmd1;
638 delete gradientCmd1;
639 delete interceptCmd1;
640 delete arbeintCmd1;
641 delete calculateCmd1;
642 delete energyspecCmd1;
643 delete diffspecCmd1;
644 delete applyEnergyWeightCmd1;
645
646 delete histDirectory;
647 delete histnameCmd1;
648 delete resethistCmd1;
649 delete histpointCmd1;
650 delete histfileCmd1;
651 delete arbintCmd1;
652
653 delete verbosityCmd;
654 delete ionCmd;
655 delete ionLvlCmd;
656 delete particleCmd;
657 delete timeCmd;
658 delete polCmd;
659 delete numberCmd;
660 delete positionCmd;
661 delete directionCmd;
662 delete energyCmd;
663 delete listCmd;
664
665 delete sourceDirectory;
666 delete addsourceCmd;
667 delete listsourceCmd;
668 delete clearsourceCmd;
669 delete getsourceCmd;
670 delete setsourceCmd;
671 delete setintensityCmd;
672 delete deletesourceCmd;
673 delete multiplevertexCmd;
674 delete flatsamplingCmd;
675
676 delete gpsDirectory;
677 theInstance = nullptr;
678}
679
680#define CHECKPG() { if (fParticleGun==nullptr) { \
681 G4ExceptionDescription msg; \
682 msg << "Command "<< command->GetCommandPath()<<"/";\
683 msg << command->GetCommandName(); \
684 msg << " used but no particle sources are set.";\
685 msg <<" Add at least a source with: /gps/source/add.";\
686 G4Exception("G4GeneralParticleSourceMessenger::SetNewValue","G4GPS003",\
687 FatalException,msg); return;\
688 } }
689
691{
692// if(command == typeCmd)
693// {
694// CHECKPG(); fParticleGun->GetPosDist()->SetPosDisType(newValues);
695// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
696// << " The command is obsolete and will be removed soon." << G4endl
697// << " Please try to use the new structured commands!" << G4endl;
698// }
699// else if(command == shapeCmd)
700// {
701// CHECKPG(); fParticleGun->GetPosDist()->SetPosDisShape(newValues);
702// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
703// << " The command is obsolete and will be removed soon." << G4endl
704// << " Please try to use the new structured commands!" << G4endl;
705// }
706// else if(command == centreCmd)
707// {
708// CHECKPG(); fParticleGun->GetPosDist()->SetCentreCoords(centreCmd->GetNew3VectorValue(newValues));
709// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
710// << " The command is obsolete and will be removed soon." << G4endl
711// << " Please try to use the new structured commands!" << G4endl;
712// }
713// else if(command == posrot1Cmd)
714// {
715// CHECKPG(); fParticleGun->GetPosDist()->SetPosRot1(posrot1Cmd->GetNew3VectorValue(newValues));
716// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
717// << " The command is obsolete and will be removed soon." << G4endl
718// << " Please try to use the new structured commands!" << G4endl;
719// }
720// else if(command == posrot2Cmd)
721// {
722// CHECKPG(); fParticleGun->GetPosDist()->SetPosRot2(posrot2Cmd->GetNew3VectorValue(newValues));
723// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
724// << " The command is obsolete and will be removed soon." << G4endl
725// << " Please try to use the new structured commands!" << G4endl;
726// }
727// else if(command == halfxCmd)
728// {
729// CHECKPG(); fParticleGun->GetPosDist()->SetHalfX(halfxCmd->GetNewDoubleValue(newValues));
730// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
731// << " The command is obsolete and will be removed soon." << G4endl
732// << " Please try to use the new structured commands!" << G4endl;
733// }
734// else if(command == halfyCmd)
735// {
736// CHECKPG(); fParticleGun->GetPosDist()->SetHalfY(halfyCmd->GetNewDoubleValue(newValues));
737// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
738// << " The command is obsolete and will be removed soon." << G4endl
739// << " Please try to use the new structured commands!" << G4endl;
740// }
741// else if(command == halfzCmd)
742// {
743// CHECKPG(); fParticleGun->GetPosDist()->SetHalfZ(halfzCmd->GetNewDoubleValue(newValues));
744// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
745// << " The command is obsolete and will be removed soon." << G4endl
746// << " Please try to use the new structured commands!" << G4endl;
747// }
748// else if(command == radiusCmd)
749// {
750// CHECKPG(); fParticleGun->GetPosDist()->SetRadius(radiusCmd->GetNewDoubleValue(newValues));
751// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
752// << " The command is obsolete and will be removed soon." << G4endl
753// << " Please try to use the new structured commands!" << G4endl;
754// }
755// else if(command == radius0Cmd)
756// {
757// CHECKPG(); fParticleGun->GetPosDist()->SetRadius0(radius0Cmd->GetNewDoubleValue(newValues));
758// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
759// << " The command is obsolete and will be removed soon." << G4endl
760// << " Please try to use the new structured commands!" << G4endl;
761// }
762// else if(command == possigmarCmd)
763// {
764// CHECKPG(); fParticleGun->GetPosDist()->SetBeamSigmaInR(possigmarCmd->GetNewDoubleValue(newValues));
765// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
766// << " The command is obsolete and will be removed soon." << G4endl
767// << " Please try to use the new structured commands!" << G4endl;
768// }
769// else if(command == possigmaxCmd)
770// {
771// CHECKPG(); fParticleGun->GetPosDist()->SetBeamSigmaInX(possigmaxCmd->GetNewDoubleValue(newValues));
772// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
773// << " The command is obsolete and will be removed soon." << G4endl
774// << " Please try to use the new structured commands!" << G4endl;
775// }
776// else if(command == possigmayCmd)
777// {
778// CHECKPG(); fParticleGun->GetPosDist()->SetBeamSigmaInY(possigmayCmd->GetNewDoubleValue(newValues));
779// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
780// << " The command is obsolete and will be removed soon." << G4endl
781// << " Please try to use the new structured commands!" << G4endl;
782// }
783// else if(command == paralpCmd)
784// {
785// CHECKPG(); fParticleGun->GetPosDist()->SetParAlpha(paralpCmd->GetNewDoubleValue(newValues));
786// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
787// << " The command is obsolete and will be removed soon." << G4endl
788// << " Please try to use the new structured commands!" << G4endl;
789// }
790// else if(command == partheCmd)
791// {
792// CHECKPG(); fParticleGun->GetPosDist()->SetParTheta(partheCmd->GetNewDoubleValue(newValues));
793// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
794// << " The command is obsolete and will be removed soon." << G4endl
795// << " Please try to use the new structured commands!" << G4endl;
796// }
797// else if(command == parphiCmd)
798// {
799// CHECKPG(); fParticleGun->GetPosDist()->SetParPhi(parphiCmd->GetNewDoubleValue(newValues));
800// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
801// << " The command is obsolete and will be removed soon." << G4endl
802// << " Please try to use the new structured commands!" << G4endl;
803// }
804// else if(command == confineCmd)
805// {
806// CHECKPG(); fParticleGun->GetPosDist()->ConfineSourceToVolume(newValues);
807// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
808// << " The command is obsolete and will be removed soon." << G4endl
809// << " Please try to use the new structured commands!" << G4endl;
810// }
811// else if(command == angtypeCmd)
812// {
813// CHECKPG(); fParticleGun->GetAngDist()->SetAngDistType(newValues);
814// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
815// << " The command is obsolete and will be removed soon." << G4endl
816// << " Please try to use the new structured commands!" << G4endl;
817// }
818// else if(command == angrot1Cmd)
819// {
820// CHECKPG();
821// G4String a = "angref1";
822// fParticleGun->GetAngDist()->DefineAngRefAxes(a,angrot1Cmd->GetNew3VectorValue(newValues));
823// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
824// << " The command is obsolete and will be removed soon." << G4endl
825// << " Please try to use the new structured commands!" << G4endl;
826// }
827// else if(command == angrot2Cmd)
828// {
829// CHECKPG();
830// G4String a = "angref2";
831// fParticleGun->GetAngDist()->DefineAngRefAxes(a,angrot2Cmd->GetNew3VectorValue(newValues));
832// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
833// << " The command is obsolete and will be removed soon." << G4endl
834// << " Please try to use the new structured commands!" << G4endl;
835// }
836// else if(command == minthetaCmd)
837// {
838// CHECKPG(); fParticleGun->GetAngDist()->SetMinTheta(minthetaCmd->GetNewDoubleValue(newValues));
839// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
840// << " The command is obsolete and will be removed soon." << G4endl
841// << " Please try to use the new structured commands!" << G4endl;
842// }
843// else if(command == minphiCmd)
844// {
845// CHECKPG(); fParticleGun->GetAngDist()->SetMinPhi(minphiCmd->GetNewDoubleValue(newValues));
846// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
847// << " The command is obsolete and will be removed soon." << G4endl
848// << " Please try to use the new structured commands!" << G4endl;
849// }
850// else if(command == maxthetaCmd)
851// {
852// CHECKPG(); fParticleGun->GetAngDist()->SetMaxTheta(maxthetaCmd->GetNewDoubleValue(newValues));
853// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
854// << " The command is obsolete and will be removed soon." << G4endl
855// << " Please try to use the new structured commands!" << G4endl;
856// }
857// else if(command == maxphiCmd)
858// {
859// CHECKPG(); fParticleGun->GetAngDist()->SetMaxPhi(maxphiCmd->GetNewDoubleValue(newValues));
860// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
861// << " The command is obsolete and will be removed soon." << G4endl
862// << " Please try to use the new structured commands!" << G4endl;
863// }
864// else if(command == angsigmarCmd)
865// {
866// CHECKPG(); fParticleGun->GetAngDist()->SetBeamSigmaInAngR(angsigmarCmd->GetNewDoubleValue(newValues));
867// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
868// << " The command is obsolete and will be removed soon." << G4endl
869// << " Please try to use the new structured commands!" << G4endl;
870// }
871// else if(command == angsigmaxCmd)
872// {
873// CHECKPG(); fParticleGun->GetAngDist()->SetBeamSigmaInAngX(angsigmaxCmd->GetNewDoubleValue(newValues));
874// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
875// << " The command is obsolete and will be removed soon." << G4endl
876// << " Please try to use the new structured commands!" << G4endl;
877// }
878// else if(command == angsigmayCmd)
879// {
880// CHECKPG(); fParticleGun->GetAngDist()->SetBeamSigmaInAngY(angsigmayCmd->GetNewDoubleValue(newValues));
881// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
882// << " The command is obsolete and will be removed soon." << G4endl
883// << " Please try to use the new structured commands!" << G4endl;
884// }
885// else if(command == useuserangaxisCmd)
886// {
887// CHECKPG(); fParticleGun->GetAngDist()->SetUseUserAngAxis(useuserangaxisCmd->GetNewBoolValue(newValues));
888// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
889// << " The command is obsolete and will be removed soon." << G4endl
890// << " Please try to use the new structured commands!" << G4endl;
891// }
892// else if(command == surfnormCmd)
893// {
894// CHECKPG(); fParticleGun->GetAngDist()->SetUserWRTSurface(surfnormCmd->GetNewBoolValue(newValues));
895// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
896// << " The command is obsolete and will be removed soon." << G4endl
897// << " Please try to use the new structured commands!" << G4endl;
898// }
899// else if(command == energytypeCmd)
900// {
901// CHECKPG(); fParticleGun->GetEneDist()->SetEnergyDisType(newValues);
902// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
903// << " The command is obsolete and will be removed soon." << G4endl
904// << " Please try to use the new structured commands!" << G4endl;
905// }
906// else if(command == eminCmd)
907// {
908// CHECKPG(); fParticleGun->GetEneDist()->SetEmin(eminCmd->GetNewDoubleValue(newValues));
909// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
910// << " The command is obsolete and will be removed soon." << G4endl
911// << " Please try to use the new structured commands!" << G4endl;
912// }
913// else if(command == emaxCmd)
914// {
915// CHECKPG(); fParticleGun->GetEneDist()->SetEmax(emaxCmd->GetNewDoubleValue(newValues));
916// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
917// << " The command is obsolete and will be removed soon." << G4endl
918// << " Please try to use the new structured commands!" << G4endl;
919// }
920// else if(command == monoenergyCmd)
921// {
922// CHECKPG(); fParticleGun->GetEneDist()->SetMonoEnergy(monoenergyCmd->GetNewDoubleValue(newValues));
923// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
924// << " The command is obsolete and will be removed soon." << G4endl
925// << " Please try to use the new structured commands!" << G4endl;
926// }
927// else if(command == engsigmaCmd)
928// {
929// CHECKPG(); fParticleGun->GetEneDist()->SetBeamSigmaInE(engsigmaCmd->GetNewDoubleValue(newValues));
930// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
931// << " The command is obsolete and will be removed soon." << G4endl
932// << " Please try to use the new structured commands!" << G4endl;
933// }
934// else if(command == alphaCmd)
935// {
936// CHECKPG(); fParticleGun->GetEneDist()->SetAlpha(alphaCmd->GetNewDoubleValue(newValues));
937// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
938// << " The command is obsolete and will be removed soon." << G4endl
939// << " Please try to use the new structured commands!" << G4endl;
940// }
941// else if(command == tempCmd)
942// {
943// CHECKPG(); fParticleGun->GetEneDist()->SetTemp(tempCmd->GetNewDoubleValue(newValues));
944// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
945// << " The command is obsolete and will be removed soon." << G4endl
946// << " Please try to use the new structured commands!" << G4endl;
947// }
948// else if(command == ezeroCmd)
949// {
950// CHECKPG(); fParticleGun->GetEneDist()->SetEzero(ezeroCmd->GetNewDoubleValue(newValues));
951// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
952// << " The command is obsolete and will be removed soon." << G4endl
953// << " Please try to use the new structured commands!" << G4endl;
954// }
955// else if(command == gradientCmd)
956// {
957// CHECKPG(); fParticleGun->GetEneDist()->SetGradient(gradientCmd->GetNewDoubleValue(newValues));
958// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
959// << " The command is obsolete and will be removed soon." << G4endl
960// << " Please try to use the new structured commands!" << G4endl;
961// }
962// else if(command == interceptCmd)
963// {
964// CHECKPG(); fParticleGun->GetEneDist()->SetInterCept(interceptCmd->GetNewDoubleValue(newValues));
965// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
966// << " The command is obsolete and will be removed soon." << G4endl
967// << " Please try to use the new structured commands!" << G4endl;
968// }
969// else if(command == calculateCmd)
970// {
971// CHECKPG(); fParticleGun->GetEneDist()->Calculate();
972// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
973// << " The command is obsolete and will be removed soon." << G4endl
974// << " Please try to use the new structured commands!" << G4endl;
975// }
976// else if(command == energyspecCmd)
977// {
978// CHECKPG(); fParticleGun->GetEneDist()->InputEnergySpectra(energyspecCmd->GetNewBoolValue(newValues));
979// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
980// << " The command is obsolete and will be removed soon." << G4endl
981// << " Please try to use the new structured commands!" << G4endl;
982// }
983// else if(command == diffspecCmd)
984// {
985// CHECKPG(); fParticleGun->GetEneDist()->InputDifferentialSpectra(diffspecCmd->GetNewBoolValue(newValues));
986// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
987// << " The command is obsolete and will be removed soon." << G4endl
988// << " Please try to use the new structured commands!" << G4endl;
989// }
990// else if(command == histnameCmd)
991// {
992// histtype = newValues;
993// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
994// << " The command is obsolete and will be removed soon." << G4endl
995// << " Please try to use the new structured commands!" << G4endl;
996// }
997// else
998// if(command == histpointCmd)
999// {
1000// CHECKPG();
1001// if(histtype == "biasx")
1002// fParticleGun->GetBiasRndm()->SetXBias(histpointCmd->GetNew3VectorValue(newValues));
1003// if(histtype == "biasy")
1004// fParticleGun->GetBiasRndm()->SetYBias(histpointCmd->GetNew3VectorValue(newValues));
1005// if(histtype == "biasz")
1006// fParticleGun->GetBiasRndm()->SetZBias(histpointCmd->GetNew3VectorValue(newValues));
1007// if(histtype == "biast")
1008// fParticleGun->GetBiasRndm()->SetThetaBias(histpointCmd->GetNew3VectorValue(newValues));
1009// if(histtype == "biasp")
1010// fParticleGun->GetBiasRndm()->SetPhiBias(histpointCmd->GetNew3VectorValue(newValues));
1011// if(histtype == "biase")
1012// fParticleGun->GetBiasRndm()->SetEnergyBias(histpointCmd->GetNew3VectorValue(newValues));
1013// if(histtype == "theta")
1014// fParticleGun->GetAngDist()->UserDefAngTheta(histpointCmd->GetNew3VectorValue(newValues));
1015// if(histtype == "phi")
1016// fParticleGun->GetAngDist()->UserDefAngPhi(histpointCmd->GetNew3VectorValue(newValues));
1017// if(histtype == "energy")
1018// fParticleGun->GetEneDist()->UserEnergyHisto(histpointCmd->GetNew3VectorValue(newValues));
1019// if(histtype == "arb")
1020// fParticleGun->GetEneDist()->ArbEnergyHisto(histpointCmd->GetNew3VectorValue(newValues));
1021// if(histtype == "epn")
1022// fParticleGun->GetEneDist()->EpnEnergyHisto(histpointCmd->GetNew3VectorValue(newValues));
1023// G4cout << " G4GeneralParticleSourceMessenger - Warning: The command is obsolete and will be removed soon. Please try to use the new structured commands!" << G4endl;
1024// }
1025// else if(command == resethistCmd)
1026// {
1027// CHECKPG();
1028// if(newValues == "theta" || newValues == "phi") {
1029// fParticleGun->GetAngDist()->ReSetHist(newValues);
1030// } else if (newValues == "energy" || newValues == "arb" || newValues == "epn") {
1031// fParticleGun->GetEneDist()->ReSetHist(newValues);
1032// } else {
1033// fParticleGun->GetBiasRndm()->ReSetHist(newValues);
1034// }
1035// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1036// << " The command is obsolete and will be removed soon." << G4endl
1037// << " Please try to use the new structured commands!" << G4endl;
1038// }
1039// else if(command == arbintCmd)
1040// {
1041// CHECKPG();
1042// fParticleGun->GetEneDist()->ArbInterpolate(newValues);
1043// G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
1044// << " The command is obsolete and will be removed soon." << G4endl
1045// << " Please try to use the new structured commands!" << G4endl;
1046// }
1047// else
1048 if( command==directionCmd )
1049 {
1050 CHECKPG();
1051 fParticleGun->GetAngDist()->SetAngDistType("planar");
1052 fParticleGun->GetAngDist()->SetParticleMomentumDirection(directionCmd->GetNew3VectorValue(newValues));
1053 }
1054 else if( command==energyCmd )
1055 {
1056 CHECKPG();
1057 fParticleGun->GetEneDist()->SetEnergyDisType("Mono");
1058 fParticleGun->GetEneDist()->SetMonoEnergy(energyCmd->GetNewDoubleValue(newValues));
1059 }
1060 else if( command==positionCmd )
1061 {
1062 CHECKPG();
1063 fParticleGun->GetPosDist()->SetPosDisType("Point");
1064 fParticleGun->GetPosDist()->SetCentreCoords(positionCmd->GetNew3VectorValue(newValues));
1065 }
1066 else if(command == verbosityCmd)
1067 {
1068 fGPS->SetVerbosity(verbosityCmd->GetNewIntValue(newValues));
1069 // CHECKPG();
1070 // fParticleGun->SetVerbosity(verbosityCmd->GetNewIntValue(newValues));
1071 }
1072 else if( command==particleCmd )
1073 {
1074 if (newValues =="ion")
1075 {
1076 fShootIon = true;
1077 }
1078 else
1079 {
1080 fShootIon = false;
1081 G4ParticleDefinition* pd = particleTable->FindParticle(newValues);
1082 if(pd != nullptr)
1083 {
1084 CHECKPG();
1085 fParticleGun->SetParticleDefinition( pd );
1086 }
1087 }
1088 }
1089 else if( command==timeCmd )
1090 {
1091 CHECKPG();
1092 fParticleGun->SetParticleTime(timeCmd->GetNewDoubleValue(newValues));
1093 }
1094 else if( command==polCmd )
1095 {
1096 CHECKPG();
1097 fParticleGun->SetParticlePolarization(polCmd->GetNew3VectorValue(newValues));
1098 }
1099 else if( command==numberCmd )
1100 {
1101 CHECKPG();
1102 fParticleGun->SetNumberOfParticles(numberCmd->GetNewIntValue(newValues));
1103 }
1104 else if( command==ionCmd )
1105 {
1106 IonCommand(newValues);
1107 }
1108 else if( command==ionLvlCmd )
1109 {
1110 IonLvlCommand(newValues);
1111 }
1112 else if( command==listCmd )
1113 {
1114 particleTable->DumpTable();
1115 }
1116 else if( command==addsourceCmd )
1117 {
1118 fGPS->AddaSource(addsourceCmd->GetNewDoubleValue(newValues));
1119 }
1120 else if( command==listsourceCmd )
1121 {
1122 fGPS->ListSource();
1123 }
1124 else if( command==clearsourceCmd )
1125 {
1126 fGPS->ClearAll();
1127 fParticleGun = nullptr;
1128 }
1129 else if( command==getsourceCmd )
1130 {
1131 G4cout << " Current source index:" << fGPS->GetCurrentSourceIndex()
1132 << " ; Intensity:" << fGPS->GetCurrentSourceIntensity() << G4endl;
1133 }
1134 else if( command==setsourceCmd )
1135 {
1136 // NOTE: This will also sets fParticleGun to the courrent source
1137 // Not very clean, the GPS::SetCurrentSourceto will call:
1138 // SetParticleSource( G4ParticleSource* )
1139 // The point is that GPS has no public API to get a source by index
1140 // TODO: Can we add this API?
1141 const G4int sn = setsourceCmd->GetNewIntValue(newValues);
1142 if ( sn >= fGPS->GetNumberofSource() )
1143 {
1145 msg << "Using command " << setsourceCmd->GetCommandPath() << "/"
1146 << setsourceCmd->GetCommandName() << " " << sn;
1147 msg << " is invalid " << fGPS->GetNumberofSource()
1148 << " source(s) are defined.";
1149 G4Exception("G4GeneralParticleSourceMessenger::SetNewValue",
1150 "G4GPS005", FatalException, msg);
1151 }
1152 fGPS->SetCurrentSourceto(setsourceCmd->GetNewIntValue(newValues));
1153 }
1154 else if( command==setintensityCmd )
1155 {
1156 fGPS->SetCurrentSourceIntensity(setintensityCmd->GetNewDoubleValue(newValues));
1157 }
1158 else if( command==deletesourceCmd )
1159 {
1160 fGPS->DeleteaSource(deletesourceCmd->GetNewIntValue(newValues));
1161 }
1162 else if(command == multiplevertexCmd)
1163 {
1164 fGPS->SetMultipleVertex(multiplevertexCmd->GetNewBoolValue(newValues));
1165 }
1166 else if(command == flatsamplingCmd)
1167 {
1168 fGPS->SetFlatSampling(flatsamplingCmd->GetNewBoolValue(newValues));
1169 }
1170 //
1171 // new implementations
1172 //
1173 else if(command == typeCmd1)
1174 {
1175 CHECKPG();
1176 fParticleGun->GetPosDist()->SetPosDisType(newValues);
1177 }
1178 else if(command == shapeCmd1)
1179 {
1180 CHECKPG();
1181 fParticleGun->GetPosDist()->SetPosDisShape(newValues);
1182 }
1183 else if(command == centreCmd1)
1184 {
1185 CHECKPG();
1186 fParticleGun->GetPosDist()->SetCentreCoords(centreCmd1->GetNew3VectorValue(newValues));
1187 }
1188 else if(command == posrot1Cmd1)
1189 {
1190 CHECKPG();
1191 fParticleGun->GetPosDist()->SetPosRot1(posrot1Cmd1->GetNew3VectorValue(newValues));
1192 }
1193 else if(command == posrot2Cmd1)
1194 {
1195 CHECKPG();
1196 fParticleGun->GetPosDist()->SetPosRot2(posrot2Cmd1->GetNew3VectorValue(newValues));
1197 }
1198 else if(command == halfxCmd1)
1199 {
1200 CHECKPG();
1201 fParticleGun->GetPosDist()->SetHalfX(halfxCmd1->GetNewDoubleValue(newValues));
1202 }
1203 else if(command == halfyCmd1)
1204 {
1205 CHECKPG();
1206 fParticleGun->GetPosDist()->SetHalfY(halfyCmd1->GetNewDoubleValue(newValues));
1207 }
1208 else if(command == halfzCmd1)
1209 {
1210 CHECKPG();
1211 fParticleGun->GetPosDist()->SetHalfZ(halfzCmd1->GetNewDoubleValue(newValues));
1212 }
1213 else if(command == radiusCmd1)
1214 {
1215 CHECKPG();
1216 fParticleGun->GetPosDist()->SetRadius(radiusCmd1->GetNewDoubleValue(newValues));
1217 }
1218 else if(command == radius0Cmd1)
1219 {
1220 CHECKPG();
1221 fParticleGun->GetPosDist()->SetRadius0(radius0Cmd1->GetNewDoubleValue(newValues));
1222 }
1223 else if(command == possigmarCmd1)
1224 {
1225 CHECKPG();
1226 fParticleGun->GetPosDist()->SetBeamSigmaInR(possigmarCmd1->GetNewDoubleValue(newValues));
1227 }
1228 else if(command == possigmaxCmd1)
1229 {
1230 CHECKPG();
1231 fParticleGun->GetPosDist()->SetBeamSigmaInX(possigmaxCmd1->GetNewDoubleValue(newValues));
1232 }
1233 else if(command == possigmayCmd1)
1234 {
1235 CHECKPG();
1236 fParticleGun->GetPosDist()->SetBeamSigmaInY(possigmayCmd1->GetNewDoubleValue(newValues));
1237 }
1238 else if(command == paralpCmd1)
1239 {
1240 CHECKPG();
1241 fParticleGun->GetPosDist()->SetParAlpha(paralpCmd1->GetNewDoubleValue(newValues));
1242 }
1243 else if(command == partheCmd1)
1244 {
1245 CHECKPG();
1246 fParticleGun->GetPosDist()->SetParTheta(partheCmd1->GetNewDoubleValue(newValues));
1247 }
1248 else if(command == parphiCmd1)
1249 {
1250 CHECKPG();
1251 fParticleGun->GetPosDist()->SetParPhi(parphiCmd1->GetNewDoubleValue(newValues));
1252 }
1253 else if(command == confineCmd1)
1254 {
1255 CHECKPG();
1256 fParticleGun->GetPosDist()->ConfineSourceToVolume(newValues);
1257 }
1258 else if(command == angtypeCmd1)
1259 {
1260 CHECKPG();
1261 fParticleGun->GetAngDist()->SetAngDistType(newValues);
1262 }
1263 else if(command == angrot1Cmd1)
1264 {
1265 CHECKPG();
1266 G4String a = "angref1";
1267 fParticleGun->GetAngDist()->DefineAngRefAxes(a,angrot1Cmd1->GetNew3VectorValue(newValues));
1268 }
1269 else if(command == angrot2Cmd1)
1270 {
1271 CHECKPG();
1272 G4String a = "angref2";
1273 fParticleGun->GetAngDist()->DefineAngRefAxes(a,angrot2Cmd1->GetNew3VectorValue(newValues));
1274 }
1275 else if(command == minthetaCmd1)
1276 {
1277 CHECKPG();
1278 fParticleGun->GetAngDist()->SetMinTheta(minthetaCmd1->GetNewDoubleValue(newValues));
1279 }
1280 else if(command == minphiCmd1)
1281 {
1282 CHECKPG();
1283 fParticleGun->GetAngDist()->SetMinPhi(minphiCmd1->GetNewDoubleValue(newValues));
1284 }
1285 else if(command == maxthetaCmd1)
1286 {
1287 CHECKPG();
1288 fParticleGun->GetAngDist()->SetMaxTheta(maxthetaCmd1->GetNewDoubleValue(newValues));
1289 }
1290 else if(command == maxphiCmd1)
1291 {
1292 CHECKPG();
1293 fParticleGun->GetAngDist()->SetMaxPhi(maxphiCmd1->GetNewDoubleValue(newValues));
1294 }
1295 else if(command == angsigmarCmd1)
1296 {
1297 CHECKPG();
1298 fParticleGun->GetAngDist()->SetBeamSigmaInAngR(angsigmarCmd1->GetNewDoubleValue(newValues));
1299 }
1300 else if(command == angsigmaxCmd1)
1301 {
1302 CHECKPG();
1303 fParticleGun->GetAngDist()->SetBeamSigmaInAngX(angsigmaxCmd1->GetNewDoubleValue(newValues));
1304 }
1305 else if(command == angsigmayCmd1)
1306 {
1307 CHECKPG();
1308 fParticleGun->GetAngDist()->SetBeamSigmaInAngY(angsigmayCmd1->GetNewDoubleValue(newValues));
1309 }
1310 else if(command == angfocusCmd)
1311 {
1312 CHECKPG();
1313 fParticleGun->GetAngDist()->SetFocusPoint(angfocusCmd->GetNew3VectorValue(newValues));
1314 }
1315 else if(command == useuserangaxisCmd1)
1316 {
1317 CHECKPG();
1318 fParticleGun->GetAngDist()->SetUseUserAngAxis(useuserangaxisCmd1->GetNewBoolValue(newValues));
1319 }
1320 else if(command == surfnormCmd1)
1321 {
1322 CHECKPG();
1323 fParticleGun->GetAngDist()->SetUserWRTSurface(surfnormCmd1->GetNewBoolValue(newValues));
1324 }
1325 else if(command == energytypeCmd1)
1326 {
1327 CHECKPG();
1328 if(newValues=="LW")
1329 {
1330 fParticleGun->GetEneDist()->SetEnergyDisType("Lin");
1331 fParticleGun->GetEneDist()->SetGradient(0.);
1332 fParticleGun->GetEneDist()->SetInterCept(1.);
1333 fParticleGun->GetEneDist()->ApplyEnergyWeight(true);
1334 }
1335 else
1336 {
1337 fParticleGun->GetEneDist()->SetEnergyDisType(newValues);
1338 fParticleGun->GetEneDist()->ApplyEnergyWeight(false);
1339 }
1340 }
1341 else if(command == eminCmd1)
1342 {
1343 CHECKPG();
1344 fParticleGun->GetEneDist()->SetEmin(eminCmd1->GetNewDoubleValue(newValues));
1345 }
1346 else if(command == emaxCmd1)
1347 {
1348 CHECKPG();
1349 fParticleGun->GetEneDist()->SetEmax(emaxCmd1->GetNewDoubleValue(newValues));
1350 }
1351 else if(command == monoenergyCmd1)
1352 {
1353 CHECKPG();
1354 fParticleGun->GetEneDist()->SetMonoEnergy(monoenergyCmd1->GetNewDoubleValue(newValues));
1355 }
1356 else if(command == engsigmaCmd1)
1357 {
1358 CHECKPG();
1359 fParticleGun->GetEneDist()->SetBeamSigmaInE(engsigmaCmd1->GetNewDoubleValue(newValues));
1360 }
1361 else if(command == alphaCmd1)
1362 {
1363 CHECKPG();
1364 fParticleGun->GetEneDist()->SetAlpha(alphaCmd1->GetNewDoubleValue(newValues));
1365 }
1366 else if(command == tempCmd1)
1367 {
1368 CHECKPG();
1369 fParticleGun->GetEneDist()->SetTemp(tempCmd1->GetNewDoubleValue(newValues));
1370 }
1371 else if(command == ezeroCmd1)
1372 {
1373 CHECKPG();
1374 fParticleGun->GetEneDist()->SetEzero(ezeroCmd1->GetNewDoubleValue(newValues));
1375 }
1376 else if(command == gradientCmd1)
1377 {
1378 CHECKPG();
1379 fParticleGun->GetEneDist()->SetGradient(gradientCmd1->GetNewDoubleValue(newValues));
1380 }
1381 else if(command == interceptCmd1)
1382 {
1383 CHECKPG();
1384 fParticleGun->GetEneDist()->SetInterCept(interceptCmd1->GetNewDoubleValue(newValues));
1385 }
1386 else if(command == arbeintCmd1)
1387 {
1388 CHECKPG();
1389 fParticleGun->GetEneDist()->SetBiasAlpha(arbeintCmd1->GetNewDoubleValue(newValues));
1390 }
1391 else if(command == calculateCmd1)
1392 {
1393 CHECKPG();
1394 fParticleGun->GetEneDist()->Calculate();
1395 }
1396 else if(command == energyspecCmd1)
1397 {
1398 CHECKPG();
1399 fParticleGun->GetEneDist()->InputEnergySpectra(energyspecCmd1->GetNewBoolValue(newValues));
1400 }
1401 else if(command == diffspecCmd1)
1402 {
1403 CHECKPG();
1404 fParticleGun->GetEneDist()->InputDifferentialSpectra(diffspecCmd1->GetNewBoolValue(newValues));
1405 }
1406 else if(command == applyEnergyWeightCmd1)
1407 {
1408 CHECKPG();
1409 auto eDisType = fParticleGun->GetEneDist()->GetEnergyDisType();
1410 if(eDisType != "Lin")
1411 {
1413 ed << "Energy distribution is defined as " << eDisType << ". /gps/ene/applyEneWeight is available only for Linear distribution.";
1414 command->CommandFailed(ed);
1415 return;
1416 }
1417 fParticleGun->GetEneDist()->ApplyEnergyWeight(applyEnergyWeightCmd1->GetNewBoolValue(newValues));
1418 }
1419 else if(command == histnameCmd1)
1420 {
1421 histtype = newValues;
1422 }
1423 else if(command == histfileCmd1)
1424 {
1425 histtype = "arb";
1426 CHECKPG();
1427 fParticleGun->GetEneDist()->ArbEnergyHistoFile(newValues);
1428 }
1429 else if(command == histpointCmd1)
1430 {
1431 CHECKPG();
1432 if(histtype == "biasx")
1433 fParticleGun->GetBiasRndm()->SetXBias(histpointCmd1->GetNew3VectorValue(newValues));
1434 if(histtype == "biasy")
1435 fParticleGun->GetBiasRndm()->SetYBias(histpointCmd1->GetNew3VectorValue(newValues));
1436 if(histtype == "biasz")
1437 fParticleGun->GetBiasRndm()->SetZBias(histpointCmd1->GetNew3VectorValue(newValues));
1438 if(histtype == "biast")
1439 fParticleGun->GetBiasRndm()->SetThetaBias(histpointCmd1->GetNew3VectorValue(newValues));
1440 if(histtype == "biasp")
1441 fParticleGun->GetBiasRndm()->SetPhiBias(histpointCmd1->GetNew3VectorValue(newValues));
1442 if(histtype == "biaspt")
1443 fParticleGun->GetBiasRndm()->SetPosThetaBias(histpointCmd1->GetNew3VectorValue(newValues));
1444 if(histtype == "biaspp")
1445 fParticleGun->GetBiasRndm()->SetPosPhiBias(histpointCmd1->GetNew3VectorValue(newValues));
1446 if(histtype == "biase")
1447 fParticleGun->GetBiasRndm()->SetEnergyBias(histpointCmd1->GetNew3VectorValue(newValues));
1448 if(histtype == "theta")
1449 fParticleGun->GetAngDist()->UserDefAngTheta(histpointCmd1->GetNew3VectorValue(newValues));
1450 if(histtype == "phi")
1451 fParticleGun->GetAngDist()->UserDefAngPhi(histpointCmd1->GetNew3VectorValue(newValues));
1452 if(histtype == "energy")
1453 fParticleGun->GetEneDist()->UserEnergyHisto(histpointCmd1->GetNew3VectorValue(newValues));
1454 if(histtype == "arb")
1455 fParticleGun->GetEneDist()->ArbEnergyHisto(histpointCmd1->GetNew3VectorValue(newValues));
1456 if(histtype == "epn")
1457 fParticleGun->GetEneDist()->EpnEnergyHisto(histpointCmd1->GetNew3VectorValue(newValues));
1458 }
1459 else if(command == resethistCmd1)
1460 {
1461 CHECKPG();
1462 if(newValues == "theta" || newValues == "phi")
1463 {
1464 fParticleGun->GetAngDist()->ReSetHist(newValues);
1465 }
1466 else if (newValues == "energy" || newValues == "arb" || newValues == "epn")
1467 {
1468 fParticleGun->GetEneDist()->ReSetHist(newValues);
1469 }
1470 else
1471 {
1472 fParticleGun->GetBiasRndm()->ReSetHist(newValues);
1473 }
1474 }
1475 else if(command == arbintCmd1)
1476 {
1477 CHECKPG(); fParticleGun->GetEneDist()->ArbInterpolate(newValues);
1478 }
1479 else
1480 {
1481 G4cout << "Error entering command" << G4endl;
1482 }
1483}
1484
1486{
1487 G4String cv;
1488
1489 // if( command==directionCmd )
1490 // { cv = directionCmd->ConvertToString(fParticleGun->GetParticleMomentumDirection()); }
1491 // else if( command==energyCmd )
1492 // { cv = energyCmd->ConvertToString(fParticleGun->GetParticleEnergy(),"GeV"); }
1493 // else if( command==positionCmd )
1494 // { cv = positionCmd->ConvertToString(fParticleGun->GetParticlePosition(),"cm"); }
1495 // else if( command==timeCmd )
1496 // { cv = timeCmd->ConvertToString(fParticleGun->GetParticleTime(),"ns"); }
1497 // else if( command==polCmd )
1498 // { cv = polCmd->ConvertToString(fParticleGun->GetParticlePolarization()); }
1499 // else if( command==numberCmd )
1500 // { cv = numberCmd->ConvertToString(fParticleGun->GetNumberOfParticles()); }
1501
1502 cv = "Not implemented yet";
1503
1504 return cv;
1505}
1506
1507void G4GeneralParticleSourceMessenger::IonCommand(G4String newValues)
1508{
1509 if (fShootIon)
1510 {
1511 G4Tokenizer next( newValues );
1512 // check argument
1513 fAtomicNumber = StoI(next());
1514 fAtomicMass = StoI(next());
1515 G4String sQ = next();
1516 if (sQ.isNull())
1517 {
1518 fIonCharge = fAtomicNumber;
1519 }
1520 else
1521 {
1522 fIonCharge = StoI(sQ);
1523 sQ = next();
1524 if (sQ.isNull())
1525 {
1526 fIonExciteEnergy = 0.0;
1527 }
1528 else
1529 {
1530 fIonExciteEnergy = StoD(sQ) * keV;
1531 }
1532 }
1534 ->GetIon(fAtomicNumber, fAtomicMass, fIonExciteEnergy);
1535 if (ion==nullptr)
1536 {
1538 ed << "Ion with Z=" << fAtomicNumber;
1539 ed << " A=" << fAtomicMass << " is not defined";
1540 ionCmd->CommandFailed(ed);
1541 }
1542 else
1543 {
1544 fParticleGun->SetParticleDefinition(ion);
1545 fParticleGun->SetParticleCharge(fIonCharge*eplus);
1546 }
1547 }
1548 else
1549 {
1551 ed << "Set /gps/particle to ion before using /gps/ion command";
1552 ionCmd->CommandFailed(ed);
1553 }
1554}
1555
1556void G4GeneralParticleSourceMessenger::IonLvlCommand(G4String newValues)
1557{
1558 if (fShootIon)
1559 {
1560 G4Tokenizer next(newValues);
1561 // check argument
1562 fAtomicNumberL = StoI(next());
1563 fAtomicMassL = StoI(next());
1564 G4String sQ = next();
1565 if (sQ.isNull())
1566 {
1567 fIonChargeL = fAtomicNumberL;
1568 }
1569 else
1570 {
1571 fIonChargeL = StoI(sQ);
1572 sQ = next();
1573 if (sQ.isNull())
1574 {
1575 fIonEnergyLevel = 0;
1576 }
1577 else
1578 {
1579 fIonEnergyLevel = StoI(sQ);
1580 }
1581 }
1582
1584 ->GetIon(fAtomicNumberL, fAtomicMassL, fIonEnergyLevel);
1585 if (ion == nullptr)
1586 {
1588 ed << "Ion with Z=" << fAtomicNumberL;
1589 ed << " A=" << fAtomicMassL << " is not defined";
1590 ionLvlCmd->CommandFailed(ed);
1591 }
1592 else
1593 {
1594 fParticleGun->SetParticleDefinition(ion);
1595 fParticleGun->SetParticleCharge(fIonChargeL*eplus);
1596 }
1597
1598 }
1599 else
1600 {
1602 ed << "Set /gps/particle to ion before using /gps/ionLvl command";
1603 ionLvlCmd->CommandFailed(ed);
1604 }
1605}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
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)
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
G4bool isNull() 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:136
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:146
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:156
void CommandFailed(G4int errCode, G4ExceptionDescription &ed)
Definition: G4UIcommand.hh:179
void SetRange(const char *rs)
Definition: G4UIcommand.hh:120
const G4String & GetCommandName() const
Definition: G4UIcommand.hh:137
G4int StoI(G4String s)
G4bool commandsShouldBeInMaster
G4double StoD(G4String s)
void SetDefaultValue(const char *theDefaultValue)