G4GeneralParticleSourceMessenger Class Reference

#include <G4GeneralParticleSourceMessenger.hh>

Inheritance diagram for G4GeneralParticleSourceMessenger:

G4UImessenger

Public Member Functions

 G4GeneralParticleSourceMessenger (G4GeneralParticleSource *)
 ~G4GeneralParticleSourceMessenger ()
void SetParticleGun (G4SingleParticleSource *fpg)
void SetNewValue (G4UIcommand *command, G4String newValues)
G4String GetCurrentValue (G4UIcommand *command)

Detailed Description

Definition at line 101 of file G4GeneralParticleSourceMessenger.hh.


Constructor & Destructor Documentation

G4GeneralParticleSourceMessenger::G4GeneralParticleSourceMessenger ( G4GeneralParticleSource  ) 

Definition at line 76 of file G4GeneralParticleSourceMessenger.cc.

References G4ParticleTable::GetParticleTable(), G4INCL::Math::pi, and G4UIparameter::SetDefaultValue().

00077     : fGPS(fPtclGun),fShootIon(false)
00078 {
00079   particleTable = G4ParticleTable::GetParticleTable();
00080   histtype = "biasx";
00081 
00082   gpsDirectory = new G4UIdirectory("/gps/");
00083   gpsDirectory->SetGuidance("General Paricle Source control commands.");
00084   //  gpsDirectory->SetGuidance(" The first 9 commands are the same as in G4ParticleGun ");
00085 
00086   // now the commands for mutiple sources
00087   sourceDirectory = new G4UIdirectory("/gps/source/");
00088   sourceDirectory->SetGuidance("Multiple source control sub-directory");
00089 
00090   addsourceCmd = new G4UIcmdWithADouble("/gps/source/add",this);
00091   addsourceCmd->SetGuidance("add a new source defintion to the particle gun with the specified intensity");
00092   addsourceCmd->SetParameterName("addsource",false,false);
00093   addsourceCmd->SetRange("addsource > 0.");
00094 
00095   listsourceCmd = new G4UIcmdWithoutParameter("/gps/source/list",this);
00096   listsourceCmd->SetGuidance("List the defined particle sources");
00097 
00098   clearsourceCmd = new G4UIcmdWithoutParameter("/gps/source/clear",this);
00099   clearsourceCmd->SetGuidance("Remove all the defined particle sources");
00100 
00101   getsourceCmd = new G4UIcmdWithoutParameter("/gps/source/show",this);
00102   getsourceCmd->SetGuidance("Show the current source index and intensity");
00103 
00104   setsourceCmd = new G4UIcmdWithAnInteger("/gps/source/set",this);
00105   setsourceCmd->SetGuidance("set the indexed source as the current one");
00106   setsourceCmd->SetGuidance("  so one can change its source definition");
00107   setsourceCmd->SetParameterName("setsource",false,false);
00108   setsourceCmd->SetRange("setsource >= 0");
00109 
00110   deletesourceCmd = new G4UIcmdWithAnInteger("/gps/source/delete",this);
00111   deletesourceCmd->SetGuidance("delete the indexed source from the list");
00112   deletesourceCmd->SetParameterName("deletesource",false,false);
00113   deletesourceCmd->SetRange("deletesource > 0");
00114 
00115   setintensityCmd = new G4UIcmdWithADouble("/gps/source/intensity",this);
00116   setintensityCmd->SetGuidance("reset the current source to the specified intensity");
00117   setintensityCmd->SetParameterName("setintensity",false,false);
00118   setintensityCmd->SetRange("setintensity > 0."); 
00119 
00120   multiplevertexCmd = new G4UIcmdWithABool("/gps/source/multiplevertex",this);
00121   multiplevertexCmd->SetGuidance("true for simulaneous generation mutiple vertex");
00122   multiplevertexCmd->SetGuidance("Default is false");
00123   multiplevertexCmd->SetParameterName("multiplevertex",true);
00124   multiplevertexCmd->SetDefaultValue(false);
00125 
00126   flatsamplingCmd = new G4UIcmdWithABool("/gps/source/flatsampling",this);
00127   flatsamplingCmd->SetGuidance("true for appling flat (biased) sampling among the sources");
00128   flatsamplingCmd->SetGuidance("Default is false");
00129   flatsamplingCmd->SetParameterName("flatsampling",true);
00130   flatsamplingCmd->SetDefaultValue(false);
00131 
00132   // below we reproduce commands awailable in G4Particle Gun
00133 
00134   listCmd = new G4UIcmdWithoutParameter("/gps/List",this);
00135   listCmd->SetGuidance("List available particles.");
00136   listCmd->SetGuidance(" Invoke G4ParticleTable.");
00137 
00138   particleCmd = new G4UIcmdWithAString("/gps/particle",this);
00139   particleCmd->SetGuidance("Set particle to be generated.");
00140   particleCmd->SetGuidance(" (geantino is default)");
00141   particleCmd->SetGuidance(" (ion can be specified for shooting ions)");
00142   particleCmd->SetParameterName("particleName",true);
00143   particleCmd->SetDefaultValue("geantino");
00144   G4String candidateList; 
00145   G4int nPtcl = particleTable->entries();
00146   for(G4int i=0;i<nPtcl;i++)
00147   {
00148     candidateList += particleTable->GetParticleName(i);
00149     candidateList += " ";
00150   }
00151   candidateList += "ion ";
00152   particleCmd->SetCandidates(candidateList);
00153 
00154   directionCmd = new G4UIcmdWith3Vector("/gps/direction",this);
00155   directionCmd->SetGuidance("Set momentum direction.");
00156   directionCmd->SetGuidance("Direction needs not to be a unit vector.");
00157   directionCmd->SetParameterName("Px","Py","Pz",false,false); 
00158   directionCmd->SetRange("Px != 0 || Py != 0 || Pz != 0");
00159   
00160   energyCmd = new G4UIcmdWithADoubleAndUnit("/gps/energy",this);
00161   energyCmd->SetGuidance("Set kinetic energy.");
00162   energyCmd->SetParameterName("Energy",false,false);
00163   energyCmd->SetDefaultUnit("GeV");
00164   //energyCmd->SetUnitCategory("Energy");
00165   //energyCmd->SetUnitCandidates("eV keV MeV GeV TeV");
00166 
00167   positionCmd = new G4UIcmdWith3VectorAndUnit("/gps/position",this);
00168   positionCmd->SetGuidance("Set starting position of the particle.");
00169   positionCmd->SetParameterName("X","Y","Z",false,false);
00170   positionCmd->SetDefaultUnit("cm");
00171   //positionCmd->SetUnitCategory("Length");
00172   //positionCmd->SetUnitCandidates("microm mm cm m km");
00173 
00174   ionCmd = new G4UIcommand("/gps/ion",this);
00175   ionCmd->SetGuidance("Set properties of ion to be generated.");
00176   ionCmd->SetGuidance("[usage] /gps/ion Z A Q E");
00177   ionCmd->SetGuidance("        Z:(int) AtomicNumber");
00178   ionCmd->SetGuidance("        A:(int) AtomicMass");
00179   ionCmd->SetGuidance("        Q:(int) Charge of Ion (in unit of e)");
00180   ionCmd->SetGuidance("        E:(double) Excitation energy (in keV)");
00181   
00182   G4UIparameter* param;
00183   param = new G4UIparameter("Z",'i',false);
00184   param->SetDefaultValue("1");
00185   ionCmd->SetParameter(param);
00186   param = new G4UIparameter("A",'i',false);
00187   param->SetDefaultValue("1");
00188   ionCmd->SetParameter(param);
00189   param = new G4UIparameter("Q",'i',true);
00190   param->SetDefaultValue("0");
00191   ionCmd->SetParameter(param);
00192   param = new G4UIparameter("E",'d',true);
00193   param->SetDefaultValue("0.0");
00194   ionCmd->SetParameter(param);
00195 
00196 
00197   timeCmd = new G4UIcmdWithADoubleAndUnit("/gps/time",this);
00198   timeCmd->SetGuidance("Set initial time of the particle.");
00199   timeCmd->SetParameterName("t0",false,false);
00200   timeCmd->SetDefaultUnit("ns");
00201   //timeCmd->SetUnitCategory("Time");
00202   //timeCmd->SetUnitCandidates("ns ms s");
00203   
00204   polCmd = new G4UIcmdWith3Vector("/gps/polarization",this);
00205   polCmd->SetGuidance("Set polarization.");
00206   polCmd->SetParameterName("Px","Py","Pz",false,false); 
00207   polCmd->SetRange("Px>=-1.&&Px<=1.&&Py>=-1.&&Py<=1.&&Pz>=-1.&&Pz<=1.");
00208 
00209   numberCmd = new G4UIcmdWithAnInteger("/gps/number",this);
00210   numberCmd->SetGuidance("Set number of particles to be generated per vertex.");
00211   numberCmd->SetParameterName("N",false,false);
00212   numberCmd->SetRange("N>0");
00213 
00214   // verbosity
00215   verbosityCmd = new G4UIcmdWithAnInteger("/gps/verbose",this);
00216   verbosityCmd->SetGuidance("Set Verbose level for GPS");
00217   verbosityCmd->SetGuidance(" 0 : Silent");
00218   verbosityCmd->SetGuidance(" 1 : Limited information");
00219   verbosityCmd->SetGuidance(" 2 : Detailed information");
00220   verbosityCmd->SetParameterName("level",false);
00221   verbosityCmd->SetRange("level>=0 && level <=2");
00222 
00223   // now extended commands
00224   // Positional ones:
00225   positionDirectory = new G4UIdirectory("/gps/pos/");
00226   positionDirectory->SetGuidance("Positional commands sub-directory");
00227 
00228   typeCmd1 = new G4UIcmdWithAString("/gps/pos/type",this);
00229   typeCmd1->SetGuidance("Sets source distribution type.");
00230   typeCmd1->SetGuidance("Either Point, Beam, Plane, Surface or Volume");
00231   typeCmd1->SetParameterName("DisType",false,false);
00232   typeCmd1->SetDefaultValue("Point");
00233   typeCmd1->SetCandidates("Point Beam Plane Surface Volume");
00234 
00235   shapeCmd1 = new G4UIcmdWithAString("/gps/pos/shape",this);
00236   shapeCmd1->SetGuidance("Sets source shape for Plan, Surface or Volume type source.");
00237   shapeCmd1->SetParameterName("Shape",false,false);
00238   shapeCmd1->SetDefaultValue("NULL");
00239   shapeCmd1->SetCandidates("Circle Annulus Ellipse Square Rectangle Sphere Ellipsoid Cylinder Para");
00240 
00241   centreCmd1 = new G4UIcmdWith3VectorAndUnit("/gps/pos/centre",this);
00242   centreCmd1->SetGuidance("Set centre coordinates of source.");
00243   centreCmd1->SetGuidance("   same effect as the /gps/position command");
00244   centreCmd1->SetParameterName("X","Y","Z",false,false);
00245   centreCmd1->SetDefaultUnit("cm");
00246   //  centreCmd1->SetUnitCandidates("micron mm cm m km");
00247 
00248   posrot1Cmd1 = new G4UIcmdWith3Vector("/gps/pos/rot1",this);
00249   posrot1Cmd1->SetGuidance("Set the 1st vector defining the rotation matrix'.");
00250   posrot1Cmd1->SetGuidance("It does not need to be a unit vector.");
00251   posrot1Cmd1->SetParameterName("R1x","R1y","R1z",false,false); 
00252   posrot1Cmd1->SetRange("R1x != 0 || R1y != 0 || R1z != 0");
00253 
00254   posrot2Cmd1 = new G4UIcmdWith3Vector("/gps/pos/rot2",this);
00255   posrot2Cmd1->SetGuidance("Set the 2nd vector defining the rotation matrix'.");
00256   posrot2Cmd1->SetGuidance("It does not need to be a unit vector.");
00257   posrot2Cmd1->SetParameterName("R2x","R2y","R2z",false,false); 
00258   posrot2Cmd1->SetRange("R2x != 0 || R2y != 0 || R2z != 0");
00259 
00260   halfxCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/halfx",this);
00261   halfxCmd1->SetGuidance("Set x half length of source.");
00262   halfxCmd1->SetParameterName("Halfx",false,false);
00263   halfxCmd1->SetDefaultUnit("cm");
00264   //  halfxCmd1->SetUnitCandidates("micron mm cm m km");
00265 
00266   halfyCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/halfy",this);
00267   halfyCmd1->SetGuidance("Set y half length of source.");
00268   halfyCmd1->SetParameterName("Halfy",false,false);
00269   halfyCmd1->SetDefaultUnit("cm");
00270   //  halfyCmd1->SetUnitCandidates("micron mm cm m km");
00271 
00272   halfzCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/halfz",this);
00273   halfzCmd1->SetGuidance("Set z half length of source.");
00274   halfzCmd1->SetParameterName("Halfz",false,false);
00275   halfzCmd1->SetDefaultUnit("cm");
00276   //  halfzCmd1->SetUnitCandidates("micron mm cm m km");
00277 
00278   radiusCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/radius",this);
00279   radiusCmd1->SetGuidance("Set radius of source.");
00280   radiusCmd1->SetParameterName("Radius",false,false);
00281   radiusCmd1->SetDefaultUnit("cm");
00282   //  radiusCmd1->SetUnitCandidates("micron mm cm m km");
00283 
00284   radius0Cmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/inner_radius",this);
00285   radius0Cmd1->SetGuidance("Set inner radius of source when required.");
00286   radius0Cmd1->SetParameterName("Radius0",false,false);
00287   radius0Cmd1->SetDefaultUnit("cm");
00288   //  radius0Cmd1->SetUnitCandidates("micron mm cm m km");
00289 
00290   possigmarCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/sigma_r",this);
00291   possigmarCmd1->SetGuidance("Set standard deviation in radial of the beam positional profile");
00292   possigmarCmd1->SetGuidance(" applicable to Beam type source only");
00293   possigmarCmd1->SetParameterName("Sigmar",false,false);
00294   possigmarCmd1->SetDefaultUnit("cm");
00295   //  possigmarCmd1->SetUnitCandidates("micron mm cm m km");
00296 
00297   possigmaxCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/sigma_x",this);
00298   possigmaxCmd1->SetGuidance("Set standard deviation of beam positional profile in x-dir");
00299   possigmaxCmd1->SetGuidance(" applicable to Beam type source only");
00300   possigmaxCmd1->SetParameterName("Sigmax",false,false);
00301   possigmaxCmd1->SetDefaultUnit("cm");
00302   //  possigmaxCmd1->SetUnitCandidates("micron mm cm m km");
00303 
00304   possigmayCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/sigma_y",this);
00305   possigmayCmd1->SetGuidance("Set standard deviation of beam positional profile in y-dir");
00306   possigmayCmd1->SetGuidance(" applicable to Beam type source only");
00307   possigmayCmd1->SetParameterName("Sigmay",false,false);
00308   possigmayCmd1->SetDefaultUnit("cm");
00309   //  possigmayCmd1->SetUnitCandidates("micron mm cm m km");
00310 
00311   paralpCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/paralp",this);
00312   paralpCmd1->SetGuidance("Angle from y-axis of y' in Para");
00313   paralpCmd1->SetParameterName("paralp",false,false);
00314   paralpCmd1->SetDefaultUnit("rad");
00315   //  paralpCmd1->SetUnitCandidates("rad deg");
00316 
00317   partheCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/parthe",this);
00318   partheCmd1->SetGuidance("Polar angle through centres of z faces");
00319   partheCmd1->SetParameterName("parthe",false,false);
00320   partheCmd1->SetDefaultUnit("rad");
00321   //  partheCmd1->SetUnitCandidates("rad deg");
00322 
00323   parphiCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/pos/parphi",this);
00324   parphiCmd1->SetGuidance("Azimuth angle through centres of z faces");
00325   parphiCmd1->SetParameterName("parphi",false,false);
00326   parphiCmd1->SetDefaultUnit("rad");
00327   //  parphiCmd1->SetUnitCandidates("rad deg");
00328 
00329   confineCmd1 = new G4UIcmdWithAString("/gps/pos/confine",this);
00330   confineCmd1->SetGuidance("Confine source to volume (NULL to unset).");
00331   confineCmd1->SetGuidance("usage: confine VolName");
00332   confineCmd1->SetParameterName("VolName",false,false);
00333   confineCmd1->SetDefaultValue("NULL");
00334 
00335   // old implementations
00336   typeCmd = new G4UIcmdWithAString("/gps/type",this);
00337   typeCmd->SetGuidance("Sets source distribution type. (obsolete!)");
00338   typeCmd->SetGuidance("Either Point, Beam, Plane, Surface or Volume");
00339   typeCmd->SetParameterName("DisType",false,false);
00340   typeCmd->SetDefaultValue("Point");
00341   typeCmd->SetCandidates("Point Beam Plane Surface Volume");
00342 
00343   shapeCmd = new G4UIcmdWithAString("/gps/shape",this);
00344   shapeCmd->SetGuidance("Sets source shape type.(obsolete!)");
00345   shapeCmd->SetParameterName("Shape",false,false);
00346   shapeCmd->SetDefaultValue("NULL");
00347   shapeCmd->SetCandidates("Circle Annulus Ellipse Square Rectangle Sphere Ellipsoid Cylinder Para");
00348 
00349   centreCmd = new G4UIcmdWith3VectorAndUnit("/gps/centre",this);
00350   centreCmd->SetGuidance("Set centre coordinates of source.(obsolete!)");
00351   centreCmd->SetParameterName("X","Y","Z",false,false);
00352   centreCmd->SetDefaultUnit("cm");
00353   //  centreCmd->SetUnitCandidates("micron mm cm m km");
00354 
00355   posrot1Cmd = new G4UIcmdWith3Vector("/gps/posrot1",this);
00356   posrot1Cmd->SetGuidance("Set rotation matrix of x'.(obsolete!)");
00357   posrot1Cmd->SetGuidance("Posrot1 does not need to be a unit vector.");
00358   posrot1Cmd->SetParameterName("R1x","R1y","R1z",false,false); 
00359   posrot1Cmd->SetRange("R1x != 0 || R1y != 0 || R1z != 0");
00360 
00361   posrot2Cmd = new G4UIcmdWith3Vector("/gps/posrot2",this);
00362   posrot2Cmd->SetGuidance("Set rotation matrix of y'.(obsolete!)");
00363   posrot2Cmd->SetGuidance("Posrot2 does not need to be a unit vector.");
00364   posrot2Cmd->SetParameterName("R2x","R2y","R2z",false,false); 
00365   posrot2Cmd->SetRange("R2x != 0 || R2y != 0 || R2z != 0");
00366 
00367   halfxCmd = new G4UIcmdWithADoubleAndUnit("/gps/halfx",this);
00368   halfxCmd->SetGuidance("Set x half length of source.(obsolete!)");
00369   halfxCmd->SetParameterName("Halfx",false,false);
00370   halfxCmd->SetDefaultUnit("cm");
00371   //  halfxCmd->SetUnitCandidates("micron mm cm m km");
00372 
00373   halfyCmd = new G4UIcmdWithADoubleAndUnit("/gps/halfy",this);
00374   halfyCmd->SetGuidance("Set y half length of source.(obsolete!)");
00375   halfyCmd->SetParameterName("Halfy",false,false);
00376   halfyCmd->SetDefaultUnit("cm");
00377   //  halfyCmd->SetUnitCandidates("micron mm cm m km");
00378 
00379   halfzCmd = new G4UIcmdWithADoubleAndUnit("/gps/halfz",this);
00380   halfzCmd->SetGuidance("Set z half length of source.(obsolete!)");
00381   halfzCmd->SetParameterName("Halfz",false,false);
00382   halfzCmd->SetDefaultUnit("cm");
00383   //  halfzCmd->SetUnitCandidates("micron mm cm m km");
00384 
00385   radiusCmd = new G4UIcmdWithADoubleAndUnit("/gps/radius",this);
00386   radiusCmd->SetGuidance("Set radius of source.(obsolete!)");
00387   radiusCmd->SetParameterName("Radius",false,false);
00388   radiusCmd->SetDefaultUnit("cm");
00389   //  radiusCmd->SetUnitCandidates("micron mm cm m km");
00390 
00391   radius0Cmd = new G4UIcmdWithADoubleAndUnit("/gps/radius0",this);
00392   radius0Cmd->SetGuidance("Set inner radius of source.(obsolete!)");
00393   radius0Cmd->SetParameterName("Radius0",false,false);
00394   radius0Cmd->SetDefaultUnit("cm");
00395   //  radius0Cmd->SetUnitCandidates("micron mm cm m km");
00396 
00397   possigmarCmd = new G4UIcmdWithADoubleAndUnit("/gps/sigmaposr",this);
00398   possigmarCmd->SetGuidance("Set standard deviation of beam position in radial(obsolete!)");
00399   possigmarCmd->SetParameterName("Sigmar",false,false);
00400   possigmarCmd->SetDefaultUnit("cm");
00401   // possigmarCmd->SetUnitCandidates("micron mm cm m km");
00402 
00403   possigmaxCmd = new G4UIcmdWithADoubleAndUnit("/gps/sigmaposx",this);
00404   possigmaxCmd->SetGuidance("Set standard deviation of beam position in x-dir(obsolete!)");
00405   possigmaxCmd->SetParameterName("Sigmax",false,false);
00406   possigmaxCmd->SetDefaultUnit("cm");
00407   //  possigmaxCmd->SetUnitCandidates("micron mm cm m km");
00408 
00409   possigmayCmd = new G4UIcmdWithADoubleAndUnit("/gps/sigmaposy",this);
00410   possigmayCmd->SetGuidance("Set standard deviation of beam position in y-dir(obsolete!)");
00411   possigmayCmd->SetParameterName("Sigmay",false,false);
00412   possigmayCmd->SetDefaultUnit("cm");
00413   //  possigmayCmd->SetUnitCandidates("micron mm cm m km");
00414 
00415   paralpCmd = new G4UIcmdWithADoubleAndUnit("/gps/paralp",this);
00416   paralpCmd->SetGuidance("Angle from y-axis of y' in Para(obsolete!)");
00417   paralpCmd->SetParameterName("paralp",false,false);
00418   paralpCmd->SetDefaultUnit("rad");
00419   //  paralpCmd->SetUnitCandidates("rad deg");
00420 
00421   partheCmd = new G4UIcmdWithADoubleAndUnit("/gps/parthe",this);
00422   partheCmd->SetGuidance("Polar angle through centres of z faces(obsolete!)");
00423   partheCmd->SetParameterName("parthe",false,false);
00424   partheCmd->SetDefaultUnit("rad");
00425   //  partheCmd->SetUnitCandidates("rad deg");
00426 
00427   parphiCmd = new G4UIcmdWithADoubleAndUnit("/gps/parphi",this);
00428   parphiCmd->SetGuidance("Azimuth angle through centres of z faces(obsolete!)");
00429   parphiCmd->SetParameterName("parphi",false,false);
00430   parphiCmd->SetDefaultUnit("rad");
00431   //  parphiCmd->SetUnitCandidates("rad deg");
00432 
00433   confineCmd = new G4UIcmdWithAString("/gps/confine",this);
00434   confineCmd->SetGuidance("Confine source to volume (NULL to unset)(obsolete!) .");
00435   confineCmd->SetGuidance("usage: confine VolName");
00436   confineCmd->SetParameterName("VolName",false,false);
00437   confineCmd->SetDefaultValue("NULL");
00438 
00439   // Angular distribution commands
00440   angularDirectory = new G4UIdirectory("/gps/ang/");
00441   angularDirectory->SetGuidance("Angular commands sub-directory");
00442 
00443   angtypeCmd1 = new G4UIcmdWithAString("/gps/ang/type",this);
00444   angtypeCmd1->SetGuidance("Sets angular source distribution type");
00445   angtypeCmd1->SetGuidance("Possible variables are: iso, cos, planar, beam1d, beam2d, focused or user");
00446   angtypeCmd1->SetParameterName("AngDis",false,false);
00447   angtypeCmd1->SetDefaultValue("iso");
00448   angtypeCmd1->SetCandidates("iso cos planar beam1d beam2d focused user");
00449 
00450   angrot1Cmd1 = new G4UIcmdWith3Vector("/gps/ang/rot1",this);
00451   angrot1Cmd1->SetGuidance("Sets the 1st vector for angular distribution rotation matrix");
00452   angrot1Cmd1->SetGuidance("Need not be a unit vector");
00453   angrot1Cmd1->SetParameterName("AR1x","AR1y","AR1z",false,false);
00454   angrot1Cmd1->SetRange("AR1x != 0 || AR1y != 0 || AR1z != 0");
00455 
00456   angrot2Cmd1 = new G4UIcmdWith3Vector("/gps/ang/rot2",this);
00457   angrot2Cmd1->SetGuidance("Sets the 2nd vector for angular distribution rotation matrix");
00458   angrot2Cmd1->SetGuidance("Need not be a unit vector");
00459   angrot2Cmd1->SetParameterName("AR2x","AR2y","AR2z",false,false);
00460   angrot2Cmd1->SetRange("AR2x != 0 || AR2y != 0 || AR2z != 0");
00461 
00462   minthetaCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ang/mintheta",this);
00463   minthetaCmd1->SetGuidance("Set minimum theta");
00464   minthetaCmd1->SetParameterName("MinTheta",false,false);
00465   minthetaCmd1->SetDefaultValue(0.);
00466   minthetaCmd1->SetDefaultUnit("rad");
00467   //  minthetaCmd1->SetUnitCandidates("rad deg");
00468 
00469   maxthetaCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ang/maxtheta",this);
00470   maxthetaCmd1->SetGuidance("Set maximum theta");
00471   maxthetaCmd1->SetParameterName("MaxTheta",false,false);
00472   maxthetaCmd1->SetDefaultValue(pi);
00473   maxthetaCmd1->SetDefaultUnit("rad");
00474   //  maxthetaCmd1->SetUnitCandidates("rad deg");
00475 
00476   minphiCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ang/minphi",this);
00477   minphiCmd1->SetGuidance("Set minimum phi");
00478   minphiCmd1->SetParameterName("MinPhi",false,false);
00479   minphiCmd1->SetDefaultUnit("rad");
00480   //  minphiCmd1->SetUnitCandidates("rad deg");
00481 
00482   maxphiCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ang/maxphi",this);
00483   maxphiCmd1->SetGuidance("Set maximum phi");
00484   maxphiCmd1->SetParameterName("MaxPhi",false,false);
00485   maxphiCmd1->SetDefaultValue(pi);
00486   maxphiCmd1->SetDefaultUnit("rad");
00487   //  maxphiCmd1->SetUnitCandidates("rad deg");
00488 
00489   angsigmarCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ang/sigma_r",this);
00490   angsigmarCmd1->SetGuidance("Set standard deviation in direction for 1D beam.");
00491   angsigmarCmd1->SetParameterName("Sigmara",false,false);
00492   angsigmarCmd1->SetDefaultUnit("rad");
00493   //  angsigmarCmd1->SetUnitCandidates("rad deg");
00494   
00495   angsigmaxCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ang/sigma_x",this);
00496   angsigmaxCmd1->SetGuidance("Set standard deviation in direction in x-direc. for 2D beam");
00497   angsigmaxCmd1->SetParameterName("Sigmaxa",false,false);
00498   angsigmaxCmd1->SetDefaultUnit("rad");
00499   //  angsigmaxCmd1->SetUnitCandidates("rad deg");
00500 
00501   angsigmayCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ang/sigma_y",this);
00502   angsigmayCmd1->SetGuidance("Set standard deviation in direction in y-direc. for 2D beam");
00503   angsigmayCmd1->SetParameterName("Sigmaya",false,false);
00504   angsigmayCmd1->SetDefaultUnit("rad");
00505   //  angsigmayCmd1->SetUnitCandidates("rad deg");
00506 
00507   angfocusCmd = new G4UIcmdWith3VectorAndUnit("/gps/ang/focuspoint",this);
00508   angfocusCmd->SetGuidance("Set the focusing point for the beam");
00509   angfocusCmd->SetParameterName("x","y","z",false,false);
00510   angfocusCmd->SetDefaultUnit("cm");
00511   //  angfocusCmd->SetUnitCandidates("micron mm cm m km");
00512 
00513   useuserangaxisCmd1 = new G4UIcmdWithABool("/gps/ang/user_coor",this);
00514   useuserangaxisCmd1->SetGuidance("true for using user defined angular co-ordinates");
00515   useuserangaxisCmd1->SetGuidance("Default is false");
00516   useuserangaxisCmd1->SetParameterName("useuserangaxis",true);
00517   useuserangaxisCmd1->SetDefaultValue(false);
00518 
00519   surfnormCmd1 = new G4UIcmdWithABool("/gps/ang/surfnorm",this);
00520   surfnormCmd1->SetGuidance("Makes a user-defined distribution with respect to surface normals rather than x,y,z axes.");
00521   surfnormCmd1->SetGuidance("Default is false");
00522   surfnormCmd1->SetParameterName("surfnorm",true);
00523   surfnormCmd1->SetDefaultValue(false);
00524 
00525   // old ones
00526   angtypeCmd = new G4UIcmdWithAString("/gps/angtype",this);
00527   angtypeCmd->SetGuidance("Sets angular source distribution type (obsolete!)");
00528   angtypeCmd->SetGuidance("Possible variables are: iso, cos planar beam1d beam2d or user");
00529   angtypeCmd->SetParameterName("AngDis",false,false);
00530   angtypeCmd->SetDefaultValue("iso");
00531   angtypeCmd->SetCandidates("iso cos planar beam1d beam2d user");
00532 
00533   angrot1Cmd = new G4UIcmdWith3Vector("/gps/angrot1",this);
00534   angrot1Cmd->SetGuidance("Sets the x' vector for angular distribution(obsolete!) ");
00535   angrot1Cmd->SetGuidance("Need not be a unit vector");
00536   angrot1Cmd->SetParameterName("AR1x","AR1y","AR1z",false,false);
00537   angrot1Cmd->SetRange("AR1x != 0 || AR1y != 0 || AR1z != 0");
00538 
00539   angrot2Cmd = new G4UIcmdWith3Vector("/gps/angrot2",this);
00540   angrot2Cmd->SetGuidance("Sets the y' vector for angular distribution (obsolete!)");
00541   angrot2Cmd->SetGuidance("Need not be a unit vector");
00542   angrot2Cmd->SetParameterName("AR2x","AR2y","AR2z",false,false);
00543   angrot2Cmd->SetRange("AR2x != 0 || AR2y != 0 || AR2z != 0");
00544 
00545   minthetaCmd = new G4UIcmdWithADoubleAndUnit("/gps/mintheta",this);
00546   minthetaCmd->SetGuidance("Set minimum theta (obsolete!)");
00547   minthetaCmd->SetParameterName("MinTheta",false,false);
00548   minthetaCmd->SetDefaultUnit("rad");
00549   //  minthetaCmd->SetUnitCandidates("rad deg");
00550 
00551   maxthetaCmd = new G4UIcmdWithADoubleAndUnit("/gps/maxtheta",this);
00552   maxthetaCmd->SetGuidance("Set maximum theta (obsolete!)");
00553   maxthetaCmd->SetParameterName("MaxTheta",false,false);
00554   maxthetaCmd->SetDefaultValue(3.1416);
00555   maxthetaCmd->SetDefaultUnit("rad");
00556   //  maxthetaCmd->SetUnitCandidates("rad deg");
00557 
00558   minphiCmd = new G4UIcmdWithADoubleAndUnit("/gps/minphi",this);
00559   minphiCmd->SetGuidance("Set minimum phi (obsolete!)");
00560   minphiCmd->SetParameterName("MinPhi",false,false);
00561   minphiCmd->SetDefaultUnit("rad");
00562   //  minphiCmd->SetUnitCandidates("rad deg");
00563 
00564   maxphiCmd = new G4UIcmdWithADoubleAndUnit("/gps/maxphi",this);
00565   maxphiCmd->SetGuidance("Set maximum phi(obsolete!)");
00566   maxphiCmd->SetParameterName("MaxPhi",false,false);
00567   maxphiCmd->SetDefaultUnit("rad");
00568   //  maxphiCmd->SetUnitCandidates("rad deg");
00569 
00570   angsigmarCmd = new G4UIcmdWithADoubleAndUnit("/gps/sigmaangr",this);
00571   angsigmarCmd->SetGuidance("Set standard deviation of beam direction in radial(obsolete!).");
00572   angsigmarCmd->SetParameterName("Sigmara",false,false);
00573   angsigmarCmd->SetDefaultUnit("rad");
00574   //  angsigmarCmd->SetUnitCandidates("rad deg");
00575 
00576   angsigmaxCmd = new G4UIcmdWithADoubleAndUnit("/gps/sigmaangx",this);
00577   angsigmaxCmd->SetGuidance("Set standard deviation of beam direction in x-direc(obsolete!).");
00578   angsigmaxCmd->SetParameterName("Sigmaxa",false,false);
00579   angsigmaxCmd->SetDefaultUnit("rad");
00580   //  angsigmaxCmd->SetUnitCandidates("rad deg");
00581 
00582   angsigmayCmd = new G4UIcmdWithADoubleAndUnit("/gps/sigmaangy",this);
00583   angsigmayCmd->SetGuidance("Set standard deviation of beam direction in y-direc.(obsolete!)");
00584   angsigmayCmd->SetParameterName("Sigmaya",false,false);
00585   angsigmayCmd->SetDefaultUnit("rad");
00586   //  angsigmayCmd->SetUnitCandidates("rad deg");
00587 
00588   useuserangaxisCmd = new G4UIcmdWithABool("/gps/useuserangaxis",this);
00589   useuserangaxisCmd->SetGuidance("true for using user defined angular co-ordinates(obsolete!)");
00590   useuserangaxisCmd->SetGuidance("Default is false");
00591   useuserangaxisCmd->SetParameterName("useuserangaxis",true);
00592   useuserangaxisCmd->SetDefaultValue(false);
00593 
00594   surfnormCmd = new G4UIcmdWithABool("/gps/surfnorm",this);
00595   surfnormCmd->SetGuidance("Makes a user-defined distribution with respect to surface normals rather than x,y,z axes (obsolete!).");
00596   surfnormCmd->SetGuidance("Default is false");
00597   surfnormCmd->SetParameterName("surfnorm",true);
00598   surfnormCmd->SetDefaultValue(false);
00599 
00600   // Energy commands
00601 
00602   energyDirectory = new G4UIdirectory("/gps/ene/");
00603   energyDirectory->SetGuidance("Spectral commands sub-directory");
00604 
00605   energytypeCmd1 = new G4UIcmdWithAString("/gps/ene/type",this);
00606   energytypeCmd1->SetGuidance("Sets energy distribution type");
00607   energytypeCmd1->SetParameterName("EnergyDis",false,false);
00608   energytypeCmd1->SetDefaultValue("Mono");
00609   energytypeCmd1->SetCandidates("Mono Lin Pow Exp Gauss Brem Bbody Cdg User Arb Epn");
00610 
00611   eminCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ene/min",this);
00612   eminCmd1->SetGuidance("Sets minimum energy");
00613   eminCmd1->SetParameterName("emin",false,false);
00614   eminCmd1->SetDefaultUnit("keV");
00615   //  eminCmd1->SetUnitCandidates("eV keV MeV GeV TeV PeV");
00616 
00617   emaxCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ene/max",this);
00618   emaxCmd1->SetGuidance("Sets maximum energy");
00619   emaxCmd1->SetParameterName("emax",false,false);
00620   emaxCmd1->SetDefaultUnit("keV");
00621   //  emaxCmd1->SetUnitCandidates("eV keV MeV GeV TeV PeV");
00622 
00623   monoenergyCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ene/mono",this);
00624   monoenergyCmd1->SetGuidance("Sets a monocromatic energy (same as  gps/energy)");
00625   monoenergyCmd1->SetParameterName("monoenergy",false,false);
00626   monoenergyCmd1->SetDefaultUnit("keV");
00627   //  monoenergyCmd1->SetUnitCandidates("eV keV MeV GeV TeV PeV");
00628 
00629   engsigmaCmd1 = new G4UIcmdWithADoubleAndUnit("/gps/ene/sigma",this);
00630   engsigmaCmd1->SetGuidance("Sets the standard deviation for Gaussian energy dist.");
00631   engsigmaCmd1->SetParameterName("Sigmae",false,false);
00632   engsigmaCmd1->SetDefaultUnit("keV");
00633   //  engsigmaCmd1->SetUnitCandidates("eV keV MeV GeV TeV PeV");
00634 
00635   alphaCmd1 = new G4UIcmdWithADouble("/gps/ene/alpha",this);
00636   alphaCmd1->SetGuidance("Sets Alpha (index) for power-law energy dist.");
00637   alphaCmd1->SetParameterName("alpha",false,false);
00638   
00639   tempCmd1 = new G4UIcmdWithADouble("/gps/ene/temp",this);
00640   tempCmd1->SetGuidance("Sets the temperature for Brem and BBody distributions (in Kelvin)");
00641   tempCmd1->SetParameterName("temp",false,false);
00642 
00643   ezeroCmd1 = new G4UIcmdWithADouble("/gps/ene/ezero",this);
00644   ezeroCmd1->SetGuidance("Sets E_0 for exponential distribution (in MeV)");
00645   ezeroCmd1->SetParameterName("ezero",false,false);
00646 
00647   gradientCmd1 = new G4UIcmdWithADouble("/gps/ene/gradient",this);
00648   gradientCmd1->SetGuidance("Sets the gradient for Lin distribution (in 1/MeV)");
00649   gradientCmd1->SetParameterName("gradient",false,false);
00650 
00651   interceptCmd1 = new G4UIcmdWithADouble("/gps/ene/intercept",this);
00652   interceptCmd1->SetGuidance("Sets the intercept for Lin distributions (in MeV)");
00653   interceptCmd1->SetParameterName("intercept",false,false);
00654 
00655   arbeintCmd1 = new G4UIcmdWithADouble("/gps/ene/biasAlpha",this);
00656   arbeintCmd1->SetGuidance("Set the power-law index for the energy sampling distri. )");
00657   arbeintCmd1->SetParameterName("arbeint",false,false);
00658 
00659   calculateCmd1 = new G4UIcmdWithoutParameter("/gps/ene/calculate",this);
00660   calculateCmd1->SetGuidance("Calculates the distributions for Cdg and BBody");
00661 
00662   energyspecCmd1 = new G4UIcmdWithABool("/gps/ene/emspec",this);
00663   energyspecCmd1->SetGuidance("True for energy and false for momentum spectra");
00664   energyspecCmd1->SetParameterName("energyspec",true);
00665   energyspecCmd1->SetDefaultValue(true);
00666 
00667   diffspecCmd1 = new G4UIcmdWithABool("/gps/ene/diffspec",this);
00668   diffspecCmd1->SetGuidance("True for differential and flase for integral spectra");
00669   diffspecCmd1->SetParameterName("diffspec",true);
00670   diffspecCmd1->SetDefaultValue(true);
00671 
00672   //old ones
00673   energytypeCmd = new G4UIcmdWithAString("/gps/energytype",this);
00674   energytypeCmd->SetGuidance("Sets energy distribution type (obsolete!)");
00675   energytypeCmd->SetParameterName("EnergyDis",false,false);
00676   energytypeCmd->SetDefaultValue("Mono");
00677   energytypeCmd->SetCandidates("Mono Lin Pow Exp Gauss Brem Bbody Cdg User Arb Epn");
00678 
00679   eminCmd = new G4UIcmdWithADoubleAndUnit("/gps/emin",this);
00680   eminCmd->SetGuidance("Sets Emin (obsolete!)");
00681   eminCmd->SetParameterName("emin",false,false);
00682   eminCmd->SetDefaultUnit("keV");
00683   //  eminCmd->SetUnitCandidates("eV keV MeV GeV TeV PeV");
00684 
00685   emaxCmd = new G4UIcmdWithADoubleAndUnit("/gps/emax",this);
00686   emaxCmd->SetGuidance("Sets Emax (obsolete!)");
00687   emaxCmd->SetParameterName("emax",false,false);
00688   emaxCmd->SetDefaultUnit("keV");
00689   //  emaxCmd->SetUnitCandidates("eV keV MeV GeV TeV PeV");
00690 
00691   monoenergyCmd = new G4UIcmdWithADoubleAndUnit("/gps/monoenergy",this);
00692   monoenergyCmd->SetGuidance("Sets Monoenergy (obsolete, use gps/energy instead!)");
00693   monoenergyCmd->SetParameterName("monoenergy",false,false);
00694   monoenergyCmd->SetDefaultUnit("keV");
00695   //  monoenergyCmd->SetUnitCandidates("eV keV MeV GeV TeV PeV");
00696 
00697   engsigmaCmd = new G4UIcmdWithADoubleAndUnit("/gps/sigmae",this);
00698   engsigmaCmd->SetGuidance("Sets the standard deviation for Gaussian energy dist.(obsolete!)");
00699   engsigmaCmd->SetParameterName("Sigmae",false,false);
00700   engsigmaCmd->SetDefaultUnit("keV");
00701   //  engsigmaCmd->SetUnitCandidates("eV keV MeV GeV TeV PeV");
00702 
00703   alphaCmd = new G4UIcmdWithADouble("/gps/alpha",this);
00704   alphaCmd->SetGuidance("Sets Alpha (index) for power-law energy dist(obsolete!).");
00705   alphaCmd->SetParameterName("alpha",false,false);
00706   
00707   tempCmd = new G4UIcmdWithADouble("/gps/temp",this);
00708   tempCmd->SetGuidance("Sets the temperature for Brem and BBody (in Kelvin)(obsolete!)");
00709   tempCmd->SetParameterName("temp",false,false);
00710 
00711   ezeroCmd = new G4UIcmdWithADouble("/gps/ezero",this);
00712   ezeroCmd->SetGuidance("Sets ezero exponential distributions (in MeV)(obsolete!)");
00713   ezeroCmd->SetParameterName("ezero",false,false);
00714 
00715   gradientCmd = new G4UIcmdWithADouble("/gps/gradient",this);
00716   gradientCmd->SetGuidance("Sets the gradient for Lin distributions (in 1/MeV)(obsolete!)");
00717   gradientCmd->SetParameterName("gradient",false,false);
00718 
00719   interceptCmd = new G4UIcmdWithADouble("/gps/intercept",this);
00720   interceptCmd->SetGuidance("Sets the intercept for Lin distributions (in MeV)(obsolete!)");
00721   interceptCmd->SetParameterName("intercept",false,false);
00722 
00723   calculateCmd = new G4UIcmdWithoutParameter("/gps/calculate",this);
00724   calculateCmd->SetGuidance("Calculates distributions for Cdg and BBody(obsolete!)");
00725 
00726   energyspecCmd = new G4UIcmdWithABool("/gps/energyspec",this);
00727   energyspecCmd->SetGuidance("True for energy and false for momentum spectra(obsolete!)");
00728   energyspecCmd->SetParameterName("energyspec",true);
00729   energyspecCmd->SetDefaultValue(true);
00730 
00731   diffspecCmd = new G4UIcmdWithABool("/gps/diffspec",this);
00732   diffspecCmd->SetGuidance("True for differential and flase for integral spectra(obsolete!)");
00733   diffspecCmd->SetParameterName("diffspec",true);
00734   diffspecCmd->SetDefaultValue(true);
00735 
00736   // Biasing + histograms in general
00737   histDirectory = new G4UIdirectory("/gps/hist/");
00738   histDirectory->SetGuidance("Histogram, biasing commands sub-directory");
00739 
00740   histnameCmd1 = new G4UIcmdWithAString("/gps/hist/type",this);
00741   histnameCmd1->SetGuidance("Sets histogram type");
00742   histnameCmd1->SetParameterName("HistType",false,false);
00743   histnameCmd1->SetDefaultValue("biasx");
00744   histnameCmd1->SetCandidates("biasx biasy biasz biast biasp biase biaspt biaspp theta phi energy arb epn");
00745 
00746   resethistCmd1 = new G4UIcmdWithAString("/gps/hist/reset",this);
00747   resethistCmd1->SetGuidance("Reset (clean) the histogram ");
00748   resethistCmd1->SetParameterName("HistType",false,false);
00749   resethistCmd1->SetDefaultValue("energy");
00750   resethistCmd1->SetCandidates("biasx biasy biasz biast biasp biase biaspt biaspp theta phi energy arb epn");
00751 
00752   histpointCmd1 = new G4UIcmdWith3Vector("/gps/hist/point",this);
00753   histpointCmd1->SetGuidance("Allows user to define a histogram");
00754   histpointCmd1->SetGuidance("Enter: Ehi Weight");
00755   histpointCmd1->SetParameterName("Ehi","Weight","Junk",true,true);
00756   histpointCmd1->SetRange("Ehi >= 0. && Weight >= 0.");
00757 
00758   histfileCmd1 = new G4UIcmdWithAString("/gps/hist/file",this);
00759   histfileCmd1->SetGuidance("import the arb energy hist in an ASCII file");
00760   histfileCmd1->SetParameterName("HistFile",false,false);
00761 
00762   arbintCmd1 = new G4UIcmdWithAString("/gps/hist/inter",this);
00763   arbintCmd1->SetGuidance("Sets the interpolation method for arbitrary distribution.");
00764   arbintCmd1->SetParameterName("int",false,false);
00765   arbintCmd1->SetDefaultValue("Lin");
00766   arbintCmd1->SetCandidates("Lin Log Exp Spline");
00767 
00768   // old ones
00769   histnameCmd = new G4UIcmdWithAString("/gps/histname",this);
00770   histnameCmd->SetGuidance("Sets histogram type (obsolete!)");
00771   histnameCmd->SetParameterName("HistType",false,false);
00772   histnameCmd->SetDefaultValue("biasx");
00773   histnameCmd->SetCandidates("biasx biasy biasz biast biasp biase biaspt biaspp theta phi energy arb epn");
00774 
00775   // re-set the histograms
00776   resethistCmd = new G4UIcmdWithAString("/gps/resethist",this);
00777   resethistCmd->SetGuidance("Re-Set the histogram (obsolete!)");
00778   resethistCmd->SetParameterName("HistType",false,false);
00779   resethistCmd->SetDefaultValue("energy");
00780   resethistCmd->SetCandidates("biasx biasy biasz biast biasp biase biaspt biaspp theta phi energy arb epn");
00781 
00782   histpointCmd = new G4UIcmdWith3Vector("/gps/histpoint",this);
00783   histpointCmd->SetGuidance("Allows user to define a histogram (obsolete!)");
00784   histpointCmd->SetGuidance("Enter: Ehi Weight");
00785   histpointCmd->SetParameterName("Ehi","Weight","Junk",false,false);
00786   histpointCmd->SetRange("Ehi >= 0. && Weight >= 0.");
00787 
00788   arbintCmd = new G4UIcmdWithAString("/gps/arbint",this);
00789   arbintCmd->SetGuidance("Sets Arbitrary Interpolation type.(obsolete!) ");
00790   arbintCmd->SetParameterName("int",false,false);
00791   arbintCmd->SetDefaultValue("NULL");
00792   arbintCmd->SetCandidates("Lin Log Exp Spline");
00793 
00794 }

