Geant4-11
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//
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// //
35
37#include "G4NuclearLevelData.hh"
38#include <sstream>
40
41
43 :theRadioactivationContainer(theRadioactivationContainer1)
44{
45 rdmDirectory = new G4UIdirectory("/process/had/rdm/");
46 rdmDirectory->SetGuidance("Controls the biased version of radioactive decay");
47
48 // Command to turn on/off variance reduction options
49 analoguemcCmd = new G4UIcmdWithABool("/process/had/rdm/analogueMC",this);
50 analoguemcCmd->SetGuidance("false: variance reduction method; true: analogue method");
51 analoguemcCmd->SetParameterName("AnalogueMC",true);
53
54 // Command to use branching ratio biasing or not
55 brbiasCmd = new G4UIcmdWithABool("/process/had/rdm/BRbias",this);
56 brbiasCmd->SetGuidance("false: no biasing; true: all branches are treated as equal");
57 brbiasCmd->SetParameterName("BRBias",true);
59
60 // Command to set the half-life thresold for isomer production
61 hlthCmd = new G4UIcmdWithADoubleAndUnit("/process/had/rdm/hlThreshold",this);
62 hlthCmd->SetGuidance("Set the h-l threshold for isomer production");
63 hlthCmd->SetParameterName("hlThreshold",false);
64 hlthCmd->SetUnitCategory("Time");
65
66 // Command to define the incident particle source time profile
67 sourcetimeprofileCmd = new G4UIcmdWithAString("/process/had/rdm/sourceTimeProfile",this);
69 ("Supply the name of the ascii file containing the source particle time profile");
70 sourcetimeprofileCmd->SetParameterName("STimeProfile",true);
72
73 // Command to define the incident particle source time profile
74 decaybiasprofileCmd = new G4UIcmdWithAString("/process/had/rdm/decayBiasProfile",this);
76 ("Supply the name of the ascii file containing the decay bias time profile");
77 decaybiasprofileCmd->SetParameterName("DBiasProfile",true);
79
80 // Command to set nuclei splitting parameter
81 splitnucleiCmd = new G4UIcmdWithAnInteger("/process/had/rdm/splitNuclei",this);
82 splitnucleiCmd->SetGuidance("Set number of splitting for the isotopes.");
83 splitnucleiCmd->SetParameterName("NSplit",true);
85 splitnucleiCmd->SetRange("NSplit>=1");
86}
87
88
90{
91 delete rdmDirectory;
92 delete analoguemcCmd;
95 delete brbiasCmd;
96 delete splitnucleiCmd;
97 delete hlthCmd;
98}
99
100
102{
103 if ( command == analoguemcCmd ) { theRadioactivationContainer->
104 SetAnalogueMonteCarlo( analoguemcCmd->GetNewBoolValue( newValues ) );
105 } else if ( command == brbiasCmd ) { theRadioactivationContainer->
106 SetBRBias( brbiasCmd->GetNewBoolValue( newValues ) );
107 } else if ( command == sourcetimeprofileCmd ) { theRadioactivationContainer->
108 SetSourceTimeProfile( newValues );
109 } else if ( command == decaybiasprofileCmd ) { theRadioactivationContainer->
110 SetDecayBias( newValues );
111 } else if ( command == splitnucleiCmd ) { theRadioactivationContainer->
112 SetSplitNuclei( splitnucleiCmd->GetNewIntValue( newValues ) );
113 } else if ( command == hlthCmd ) { theRadioactivationContainer->
114 SetHLThreshold( hlthCmd->GetNewDoubleValue( newValues ) );
115 }
116}
117
G4Radioactivation * theRadioactivationContainer
G4RadioactivationMessenger(G4Radioactivation *theRadioactivationContainer)
void SetNewValue(G4UIcommand *command, G4String newValues)
G4UIcmdWithADoubleAndUnit * hlthCmd
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