Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4RadioactiveDecaymessenger.cc
이 파일의 문서화 페이지로 가기
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 
28 #include "G4NuclearLevelData.hh"
29 
30 #include <sstream>
31 
32 
34 (G4RadioactiveDecay* theRadioactiveDecayContainer1)
35 :theRadioactiveDecayContainer(theRadioactiveDecayContainer1)
36 {
37  // main directory for control of the RDM
38  grdmDirectory = new G4UIdirectory("/grdm/");
39  grdmDirectory->SetGuidance("Controls for the Radioactive Decay Module.");
40 
41  // Command to define the limits on nucleus the RDM will treat.
42  nucleuslimitsCmd = new
43  G4UIcmdWithNucleusLimits("/grdm/nucleusLimits",this);
44  nucleuslimitsCmd->SetGuidance
45  ("Set the atomic weight and number limits for the RDM.");
46  nucleuslimitsCmd->SetParameterName("aMin","aMax","zMin","zMax",true);
47 
48  // The next command contols whether the decay will be treated analoguely or
49  // with variance reduction
50  analoguemcCmd = new G4UIcmdWithABool ("/grdm/analogueMC",this);
51  analoguemcCmd->SetGuidance("false: variance reduction method; true: analogue method");
52  analoguemcCmd->SetParameterName("AnalogueMC",true);
53  analoguemcCmd->SetDefaultValue(true);
54  //
55  // The next command contols whether beta decay will be treated faithfully or
56  // in fast mode
57  //
58  fbetaCmd = new G4UIcmdWithABool ("/grdm/fBeta",this);
59  fbetaCmd->SetGuidance("false: use 3-body decay, true: use histogram method");
60  fbetaCmd->SetParameterName("fBeta",true);
61  fbetaCmd->SetDefaultValue(false);
62 
63  //
64  //
65  // Command to selete a logical volume for RDM.
66  //
67  avolumeCmd = new
68  G4UIcmdWithAString("/grdm/selectVolume",this);
69  avolumeCmd->SetGuidance
70  ("Suppply a logical volumes name to add it to the RDM apply list");
71  avolumeCmd->SetParameterName("aVolume",false);
72  //
73  //
74  //
75  // Command to de-selete a logical volume for RDM.
76  //
77  deavolumeCmd = new
78  G4UIcmdWithAString("/grdm/deselectVolume",this);
79  deavolumeCmd->SetGuidance
80  ("Suppply a logical volumes name to remove it from the RDM apply list");
81  deavolumeCmd->SetParameterName("aVolume",false);
82  //
83  //
84  // Command to selete all logical volumes for RDM.
85  //
86  allvolumesCmd = new
87  G4UIcmdWithoutParameter("/grdm/allVolumes",this);
88  allvolumesCmd->SetGuidance
89  (" apply RDM to all logical volumes. No parameter required.");
90  // allvolumeCmd->SetParameterName("AddAVolume",true);
91 
92  //
93  // Command to de-selete a logical volume for RDM.
94  //
95  deallvolumesCmd = new
96  G4UIcmdWithoutParameter("/grdm/noVolumes",this);
97  deallvolumesCmd->SetGuidance
98  (" RDM is not applied to any logical volumes");
99 
100  // deallvolumesCmd->SetParameterName("RemoveAVolume",true);
101  //
102  // The next command contols whether the branching ratio biasing will be applied or not
103  //
104  brbiasCmd = new G4UIcmdWithABool ("/grdm/BRbias",this);
105  brbiasCmd->SetGuidance("false: no biasing; true: all branches are treated as equal");
106  brbiasCmd->SetParameterName("BRBias",true);
107  brbiasCmd->SetDefaultValue(true);
108  //
109  // Command contols whether ICM will be applied or not
110  //
111  icmCmd = new G4UIcmdWithABool ("/grdm/applyICM",this);
112  icmCmd->SetGuidance("True: ICM is applied; false: no");
113  icmCmd->SetParameterName("applyICM",true);
114  icmCmd->SetDefaultValue(true);
115  //
116  // Command contols whether ARM will be applied or not
117  //
118  armCmd = new G4UIcmdWithABool ("/grdm/applyARM",this);
119  armCmd->SetGuidance("True: ARM is applied; false: no");
120  armCmd->SetParameterName("applyARM",true);
121  armCmd->SetDefaultValue(true);
122  //armCmd->AvailableForStates(G4State_PreInit);
123  //
124  // Command to set the h-l thresold for isomer production
125  //
126  hlthCmd = new G4UIcmdWithADoubleAndUnit("/grdm/hlThreshold",this);
127  hlthCmd->SetGuidance("Set the h-l threshold for isomer production");
128  hlthCmd->SetParameterName("hlThreshold",false);
129  // hlthCmd->SetRange("hlThreshold>0.");
130  hlthCmd->SetUnitCategory("Time");
131  // hlthCmd->AvailableForStates(G4State_PreInit);
132  //
133  // Command to define the incident particle source time profile.
134  //
135  sourcetimeprofileCmd = new
136  G4UIcmdWithAString("/grdm/sourceTimeProfile",this);
137  sourcetimeprofileCmd->SetGuidance
138  ("Supply the name of the ascii file containing the source particle time profile");
139  sourcetimeprofileCmd->SetParameterName("STimeProfile",true);
140  sourcetimeprofileCmd->SetDefaultValue("source.data");
141  //
142  //
143  // Command to define the incident particle source time profile.
144  //
145  decaybiasprofileCmd = new
146  G4UIcmdWithAString("/grdm/decayBiasProfile",this);
147  decaybiasprofileCmd->SetGuidance
148  ("Supply the name of the ascii file containing the decay bias time profile");
149  decaybiasprofileCmd->SetParameterName("DBiasProfile",true);
150  decaybiasprofileCmd->SetDefaultValue("bias.data");
151 
152  //
153  // Command to set the directional bias (collimation) vector
154  //
155  colldirCmd = new G4UIcmdWith3Vector("/grdm/decayDirection",this);
156  colldirCmd->SetGuidance("Supply the direction vector for decay products");
157  colldirCmd->SetParameterName("X","Y","Z",false);
158 
159  //
160  // Command to set the directional bias (collimation) half angle ("cone")
161  //
162  collangleCmd = new G4UIcmdWithADoubleAndUnit("/grdm/decayHalfAngle",this);
163  collangleCmd->SetGuidance
164  ("Supply maximum angle from direction vector for decay products");
165  collangleCmd->SetParameterName("halfAngle",false);
166  collangleCmd->SetUnitCategory("Angle");
167 
168  //
169  // This command setup the nuclei spliting parameter
170  //
171  splitnucleiCmd = new G4UIcmdWithAnInteger("/grdm/splitNuclei",this);
172  splitnucleiCmd->SetGuidance("Set number of spliting for the isotopes.");
173  splitnucleiCmd->SetParameterName("NSplit",true);
174  splitnucleiCmd->SetDefaultValue(1);
175  splitnucleiCmd->SetRange("NSplit>=1");
176 
177  //
178  // This command setup the verbose level of radioactive decay
179  //
180  verboseCmd = new G4UIcmdWithAnInteger("/grdm/verbose",this);
181  verboseCmd->SetGuidance("Set verbose level: 0, 1, 2 or 3");
182  verboseCmd->SetParameterName("VerboseLevel",true);
183  verboseCmd->SetDefaultValue(1);
184  verboseCmd->SetRange("VerboseLevel>=0");
185 
186  //
187  //This commansd allows the user to define its own decay datafile for
188  // a given isotope
189  //
190  userDecayDataCmd = new G4UIcommand("/grdm/setRadioactiveDecayFile",this);
191  G4UIparameter* Z_para= new G4UIparameter("Z_isotope",'i',true);
192  Z_para->SetParameterRange("Z_isotope > 0");
193  Z_para->SetGuidance("Z: Charge number of isotope");
194 
195 
196  G4UIparameter* A_para= new G4UIparameter("A_isotope",'i',true);
197  A_para->SetParameterRange("A_isotope > 1");
198  A_para->SetGuidance("A: mass number of isotope");
199 
200  G4UIparameter* FileName_para= new G4UIparameter("file_name",'s',true);
201  FileName_para->SetGuidance("Name of the user data file");
202  userDecayDataCmd->SetParameter(Z_para);
203  userDecayDataCmd->SetParameter(A_para);
204  userDecayDataCmd->SetParameter(FileName_para);
205 
206  //
207  //This commands allows the user to define its own evaporation data file for
208  // a given isotope
209  //
210  userEvaporationDataCmd = new G4UIcommand("/grdm/setPhotoEvaporationFile",this);
211  userEvaporationDataCmd->SetParameter(Z_para);
212  userEvaporationDataCmd->SetParameter(A_para);
213  userEvaporationDataCmd->SetParameter(FileName_para);
214 
215 
216 }
218 //
220 {
221  delete grdmDirectory;
222  delete nucleuslimitsCmd;
223  delete sourcetimeprofileCmd;
224  delete decaybiasprofileCmd;
225  delete analoguemcCmd;
226  delete fbetaCmd;
227  delete brbiasCmd;
228  delete splitnucleiCmd;
229  delete verboseCmd;
230  delete avolumeCmd;
231  delete deavolumeCmd;
232  delete allvolumesCmd;
233  delete deallvolumesCmd;
234  delete icmCmd;
235  delete armCmd;
236  delete hlthCmd;
237  delete userDecayDataCmd;
238  delete userEvaporationDataCmd;
239  delete colldirCmd;
240  delete collangleCmd;
241 
242 }
244 //
246 {
248  SetNucleusLimits(nucleuslimitsCmd->GetNewNucleusLimitsValue(newValues));
249 
250  } else if (command==analoguemcCmd) {theRadioactiveDecayContainer->
251  SetAnalogueMonteCarlo(analoguemcCmd->GetNewBoolValue(newValues));
252 
253  } else if (command==fbetaCmd) {theRadioactiveDecayContainer->
254  SetFBeta(fbetaCmd->GetNewBoolValue(newValues));
255 
256  } else if (command==avolumeCmd) {theRadioactiveDecayContainer->
257  SelectAVolume(newValues);
258 
259  } else if (command==deavolumeCmd) {theRadioactiveDecayContainer->
260  DeselectAVolume(newValues);}
261  else if (command==allvolumesCmd) {theRadioactiveDecayContainer->
262  SelectAllVolumes();}
263  else if (command==deallvolumesCmd) {theRadioactiveDecayContainer->
264  DeselectAllVolumes();}
265  else if (command==brbiasCmd) {theRadioactiveDecayContainer->
266  SetBRBias(brbiasCmd->GetNewBoolValue(newValues));}
268  SetSourceTimeProfile(newValues);}
270  SetDecayBias(newValues);}
271  else if (command==splitnucleiCmd) {theRadioactiveDecayContainer->
272  SetSplitNuclei(splitnucleiCmd->GetNewIntValue(newValues));}
273  else if (command==verboseCmd) {theRadioactiveDecayContainer->
274  SetVerboseLevel(verboseCmd->GetNewIntValue(newValues));}
275  else if (command==icmCmd ) {theRadioactiveDecayContainer->
276  SetICM(icmCmd->GetNewBoolValue(newValues));}
277  else if (command==armCmd ) {theRadioactiveDecayContainer->
278  SetARM(armCmd->GetNewBoolValue(newValues));}
279  else if (command==hlthCmd ) {theRadioactiveDecayContainer->
280  SetHLThreshold(hlthCmd->GetNewDoubleValue(newValues));
281 
282  } else if (command ==userDecayDataCmd){
283  G4int Z,A;
284  G4String file_name;
285  const char* nv = (const char*)newValues;
286  std::istringstream is(nv);
287  is >> Z >> A >> file_name;
289 
290  } else if (command ==userEvaporationDataCmd){
291  G4int Z,A;
292  G4String file_name;
293  const char* nv = (const char*)newValues;
294  std::istringstream is(nv);
295  is >> Z >> A >> file_name;
297 
298  } else if (command==colldirCmd) {theRadioactiveDecayContainer->
299  SetDecayDirection(colldirCmd->GetNew3VectorValue(newValues));
300 
301  } else if (command==collangleCmd) {theRadioactiveDecayContainer->
302  SetDecayHalfAngle(collangleCmd->GetNewDoubleValue(newValues));
303  }
304 }
305 
306 
307 
308 
309 
310 
G4NucleusLimits GetNewNucleusLimitsValue(G4String paramString)
G4UIcmdWithADoubleAndUnit * hlthCmd
G4UIcmdWithoutParameter * allvolumesCmd
G4UIcmdWithNucleusLimits * nucleuslimitsCmd
static G4NuclearLevelData * GetInstance()
G4bool AddPrivateData(G4int Z, G4int A, const G4String &filename)
void SetGuidance(const char *theGuidance)
static G4bool GetNewBoolValue(const char *paramString)
void SetParameterRange(const char *theRange)
Float_t Z
G4RadioactiveDecaymessenger(G4RadioactiveDecay *theRadioactiveDecayContainer)
static G4double GetNewDoubleValue(const char *paramString)
double A(double temperature)
G4RadioactiveDecay * theRadioactiveDecayContainer
static G4ThreeVector GetNew3VectorValue(const char *paramString)
G4UIcmdWithoutParameter * deallvolumesCmd
G4UIcmdWithADoubleAndUnit * collangleCmd
void SetNewValue(G4UIcommand *command, G4String newValues)
int G4int
Definition: G4Types.hh:78
static G4int GetNewIntValue(const char *paramString)
void AddUserDecayDataFile(G4int Z, G4int A, G4String filename)