G4GeneralParticleSourceMessenger::~G4GeneralParticleSourceMessenger (  ) 

Definition at line 796 of file G4GeneralParticleSourceMessenger.cc.

00797 {
00798   delete positionDirectory;
00799   delete typeCmd;
00800   delete shapeCmd;
00801   delete centreCmd;
00802   delete posrot1Cmd;
00803   delete posrot2Cmd;
00804   delete halfxCmd;
00805   delete halfyCmd;
00806   delete halfzCmd;
00807   delete radiusCmd;
00808   delete radius0Cmd;
00809   delete possigmarCmd;
00810   delete possigmaxCmd;
00811   delete possigmayCmd;
00812   delete paralpCmd;
00813   delete partheCmd;
00814   delete parphiCmd;
00815   delete confineCmd;
00816   delete typeCmd1;
00817   delete shapeCmd1;
00818   delete centreCmd1;
00819   delete posrot1Cmd1;
00820   delete posrot2Cmd1;
00821   delete halfxCmd1;
00822   delete halfyCmd1;
00823   delete halfzCmd1;
00824   delete radiusCmd1;
00825   delete radius0Cmd1;
00826   delete possigmarCmd1;
00827   delete possigmaxCmd1;
00828   delete possigmayCmd1;
00829   delete paralpCmd1;
00830   delete partheCmd1;
00831   delete parphiCmd1;
00832   delete confineCmd1;
00833 
00834   delete angularDirectory;
00835   delete angtypeCmd;
00836   delete angrot1Cmd;
00837   delete angrot2Cmd;
00838   delete minthetaCmd;
00839   delete maxthetaCmd;
00840   delete minphiCmd;
00841   delete maxphiCmd;
00842   delete angsigmarCmd;
00843   delete angsigmaxCmd;
00844   delete angsigmayCmd;
00845   delete useuserangaxisCmd;
00846   delete surfnormCmd;
00847   delete angtypeCmd1;
00848   delete angrot1Cmd1;
00849   delete angrot2Cmd1;
00850   delete minthetaCmd1;
00851   delete maxthetaCmd1;
00852   delete minphiCmd1;
00853   delete maxphiCmd1;
00854   delete angsigmarCmd1;
00855   delete angsigmaxCmd1;
00856   delete angsigmayCmd1;
00857   delete angfocusCmd;
00858   delete useuserangaxisCmd1;
00859   delete surfnormCmd1;
00860 
00861   delete energyDirectory;
00862   delete energytypeCmd;
00863   delete eminCmd;
00864   delete emaxCmd;
00865   delete monoenergyCmd;
00866   delete engsigmaCmd;
00867   delete alphaCmd;
00868   delete tempCmd;
00869   delete ezeroCmd;
00870   delete gradientCmd;
00871   delete interceptCmd;
00872   delete calculateCmd;
00873   delete energyspecCmd;
00874   delete diffspecCmd;
00875   delete energytypeCmd1;
00876   delete eminCmd1;
00877   delete emaxCmd1;
00878   delete monoenergyCmd1;
00879   delete engsigmaCmd1;
00880   delete alphaCmd1;
00881   delete tempCmd1;
00882   delete ezeroCmd1;
00883   delete gradientCmd1;
00884   delete interceptCmd1;
00885   delete arbeintCmd1;
00886   delete calculateCmd1;
00887   delete energyspecCmd1;
00888   delete diffspecCmd1;
00889 
00890   delete histDirectory;
00891   delete histnameCmd;
00892   delete resethistCmd;
00893   delete histpointCmd;
00894   delete arbintCmd;
00895   delete histnameCmd1;
00896   delete resethistCmd1;
00897   delete histpointCmd1;
00898   delete histfileCmd1;
00899   delete arbintCmd1;
00900 
00901   delete verbosityCmd;
00902   delete ionCmd;
00903   delete particleCmd;
00904   delete timeCmd;
00905   delete polCmd;
00906   delete numberCmd;
00907   delete positionCmd;
00908   delete directionCmd;
00909   delete energyCmd;
00910   delete listCmd;
00911 
00912   delete sourceDirectory;
00913   delete addsourceCmd;
00914   delete listsourceCmd;
00915   delete clearsourceCmd;
00916   delete getsourceCmd;
00917   delete setsourceCmd;
00918   delete setintensityCmd;
00919   delete deletesourceCmd;
00920   delete multiplevertexCmd;
00921   delete flatsamplingCmd;
00922 
00923   delete gpsDirectory;
00924   
00925 }


Member Function Documentation

G4String G4GeneralParticleSourceMessenger::GetCurrentValue ( G4UIcommand command  )  [virtual]

Reimplemented from G4UImessenger.

Definition at line 1596 of file G4GeneralParticleSourceMessenger.cc.

01597 {
01598   G4String cv;
01599   
01600   //  if( command==directionCmd )
01601   //  { cv = directionCmd->ConvertToString(fParticleGun->GetParticleMomentumDirection()); }
01602   //  else if( command==energyCmd )
01603   //  { cv = energyCmd->ConvertToString(fParticleGun->GetParticleEnergy(),"GeV"); }
01604   //  else if( command==positionCmd )
01605   //  { cv = positionCmd->ConvertToString(fParticleGun->GetParticlePosition(),"cm"); }
01606   //  else if( command==timeCmd )
01607   //  { cv = timeCmd->ConvertToString(fParticleGun->GetParticleTime(),"ns"); }
01608   //  else if( command==polCmd )
01609   //  { cv = polCmd->ConvertToString(fParticleGun->GetParticlePolarization()); }
01610   //  else if( command==numberCmd )
01611   //  { cv = numberCmd->ConvertToString(fParticleGun->GetNumberOfParticles()); }
01612   
01613   cv = "Not implemented yet";
01614 
01615   return cv;
01616 }

