Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4RadioactivationMessenger.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// File: G4RadioactivationMessenger.cc //
29// Author: D.H. Wright (SLAC) //
30// Date: 29 August 2017 //
31// Description: messenger class for biased version of G4RadioactiveDecay. //
32// Based on the code of F. Lei and P.R. Truscott. //
33// //
34////////////////////////////////////////////////////////////////////////////////
35
37#include "G4NuclearLevelData.hh"
38#include <sstream>
40
41
43 :theRadioactivationContainer(theRadioactivationContainer1)
44{
45 old_grdmDirectory = new G4UIdirectory("/grdm/"); // To be removed in G4 11.0
46 old_grdmDirectory->SetGuidance("Controls the biased version of radioactive decay");
47
48 rdmDirectory = new G4UIdirectory("/process/had/rdm/");
49 rdmDirectory->SetGuidance("Controls the biased version of radioactive decay");
50
51 // Command to turn on/off variance reduction options
52 old_analoguemcCmd = new G4UIcmdWithABool("/grdm/analogueMC",this); // To be removed in G4 11.0
53 old_analoguemcCmd->SetGuidance("false: variance reduction method; true: analogue method");
54 old_analoguemcCmd->SetParameterName("OldAnalogueMC",true);
55 old_analoguemcCmd->SetDefaultValue(true);
56
57 analoguemcCmd = new G4UIcmdWithABool("/process/had/rdm/analogueMC",this);
58 analoguemcCmd->SetGuidance("false: variance reduction method; true: analogue method");
59 analoguemcCmd->SetParameterName("AnalogueMC",true);
60 analoguemcCmd->SetDefaultValue(true);
61
62 // Command to use branching ratio biasing or not
63 old_brbiasCmd = new G4UIcmdWithABool("/grdm/BRbias",this); // To be removed in G4 11.0
64 old_brbiasCmd->SetGuidance("false: no biasing; true: all branches are treated as equal");
65 old_brbiasCmd->SetParameterName("OldBRBias",true);
66 old_brbiasCmd->SetDefaultValue(true);
67
68 brbiasCmd = new G4UIcmdWithABool("/process/had/rdm/BRbias",this);
69 brbiasCmd->SetGuidance("false: no biasing; true: all branches are treated as equal");
70 brbiasCmd->SetParameterName("BRBias",true);
71 brbiasCmd->SetDefaultValue(true);
72
73 // Command to set the half-life thresold for isomer production
74 old_hlthCmd = new G4UIcmdWithADoubleAndUnit("/grdm/hlThreshold",this); // To be removed in G4 11.0
75 old_hlthCmd->SetGuidance("Set the h-l threshold for isomer production");
76 old_hlthCmd->SetParameterName("OldhlThreshold",false);
77 old_hlthCmd->SetUnitCategory("Time");
78
79 hlthCmd = new G4UIcmdWithADoubleAndUnit("/process/had/rdm/hlThreshold",this);
80 hlthCmd->SetGuidance("Set the h-l threshold for isomer production");
81 hlthCmd->SetParameterName("hlThreshold",false);
82 hlthCmd->SetUnitCategory("Time");
83
84 // Command to define the incident particle source time profile
85 old_sourcetimeprofileCmd = new G4UIcmdWithAString("/grdm/sourceTimeProfile",this); // To be removed in G4 11.0
86 old_sourcetimeprofileCmd->SetGuidance
87 ("Supply the name of the ascii file containing the source particle time profile");
88 old_sourcetimeprofileCmd->SetParameterName("OldSTimeProfile",true);
89 old_sourcetimeprofileCmd->SetDefaultValue("source.data");
90
91 sourcetimeprofileCmd = new G4UIcmdWithAString("/process/had/rdm/sourceTimeProfile",this);
92 sourcetimeprofileCmd->SetGuidance
93 ("Supply the name of the ascii file containing the source particle time profile");
94 sourcetimeprofileCmd->SetParameterName("STimeProfile",true);
95 sourcetimeprofileCmd->SetDefaultValue("source.data");
96
97 // Command to define the incident particle source time profile
98 old_decaybiasprofileCmd = new G4UIcmdWithAString("/grdm/decayBiasProfile",this); // To be removed in G4 11.0
99 old_decaybiasprofileCmd->SetGuidance
100 ("Supply the name of the ascii file containing the decay bias time profile");
101 old_decaybiasprofileCmd->SetParameterName("OldDBiasProfile",true);
102 old_decaybiasprofileCmd->SetDefaultValue("bias.data");
103
104 decaybiasprofileCmd = new G4UIcmdWithAString("/process/had/rdm/decayBiasProfile",this);
105 decaybiasprofileCmd->SetGuidance
106 ("Supply the name of the ascii file containing the decay bias time profile");
107 decaybiasprofileCmd->SetParameterName("DBiasProfile",true);
108 decaybiasprofileCmd->SetDefaultValue("bias.data");
109
110 // Command to set nuclei splitting parameter
111 old_splitnucleiCmd = new G4UIcmdWithAnInteger("/grdm/splitNuclei",this); // To be removed in G4 11.0
112 old_splitnucleiCmd->SetGuidance("Set number of splitting for the isotopes.");
113 old_splitnucleiCmd->SetParameterName("OldNSplit",true);
114 old_splitnucleiCmd->SetDefaultValue(1);
115 old_splitnucleiCmd->SetRange("OldNSplit>=1");
116
117 splitnucleiCmd = new G4UIcmdWithAnInteger("/process/had/rdm/splitNuclei",this);
118 splitnucleiCmd->SetGuidance("Set number of splitting for the isotopes.");
119 splitnucleiCmd->SetParameterName("NSplit",true);
120 splitnucleiCmd->SetDefaultValue(1);
121 splitnucleiCmd->SetRange("NSplit>=1");
122}
123
124
126{
127 delete old_grdmDirectory; // To be removed in G4 11.0
128 delete rdmDirectory;
129 delete old_analoguemcCmd; // To be removed in G4 11.0
130 delete analoguemcCmd;
131 delete old_sourcetimeprofileCmd; // To be removed in G4 11.0
132 delete sourcetimeprofileCmd;
133 delete old_decaybiasprofileCmd; // To be removed in G4 11.0
134 delete decaybiasprofileCmd;
135 delete old_brbiasCmd; // To be removed in G4 11.0
136 delete brbiasCmd;
137 delete old_splitnucleiCmd; // To be removed in G4 11.0
138 delete splitnucleiCmd;
139 delete old_hlthCmd; // To be removed in G4 11.0
140 delete hlthCmd;
141}
142
143
145{
146 // Old commands to be removed in G4 11.0
147 if ( command == old_analoguemcCmd ) { theRadioactivationContainer->
148 SetAnalogueMonteCarlo( old_analoguemcCmd->GetNewBoolValue( newValues ) );
150 ed << "This command is valid but deprecated and will be replaced with the command:\n"
151 << "/process/had/rdm/analogueMC in the next major release, Geant4 version 11.0";
152 G4Exception( "G4RadioactivationMessenger", "HAD_RDM_871", JustWarning, ed );
153 } else if ( command == old_brbiasCmd ) { theRadioactivationContainer->
154 SetBRBias( old_brbiasCmd->GetNewBoolValue( newValues ) );
156 ed << "This command is valid but deprecated and will be replaced with the command:\n"
157 << "/process/had/rdm/BRbias in the next major release, Geant4 version 11.0";
158 G4Exception( "G4RadioactivationMessenger", "HAD_RDM_872", JustWarning, ed );
159 } else if ( command == old_sourcetimeprofileCmd ) { theRadioactivationContainer->
160 SetSourceTimeProfile( newValues );
162 ed << "This command is valid but deprecated and will be replaced with the command:\n"
163 << "/process/had/rdm/sourceTimeProfile in the next major release, Geant4 version 11.0";
164 G4Exception( "G4RadioactivationMessenger", "HAD_RDM_873", JustWarning, ed );
165 } else if ( command == old_decaybiasprofileCmd ) { theRadioactivationContainer->
166 SetDecayBias( newValues );
168 ed << "This command is valid but deprecated and will be replaced with the command:\n"
169 << "/process/had/rdm/decayBiasProfile in the next major release, Geant4 version 11.0";
170 G4Exception( "G4RadioactivationMessenger", "HAD_RDM_874", JustWarning, ed );
171 } else if ( command == old_splitnucleiCmd ) { theRadioactivationContainer->
172 SetSplitNuclei( old_splitnucleiCmd->GetNewIntValue( newValues ) );
174 ed << "This command is valid but deprecated and will be replaced with the command:\n"
175 << "/process/had/rdm/splitNuclei in the next major release, Geant4 version 11.0";
176 G4Exception( "G4RadioactivationMessenger", "HAD_RDM_875", JustWarning, ed );
177 } else if ( command == old_hlthCmd ) { theRadioactivationContainer->
178 SetHLThreshold( old_hlthCmd->GetNewDoubleValue( newValues ) );
180 ed << "This command is valid but deprecated and will be replaced with the command:\n"
181 << "/process/had/rdm/hlThreshold in the next major release, Geant4 version 11.0";
182 G4Exception( "G4RadioactivationMessenger", "HAD_RDM_876", JustWarning, ed );
183 }
184
185 // New commands
186 if ( command == analoguemcCmd ) { theRadioactivationContainer->
187 SetAnalogueMonteCarlo( analoguemcCmd->GetNewBoolValue( newValues ) );
188 } else if ( command == brbiasCmd ) { theRadioactivationContainer->
189 SetBRBias( brbiasCmd->GetNewBoolValue( newValues ) );
190 } else if ( command == sourcetimeprofileCmd ) { theRadioactivationContainer->
191 SetSourceTimeProfile( newValues );
192 } else if ( command == decaybiasprofileCmd ) { theRadioactivationContainer->
193 SetDecayBias( newValues );
194 } else if ( command == splitnucleiCmd ) { theRadioactivationContainer->
195 SetSplitNuclei( splitnucleiCmd->GetNewIntValue( newValues ) );
196 } else if ( command == hlthCmd ) { theRadioactivationContainer->
197 SetHLThreshold( hlthCmd->GetNewDoubleValue( newValues ) );
198 }
199}
200
@ 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
G4RadioactivationMessenger(G4Radioactivation *theRadioactivationContainer)
void SetNewValue(G4UIcommand *command, G4String newValues)
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 SetDefaultValue(const char *defVal)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
static G4int GetNewIntValue(const char *paramString)
void SetDefaultValue(G4int defVal)
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:156
void SetRange(const char *rs)
Definition: G4UIcommand.hh:120