Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4RadioactiveDecayBaseMessenger.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26
27////////////////////////////////////////////////////////////////////////////////
28// //
29// File: G4RadioactiveDecayBaseMessenger.cc //
30// Author: D.H. Wright (SLAC) //
31// Date: 29 August 2017 //
32// Description: messenger class for the non-biased version of //
33// G4RadioactiveDecay. Based on the code of F. Lei and //
34// P.R. Truscott. //
35// //
36////////////////////////////////////////////////////////////////////////////////
37
39#include "G4NuclearLevelData.hh"
40#include <sstream>
42
43
45(G4RadioactiveDecayBase* theRadioactiveDecayContainer1)
46:theRadioactiveDecayContainer(theRadioactiveDecayContainer1)
47{
48 // main directory for control of the RDM
49 old_grdmDirectory = new G4UIdirectory("/grdm/"); // To be removed in G4 11.0
50 old_grdmDirectory->SetGuidance("Controls for the Radioactive Decay Module.");
51
52 rdmDirectory = new G4UIdirectory("/process/had/rdm/");
53 rdmDirectory->SetGuidance("Controls for the Radioactive Decay Module.");
54
55 // Command to define the limits on nucleus the RDM will treat.
56 old_nucleuslimitsCmd = new G4UIcmdWithNucleusLimits("/grdm/nucleusLimits",this); // To be removed in G4 11.0
57 old_nucleuslimitsCmd->SetGuidance("Set the atomic weight and number limits for the RDM.");
58 old_nucleuslimitsCmd->SetParameterName("OldAMin","OldAMax","OldZMin","OldZMax",true);
59
60 nucleuslimitsCmd = new G4UIcmdWithNucleusLimits("/process/had/rdm/nucleusLimits",this);
61 nucleuslimitsCmd->SetGuidance("Set the atomic weight and number limits for the RDM.");
62 nucleuslimitsCmd->SetParameterName("AMin","AMax","ZMin","ZMax",true);
63
64 // Select a logical volume for RDM
65 old_avolumeCmd = new G4UIcmdWithAString("/grdm/selectVolume",this); // To be removed in G4 11.0
66 old_avolumeCmd->SetGuidance("Supply a logical volumes name to add it to the RDM apply list");
67 old_avolumeCmd->SetParameterName("OldAVolume",false);
68
69 avolumeCmd = new G4UIcmdWithAString("/process/had/rdm/selectVolume",this);
70 avolumeCmd->SetGuidance("Supply a logical volumes name to add it to the RDM apply list");
71 avolumeCmd->SetParameterName("AVolume",false);
72
73 // De-select a logical volume for RDM
74 old_deavolumeCmd = new G4UIcmdWithAString("/grdm/deselectVolume",this); // To be removed in G4 11.0
75 old_deavolumeCmd->SetGuidance("Supply a logical volumes name to remove it from the RDM apply list");
76 old_deavolumeCmd->SetParameterName("OldAVolume",false);
77
78 deavolumeCmd = new G4UIcmdWithAString("/process/had/rdm/deselectVolume",this);
79 deavolumeCmd->SetGuidance("Supply a logical volumes name to remove it from the RDM apply list");
80 deavolumeCmd->SetParameterName("AVolume",false);
81
82 // Select all logical volumes for RDM
83 old_allvolumesCmd = new G4UIcmdWithoutParameter("/grdm/allVolumes",this); // To be removed in G4 11.0
84 old_allvolumesCmd->SetGuidance(" apply RDM to all logical volumes. No parameter required.");
85 // old_allvolumeCmd->SetParameterName("OldAddAVolume",true);
86
87 allvolumesCmd = new G4UIcmdWithoutParameter("/process/had/rdm/allVolumes",this);
88 allvolumesCmd->SetGuidance(" apply RDM to all logical volumes. No parameter required.");
89 // allvolumeCmd->SetParameterName("AddAVolume",true);
90
91 // De-select all logical volumes for RDM
92 old_deallvolumesCmd = new G4UIcmdWithoutParameter("/grdm/noVolumes",this); // To be removed in G4 11.0
93 old_deallvolumesCmd->SetGuidance(" RDM is not applied to any logical volumes");
94 // old_deallvolumesCmd->SetParameterName("OldRemoveAVolume",true);
95
96 deallvolumesCmd = new G4UIcmdWithoutParameter("/process/had/rdm/noVolumes",this);
97 deallvolumesCmd->SetGuidance(" RDM is not applied to any logical volumes");
98 // deallvolumesCmd->SetParameterName("RemoveAVolume",true);
99
100 // Command to invoke internal conversion or not
101 old_icmCmd = new G4UIcmdWithABool("/grdm/applyICM",this); // To be removed in G4 11.0
102 old_icmCmd->SetGuidance("Command not active; kept for backward compatibility.");
103 old_icmCmd->SetGuidance("Internal conversion is always turned on.");
104 old_icmCmd->SetParameterName("OldApplyICM",true);
105 old_icmCmd->SetDefaultValue(true);
106
107 icmCmd = new G4UIcmdWithABool("/process/had/rdm/applyICM",this);
108 icmCmd->SetGuidance("Command not active; kept for backward compatibility.");
109 icmCmd->SetGuidance("Internal conversion is always turned on.");
110 icmCmd->SetParameterName("ApplyICM",true);
111 icmCmd->SetDefaultValue(true);
112
113 // Command to invoke atomic relaxation or not
114 old_armCmd = new G4UIcmdWithABool("/grdm/applyARM",this); // To be removed in G4 11.0
115 old_armCmd->SetGuidance("True: ARM is applied; false: no");
116 old_armCmd->SetParameterName("OldApplyARM",true);
117 old_armCmd->SetDefaultValue(true);
118 //old_armCmd->AvailableForStates(G4State_PreInit);
119
120 armCmd = new G4UIcmdWithABool("/process/had/rdm/applyARM",this);
121 armCmd->SetGuidance("True: ARM is applied; false: no");
122 armCmd->SetParameterName("ApplyARM",true);
123 armCmd->SetDefaultValue(true);
124 //armCmd->AvailableForStates(G4State_PreInit);
125
126 // Command to set the directional bias (collimation) vector
127 old_colldirCmd = new G4UIcmdWith3Vector("/grdm/decayDirection",this); // To be removed in G4 11.0
128 old_colldirCmd->SetGuidance("Supply the direction vector for decay products");
129 old_colldirCmd->SetParameterName("OldX","OldY","OldZ",false);
130
131 colldirCmd = new G4UIcmdWith3Vector("/process/had/rdm/decayDirection",this);
132 colldirCmd->SetGuidance("Supply the direction vector for decay products");
133 colldirCmd->SetParameterName("X","Y","Z",false);
134
135 // Command to set the directional bias (collimation) half angle ("cone")
136 old_collangleCmd = new G4UIcmdWithADoubleAndUnit("/grdm/decayHalfAngle",this); // To be removed in G4 11.0
137 old_collangleCmd->SetGuidance("Supply maximum angle from direction vector for decay products");
138 old_collangleCmd->SetParameterName("OldHalfAngle",false);
139 old_collangleCmd->SetUnitCategory("Angle");
140
141 collangleCmd = new G4UIcmdWithADoubleAndUnit("/process/had/rdm/decayHalfAngle",this);
142 collangleCmd->SetGuidance("Supply maximum angle from direction vector for decay products");
143 collangleCmd->SetParameterName("HalfAngle",false);
144 collangleCmd->SetUnitCategory("Angle");
145
146 // This command setup the verbose level of radioactive decay
147 old_verboseCmd = new G4UIcmdWithAnInteger("/grdm/verbose",this); // To be removed in G4 11.0
148 old_verboseCmd->SetGuidance("Set verbose level: 0, 1, 2 or 3");
149 old_verboseCmd->SetParameterName("OldVerboseLevel",true);
150 old_verboseCmd->SetDefaultValue(1);
151 old_verboseCmd->SetRange("OldVerboseLevel>=0");
152
153 verboseCmd = new G4UIcmdWithAnInteger("/process/had/rdm/verbose",this);
154 verboseCmd->SetGuidance("Set verbose level: 0, 1, 2 or 3");
155 verboseCmd->SetParameterName("VerboseLevel",true);
156 verboseCmd->SetDefaultValue(1);
157 verboseCmd->SetRange("VerboseLevel>=0");
158
159 // Use a user-defined decay datafile for a given isotope
160 old_userDecayDataCmd = new G4UIcommand("/grdm/setRadioactiveDecayFile",this); // To be removed in G4 11.0
161 old_userDecayDataCmd->SetGuidance("Supply user-defined radioactive decay data file");
162 G4UIparameter* old_Z_para= new G4UIparameter("Z_isotope",'i',true);
163 old_Z_para->SetParameterRange("Z_isotope > 0");
164 old_Z_para->SetGuidance("Z: Charge number of isotope");
165 G4UIparameter* old_A_para= new G4UIparameter("A_isotope",'i',true);
166 old_A_para->SetParameterRange("A_isotope > 1");
167 old_A_para->SetGuidance("A: mass number of isotope");
168 G4UIparameter* old_FileName_para= new G4UIparameter("file_name",'s',true);
169 old_FileName_para->SetGuidance("Name of the user data file");
170 old_userDecayDataCmd->SetParameter(old_Z_para);
171 old_userDecayDataCmd->SetParameter(old_A_para);
172 old_userDecayDataCmd->SetParameter(old_FileName_para);
173
174 userDecayDataCmd = new G4UIcommand("/process/had/rdm/setRadioactiveDecayFile",this);
175 userDecayDataCmd->SetGuidance("Supply user-defined radioactive decay data file");
176 G4UIparameter* Z_para= new G4UIparameter("Z_isotope",'i',true);
177 Z_para->SetParameterRange("Z_isotope > 0");
178 Z_para->SetGuidance("Z: Charge number of isotope");
179 G4UIparameter* A_para= new G4UIparameter("A_isotope",'i',true);
180 A_para->SetParameterRange("A_isotope > 1");
181 A_para->SetGuidance("A: mass number of isotope");
182 G4UIparameter* FileName_para= new G4UIparameter("file_name",'s',true);
183 FileName_para->SetGuidance("Name of the user data file");
184 userDecayDataCmd->SetParameter(Z_para);
185 userDecayDataCmd->SetParameter(A_para);
186 userDecayDataCmd->SetParameter(FileName_para);
187
188 // Use a user-defined evaporation data file for a given isotope
189 old_userEvaporationDataCmd = new G4UIcommand("/grdm/setPhotoEvaporationFile",this); // To be removed in G4 11.0
190 old_userEvaporationDataCmd->SetGuidance("Supply user-defined photon evaporation data file");
191 G4UIparameter* old_Zpara= new G4UIparameter("Z_isotope",'i',true);
192 old_Zpara->SetParameterRange("Z_isotope > 0");
193 old_Zpara->SetGuidance("Z: Charge number of isotope");
194 G4UIparameter* old_Apara= new G4UIparameter("A_isotope",'i',true);
195 old_Apara->SetParameterRange("A_isotope > 1");
196 old_Apara->SetGuidance("A: mass number of isotope");
197 G4UIparameter* old_FileNamepara= new G4UIparameter("file_name",'s',true);
198 old_FileNamepara->SetGuidance("Name of the user data file");
199 old_userEvaporationDataCmd->SetParameter(old_Zpara);
200 old_userEvaporationDataCmd->SetParameter(old_Apara);
201 old_userEvaporationDataCmd->SetParameter(old_FileNamepara);
202
203 userEvaporationDataCmd = new G4UIcommand("/process/had/rdm/setPhotoEvaporationFile",this);
204 userEvaporationDataCmd->SetGuidance("Supply user-defined photon evaporation data file");
205 G4UIparameter* Zpara= new G4UIparameter("Z_isotope",'i',true);
206 Zpara->SetParameterRange("Z_isotope > 0");
207 Zpara->SetGuidance("Z: Charge number of isotope");
208 G4UIparameter* Apara= new G4UIparameter("A_isotope",'i',true);
209 Apara->SetParameterRange("A_isotope > 1");
210 Apara->SetGuidance("A: mass number of isotope");
211 G4UIparameter* FileNamepara= new G4UIparameter("file_name",'s',true);
212 FileNamepara->SetGuidance("Name of the user data file");
213 userEvaporationDataCmd->SetParameter(Zpara);
214 userEvaporationDataCmd->SetParameter(Apara);
215 userEvaporationDataCmd->SetParameter(FileNamepara);
216}
217
218
220{
221 delete old_grdmDirectory; // To be removed in G4 11.0
222 delete rdmDirectory;
223 delete old_nucleuslimitsCmd; // To be removed in G4 11.0
224 delete nucleuslimitsCmd;
225 delete old_verboseCmd; // To be removed in G4 11.0
226 delete verboseCmd;
227 delete old_avolumeCmd; // To be removed in G4 11.0
228 delete avolumeCmd;
229 delete old_deavolumeCmd; // To be removed in G4 11.0
230 delete deavolumeCmd;
231 delete old_allvolumesCmd; // To be removed in G4 11.0
232 delete allvolumesCmd;
233 delete old_deallvolumesCmd; // To be removed in G4 11.0
234 delete deallvolumesCmd;
235 delete old_icmCmd; // To be removed in G4 11.0
236 delete icmCmd;
237 delete old_armCmd; // To be removed in G4 11.0
238 delete armCmd;
239 delete old_userDecayDataCmd; // To be removed in G4 11.0
240 delete userDecayDataCmd; //<---
241 delete old_userEvaporationDataCmd; // To be removed in G4 11.0
242 delete userEvaporationDataCmd;
243 delete old_colldirCmd; // To be removed in G4 11.0
244 delete colldirCmd;
245 delete old_collangleCmd; // To be removed in G4 11.0
246 delete collangleCmd;
247}
248
249
250void
252{
253 // Old commands to be removed in G4 11.0
254 if ( command == old_nucleuslimitsCmd ) {
255 theRadioactiveDecayContainer->
256 SetNucleusLimits( old_nucleuslimitsCmd->GetNewNucleusLimitsValue( newValues ) );
258 ed << "This command is valid but deprecated and will be replaced with the command:\n"
259 << "/process/had/rdm/nucleusLimits in the next major release, Geant4 version 11.0";
260 G4Exception( "G4RadioactiveDecayBaseMessenger", "HAD_RDM_671", JustWarning, ed );
261 } else if ( command == old_avolumeCmd ) {
262 theRadioactiveDecayContainer->SelectAVolume( newValues );
264 ed << "This command is valid but deprecated and will be replaced with the command:\n"
265 << "/process/had/rdm/selectVolume in the next major release, Geant4 version 11.0";
266 G4Exception( "G4RadioactiveDecayBaseMessenger", "HAD_RDM_672", JustWarning, ed );
267 } else if ( command == old_deavolumeCmd ) {
268 theRadioactiveDecayContainer->DeselectAVolume( newValues );
270 ed << "This command is valid but deprecated and will be replaced with the command:\n"
271 << "/process/had/rdm/deselectVolume in the next major release, Geant4 version 11.0";
272 G4Exception( "G4RadioactiveDecayBaseMessenger", "HAD_RDM_673", JustWarning, ed );
273 } else if ( command == old_allvolumesCmd ) {
274 theRadioactiveDecayContainer->SelectAllVolumes();
276 ed << "This command is valid but deprecated and will be replaced with the command:\n"
277 << "/process/had/rdm/allVolumes in the next major release, Geant4 version 11.0";
278 G4Exception( "G4RadioactiveDecayBaseMessenger", "HAD_RDM_674", JustWarning, ed );
279 } else if ( command == old_deallvolumesCmd ) {
280 theRadioactiveDecayContainer->DeselectAllVolumes();
282 ed << "This command is valid but deprecated and will be replaced with the command:\n"
283 << "/process/had/rdm/noVolumes in the next major release, Geant4 version 11.0";
284 G4Exception( "G4RadioactiveDecayBaseMessenger", "HAD_RDM_675", JustWarning, ed );
285 } else if ( command == old_verboseCmd ) {
286 theRadioactiveDecayContainer->SetVerboseLevel( old_verboseCmd->GetNewIntValue( newValues ) );
288 ed << "This command is valid but deprecated and will be replaced with the command:\n"
289 << "/process/had/rdm/verbose in the next major release, Geant4 version 11.0";
290 G4Exception( "G4RadioactiveDecayBaseMessenger", "HAD_RDM_676", JustWarning, ed );
291 } else if ( command == old_icmCmd ) {
292 theRadioactiveDecayContainer->SetICM( old_icmCmd->GetNewBoolValue( newValues ) );
294 ed << "This command is valid but deprecated and will be replaced with the command:\n"
295 << "/process/had/rdm/applyICM in the next major release, Geant4 version 11.0";
296 G4Exception( "G4RadioactiveDecayBaseMessenger", "HAD_RDM_677", JustWarning, ed );
297 } else if ( command == old_armCmd ) {
298 theRadioactiveDecayContainer->SetARM( old_armCmd->GetNewBoolValue( newValues ) );
300 ed << "This command is valid but deprecated and will be replaced with the command:\n"
301 << "/process/had/rdm/applyARM in the next major release, Geant4 version 11.0";
302 G4Exception( "G4RadioactiveDecayBaseMessenger", "HAD_RDM_678", JustWarning, ed );
303 } else if ( command == old_userDecayDataCmd ) {
304 G4int Z,A;
305 G4String file_name;
306 const char* nv = (const char*)newValues;
307 std::istringstream is(nv);
308 is >> Z >> A >> file_name;
309 theRadioactiveDecayContainer->AddUserDecayDataFile(Z,A,file_name);
311 ed << "This command is valid but deprecated and will be replaced with the command:\n"
312 << "/process/had/rdm/setRadioactiveDecayFile in the next major release, Geant4 version 11.0";
313 G4Exception( "G4RadioactiveDecayBaseMessenger", "HAD_RDM_679", JustWarning, ed );
314 } else if ( command == old_userEvaporationDataCmd ) {
315 G4int Z,A;
316 G4String file_name;
317 const char* nv = (const char*)newValues;
318 std::istringstream is(nv);
319 is >> Z >> A >> file_name;
322 ed << "This command is valid but deprecated and will be replaced with the command:\n"
323 << "/process/had/rdm/setPhotoEvaporationFile in the next major release, Geant4 version 11.0";
324 G4Exception( "G4RadioactiveDecayBaseMessenger", "HAD_RDM_680", JustWarning, ed );
325 } else if ( command == old_colldirCmd ) {
326 theRadioactiveDecayContainer->SetDecayDirection( old_colldirCmd->GetNew3VectorValue( newValues ) );
328 ed << "This command is valid but deprecated and will be replaced with the command:\n"
329 << "/process/had/rdm/decayDirection in the next major release, Geant4 version 11.0";
330 G4Exception( "G4RadioactiveDecayBaseMessenger", "HAD_RDM_681", JustWarning, ed );
331 } else if ( command == old_collangleCmd ) {
332 theRadioactiveDecayContainer->SetDecayHalfAngle( old_collangleCmd->GetNewDoubleValue( newValues ) );
334 ed << "This command is valid but deprecated and will be replaced with the command:\n"
335 << "/process/had/rdm/decayHalfAngle in the next major release, Geant4 version 11.0";
336 G4Exception( "G4RadioactiveDecayBaseMessenger", "HAD_RDM_682", JustWarning, ed );
337 }
338
339 // New commands
340 if ( command == nucleuslimitsCmd ) {
341 theRadioactiveDecayContainer->
342 SetNucleusLimits( nucleuslimitsCmd->GetNewNucleusLimitsValue( newValues ) );
343 } else if ( command == avolumeCmd ) {
344 theRadioactiveDecayContainer->SelectAVolume( newValues );
345 } else if ( command == deavolumeCmd ) {
346 theRadioactiveDecayContainer->DeselectAVolume( newValues );
347 } else if ( command == allvolumesCmd ) {
348 theRadioactiveDecayContainer->SelectAllVolumes();
349 } else if ( command == deallvolumesCmd ) {
350 theRadioactiveDecayContainer->DeselectAllVolumes();
351 } else if ( command == verboseCmd ) {
352 theRadioactiveDecayContainer->SetVerboseLevel( verboseCmd->GetNewIntValue( newValues ) );
353 } else if ( command == icmCmd ) {
354 theRadioactiveDecayContainer->SetICM( icmCmd->GetNewBoolValue( newValues ) );
355 } else if ( command == armCmd ) {
356 theRadioactiveDecayContainer->SetARM( armCmd->GetNewBoolValue( newValues ) );
357 } else if ( command == userDecayDataCmd ) {
358 G4int Z,A;
359 G4String file_name;
360 const char* nv = (const char*)newValues;
361 std::istringstream is(nv);
362 is >> Z >> A >> file_name;
363 theRadioactiveDecayContainer->AddUserDecayDataFile(Z,A,file_name);
364 } else if ( command == userEvaporationDataCmd ) {
365 G4int Z,A;
366 G4String file_name;
367 const char* nv = (const char*)newValues;
368 std::istringstream is(nv);
369 is >> Z >> A >> file_name;
371 } else if ( command == colldirCmd ) {
372 theRadioactiveDecayContainer->SetDecayDirection( colldirCmd->GetNew3VectorValue( newValues ) );
373 } else if ( command == collangleCmd ) {
374 theRadioactiveDecayContainer->SetDecayHalfAngle( collangleCmd->GetNewDoubleValue( newValues ) );
375 }
376}
377
double A(double temperature)
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
int G4int
Definition: G4Types.hh:85
G4bool AddPrivateData(G4int Z, G4int A, const G4String &filename)
static G4NuclearLevelData * GetInstance()
G4RadioactiveDecayBaseMessenger(G4RadioactiveDecayBase *theRadioactiveDecayContainer)
void SetNewValue(G4UIcommand *command, G4String newValues)
void AddUserDecayDataFile(G4int Z, G4int A, G4String filename)
void DeselectAVolume(const G4String aVolume)
void SelectAVolume(const G4String aVolume)
void SetDecayHalfAngle(G4double halfAngle=0.*CLHEP::deg)
void SetDecayDirection(const G4ThreeVector &theDir)
static G4ThreeVector GetNew3VectorValue(const char *paramString)
void SetParameterName(const char *theNameX, const char *theNameY, const char *theNameZ, G4bool omittable, G4bool currentAsDefault=false)
static G4bool GetNewBoolValue(const char *paramString)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetDefaultValue(G4bool defVal)
void SetUnitCategory(const char *unitCategory)
static G4double GetNewDoubleValue(const char *paramString)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
static G4int GetNewIntValue(const char *paramString)
void SetDefaultValue(G4int defVal)
G4NucleusLimits GetNewNucleusLimitsValue(G4String paramString)
void SetParameterName(const char *theNameAMin, const char *theNameAMax, const char *theNameZMin, const char *theNameZMax, G4bool omittable, G4bool currentAsDefault=true)
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:146
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:156
void SetRange(const char *rs)
Definition: G4UIcommand.hh:120
void SetGuidance(const char *theGuidance)
void SetParameterRange(const char *theRange)