void G4GeneralParticleSourceMessenger::SetNewValue ( G4UIcommand command,
G4String  newValues 
) [virtual]

Reimplemented from G4UImessenger.

Definition at line 927 of file G4GeneralParticleSourceMessenger.cc.

References G4GeneralParticleSource::AddaSource(), G4SPSEneDistribution::ArbEnergyHisto(), G4SPSEneDistribution::ArbEnergyHistoFile(), G4SPSEneDistribution::ArbInterpolate(), G4SPSEneDistribution::Calculate(), G4GeneralParticleSource::ClearAll(), G4SPSPosDistribution::ConfineSourceToVolume(), G4SPSAngDistribution::DefineAngRefAxes(), G4GeneralParticleSource::DeleteaSource(), G4ParticleTable::DumpTable(), G4SPSEneDistribution::EpnEnergyHisto(), G4ParticleTable::FindParticle(), G4cout, G4endl, G4SingleParticleSource::GetAngDist(), G4SingleParticleSource::GetBiasRndm(), G4GeneralParticleSource::GetCurrentSourceIndex(), G4GeneralParticleSource::GetCurrentSourceIntensity(), G4SingleParticleSource::GetEneDist(), G4UIcmdWith3Vector::GetNew3VectorValue(), G4UIcmdWith3VectorAndUnit::GetNew3VectorValue(), G4UIcmdWithABool::GetNewBoolValue(), G4UIcmdWithADouble::GetNewDoubleValue(), G4UIcmdWithADoubleAndUnit::GetNewDoubleValue(), G4UIcmdWithAnInteger::GetNewIntValue(), G4SingleParticleSource::GetPosDist(), G4SPSEneDistribution::InputDifferentialSpectra(), G4SPSEneDistribution::InputEnergySpectra(), G4GeneralParticleSource::ListSource(), G4SPSRandomGenerator::ReSetHist(), G4SPSEneDistribution::ReSetHist(), G4SPSAngDistribution::ReSetHist(), G4SPSEneDistribution::SetAlpha(), G4SPSAngDistribution::SetAngDistType(), G4SPSAngDistribution::SetBeamSigmaInAngR(), G4SPSAngDistribution::SetBeamSigmaInAngX(), G4SPSAngDistribution::SetBeamSigmaInAngY(), G4SPSEneDistribution::SetBeamSigmaInE(), G4SPSPosDistribution::SetBeamSigmaInR(), G4SPSPosDistribution::SetBeamSigmaInX(), G4SPSPosDistribution::SetBeamSigmaInY(), G4SPSEneDistribution::SetBiasAlpha(), G4SPSPosDistribution::SetCentreCoords(), G4GeneralParticleSource::SetCurrentSourceIntensity(), G4GeneralParticleSource::SetCurrentSourceto(), G4SPSEneDistribution::SetEmax(), G4SPSEneDistribution::SetEmin(), G4SPSRandomGenerator::SetEnergyBias(), G4SPSEneDistribution::SetEnergyDisType(), G4SPSEneDistribution::SetEzero(), G4GeneralParticleSource::SetFlatSampling(), G4SPSAngDistribution::SetFocusPoint(), G4SPSEneDistribution::SetGradient(), G4SPSPosDistribution::SetHalfX(), G4SPSPosDistribution::SetHalfY(), G4SPSPosDistribution::SetHalfZ(), G4SPSEneDistribution::SetInterCept(), G4SPSAngDistribution::SetMaxPhi(), G4SPSAngDistribution::SetMaxTheta(), G4SPSAngDistribution::SetMinPhi(), G4SPSAngDistribution::SetMinTheta(), G4SPSEneDistribution::SetMonoEnergy(), G4GeneralParticleSource::SetMultipleVertex(), G4SingleParticleSource::SetNumberOfParticles(), G4SPSPosDistribution::SetParAlpha(), G4SPSPosDistribution::SetParPhi(), G4SPSPosDistribution::SetParTheta(), G4SingleParticleSource::SetParticleDefinition(), G4SPSAngDistribution::SetParticleMomentumDirection(), G4SingleParticleSource::SetParticlePolarization(), G4SingleParticleSource::SetParticleTime(), G4SPSRandomGenerator::SetPhiBias(), G4SPSPosDistribution::SetPosDisShape(), G4SPSPosDistribution::SetPosDisType(), G4SPSRandomGenerator::SetPosPhiBias(), G4SPSPosDistribution::SetPosRot1(), G4SPSPosDistribution::SetPosRot2(), G4SPSRandomGenerator::SetPosThetaBias(), G4SPSPosDistribution::SetRadius(), G4SPSPosDistribution::SetRadius0(), G4SPSEneDistribution::SetTemp(), G4SPSRandomGenerator::SetThetaBias(), G4SPSAngDistribution::SetUserWRTSurface(), G4SPSAngDistribution::SetUseUserAngAxis(), G4SingleParticleSource::SetVerbosity(), G4SPSRandomGenerator::SetXBias(), G4SPSRandomGenerator::SetYBias(), G4SPSRandomGenerator::SetZBias(), G4SPSAngDistribution::UserDefAngPhi(), G4SPSAngDistribution::UserDefAngTheta(), and G4SPSEneDistribution::UserEnergyHisto().

