Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4RadioactiveDecayMessenger Class Reference

#include <G4RadioactiveDecayMessenger.hh>

+ Inheritance diagram for G4RadioactiveDecayMessenger:

Public Member Functions

 G4RadioactiveDecayMessenger (G4RadioactiveDecay *theRadioactiveDecayContainer)
 
 ~G4RadioactiveDecayMessenger ()
 
void SetNewValue (G4UIcommand *command, G4String newValues)
 
- Public Member Functions inherited from G4UImessenger
 G4UImessenger ()=default
 
 G4UImessenger (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
virtual ~G4UImessenger ()
 
virtual G4String GetCurrentValue (G4UIcommand *command)
 
G4bool CommandsShouldBeInMaster () const
 

Additional Inherited Members

- Protected Member Functions inherited from G4UImessenger
G4String ItoS (G4int i)
 
G4String LtoS (G4long l)
 
G4String DtoS (G4double a)
 
G4String BtoS (G4bool b)
 
G4int StoI (const G4String &s)
 
G4long StoL (const G4String &s)
 
G4double StoD (const G4String &s)
 
G4bool StoB (G4String s)
 
void AddUIcommand (G4UIcommand *newCommand)
 
void CreateDirectory (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
template<typename T >
T * CreateCommand (const G4String &cname, const G4String &dsc)
 
- Protected Attributes inherited from G4UImessenger
G4UIdirectorybaseDir = nullptr
 
G4String baseDirName = ""
 
G4bool commandsShouldBeInMaster = false
 

Detailed Description

Definition at line 57 of file G4RadioactiveDecayMessenger.hh.

Constructor & Destructor Documentation

◆ G4RadioactiveDecayMessenger()

G4RadioactiveDecayMessenger::G4RadioactiveDecayMessenger ( G4RadioactiveDecay * theRadioactiveDecayContainer)

Definition at line 44 of file G4RadioactiveDecayMessenger.cc.

46:theRadioactiveDecayContainer(theRadioactiveDecayContainer1)
47{
48 rdmDirectory = new G4UIdirectory("/process/had/rdm/");
49 rdmDirectory->SetGuidance("Controls for the Radioactive Decay Module.");
50
51 // Command to define the limits on nucleus the RDM will treat.
52 nucleuslimitsCmd = new G4UIcmdWithNucleusLimits("/process/had/rdm/nucleusLimits",this);
53 nucleuslimitsCmd->SetGuidance("Set the atomic weight and number limits for the RDM.");
54 nucleuslimitsCmd->SetParameterName("AMin","AMax","ZMin","ZMax",true);
55
56 // Select a logical volume for RDM
57 avolumeCmd = new G4UIcmdWithAString("/process/had/rdm/selectVolume",this);
58 avolumeCmd->SetGuidance("Supply a logical volumes name to add it to the RDM apply list");
59 avolumeCmd->SetParameterName("AVolume",false);
60
61 // De-select a logical volume for RDM
62 deavolumeCmd = new G4UIcmdWithAString("/process/had/rdm/deselectVolume",this);
63 deavolumeCmd->SetGuidance("Supply a logical volumes name to remove it from the RDM apply list");
64 deavolumeCmd->SetParameterName("AVolume",false);
65
66 // Select all logical volumes for RDM
67 allvolumesCmd = new G4UIcmdWithoutParameter("/process/had/rdm/allVolumes",this);
68 allvolumesCmd->SetGuidance(" apply RDM to all logical volumes. No parameter required.");
69 // allvolumeCmd->SetParameterName("AddAVolume",true);
70
71 // De-select all logical volumes for RDM
72 deallvolumesCmd = new G4UIcmdWithoutParameter("/process/had/rdm/noVolumes",this);
73 deallvolumesCmd->SetGuidance(" RDM is not applied to any logical volumes");
74 // deallvolumesCmd->SetParameterName("RemoveAVolume",true);
75
76 // Command to invoke atomic relaxation or not
77 armCmd = new G4UIcmdWithABool("/process/had/rdm/applyARM",this);
78 armCmd->SetGuidance("True: ARM is applied; false: no");
79 armCmd->SetParameterName("ApplyARM",true);
80 armCmd->SetDefaultValue(true);
81 //armCmd->AvailableForStates(G4State_PreInit);
82
83 // Command to set the directional bias (collimation) vector
84 colldirCmd = new G4UIcmdWith3Vector("/process/had/rdm/decayDirection",this);
85 colldirCmd->SetGuidance("Supply the direction vector for decay products");
86 colldirCmd->SetParameterName("X","Y","Z",false);
87
88 // Command to set the directional bias (collimation) half angle ("cone")
89 collangleCmd = new G4UIcmdWithADoubleAndUnit("/process/had/rdm/decayHalfAngle",this);
90 collangleCmd->SetGuidance("Supply maximum angle from direction vector for decay products");
91 collangleCmd->SetParameterName("HalfAngle",false);
92 collangleCmd->SetUnitCategory("Angle");
93
94 // This command setup the verbose level of radioactive decay
95 verboseCmd = new G4UIcmdWithAnInteger("/process/had/rdm/verbose",this);
96 verboseCmd->SetGuidance("Set verbose level: 0, 1, 2 or 3");
97 verboseCmd->SetParameterName("VerboseLevel",true);
98 verboseCmd->SetDefaultValue(1);
99 verboseCmd->SetRange("VerboseLevel>=0");
100
101 // Use a user-defined decay datafile for a given isotope
102 userDecayDataCmd = new G4UIcommand("/process/had/rdm/setRadioactiveDecayFile",this);
103 userDecayDataCmd->SetGuidance("Supply user-defined radioactive decay data file");
104 G4UIparameter* Z_para= new G4UIparameter("Z_isotope",'i',true);
105 Z_para->SetParameterRange("Z_isotope > 0");
106 Z_para->SetGuidance("Z: Charge number of isotope");
107 G4UIparameter* A_para= new G4UIparameter("A_isotope",'i',true);
108 A_para->SetParameterRange("A_isotope > 1");
109 A_para->SetGuidance("A: mass number of isotope");
110 G4UIparameter* FileName_para= new G4UIparameter("file_name",'s',true);
111 FileName_para->SetGuidance("Name of the user data file");
112 userDecayDataCmd->SetParameter(Z_para);
113 userDecayDataCmd->SetParameter(A_para);
114 userDecayDataCmd->SetParameter(FileName_para);
115
116 // Use a user-defined evaporation data file for a given isotope
117 userEvaporationDataCmd = new G4UIcommand("/process/had/rdm/setPhotoEvaporationFile",this);
118 userEvaporationDataCmd->SetGuidance("Supply user-defined photon evaporation data file");
119 G4UIparameter* Zpara= new G4UIparameter("Z_isotope",'i',true);
120 Zpara->SetParameterRange("Z_isotope > 0");
121 Zpara->SetGuidance("Z: Charge number of isotope");
122 G4UIparameter* Apara= new G4UIparameter("A_isotope",'i',true);
123 Apara->SetParameterRange("A_isotope > 1");
124 Apara->SetGuidance("A: mass number of isotope");
125 G4UIparameter* FileNamepara= new G4UIparameter("file_name",'s',true);
126 FileNamepara->SetGuidance("Name of the user data file");
127 userEvaporationDataCmd->SetParameter(Zpara);
128 userEvaporationDataCmd->SetParameter(Apara);
129 userEvaporationDataCmd->SetParameter(FileNamepara);
130
131 // Command to set the threshold for very long decay time (i.e. radioactive
132 // decays of nuclides at rest happening later than this threshold are ignored)
133 thresholdForVeryLongDecayTimeCmd = new G4UIcmdWithADoubleAndUnit("/process/had/rdm/thresholdForVeryLongDecayTime",this);
134 thresholdForVeryLongDecayTimeCmd->SetGuidance("Ignore decays at rest of nuclides happening after this time threshold");
135 thresholdForVeryLongDecayTimeCmd->SetParameterName("ThresholdForVeryLongDecayTime",false);
136 thresholdForVeryLongDecayTimeCmd->SetUnitCategory("Time");
137}
void SetParameterName(const char *theNameX, const char *theNameY, const char *theNameZ, G4bool omittable, G4bool currentAsDefault=false)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetDefaultValue(G4bool defVal)
void SetUnitCategory(const char *unitCategory)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetDefaultValue(G4int defVal)
void SetParameterName(const char *theNameAMin, const char *theNameAMax, const char *theNameZMin, const char *theNameZMax, G4bool omittable, G4bool currentAsDefault=true)
void SetParameter(G4UIparameter *const newParameter)
void SetGuidance(const char *aGuidance)
void SetRange(const char *rs)
void SetGuidance(const char *theGuidance)
void SetParameterRange(const char *theRange)

◆ ~G4RadioactiveDecayMessenger()

G4RadioactiveDecayMessenger::~G4RadioactiveDecayMessenger ( )

Definition at line 140 of file G4RadioactiveDecayMessenger.cc.

141{
142 delete rdmDirectory;
143 delete nucleuslimitsCmd;
144 delete verboseCmd;
145 delete avolumeCmd;
146 delete deavolumeCmd;
147 delete allvolumesCmd;
148 delete deallvolumesCmd;
149 delete armCmd;
150 delete userDecayDataCmd;
151 delete userEvaporationDataCmd;
152 delete colldirCmd;
153 delete collangleCmd;
154 delete thresholdForVeryLongDecayTimeCmd;
155}

Member Function Documentation

◆ SetNewValue()

void G4RadioactiveDecayMessenger::SetNewValue ( G4UIcommand * command,
G4String newValues )
virtual

Reimplemented from G4UImessenger.

Definition at line 159 of file G4RadioactiveDecayMessenger.cc.

160{
161 if ( command == nucleuslimitsCmd ) {
162 theRadioactiveDecayContainer->
163 SetNucleusLimits( nucleuslimitsCmd->GetNewNucleusLimitsValue( newValues ) );
164 } else if ( command == avolumeCmd ) {
165 theRadioactiveDecayContainer->SelectAVolume( newValues );
166 } else if ( command == deavolumeCmd ) {
167 theRadioactiveDecayContainer->DeselectAVolume( newValues );
168 } else if ( command == allvolumesCmd ) {
169 theRadioactiveDecayContainer->SelectAllVolumes();
170 } else if ( command == deallvolumesCmd ) {
171 theRadioactiveDecayContainer->DeselectAllVolumes();
172 } else if (command == verboseCmd) {
173 theRadioactiveDecayContainer->SetVerboseLevel(verboseCmd->GetNewIntValue(newValues) );
174 } else if (command == armCmd) {
175 theRadioactiveDecayContainer->SetARM(armCmd->GetNewBoolValue(newValues) );
176 } else if ( command == userDecayDataCmd ) {
177 G4int Z,A;
178 G4String file_name;
179 const char* nv = (const char*)newValues;
180 std::istringstream is(nv);
181 is >> Z >> A >> file_name;
182 theRadioactiveDecayContainer->AddUserDecayDataFile(Z,A,file_name);
183 } else if ( command == userEvaporationDataCmd ) {
184 G4int Z,A;
185 G4String file_name;
186 const char* nv = (const char*)newValues;
187 std::istringstream is(nv);
188 is >> Z >> A >> file_name;
190 } else if ( command == colldirCmd ) {
191 theRadioactiveDecayContainer->SetDecayDirection( colldirCmd->GetNew3VectorValue( newValues ) );
192 } else if ( command == collangleCmd ) {
193 theRadioactiveDecayContainer->SetDecayHalfAngle( collangleCmd->GetNewDoubleValue( newValues ) );
194 } else if (command == thresholdForVeryLongDecayTimeCmd) {
195 theRadioactiveDecayContainer->SetThresholdForVeryLongDecayTime(thresholdForVeryLongDecayTimeCmd->GetNewDoubleValue(newValues) );
196 }
197}
int G4int
Definition G4Types.hh:85
const G4double A[17]
G4bool AddPrivateData(G4int Z, G4int A, const G4String &filename)
static G4NuclearLevelData * GetInstance()
void AddUserDecayDataFile(G4int Z, G4int A, const G4String &filename)
void SelectAVolume(const G4String &aVolume)
void SetDecayHalfAngle(G4double halfAngle=0.*CLHEP::deg)
void SetDecayDirection(const G4ThreeVector &theDir)
void SetARM(G4bool arm)
void SetThresholdForVeryLongDecayTime(const G4double inputThreshold)
void DeselectAVolume(const G4String &aVolume)
static G4ThreeVector GetNew3VectorValue(const char *paramString)
static G4bool GetNewBoolValue(const char *paramString)
static G4double GetNewDoubleValue(const char *paramString)
static G4int GetNewIntValue(const char *paramString)
G4NucleusLimits GetNewNucleusLimitsValue(G4String paramString)
void SetVerboseLevel(G4int value)

The documentation for this class was generated from the following files: