Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4RadioactiveDecayBaseMessenger.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 // //
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 // //
37 
39 #include "G4NuclearLevelData.hh"
40 #include <sstream>
41 
42 
44 (G4RadioactiveDecayBase* theRadioactiveDecayContainer1)
45 :theRadioactiveDecayContainer(theRadioactiveDecayContainer1)
46 {
47  // main directory for control of the RDM
48  grdmDirectory = new G4UIdirectory("/grdm/");
49  grdmDirectory->SetGuidance("Controls for the Radioactive Decay Module.");
50 
51  // Command to define the limits on nucleus the RDM will treat.
52  nucleuslimitsCmd = new
53  G4UIcmdWithNucleusLimits("/grdm/nucleusLimits",this);
54  nucleuslimitsCmd->SetGuidance
55  ("Set the atomic weight and number limits for the RDM.");
56  nucleuslimitsCmd->SetParameterName("aMin","aMax","zMin","zMax",true);
57 
58  // Command controlling whether beta decay will be treated faithfully or
59  // in fast mode (obsolete)
60  fbetaCmd = new G4UIcmdWithABool ("/grdm/fBeta",this);
61 // fbetaCmd->SetGuidance("false: use 3-body decay, true: use histogram method");
62  fbetaCmd->SetGuidance("Fast (approximate) beta decay no longer used.");
63  fbetaCmd->SetGuidance("Kept for backward compatibility.");
64  fbetaCmd->SetParameterName("fBeta",true);
65  fbetaCmd->SetDefaultValue(false);
66 
67  // Select a logical volume for RDM
68  avolumeCmd = new
69  G4UIcmdWithAString("/grdm/selectVolume",this);
70  avolumeCmd->SetGuidance
71  ("Suppply a logical volumes name to add it to the RDM apply list");
72  avolumeCmd->SetParameterName("aVolume",false);
73 
74  // De-select a logical volume for RDM
75  deavolumeCmd = new
76  G4UIcmdWithAString("/grdm/deselectVolume",this);
77  deavolumeCmd->SetGuidance
78  ("Suppply a logical volumes name to remove it from the RDM apply list");
79  deavolumeCmd->SetParameterName("aVolume",false);
80 
81  // Select all logical volumes for RDM
82  allvolumesCmd = new
83  G4UIcmdWithoutParameter("/grdm/allVolumes",this);
84  allvolumesCmd->SetGuidance
85  (" apply RDM to all logical volumes. No parameter required.");
86  // allvolumeCmd->SetParameterName("AddAVolume",true);
87 
88  // De-select all logical volumes for RDM
89  deallvolumesCmd = new
90  G4UIcmdWithoutParameter("/grdm/noVolumes",this);
91  deallvolumesCmd->SetGuidance
92  (" RDM is not applied to any logical volumes");
93  // deallvolumesCmd->SetParameterName("RemoveAVolume",true);z
94 
95  // Command to invoke internal conversion or not
96  icmCmd = new G4UIcmdWithABool ("/grdm/applyICM",this);
97  icmCmd->SetGuidance("Command not active; kept for backward compatibility.");
98  icmCmd->SetGuidance("Internal conversion is always turned on.");
99  icmCmd->SetParameterName("applyICM",true);
100  icmCmd->SetDefaultValue(true);
101 
102  // Command to invoke atomic relaxation or not
103  armCmd = new G4UIcmdWithABool ("/grdm/applyARM",this);
104  armCmd->SetGuidance("True: ARM is applied; false: no");
105  armCmd->SetParameterName("applyARM",true);
106  armCmd->SetDefaultValue(true);
107  //armCmd->AvailableForStates(G4State_PreInit);
108 
109  // Command to set the directional bias (collimation) vector
110  colldirCmd = new G4UIcmdWith3Vector("/grdm/decayDirection",this);
111  colldirCmd->SetGuidance("Supply the direction vector for decay products");
112  colldirCmd->SetParameterName("X","Y","Z",false);
113 
114  // Command to set the directional bias (collimation) half angle ("cone")
115  collangleCmd = new G4UIcmdWithADoubleAndUnit("/grdm/decayHalfAngle",this);
116  collangleCmd->SetGuidance
117  ("Supply maximum angle from direction vector for decay products");
118  collangleCmd->SetParameterName("halfAngle",false);
119  collangleCmd->SetUnitCategory("Angle");
120 
121  // This command setup the verbose level of radioactive decay
122  verboseCmd = new G4UIcmdWithAnInteger("/grdm/verbose",this);
123  verboseCmd->SetGuidance("Set verbose level: 0, 1, 2 or 3");
124  verboseCmd->SetParameterName("VerboseLevel",true);
125  verboseCmd->SetDefaultValue(1);
126  verboseCmd->SetRange("VerboseLevel>=0");
127 
128  // Use a user-defined decay datafile for a given isotope
129  userDecayDataCmd = new G4UIcommand("/grdm/setRadioactiveDecayFile",this);
130  userDecayDataCmd->SetGuidance("Supply user-defined radioactive decay data file");
131 
132  G4UIparameter* Z_para= new G4UIparameter("Z_isotope",'i',true);
133  Z_para->SetParameterRange("Z_isotope > 0");
134  Z_para->SetGuidance("Z: Charge number of isotope");
135 
136  G4UIparameter* A_para= new G4UIparameter("A_isotope",'i',true);
137  A_para->SetParameterRange("A_isotope > 1");
138  A_para->SetGuidance("A: mass number of isotope");
139 
140  G4UIparameter* FileName_para= new G4UIparameter("file_name",'s',true);
141  FileName_para->SetGuidance("Name of the user data file");
142 
143  userDecayDataCmd->SetParameter(Z_para);
144  userDecayDataCmd->SetParameter(A_para);
145  userDecayDataCmd->SetParameter(FileName_para);
146 
147  // Use a user-defined evaporation data file for a given isotope
148  userEvaporationDataCmd = new G4UIcommand("/grdm/setPhotoEvaporationFile",this);
149  userEvaporationDataCmd->SetGuidance("Supply user-defined photon evaporation data file");
150 
151  G4UIparameter* Zpara= new G4UIparameter("Z_isotope",'i',true);
152  Zpara->SetParameterRange("Z_isotope > 0");
153  Zpara->SetGuidance("Z: Charge number of isotope");
154 
155  G4UIparameter* Apara= new G4UIparameter("A_isotope",'i',true);
156  Apara->SetParameterRange("A_isotope > 1");
157  Apara->SetGuidance("A: mass number of isotope");
158 
159  G4UIparameter* FileNamepara= new G4UIparameter("file_name",'s',true);
160  FileNamepara->SetGuidance("Name of the user data file");
161 
162  userEvaporationDataCmd->SetParameter(Zpara);
163  userEvaporationDataCmd->SetParameter(Apara);
164  userEvaporationDataCmd->SetParameter(FileNamepara);
165 }
166 
167 
169 {
170  delete grdmDirectory;
171  delete nucleuslimitsCmd;
172  delete fbetaCmd;
173  delete verboseCmd;
174  delete avolumeCmd;
175  delete deavolumeCmd;
176  delete allvolumesCmd;
177  delete deallvolumesCmd;
178  delete icmCmd;
179  delete armCmd;
180  delete userDecayDataCmd;
181  delete userEvaporationDataCmd;
182  delete colldirCmd;
183  delete collangleCmd;
184 }
185 
186 
187 void
189 {
190  if (command==nucleuslimitsCmd) {
192  SetNucleusLimits(nucleuslimitsCmd->GetNewNucleusLimitsValue(newValues));
193 
194  } else if (command==fbetaCmd) {
196  SetFBeta(fbetaCmd->GetNewBoolValue(newValues));
197 
198  } else if (command==avolumeCmd) {
200 
201  } else if (command==deavolumeCmd) {
203 
204  } else if (command==allvolumesCmd) {
206 
207  } else if (command==deallvolumesCmd) {
209 
210  } else if (command==verboseCmd) {
212  SetVerboseLevel(verboseCmd->GetNewIntValue(newValues));
213  } else if (command==icmCmd) {
215  SetICM(icmCmd->GetNewBoolValue(newValues));
216 
217  } else if (command==armCmd) {
219  SetARM(armCmd->GetNewBoolValue(newValues));
220 
221  } else if (command ==userDecayDataCmd) {
222  G4int Z,A;
223  G4String file_name;
224  const char* nv = (const char*)newValues;
225  std::istringstream is(nv);
226  is >> Z >> A >> file_name;
228 
229  } else if (command ==userEvaporationDataCmd) {
230  G4int Z,A;
231  G4String file_name;
232  const char* nv = (const char*)newValues;
233  std::istringstream is(nv);
234  is >> Z >> A >> file_name;
236 
237  } else if (command==colldirCmd) {
239  SetDecayDirection(colldirCmd->GetNew3VectorValue(newValues));
240 
241  } else if (command==collangleCmd) {
243  SetDecayHalfAngle(collangleCmd->GetNewDoubleValue(newValues));
244  }
245 }
246 
G4NucleusLimits GetNewNucleusLimitsValue(G4String paramString)
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
static G4double GetNewDoubleValue(const char *paramString)
void SelectAVolume(const G4String aVolume)
void SetNewValue(G4UIcommand *command, G4String newValues)
double A(double temperature)
G4RadioactiveDecayBaseMessenger(G4RadioactiveDecayBase *theRadioactiveDecayContainer)
void DeselectAVolume(const G4String aVolume)
static G4ThreeVector GetNew3VectorValue(const char *paramString)
int G4int
Definition: G4Types.hh:78
void AddUserDecayDataFile(G4int Z, G4int A, G4String filename)
static G4int GetNewIntValue(const char *paramString)