00928 {
00929   if(command == typeCmd)
00930     {
00931       fParticleGun->GetPosDist()->SetPosDisType(newValues);
00932       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
00933              << " The command is obsolete and will be removed soon." << G4endl
00934              << " Please try to use the new structured commands!" << G4endl;
00935     }
00936   else if(command == shapeCmd)
00937     {
00938       fParticleGun->GetPosDist()->SetPosDisShape(newValues);
00939       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
00940              << " The command is obsolete and will be removed soon." << G4endl
00941              << " Please try to use the new structured commands!" << G4endl;
00942     }
00943   else if(command == centreCmd)
00944     {
00945       fParticleGun->GetPosDist()->SetCentreCoords(centreCmd->GetNew3VectorValue(newValues));
00946       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
00947              << " The command is obsolete and will be removed soon." << G4endl
00948              << " Please try to use the new structured commands!" << G4endl;
00949     }
00950   else if(command == posrot1Cmd)
00951     {
00952       fParticleGun->GetPosDist()->SetPosRot1(posrot1Cmd->GetNew3VectorValue(newValues));
00953       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
00954              << " The command is obsolete and will be removed soon." << G4endl
00955              << " Please try to use the new structured commands!" << G4endl;
00956     }
00957   else if(command == posrot2Cmd)
00958     {
00959       fParticleGun->GetPosDist()->SetPosRot2(posrot2Cmd->GetNew3VectorValue(newValues));
00960       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
00961              << " The command is obsolete and will be removed soon." << G4endl
00962              << " Please try to use the new structured commands!" << G4endl;
00963     }
00964   else if(command == halfxCmd)
00965     {
00966       fParticleGun->GetPosDist()->SetHalfX(halfxCmd->GetNewDoubleValue(newValues));
00967       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
00968              << " The command is obsolete and will be removed soon." << G4endl
00969              << " Please try to use the new structured commands!" << G4endl;
00970     }
00971   else if(command == halfyCmd)
00972     {
00973       fParticleGun->GetPosDist()->SetHalfY(halfyCmd->GetNewDoubleValue(newValues));
00974       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
00975              << " The command is obsolete and will be removed soon." << G4endl
00976              << " Please try to use the new structured commands!" << G4endl;
00977     }
00978   else if(command == halfzCmd)
00979     {
00980       fParticleGun->GetPosDist()->SetHalfZ(halfzCmd->GetNewDoubleValue(newValues));
00981       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
00982              << " The command is obsolete and will be removed soon." << G4endl
00983              << " Please try to use the new structured commands!" << G4endl;
00984     }
00985   else if(command == radiusCmd)
00986     {
00987       fParticleGun->GetPosDist()->SetRadius(radiusCmd->GetNewDoubleValue(newValues));
00988       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
00989              << " The command is obsolete and will be removed soon." << G4endl
00990              << " Please try to use the new structured commands!" << G4endl;
00991     }
00992   else if(command == radius0Cmd)
00993     {
00994       fParticleGun->GetPosDist()->SetRadius0(radius0Cmd->GetNewDoubleValue(newValues));
00995       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
00996              << " The command is obsolete and will be removed soon." << G4endl
00997              << " Please try to use the new structured commands!" << G4endl;
00998     }
00999   else if(command == possigmarCmd)
01000     {
01001       fParticleGun->GetPosDist()->SetBeamSigmaInR(possigmarCmd->GetNewDoubleValue(newValues));
01002       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01003              << " The command is obsolete and will be removed soon." << G4endl
01004              << " Please try to use the new structured commands!" << G4endl;
01005     }
01006   else if(command == possigmaxCmd)
01007     {
01008       fParticleGun->GetPosDist()->SetBeamSigmaInX(possigmaxCmd->GetNewDoubleValue(newValues));
01009       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01010              << " The command is obsolete and will be removed soon." << G4endl
01011              << " Please try to use the new structured commands!" << G4endl;
01012     }
01013   else if(command == possigmayCmd)
01014     {
01015       fParticleGun->GetPosDist()->SetBeamSigmaInY(possigmayCmd->GetNewDoubleValue(newValues));
01016       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01017              << " The command is obsolete and will be removed soon." << G4endl
01018              << " Please try to use the new structured commands!" << G4endl;
01019     }
01020   else if(command == paralpCmd)
01021     {
01022       fParticleGun->GetPosDist()->SetParAlpha(paralpCmd->GetNewDoubleValue(newValues));
01023       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01024              << " The command is obsolete and will be removed soon." << G4endl
01025              << " Please try to use the new structured commands!" << G4endl;
01026     }
01027   else if(command == partheCmd)
01028     {
01029       fParticleGun->GetPosDist()->SetParTheta(partheCmd->GetNewDoubleValue(newValues));
01030       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01031              << " The command is obsolete and will be removed soon." << G4endl
01032              << " Please try to use the new structured commands!" << G4endl;
01033     }
01034   else if(command == parphiCmd)
01035     {
01036       fParticleGun->GetPosDist()->SetParPhi(parphiCmd->GetNewDoubleValue(newValues));
01037       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01038              << " The command is obsolete and will be removed soon." << G4endl
01039              << " Please try to use the new structured commands!" << G4endl;
01040     }
01041   else if(command == confineCmd)
01042     {
01043       fParticleGun->GetPosDist()->ConfineSourceToVolume(newValues);
01044       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01045              << " The command is obsolete and will be removed soon." << G4endl
01046              << " Please try to use the new structured commands!" << G4endl;
01047     }
01048   else if(command == angtypeCmd)
01049     {
01050       fParticleGun->GetAngDist()->SetAngDistType(newValues);
01051       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01052              << " The command is obsolete and will be removed soon." << G4endl
01053              << " Please try to use the new structured commands!" << G4endl;
01054     }
01055   else if(command == angrot1Cmd)
01056     {
01057       G4String a = "angref1";
01058       fParticleGun->GetAngDist()->DefineAngRefAxes(a,angrot1Cmd->GetNew3VectorValue(newValues));
01059       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01060              << " The command is obsolete and will be removed soon." << G4endl
01061              << " Please try to use the new structured commands!" << G4endl;
01062     }
01063   else if(command == angrot2Cmd)
01064     {
01065       G4String a = "angref2";
01066       fParticleGun->GetAngDist()->DefineAngRefAxes(a,angrot2Cmd->GetNew3VectorValue(newValues));
01067       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01068              << " The command is obsolete and will be removed soon." << G4endl
01069              << " Please try to use the new structured commands!" << G4endl;
01070     }
01071   else if(command == minthetaCmd)
01072     {
01073       fParticleGun->GetAngDist()->SetMinTheta(minthetaCmd->GetNewDoubleValue(newValues));
01074       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01075              << " The command is obsolete and will be removed soon." << G4endl
01076              << " Please try to use the new structured commands!" << G4endl;
01077     }
01078   else if(command == minphiCmd)
01079     {
01080       fParticleGun->GetAngDist()->SetMinPhi(minphiCmd->GetNewDoubleValue(newValues));
01081       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01082              << " The command is obsolete and will be removed soon." << G4endl
01083              << " Please try to use the new structured commands!" << G4endl;
01084     }
01085   else if(command == maxthetaCmd)
01086     {
01087       fParticleGun->GetAngDist()->SetMaxTheta(maxthetaCmd->GetNewDoubleValue(newValues));
01088       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01089              << " The command is obsolete and will be removed soon." << G4endl
01090              << " Please try to use the new structured commands!" << G4endl;
01091     }
01092   else if(command == maxphiCmd)
01093     {
01094       fParticleGun->GetAngDist()->SetMaxPhi(maxphiCmd->GetNewDoubleValue(newValues));
01095       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01096              << " The command is obsolete and will be removed soon." << G4endl
01097              << " Please try to use the new structured commands!" << G4endl;
01098     }
01099   else if(command == angsigmarCmd)
01100     {
01101       fParticleGun->GetAngDist()->SetBeamSigmaInAngR(angsigmarCmd->GetNewDoubleValue(newValues));
01102       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01103              << " The command is obsolete and will be removed soon." << G4endl
01104              << " Please try to use the new structured commands!" << G4endl;
01105     }
01106   else if(command == angsigmaxCmd)
01107     {
01108       fParticleGun->GetAngDist()->SetBeamSigmaInAngX(angsigmaxCmd->GetNewDoubleValue(newValues));
01109       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01110              << " The command is obsolete and will be removed soon." << G4endl
01111              << " Please try to use the new structured commands!" << G4endl;
01112     }
01113   else if(command == angsigmayCmd)
01114     {
01115       fParticleGun->GetAngDist()->SetBeamSigmaInAngY(angsigmayCmd->GetNewDoubleValue(newValues));
01116       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01117              << " The command is obsolete and will be removed soon." << G4endl
01118              << " Please try to use the new structured commands!" << G4endl;
01119     }
01120   else if(command == useuserangaxisCmd)
01121     {
01122       fParticleGun->GetAngDist()->SetUseUserAngAxis(useuserangaxisCmd->GetNewBoolValue(newValues));
01123       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01124              << " The command is obsolete and will be removed soon." << G4endl
01125              << " Please try to use the new structured commands!" << G4endl;
01126     }
01127   else if(command == surfnormCmd)
01128     {
01129       fParticleGun->GetAngDist()->SetUserWRTSurface(surfnormCmd->GetNewBoolValue(newValues));
01130       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01131              << " The command is obsolete and will be removed soon." << G4endl
01132              << " Please try to use the new structured commands!" << G4endl;
01133     }
01134   else if(command == energytypeCmd)
01135     {
01136       fParticleGun->GetEneDist()->SetEnergyDisType(newValues);
01137       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01138              << " The command is obsolete and will be removed soon." << G4endl
01139              << " Please try to use the new structured commands!" << G4endl;
01140     }
01141   else if(command == eminCmd)
01142     {
01143       fParticleGun->GetEneDist()->SetEmin(eminCmd->GetNewDoubleValue(newValues));
01144       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01145              << " The command is obsolete and will be removed soon." << G4endl
01146              << " Please try to use the new structured commands!" << G4endl;
01147     }
01148   else if(command == emaxCmd)
01149     {
01150       fParticleGun->GetEneDist()->SetEmax(emaxCmd->GetNewDoubleValue(newValues));
01151       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01152              << " The command is obsolete and will be removed soon." << G4endl
01153              << " Please try to use the new structured commands!" << G4endl;
01154     }
01155   else if(command == monoenergyCmd)
01156     {
01157       fParticleGun->GetEneDist()->SetMonoEnergy(monoenergyCmd->GetNewDoubleValue(newValues));
01158       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01159              << " The command is obsolete and will be removed soon." << G4endl
01160              << " Please try to use the new structured commands!" << G4endl;
01161     }
01162   else if(command == engsigmaCmd)
01163     {
01164       fParticleGun->GetEneDist()->SetBeamSigmaInE(engsigmaCmd->GetNewDoubleValue(newValues));
01165       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01166              << " The command is obsolete and will be removed soon." << G4endl
01167              << " Please try to use the new structured commands!" << G4endl;
01168     }
01169   else if(command == alphaCmd)
01170     {
01171       fParticleGun->GetEneDist()->SetAlpha(alphaCmd->GetNewDoubleValue(newValues));
01172       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01173              << " The command is obsolete and will be removed soon." << G4endl
01174              << " Please try to use the new structured commands!" << G4endl;
01175     }
01176   else if(command == tempCmd)
01177     {
01178       fParticleGun->GetEneDist()->SetTemp(tempCmd->GetNewDoubleValue(newValues));
01179       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01180              << " The command is obsolete and will be removed soon." << G4endl
01181              << " Please try to use the new structured commands!" << G4endl;
01182     }
01183   else if(command == ezeroCmd)
01184     {
01185       fParticleGun->GetEneDist()->SetEzero(ezeroCmd->GetNewDoubleValue(newValues));
01186       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01187              << " The command is obsolete and will be removed soon." << G4endl
01188              << " Please try to use the new structured commands!" << G4endl;
01189     }
01190   else if(command == gradientCmd)
01191     {
01192       fParticleGun->GetEneDist()->SetGradient(gradientCmd->GetNewDoubleValue(newValues));
01193       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01194              << " The command is obsolete and will be removed soon." << G4endl
01195              << " Please try to use the new structured commands!" << G4endl;
01196     }
01197   else if(command == interceptCmd)
01198     {
01199       fParticleGun->GetEneDist()->SetInterCept(interceptCmd->GetNewDoubleValue(newValues));
01200       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01201              << " The command is obsolete and will be removed soon." << G4endl
01202              << " Please try to use the new structured commands!" << G4endl;
01203     }
01204   else if(command == calculateCmd)
01205     {
01206       fParticleGun->GetEneDist()->Calculate();
01207       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01208              << " The command is obsolete and will be removed soon." << G4endl
01209              << " Please try to use the new structured commands!" << G4endl;
01210     }
01211   else if(command == energyspecCmd)
01212     {
01213       fParticleGun->GetEneDist()->InputEnergySpectra(energyspecCmd->GetNewBoolValue(newValues));
01214       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01215              << " The command is obsolete and will be removed soon." << G4endl
01216              << " Please try to use the new structured commands!" << G4endl;
01217     }
01218   else if(command == diffspecCmd)
01219     {
01220       fParticleGun->GetEneDist()->InputDifferentialSpectra(diffspecCmd->GetNewBoolValue(newValues));
01221       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01222              << " The command is obsolete and will be removed soon." << G4endl
01223              << " Please try to use the new structured commands!" << G4endl;
01224     }
01225   else if(command == histnameCmd)
01226     {
01227       histtype = newValues;
01228       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01229              << " The command is obsolete and will be removed soon." << G4endl
01230              << " Please try to use the new structured commands!" << G4endl;
01231     }
01232   else if(command == histpointCmd)
01233     {
01234       if(histtype == "biasx")
01235         fParticleGun->GetBiasRndm()->SetXBias(histpointCmd->GetNew3VectorValue(newValues));
01236       if(histtype == "biasy")
01237         fParticleGun->GetBiasRndm()->SetYBias(histpointCmd->GetNew3VectorValue(newValues));
01238       if(histtype == "biasz")
01239         fParticleGun->GetBiasRndm()->SetZBias(histpointCmd->GetNew3VectorValue(newValues));
01240       if(histtype == "biast")
01241         fParticleGun->GetBiasRndm()->SetThetaBias(histpointCmd->GetNew3VectorValue(newValues));
01242       if(histtype == "biasp")
01243         fParticleGun->GetBiasRndm()->SetPhiBias(histpointCmd->GetNew3VectorValue(newValues));
01244       if(histtype == "biase")
01245         fParticleGun->GetBiasRndm()->SetEnergyBias(histpointCmd->GetNew3VectorValue(newValues));
01246       if(histtype == "theta")
01247         fParticleGun->GetAngDist()->UserDefAngTheta(histpointCmd->GetNew3VectorValue(newValues));
01248       if(histtype == "phi")
01249         fParticleGun->GetAngDist()->UserDefAngPhi(histpointCmd->GetNew3VectorValue(newValues));
01250       if(histtype == "energy")
01251         fParticleGun->GetEneDist()->UserEnergyHisto(histpointCmd->GetNew3VectorValue(newValues));
01252       if(histtype == "arb")
01253         fParticleGun->GetEneDist()->ArbEnergyHisto(histpointCmd->GetNew3VectorValue(newValues));
01254       if(histtype == "epn")
01255         fParticleGun->GetEneDist()->EpnEnergyHisto(histpointCmd->GetNew3VectorValue(newValues));
01256       G4cout << " G4GeneralParticleSourceMessenger - Warning: The command is obsolete and will be removed soon. Please try to use the new structured commands!" << G4endl;
01257     }
01258   else if(command == resethistCmd)
01259     {
01260       if(newValues == "theta" || newValues == "phi") {
01261         fParticleGun->GetAngDist()->ReSetHist(newValues);
01262       } else if (newValues == "energy" || newValues == "arb" || newValues == "epn") {
01263         fParticleGun->GetEneDist()->ReSetHist(newValues);
01264       } else {
01265         fParticleGun->GetBiasRndm()->ReSetHist(newValues);
01266       }
01267       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01268              << " The command is obsolete and will be removed soon." << G4endl
01269              << " Please try to use the new structured commands!" << G4endl;
01270     }
01271   else if(command == arbintCmd)
01272     {
01273       fParticleGun->GetEneDist()->ArbInterpolate(newValues);
01274       G4cout << " G4GeneralParticleSourceMessenger - Warning:" << G4endl
01275              << " The command is obsolete and will be removed soon." << G4endl
01276              << " Please try to use the new structured commands!" << G4endl;
01277     }
01278   else if( command==directionCmd )
01279     { 
01280       fParticleGun->GetAngDist()->SetAngDistType("planar");
01281       fParticleGun->GetAngDist()->SetParticleMomentumDirection(directionCmd->GetNew3VectorValue(newValues));
01282     }
01283   else if( command==energyCmd )
01284     {    
01285       fParticleGun->GetEneDist()->SetEnergyDisType("Mono");
01286       fParticleGun->GetEneDist()->SetMonoEnergy(energyCmd->GetNewDoubleValue(newValues));
01287     }
01288   else if( command==positionCmd )
01289     { 
01290       fParticleGun->GetPosDist()->SetPosDisType("Point");    
01291       fParticleGun->GetPosDist()->SetCentreCoords(positionCmd->GetNew3VectorValue(newValues));
01292     }
01293   else if(command == verbosityCmd)
01294     {
01295       fParticleGun->SetVerbosity(verbosityCmd->GetNewIntValue(newValues));
01296     }
01297   else if( command==particleCmd )
01298     {
01299       if (newValues =="ion") {
01300         fShootIon = true;
01301       } else {
01302         fShootIon = false;
01303         G4ParticleDefinition* pd = particleTable->FindParticle(newValues);
01304         if(pd != NULL)
01305           { fParticleGun->SetParticleDefinition( pd ); }
01306       }
01307     }
01308   else if( command==timeCmd )
01309     { fParticleGun->SetParticleTime(timeCmd->GetNewDoubleValue(newValues)); }
01310   else if( command==polCmd )
01311     { fParticleGun->SetParticlePolarization(polCmd->GetNew3VectorValue(newValues)); }
01312   else if( command==numberCmd )
01313     { fParticleGun->SetNumberOfParticles(numberCmd->GetNewIntValue(newValues)); }
01314   else if( command==ionCmd )
01315     { IonCommand(newValues); }
01316   else if( command==listCmd ){ 
01317     particleTable->DumpTable(); 
01318   }
01319   else if( command==addsourceCmd )
01320     { 
01321       fGPS->AddaSource(addsourceCmd->GetNewDoubleValue(newValues));
01322     }
01323   else if( command==listsourceCmd )
01324     { 
01325       fGPS->ListSource();
01326     }
01327   else if( command==clearsourceCmd )
01328     { 
01329       fGPS->ClearAll();
01330     }
01331   else if( command==getsourceCmd )
01332     { 
01333       G4cout << " Current source index:" << fGPS->GetCurrentSourceIndex() 
01334              << " ; Intensity:" << fGPS->GetCurrentSourceIntensity() << G4endl;
01335     }
01336   else if( command==setsourceCmd )
01337     { 
01338       fGPS->SetCurrentSourceto(setsourceCmd->GetNewIntValue(newValues));
01339     }
01340   else if( command==setintensityCmd )
01341     { 
01342       fGPS->SetCurrentSourceIntensity(setintensityCmd->GetNewDoubleValue(newValues));
01343     }
01344   else if( command==deletesourceCmd )
01345     { 
01346       fGPS->DeleteaSource(deletesourceCmd->GetNewIntValue(newValues));
01347     }
01348   else if(command == multiplevertexCmd)
01349     {
01350       fGPS->SetMultipleVertex(multiplevertexCmd->GetNewBoolValue(newValues));
01351     }
01352   else if(command == flatsamplingCmd)
01353     {
01354       fGPS->SetFlatSampling(flatsamplingCmd->GetNewBoolValue(newValues));
01355     }
01356   //
01357   // new implementations
01358   //
01359   //
01360   else if(command == typeCmd1)
01361     {
01362       fParticleGun->GetPosDist()->SetPosDisType(newValues);
01363     }
01364   else if(command == shapeCmd1)
01365     {
01366       fParticleGun->GetPosDist()->SetPosDisShape(newValues);
01367     }
01368   else if(command == centreCmd1)
01369     {
01370       fParticleGun->GetPosDist()->SetCentreCoords(centreCmd1->GetNew3VectorValue(newValues));
01371     }
01372   else if(command == posrot1Cmd1)
01373     {
01374       fParticleGun->GetPosDist()->SetPosRot1(posrot1Cmd1->GetNew3VectorValue(newValues));
01375     }
01376   else if(command == posrot2Cmd1)
01377     {
01378       fParticleGun->GetPosDist()->SetPosRot2(posrot2Cmd1->GetNew3VectorValue(newValues));
01379     }
01380   else if(command == halfxCmd1)
01381     {
01382       fParticleGun->GetPosDist()->SetHalfX(halfxCmd1->GetNewDoubleValue(newValues));
01383     }
01384   else if(command == halfyCmd1)
01385     {
01386       fParticleGun->GetPosDist()->SetHalfY(halfyCmd1->GetNewDoubleValue(newValues));
01387     }
01388   else if(command == halfzCmd1)
01389     {
01390       fParticleGun->GetPosDist()->SetHalfZ(halfzCmd1->GetNewDoubleValue(newValues));
01391     }
01392   else if(command == radiusCmd1)
01393     {
01394       fParticleGun->GetPosDist()->SetRadius(radiusCmd1->GetNewDoubleValue(newValues));
01395     }
01396   else if(command == radius0Cmd1)
01397     {
01398       fParticleGun->GetPosDist()->SetRadius0(radius0Cmd1->GetNewDoubleValue(newValues));
01399     }
01400   else if(command == possigmarCmd1)
01401     {
01402       fParticleGun->GetPosDist()->SetBeamSigmaInR(possigmarCmd1->GetNewDoubleValue(newValues));
01403     }
01404   else if(command == possigmaxCmd1)
01405     {
01406       fParticleGun->GetPosDist()->SetBeamSigmaInX(possigmaxCmd1->GetNewDoubleValue(newValues));
01407     }
01408   else if(command == possigmayCmd1)
01409     {
01410       fParticleGun->GetPosDist()->SetBeamSigmaInY(possigmayCmd1->GetNewDoubleValue(newValues));
01411     }
01412   else if(command == paralpCmd1)
01413     {
01414       fParticleGun->GetPosDist()->SetParAlpha(paralpCmd1->GetNewDoubleValue(newValues));
01415     }
01416   else if(command == partheCmd1)
01417     {
01418       fParticleGun->GetPosDist()->SetParTheta(partheCmd1->GetNewDoubleValue(newValues));
01419     }
01420   else if(command == parphiCmd1)
01421     {
01422       fParticleGun->GetPosDist()->SetParPhi(parphiCmd1->GetNewDoubleValue(newValues));
01423     }
01424   else if(command == confineCmd1)
01425     {
01426       fParticleGun->GetPosDist()->ConfineSourceToVolume(newValues);
01427     }
01428   else if(command == angtypeCmd1)
01429     {
01430       fParticleGun->GetAngDist()->SetAngDistType(newValues);
01431     }
01432   else if(command == angrot1Cmd1)
01433     {
01434       G4String a = "angref1";
01435       fParticleGun->GetAngDist()->DefineAngRefAxes(a,angrot1Cmd1->GetNew3VectorValue(newValues));
01436     }
01437   else if(command == angrot2Cmd1)
01438     {
01439       G4String a = "angref2";
01440       fParticleGun->GetAngDist()->DefineAngRefAxes(a,angrot2Cmd1->GetNew3VectorValue(newValues));
01441     }
01442   else if(command == minthetaCmd1)
01443     {
01444       fParticleGun->GetAngDist()->SetMinTheta(minthetaCmd1->GetNewDoubleValue(newValues));
01445     }
01446   else if(command == minphiCmd1)
01447     {
01448       fParticleGun->GetAngDist()->SetMinPhi(minphiCmd1->GetNewDoubleValue(newValues));
01449     }
01450   else if(command == maxthetaCmd1)
01451     {
01452       fParticleGun->GetAngDist()->SetMaxTheta(maxthetaCmd1->GetNewDoubleValue(newValues));
01453     }
01454   else if(command == maxphiCmd1)
01455     {
01456       fParticleGun->GetAngDist()->SetMaxPhi(maxphiCmd1->GetNewDoubleValue(newValues));
01457     }
01458   else if(command == angsigmarCmd1)
01459     {
01460       fParticleGun->GetAngDist()->SetBeamSigmaInAngR(angsigmarCmd1->GetNewDoubleValue(newValues));
01461     }
01462   else if(command == angsigmaxCmd1)
01463     {
01464       fParticleGun->GetAngDist()->SetBeamSigmaInAngX(angsigmaxCmd1->GetNewDoubleValue(newValues));
01465     }
01466   else if(command == angsigmayCmd1)
01467     {
01468       fParticleGun->GetAngDist()->SetBeamSigmaInAngY(angsigmayCmd1->GetNewDoubleValue(newValues));
01469     }
01470   else if(command == angfocusCmd)
01471     {
01472       fParticleGun->GetAngDist()->SetFocusPoint(angfocusCmd->GetNew3VectorValue(newValues));
01473     }
01474   else if(command == useuserangaxisCmd1)
01475     {
01476       fParticleGun->GetAngDist()->SetUseUserAngAxis(useuserangaxisCmd1->GetNewBoolValue(newValues));
01477     }
01478   else if(command == surfnormCmd1)
01479     {
01480       fParticleGun->GetAngDist()->SetUserWRTSurface(surfnormCmd1->GetNewBoolValue(newValues));
01481     }
01482   else if(command == energytypeCmd1)
01483     {
01484       fParticleGun->GetEneDist()->SetEnergyDisType(newValues);
01485     }
01486   else if(command == eminCmd1)
01487     {
01488       fParticleGun->GetEneDist()->SetEmin(eminCmd1->GetNewDoubleValue(newValues));
01489     }
01490   else if(command == emaxCmd1)
01491     {
01492       fParticleGun->GetEneDist()->SetEmax(emaxCmd1->GetNewDoubleValue(newValues));
01493     }
01494   else if(command == monoenergyCmd1)
01495     {
01496       fParticleGun->GetEneDist()->SetMonoEnergy(monoenergyCmd1->GetNewDoubleValue(newValues));
01497     }
01498   else if(command == engsigmaCmd1)
01499     {
01500       fParticleGun->GetEneDist()->SetBeamSigmaInE(engsigmaCmd1->GetNewDoubleValue(newValues));
01501     }
01502   else if(command == alphaCmd1)
01503     {
01504       fParticleGun->GetEneDist()->SetAlpha(alphaCmd1->GetNewDoubleValue(newValues));
01505     }
01506   else if(command == tempCmd1)
01507     {
01508       fParticleGun->GetEneDist()->SetTemp(tempCmd1->GetNewDoubleValue(newValues));
01509     }
01510   else if(command == ezeroCmd1)
01511     {
01512       fParticleGun->GetEneDist()->SetEzero(ezeroCmd1->GetNewDoubleValue(newValues));
01513     }
01514   else if(command == gradientCmd1)
01515     {
01516       fParticleGun->GetEneDist()->SetGradient(gradientCmd1->GetNewDoubleValue(newValues));
01517     }
01518   else if(command == interceptCmd1)
01519     {
01520       fParticleGun->GetEneDist()->SetInterCept(interceptCmd1->GetNewDoubleValue(newValues));
01521     }
01522   else if(command == arbeintCmd1)
01523     {
01524       fParticleGun->GetEneDist()->SetBiasAlpha(arbeintCmd1->GetNewDoubleValue(newValues));
01525     }
01526   else if(command == calculateCmd1)
01527     {
01528       fParticleGun->GetEneDist()->Calculate();
01529     }
01530   else if(command == energyspecCmd1)
01531     {
01532       fParticleGun->GetEneDist()->InputEnergySpectra(energyspecCmd1->GetNewBoolValue(newValues));
01533     }
01534   else if(command == diffspecCmd1)
01535     {
01536       fParticleGun->GetEneDist()->InputDifferentialSpectra(diffspecCmd1->GetNewBoolValue(newValues));
01537     }
01538   else if(command == histnameCmd1)
01539     {
01540       histtype = newValues;
01541     }
01542   else if(command == histfileCmd1)
01543     {
01544       histtype = "arb";
01545       fParticleGun->GetEneDist()->ArbEnergyHistoFile(newValues);
01546     }
01547   else if(command == histpointCmd1)
01548     {
01549       if(histtype == "biasx")
01550         fParticleGun->GetBiasRndm()->SetXBias(histpointCmd1->GetNew3VectorValue(newValues));
01551       if(histtype == "biasy")
01552         fParticleGun->GetBiasRndm()->SetYBias(histpointCmd1->GetNew3VectorValue(newValues));
01553       if(histtype == "biasz")
01554         fParticleGun->GetBiasRndm()->SetZBias(histpointCmd1->GetNew3VectorValue(newValues));
01555       if(histtype == "biast")
01556         fParticleGun->GetBiasRndm()->SetThetaBias(histpointCmd1->GetNew3VectorValue(newValues));
01557       if(histtype == "biasp")
01558         fParticleGun->GetBiasRndm()->SetPhiBias(histpointCmd1->GetNew3VectorValue(newValues));
01559       if(histtype == "biaspt")
01560         fParticleGun->GetBiasRndm()->SetPosThetaBias(histpointCmd1->GetNew3VectorValue(newValues));
01561       if(histtype == "biaspp")
01562         fParticleGun->GetBiasRndm()->SetPosPhiBias(histpointCmd1->GetNew3VectorValue(newValues));
01563       if(histtype == "biase")
01564         fParticleGun->GetBiasRndm()->SetEnergyBias(histpointCmd1->GetNew3VectorValue(newValues));
01565       if(histtype == "theta")
01566         fParticleGun->GetAngDist()->UserDefAngTheta(histpointCmd1->GetNew3VectorValue(newValues));
01567       if(histtype == "phi")
01568         fParticleGun->GetAngDist()->UserDefAngPhi(histpointCmd1->GetNew3VectorValue(newValues));
01569       if(histtype == "energy")
01570         fParticleGun->GetEneDist()->UserEnergyHisto(histpointCmd1->GetNew3VectorValue(newValues));
01571       if(histtype == "arb")
01572         fParticleGun->GetEneDist()->ArbEnergyHisto(histpointCmd1->GetNew3VectorValue(newValues));
01573       if(histtype == "epn")
01574         fParticleGun->GetEneDist()->EpnEnergyHisto(histpointCmd1->GetNew3VectorValue(newValues));
01575     }
01576   else if(command == resethistCmd1)
01577     {
01578       if(newValues == "theta" || newValues == "phi") {
01579         fParticleGun->GetAngDist()->ReSetHist(newValues);
01580       } else if (newValues == "energy" || newValues == "arb" || newValues == "epn") {
01581         fParticleGun->GetEneDist()->ReSetHist(newValues);
01582       } else {
01583         fParticleGun->GetBiasRndm()->ReSetHist(newValues);
01584       }
01585     }
01586   else if(command == arbintCmd1)
01587     {
01588       fParticleGun->GetEneDist()->ArbInterpolate(newValues);
01589     }
01590   else 
01591     {
01592       G4cout << "Error entering command" << G4endl;
01593     }
01594 }

void G4GeneralParticleSourceMessenger::SetParticleGun ( G4SingleParticleSource fpg  )  [inline]

Definition at line 107 of file G4GeneralParticleSourceMessenger.hh.

Referenced by G4GeneralParticleSource::AddaSource(), and G4GeneralParticleSource::SetCurrentSourceto().

00107 { fParticleGun = fpg; } ;


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:52:04 2013 for Geant4 by  doxygen 1.4.7